Support more complex fixed position selections
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 28 Jan 2015 20:06:26 +0000 (22:06 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Sat, 7 Feb 2015 00:37:21 +0000 (01:37 +0100)
Make it possible to create selections like
  [0,0,0] plus [0,1,0]
to specify a fixed set of positions.  This can be useful at least for
'gmx gangle' to calculate angles from a fixed reference vector.
Such selections were already understood properly by the parser, but
caused various crashes elsewhere in the code.  Fixed those crashes by
consistently managing the memory for the involved t_block structures
(always allocate one value in the index array, even if the block is
empty), and by not assuming that a plus-like keyword always has a child
element in the evaluation tree.

Fixes #1619.

Change-Id: I513cddbe882f269ad867c726a54583ee48b41b4d

src/gromacs/selection/evaluate.cpp
src/gromacs/selection/indexutil.cpp
src/gromacs/selection/poscalc.cpp
src/gromacs/selection/position.cpp
src/gromacs/selection/selection.cpp
src/gromacs/selection/sm_permute.cpp
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesConstantPositions.xml
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesConstantPositionsWithModifiers.xml [new file with mode: 0644]
src/gromacs/selection/tests/selectioncollection.cpp

index db2fd2cc666c7aeca092e7c65be588211f1e2673..e1f4df2f90be4d300eb740e54ab66ede04508446 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -1021,15 +1021,12 @@ _gmx_sel_evaluate_modifier(gmx_sel_evaluate_t                     *data,
         sel->u.expr.method->init_frame(data->top, data->fr, data->pbc,
                                        sel->u.expr.mdata);
     }
-    GMX_RELEASE_ASSERT(sel->child != NULL,
-                       "Modifier element with a value must have a child");
-    if (sel->child->v.type != POS_VALUE)
+    if (sel->child && sel->child->v.type != POS_VALUE)
     {
         GMX_THROW(gmx::NotImplementedError("Non-position valued modifiers not implemented"));
     }
     sel->u.expr.method->pupdate(data->top, data->fr, data->pbc,
-                                sel->child->v.u.p,
-                                &sel->v, sel->u.expr.mdata);
+                                NULL, &sel->v, sel->u.expr.mdata);
 }
 
 
index c448bc2274153490ec2553ac8e0d5369a827b75c..2e85724ac4d8e43ee679f2e55269a643fb5c03eb 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -801,8 +801,7 @@ next_group_index(int atomIndex, t_topology *top, e_index_t type, int *id)
  *   (depending on \p type) that is partially contained in the group.
  *   If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
  *
- * \p m should have been initialized somehow (calloc() is enough) unless
- * \p type is INDEX_UNKNOWN.
+ * \p m should have been initialized somehow (calloc() is enough).
  * \p g should be sorted.
  */
 void
@@ -811,8 +810,9 @@ gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
 {
     if (type == INDEX_UNKNOWN)
     {
+        sfree(t->a);
+        srenew(t->index, 2);
         t->nr           = 1;
-        snew(t->index, 2);
         t->nalloc_index = 2;
         t->index[0]     = 0;
         t->index[1]     = 0;
index f193abd42c577b1445fc4d4707a106080528671e..616ca1a0584abc48182b27ca8a33c44b54a3660e 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -1066,7 +1066,7 @@ gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p)
     {
         gmx_ana_indexmap_set_static(&p->m, &pc->b);
     }
-    gmx_ana_pos_reserve(p, p->m.mapb.nr, 0);
+    gmx_ana_pos_reserve(p, p->m.mapb.nr, -1);
     if (pc->flags & POS_VELOCITIES)
     {
         gmx_ana_pos_reserve_velocities(p);
index e5c6346f94e149531eb106c53f30752eacf96159..1315d94aaea8450357fb3bf72a78707990e0b814 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -99,7 +99,7 @@ gmx_ana_pos_reserve(gmx_ana_pos_t *pos, int n, int isize)
             srenew(pos->f, n);
         }
     }
-    if (isize > 0)
+    if (isize >= 0)
     {
         gmx_ana_indexmap_reserve(&pos->m, n, isize);
     }
@@ -200,7 +200,7 @@ gmx_ana_pos_copy(gmx_ana_pos_t *dest, gmx_ana_pos_t *src, bool bFirst)
 {
     if (bFirst)
     {
-        gmx_ana_pos_reserve(dest, src->count(), 0);
+        gmx_ana_pos_reserve(dest, src->count(), -1);
         if (src->v)
         {
             gmx_ana_pos_reserve_velocities(dest);
@@ -247,7 +247,8 @@ gmx_ana_pos_empty_init(gmx_ana_pos_t *pos)
     pos->m.mapb.nra = 0;
     pos->m.b.nr     = 0;
     pos->m.b.nra    = 0;
-    /* This should not really be necessary, but do it for safety... */
+    /* Initializing these should not really be necessary, but do it for
+     * safety... */
     pos->m.mapb.index[0] = 0;
     pos->m.b.index[0]    = 0;
     /* This function should only be used to construct all the possible
index fae8d6d64c013625b1c77d0ddc8441a496f36f4d..7b77d5c5ac4845a31d286e680517d69395f48180 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -87,6 +87,10 @@ SelectionData::SelectionData(SelectionTreeElement *elem,
         while (child->type == SEL_MODIFIER)
         {
             child = child->child;
+            if (!child)
+            {
+                break;
+            }
             if (child->type == SEL_SUBEXPRREF)
             {
                 child = child->child;
@@ -100,13 +104,16 @@ SelectionData::SelectionData(SelectionTreeElement *elem,
                 }
             }
         }
-        /* For variable references, we should skip the
-         * SEL_SUBEXPRREF and SEL_SUBEXPR elements. */
-        if (child->type == SEL_SUBEXPRREF)
+        if (child)
         {
-            child = child->child->child;
+            /* For variable references, we should skip the
+             * SEL_SUBEXPRREF and SEL_SUBEXPR elements. */
+            if (child->type == SEL_SUBEXPRREF)
+            {
+                child = child->child->child;
+            }
+            bDynamic_ = (child->child->flags & SEL_DYNAMIC);
         }
-        bDynamic_ = (child->child->flags & SEL_DYNAMIC);
     }
     initCoveredFraction(CFRAC_NONE);
 }
index a4f52fe150e3e97b472a5e59996b4446153778f5..eb7467149f342cc362c350aaa7f4dbb50ba39e23 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -238,7 +238,7 @@ free_data_permute(void *data)
 
 static void
 evaluate_permute(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
-                 gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data)
+                 gmx_ana_pos_t * /*p*/, gmx_ana_selvalue_t *out, void *data)
 {
     t_methoddata_permute *d = (t_methoddata_permute *)data;
     int                   i, j, b;
@@ -261,7 +261,7 @@ evaluate_permute(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc *
                 /* De-permute the reference ID */
                 refid = refid - (refid % d->n) + d->perm[refid % d->n];
             }
-            gmx_ana_pos_append(out->u.p, p, b, refid);
+            gmx_ana_pos_append(out->u.p, &d->p, b, refid);
         }
     }
     gmx_ana_pos_append_finish(out->u.p);
index e0e5a09a4269d598b055ef049a34cbe1227e91dc..386a27253fdeab3ee549a577d2232cceea75f205 100644 (file)
       <Sequence Name="Atoms">
         <Int Name="Length">0</Int>
       </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">1</Int>
+        <Position>
+          <Int Name="RefId">0</Int>
+          <Int Name="MappedId">0</Int>
+        </Position>
+      </Sequence>
     </Selection>
   </CompiledSelections>
   <EvaluatedSelections Name="Frame1">
@@ -28,6 +35,8 @@
             <Real Name="Y">-2</Real>
             <Real Name="Z">3.5</Real>
           </Vector>
+          <Int Name="RefId">0</Int>
+          <Int Name="MappedId">0</Int>
         </Position>
       </Sequence>
     </Selection>
diff --git a/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesConstantPositionsWithModifiers.xml b/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesConstantPositionsWithModifiers.xml
new file mode 100644 (file)
index 0000000..24f59c6
--- /dev/null
@@ -0,0 +1,57 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <ParsedSelections Name="Parsed">
+    <ParsedSelection Name="Selection1">
+      <String Name="Input">[0, 0, 0] plus [0, 1, 0]</String>
+      <String Name="Text">[0, 0, 0] plus [0, 1, 0]</String>
+      <Bool Name="Dynamic">false</Bool>
+    </ParsedSelection>
+  </ParsedSelections>
+  <CompiledSelections Name="Compiled">
+    <Selection Name="Selection1">
+      <Sequence Name="Atoms">
+        <Int Name="Length">0</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">2</Int>
+        <Position>
+          <Int Name="RefId">0</Int>
+          <Int Name="MappedId">0</Int>
+        </Position>
+        <Position>
+          <Int Name="RefId">1</Int>
+          <Int Name="MappedId">0</Int>
+        </Position>
+      </Sequence>
+    </Selection>
+  </CompiledSelections>
+  <EvaluatedSelections Name="Frame1">
+    <Selection Name="Selection1">
+      <Sequence Name="Atoms">
+        <Int Name="Length">0</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">2</Int>
+        <Position>
+          <Vector Name="Coordinates">
+            <Real Name="X">0</Real>
+            <Real Name="Y">0</Real>
+            <Real Name="Z">0</Real>
+          </Vector>
+          <Int Name="RefId">0</Int>
+          <Int Name="MappedId">0</Int>
+        </Position>
+        <Position>
+          <Vector Name="Coordinates">
+            <Real Name="X">0</Real>
+            <Real Name="Y">1</Real>
+            <Real Name="Z">0</Real>
+          </Vector>
+          <Int Name="RefId">1</Int>
+          <Int Name="MappedId">0</Int>
+        </Position>
+      </Sequence>
+    </Selection>
+  </EvaluatedSelections>
+</ReferenceData>
index 18fddac4c85c45ca98370208a74ef76c16584b41..f17c17d2b3143ebe4b6bb9af1a9741ff65d454bc 100644 (file)
@@ -1028,12 +1028,25 @@ TEST_F(SelectionCollectionDataTest, HandlesUnsortedIndexGroupsInSelectionsDelaye
     ASSERT_NO_FATAL_FAILURE(runCompiler());
 }
 
+
 TEST_F(SelectionCollectionDataTest, HandlesConstantPositions)
 {
     static const char * const selections[] = {
         "[1, -2, 3.5]"
     };
-    setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
+    setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
+             | efTestPositionMapping);
+    runTest("simple.gro", selections);
+}
+
+
+TEST_F(SelectionCollectionDataTest, HandlesConstantPositionsWithModifiers)
+{
+    static const char * const selections[] = {
+        "[0, 0, 0] plus [0, 1, 0]"
+    };
+    setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
+             | efTestPositionMapping);
     runTest("simple.gro", selections);
 }