Remove gmx_ana_index_t::name and unused code.
authorTeemu Murtola <teemu.murtola@gmail.com>
Fri, 12 Apr 2013 11:18:41 +0000 (14:18 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Mon, 29 Apr 2013 19:23:31 +0000 (21:23 +0200)
Remove the name field from gmx_ana_index_t, which was in most cases
NULL, and was a source of multiple ugly hacks in the selection code to
avoid leaking the memory allocated for the name.
Now that SelectionTreeElement has a proper std::string name member,
there is no need for this name field anywhere in the selection code
where the atom groups are handled as selection values.

Also remove unused code from indexutil.* (that would have required extra
effort to not break it).

Change-Id: I14698c17da3b479de4016a6188fcd235514b9e2a

src/gromacs/selection/compiler.cpp
src/gromacs/selection/evaluate.cpp
src/gromacs/selection/indexutil.cpp
src/gromacs/selection/indexutil.h
src/gromacs/selection/parsetree.cpp
src/gromacs/selection/poscalc.cpp
src/gromacs/selection/selection.cpp
src/gromacs/selection/selectioncollection.cpp
src/gromacs/trajectoryanalysis/runnercommon.cpp

index 2b198db6de5b6322d555682f98ec9e1d916b5b67..60a87a179579ab6568dd1ee9490d4d3c8a8f2b24 100644 (file)
@@ -458,8 +458,6 @@ void SelectionTreeElement::freeCompilerData()
         evaluate = cdata->evaluate;
         if (cdata->flags & SEL_CDATA_MINMAXALLOC)
         {
-            cdata->gmin->name = NULL;
-            cdata->gmax->name = NULL;
             gmx_ana_index_deinit(cdata->gmin);
             gmx_ana_index_deinit(cdata->gmax);
             sfree(cdata->gmin);
@@ -748,7 +746,6 @@ create_subexpression_name(const SelectionTreeElementPointer &sel, int i)
 {
     std::string name(gmx::formatString("SubExpr %d", i));
     sel->setName(name);
-    sel->u.cgrp.name = strdup(name.c_str());
 }
 
 /*! \brief
@@ -1621,12 +1618,11 @@ initialize_evalgrps(gmx_ana_selcollection_t *sc)
         if (root->child->type != SEL_SUBEXPR
             || (root->child->v.type != GROUP_VALUE && !(root->flags & SEL_ATOMVAL)))
         {
-            gmx_ana_index_set(&root->u.cgrp, -1, 0, root->u.cgrp.name, 0);
+            gmx_ana_index_set(&root->u.cgrp, -1, 0, 0);
         }
         else if (root->child->cdata->flags & SEL_CDATA_FULLEVAL)
         {
-            gmx_ana_index_set(&root->u.cgrp, sc->gall.isize, sc->gall.index,
-                              root->u.cgrp.name, 0);
+            gmx_ana_index_set(&root->u.cgrp, sc->gall.isize, sc->gall.index, 0);
         }
         root = root->next;
     }
@@ -1755,7 +1751,7 @@ make_static(const SelectionTreeElementPointer &sel)
      * */
     if (sel->v.type == GROUP_VALUE)
     {
-        gmx_ana_index_set(&sel->u.cgrp, sel->v.u.g->isize, sel->v.u.g->index, NULL, 0);
+        gmx_ana_index_set(&sel->u.cgrp, sel->v.u.g->isize, sel->v.u.g->index, 0);
     }
 }
 
@@ -1990,10 +1986,8 @@ evaluate_boolean_static_part(gmx_sel_evaluate_t                *data,
          * and hence we can overwrite it safely. */
         if (child->u.cgrp.nalloc_index > 0)
         {
-            char *name = child->u.cgrp.name;
             gmx_ana_index_copy(&child->u.cgrp, child->v.u.g, false);
             gmx_ana_index_squeeze(&child->u.cgrp);
-            child->u.cgrp.name = name;
         }
         else
         {
@@ -2088,11 +2082,8 @@ evaluate_boolean_minmax_grps(const SelectionTreeElementPointer &sel,
                 gmx_ana_index_copy(sel->child->v.u.g, gmin, false);
                 if (sel->child->u.cgrp.nalloc_index > 0)
                 {
-                    /* Keep the name as in evaluate_boolean_static_part(). */
-                    char *name = sel->child->u.cgrp.name;
                     gmx_ana_index_reserve(&sel->child->u.cgrp, gmin->isize);
                     gmx_ana_index_copy(&sel->child->u.cgrp, gmin, false);
-                    sel->child->u.cgrp.name = name;
                 }
                 else
                 {
@@ -2339,10 +2330,6 @@ analyze_static(gmx_sel_evaluate_t                *data,
     {
         gmx_ana_index_squeeze(sel->cdata->gmin);
         gmx_ana_index_squeeze(sel->cdata->gmax);
-        sfree(sel->cdata->gmin->name);
-        sfree(sel->cdata->gmax->name);
-        sel->cdata->gmin->name = NULL;
-        sel->cdata->gmax->name = NULL;
     }
 
     /* Replace the result of the evaluation */
@@ -2399,7 +2386,6 @@ init_root_item(const SelectionTreeElementPointer &root,
     }
 
     /* Set the evaluation group */
-    char *name = root->u.cgrp.name;
     if (root->evaluate)
     {
         /* Non-atom-valued non-group expressions don't care about the group, so
@@ -2407,12 +2393,12 @@ init_root_item(const SelectionTreeElementPointer &root,
         if ((expr->flags & SEL_VARNUMVAL)
             || ((expr->flags & SEL_SINGLEVAL) && expr->v.type != GROUP_VALUE))
         {
-            gmx_ana_index_set(&root->u.cgrp, -1, NULL, NULL, 0);
+            gmx_ana_index_set(&root->u.cgrp, -1, NULL, 0);
         }
         else if (expr->cdata->gmax->isize == gall->isize)
         {
             /* Save some memory by only referring to the global group. */
-            gmx_ana_index_set(&root->u.cgrp, gall->isize, gall->index, NULL, 0);
+            gmx_ana_index_set(&root->u.cgrp, gall->isize, gall->index, 0);
         }
         else
         {
@@ -2452,7 +2438,6 @@ init_root_item(const SelectionTreeElementPointer &root,
     {
         gmx_ana_index_clear(&root->u.cgrp);
     }
-    root->u.cgrp.name = name;
 }
 
 
@@ -2492,14 +2477,10 @@ postprocess_item_subexpressions(const SelectionTreeElementPointer &sel)
         && (sel->cdata->flags & SEL_CDATA_STATICEVAL)
         && !(sel->cdata->flags & SEL_CDATA_FULLEVAL))
     {
-        char *name;
-
         /* We need to free memory allocated for the group, because it is no
          * longer needed (and would be lost on next call to the evaluation
-         * function). But we need to preserve the name. */
-        name = sel->u.cgrp.name;
+         * function). */
         gmx_ana_index_deinit(&sel->u.cgrp);
-        sel->u.cgrp.name = name;
 
         sel->evaluate        = &_gmx_sel_evaluate_subexpr_staticeval;
         sel->cdata->evaluate = sel->evaluate;
index cc1751b6c5b7f225daeacde77a8183092ea85408..9a4f0744166fce60c963f014e6bc837c80235d32 100644 (file)
@@ -590,7 +590,7 @@ _gmx_sel_evaluate_subexpr_staticeval(gmx_sel_evaluate_t                *data,
         }
         else
         {
-            gmx_ana_index_set(&sel->u.cgrp, g->isize, g->index, sel->u.cgrp.name, 0);
+            gmx_ana_index_set(&sel->u.cgrp, g->isize, g->index, 0);
         }
     }
 }
@@ -629,18 +629,13 @@ _gmx_sel_evaluate_subexpr(gmx_sel_evaluate_t                *data,
             SelelemTemporaryValueAssigner assigner(sel->child, *sel);
             sel->child->evaluate(data, sel->child, g);
         }
-        /* We need to keep the name for the cgrp across the copy to avoid
-         * problems if g has a name set. */
-        char *name = sel->u.cgrp.name;
         gmx_ana_index_copy(&sel->u.cgrp, g, false);
-        sel->u.cgrp.name = name;
         gmiss.isize      = 0;
     }
     else
     {
         gmissreserver.reserve(&gmiss, g->isize);
         gmx_ana_index_difference(&gmiss, g, &sel->u.cgrp);
-        gmiss.name = NULL;
     }
     if (gmiss.isize > 0)
     {
@@ -1130,7 +1125,6 @@ _gmx_sel_evaluate_or(gmx_sel_evaluate_t                *data,
     child = child->next;
     while (child && tmp.isize > 0)
     {
-        tmp.name = NULL;
         {
             MempoolSelelemReserver reserver(child, tmp.isize);
             child->evaluate(data, child, &tmp);
index 401b737fdc333ec9d4d7d89eaab459f347c73af5..e9eda49f16944846addcccb17c3963cf342d5e04 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
  * \author Teemu Murtola <teemu.murtola@gmail.com>
  * \ingroup module_selection
  */
+#include "gromacs/selection/indexutil.h"
+
+#include <cstring>
+
+#include <string>
+#include <vector>
+
 #include "gromacs/legacyheaders/index.h"
 #include "gromacs/legacyheaders/gmx_fatal.h"
 #include "gromacs/legacyheaders/smalloc.h"
-#include "gromacs/legacyheaders/string2.h"
 #include "gromacs/legacyheaders/typedefs.h"
 
-#include "gromacs/selection/indexutil.h"
 
 /********************************************************************
  * gmx_ana_indexgrps_t functions
  */
 struct gmx_ana_indexgrps_t
 {
+    //! Initializes an empty set of groups.
+    explicit gmx_ana_indexgrps_t(int nr) : nr(nr), g(NULL)
+    {
+        names.reserve(nr);
+        snew(g, nr);
+    }
+    ~gmx_ana_indexgrps_t()
+    {
+        for (int i = 0; i < nr; ++i)
+        {
+            gmx_ana_index_deinit(&g[i]);
+        }
+        sfree(g);
+    }
+
     /** Number of index groups. */
-    int                 nr;
+    int                       nr;
     /** Array of index groups. */
-    gmx_ana_index_t    *g;
+    gmx_ana_index_t          *g;
+    /** Group names. */
+    std::vector<std::string>  names;
 };
 
-/*!
- * \param[out] g     Index group structure.
- * \param[in]  ngrps Number of groups for which memory is allocated.
- */
-void
-gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps)
-{
-    snew(*g, 1);
-    (*g)->nr = ngrps;
-    snew((*g)->g,    ngrps);
-}
-
 /*!
  * \param[out] g     Index group structure.
  * \param[in]  top   Topology structure.
@@ -93,7 +103,6 @@ gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
 {
     t_blocka *block = NULL;
     char    **names = NULL;
-    int       i, j;
 
     if (fnm)
     {
@@ -106,29 +115,44 @@ gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
     }
     else
     {
-        snew(*g, 1);
-        (*g)->nr    = 0;
-        (*g)->g     = NULL;
+        *g = new gmx_ana_indexgrps_t(0);
         return;
     }
 
-    gmx_ana_indexgrps_alloc(g, block->nr);
-    for (i = 0; i < block->nr; ++i)
+    try
     {
-        gmx_ana_index_t *grp = &(*g)->g[i];
+        *g = new gmx_ana_indexgrps_t(block->nr);
+        for (int i = 0; i < block->nr; ++i)
+        {
+            gmx_ana_index_t *grp = &(*g)->g[i];
 
-        grp->isize = block->index[i+1] - block->index[i];
-        snew(grp->index, grp->isize);
-        for (j = 0; j < grp->isize; ++j)
+            grp->isize = block->index[i+1] - block->index[i];
+            snew(grp->index, grp->isize);
+            for (int j = 0; j < grp->isize; ++j)
+            {
+                grp->index[j] = block->a[block->index[i]+j];
+            }
+            grp->nalloc_index = grp->isize;
+            (*g)->names.push_back(names[i]);
+        }
+    }
+    catch (...)
+    {
+        done_blocka(block);
+        sfree(block);
+        for (int i = 0; i < block->nr; ++i)
         {
-            grp->index[j] = block->a[block->index[i]+j];
+            sfree(names[i]);
         }
-        grp->name         = names[i];
-        grp->nalloc_index = grp->isize;
+        sfree(names);
+        throw;
     }
-
     done_blocka(block);
     sfree(block);
+    for (int i = 0; i < block->nr; ++i)
+    {
+        sfree(names[i]);
+    }
     sfree(names);
 }
 
@@ -140,39 +164,7 @@ gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
 void
 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
 {
-    int  i;
-
-    if (g->nr == 0)
-    {
-        sfree(g);
-        return;
-    }
-    for (i = 0; i < g->nr; ++i)
-    {
-        gmx_ana_index_deinit(&g->g[i]);
-    }
-    sfree(g->g);
-    g->nr    = 0;
-    g->g     = NULL;
-    sfree(g);
-}
-
-/*!
- * \param[out] dest Destination index groups.
- * \param[in]  src  Source index groups.
- *
- * A deep copy is made for all fields, including the group names.
- */
-void
-gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t **dest, gmx_ana_indexgrps_t *src)
-{
-    int g;
-
-    gmx_ana_indexgrps_alloc(dest, src->nr);
-    for (g = 0; g < src->nr; ++g)
-    {
-        gmx_ana_index_copy(&(*dest)->g[g], &src->g[g], true);
-    }
+    delete g;
 }
 
 /*!
@@ -203,53 +195,64 @@ gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
 }
 
 /*!
- * \param[out] dest Output structure.
- * \param[in]  src  Input index groups.
- * \param[in]  n    Number of the group to extract.
+ * \param[out] dest     Output structure.
+ * \param[out] destName Receives the name of the group if found.
+ * \param[in]  src      Input index groups.
+ * \param[in]  n        Number of the group to extract.
  * \returns true if \p n is a valid group in \p src, false otherwise.
  */
 bool
-gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, int n)
+gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
+                          gmx_ana_indexgrps_t *src, int n)
 {
+    destName->clear();
     if (n < 0 || n >= src->nr)
     {
         dest->isize = 0;
         return false;
     }
 
+    if (destName != NULL)
+    {
+        *destName = src->names[n];
+    }
     gmx_ana_index_copy(dest, &src->g[n], true);
     return true;
 }
 
 /*!
- * \param[out] dest Output structure.
- * \param[in]  src  Input index groups.
- * \param[in]  name Name (or part of the name) of the group to extract.
+ * \param[out] dest     Output structure.
+ * \param[out] destName Receives the name of the group if found.
+ * \param[in]  src      Input index groups.
+ * \param[in]  name     Name (or part of the name) of the group to extract.
  * \returns true if \p name is a valid group in \p src, false otherwise.
  *
  * Uses the Gromacs routine find_group() to find the actual group;
  * the comparison is case-insensitive.
  */
 bool
-gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *name)
+gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
+                       gmx_ana_indexgrps_t *src,
+                       const char *name)
 {
-    int    i;
-    char **names;
+    const char **names;
 
+    destName->clear();
     snew(names, src->nr);
-    for (i = 0; i < src->nr; ++i)
+    for (int i = 0; i < src->nr; ++i)
     {
-        names[i] = src->g[i].name;
+        names[i] = src->names[i].c_str();
     }
-    i = find_group(name, src->nr, names);
+    int n = find_group(const_cast<char *>(name), src->nr,
+                       const_cast<char **>(names));
     sfree(names);
-    if (i == NOTSET)
+    if (n == NOTSET)
     {
         dest->isize = 0;
         return false;
     }
 
-    return gmx_ana_indexgrps_extract(dest, src, i);
+    return gmx_ana_indexgrps_extract(dest, destName, src, n);
 }
 
 /*!
@@ -261,12 +264,10 @@ gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *na
 void
 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn)
 {
-    int  i;
-
-    for (i = 0; i < g->nr; ++i)
+    for (int i = 0; i < g->nr; ++i)
     {
-        fprintf(fp, " %2d: ", i);
-        gmx_ana_index_dump(fp, &g->g[i], i, maxn);
+        fprintf(fp, " Group %2d \"%s\" ", i + 1, g->names[i].c_str());
+        gmx_ana_index_dump(fp, &g->g[i], maxn);
     }
 }
 
@@ -311,7 +312,6 @@ gmx_ana_index_clear(gmx_ana_index_t *g)
 {
     g->isize        = 0;
     g->index        = NULL;
-    g->name         = NULL;
     g->nalloc_index = 0;
 }
 
@@ -319,29 +319,25 @@ gmx_ana_index_clear(gmx_ana_index_t *g)
  * \param[out] g      Output structure.
  * \param[in]  isize  Number of atoms in the new group.
  * \param[in]  index  Array of \p isize atoms (can be NULL if \p isize is 0).
- * \param[in]  name   Name for the new group (can be NULL).
  * \param[in]  nalloc Number of elements allocated for \p index
  *   (if 0, \p index is not freed in gmx_ana_index_deinit())
  *
  * No copy if \p index is made.
  */
 void
-gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, char *name,
-                  int nalloc)
+gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, int nalloc)
 {
     g->isize        = isize;
     g->index        = index;
-    g->name         = name;
     g->nalloc_index = nalloc;
 }
 
 /*!
  * \param[out] g      Output structure.
  * \param[in]  natoms Number of atoms.
- * \param[in]  name   Name for the new group (can be NULL).
  */
 void
-gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name)
+gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms)
 {
     int  i;
 
@@ -351,7 +347,6 @@ gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name)
     {
         g->index[i] = i;
     }
-    g->name         = name;
     g->nalloc_index = natoms;
 }
 
@@ -367,7 +362,6 @@ gmx_ana_index_deinit(gmx_ana_index_t *g)
     {
         sfree(g->index);
     }
-    sfree(g->name);
     gmx_ana_index_clear(g);
 }
 
@@ -390,38 +384,21 @@ gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
             snew(dest->index, dest->isize);
             dest->nalloc_index = dest->isize;
         }
-        memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
-    }
-    if (bAlloc && src->name)
-    {
-        dest->name = strdup(src->name);
-    }
-    else if (bAlloc || src->name)
-    {
-        dest->name = src->name;
+        std::memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
     }
 }
 
 /*!
  * \param[in]  fp     Where to print the output.
  * \param[in]  g      Index group to print.
- * \param[in]  i      Group number to use if the name is NULL.
  * \param[in]  maxn   Maximum number of indices to print (-1 = print all).
  */
 void
-gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int i, int maxn)
+gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int maxn)
 {
     int  j, n;
 
-    if (g->name)
-    {
-        fprintf(fp, "\"%s\"", g->name);
-    }
-    else
-    {
-        fprintf(fp, "Group %d", i+1);
-    }
-    fprintf(fp, " (%d atoms)", g->isize);
+    fprintf(fp, "(%d atoms)", g->isize);
     if (maxn != 0)
     {
         fprintf(fp, ":");
@@ -442,35 +419,6 @@ gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int i, int maxn)
     fprintf(fp, "\n");
 }
 
-/*!
- * \param[in]  g      Input index group.
- * \param[in]  natoms Number of atoms to check against.
- *
- * If any atom index in the index group is less than zero or >= \p natoms,
- * gmx_fatal() is called.
- */
-void
-gmx_ana_index_check(gmx_ana_index_t *g, int natoms)
-{
-    int  j;
-
-    for (j = 0; j < g->isize; ++j)
-    {
-        if (g->index[j] >= natoms)
-        {
-            gmx_fatal(FARGS, "Atom index (%d) in index group %s (%d atoms) "
-                      "larger than number of atoms in trajectory (%d atoms)",
-                      g->index[j], g->name, g->isize, natoms);
-        }
-        else if (g->index[j] < 0)
-        {
-            gmx_fatal(FARGS, "Atom index (%d) in index group %s (%d atoms) "
-                      "is less than zero",
-                      g->index[j], g->name, g->isize);
-        }
-    }
-}
-
 /*!
  * \param[in]  g      Index group to check.
  * \returns    true if the index group is sorted and has no duplicates,
@@ -827,7 +775,7 @@ gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
             srenew(t->a, g->isize);
             t->nalloc_a = g->isize;
         }
-        memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
+        std::memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
     }
 
     /* Allocate memory for the block index. We don't know in advance
@@ -1169,7 +1117,7 @@ gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
         m->mapid[i] = m->orgid[i];
     }
     m->mapb.nr = m->nr;
-    memcpy(m->mapb.index, m->b.index, (m->nr+1)*sizeof(*(m->mapb.index)));
+    std::memcpy(m->mapb.index, m->b.index, (m->nr+1)*sizeof(*(m->mapb.index)));
     m->bStatic    = true;
     m->bMapStatic = true;
 }
@@ -1223,15 +1171,15 @@ gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bF
         dest->type       = src->type;
         dest->b.nr       = src->b.nr;
         dest->b.nra      = src->b.nra;
-        memcpy(dest->orgid,      src->orgid,      dest->b.nr*sizeof(*dest->orgid));
-        memcpy(dest->b.index,    src->b.index,   (dest->b.nr+1)*sizeof(*dest->b.index));
-        memcpy(dest->b.a,        src->b.a,        dest->b.nra*sizeof(*dest->b.a));
+        std::memcpy(dest->orgid,      src->orgid,      dest->b.nr*sizeof(*dest->orgid));
+        std::memcpy(dest->b.index,    src->b.index,   (dest->b.nr+1)*sizeof(*dest->b.index));
+        std::memcpy(dest->b.a,        src->b.a,        dest->b.nra*sizeof(*dest->b.a));
     }
     dest->nr         = src->nr;
     dest->mapb.nr    = src->mapb.nr;
-    memcpy(dest->refid,      src->refid,      dest->nr*sizeof(*dest->refid));
-    memcpy(dest->mapid,      src->mapid,      dest->nr*sizeof(*dest->mapid));
-    memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
+    std::memcpy(dest->refid,      src->refid,      dest->nr*sizeof(*dest->refid));
+    std::memcpy(dest->mapid,      src->mapid,      dest->nr*sizeof(*dest->mapid));
+    std::memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
     dest->bStatic    = src->bStatic;
     dest->bMapStatic = src->bMapStatic;
 }
index f7f0ddbdd125220cc4721c6f56f9ac56b6ccc446..da31448525535e98285e0f67c79a180b03aaf9b0 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -61,6 +61,8 @@
 #ifndef GMX_SELECTION_INDEXUTIL_H
 #define GMX_SELECTION_INDEXUTIL_H
 
+#include <string>
+
 #include "../legacyheaders/typedefs.h"
 
 /** Stores a set of index groups. */
@@ -89,8 +91,6 @@ typedef struct gmx_ana_index_t
     int                 isize;
     /** List of atoms. */
     atom_id            *index;
-    /** Group name. */
-    char               *name;
     /** Number of items allocated for \p index. */
     int                 nalloc_index;
 } gmx_ana_index_t;
@@ -191,9 +191,6 @@ typedef struct gmx_ana_indexmap_t
 /*! \name Functions for handling gmx_ana_indexgrps_t
  */
 /*@{*/
-/** Allocate memory for index groups. */
-void
-gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps);
 /** Reads index groups from a file or constructs them from topology. */
 void
 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
@@ -201,9 +198,6 @@ gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
 /** Frees memory allocated for index groups. */
 void
 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g);
-/** Create a deep copy of \c gmx_ana_indexgrps_t. */
-void
-gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t **dest, gmx_ana_indexgrps_t *src);
 /** Returns true if the index group structure is emtpy. */
 bool
 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g);
@@ -213,10 +207,12 @@ gmx_ana_index_t *
 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n);
 /** Extracts a single index group. */
 bool
-gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, int n);
+gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
+                          gmx_ana_indexgrps_t *src, int n);
 /** Finds and extracts a single index group by name. */
 bool
-gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *name);
+gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
+                       gmx_ana_indexgrps_t *src, const char *name);
 
 /** Writes out a list of index groups. */
 void
@@ -237,11 +233,10 @@ void
 gmx_ana_index_clear(gmx_ana_index_t *g);
 /** Constructs a \c gmx_ana_index_t from given values. */
 void
-gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, char *name,
-                  int nalloc);
+gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, int nalloc);
 /** Creates a simple index group from the first to the \p natoms'th atom. */
 void
-gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name);
+gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms);
 /** Frees memory allocated for an index group. */
 void
 gmx_ana_index_deinit(gmx_ana_index_t *g);
@@ -251,11 +246,8 @@ gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc);
 
 /** Writes out the contents of a index group. */
 void
-gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int i, int maxn);
+gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int maxn);
 
-/** Checks whether all indices are between 0 and \p natoms. */
-void
-gmx_ana_index_check(gmx_ana_index_t *g, int natoms);
 /** Checks whether an index group is sorted. */
 bool
 gmx_ana_index_check_sorted(gmx_ana_index_t *g);
index c1e6614eb64419828283290bf5e45b52d62d88eb..6790798160d3bd6c0db71558c729a003841c4573 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -930,13 +930,13 @@ _gmx_sel_init_group_by_name(const char *name, yyscan_t scanner)
     }
     SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_CONST));
     _gmx_selelem_set_vtype(sel, GROUP_VALUE);
-    /* FIXME: The constness should not be cast away */
-    if (!gmx_ana_indexgrps_find(&sel->u.cgrp, grps, (char *)name))
+    std::string                 foundName;
+    if (!gmx_ana_indexgrps_find(&sel->u.cgrp, &foundName, grps, name))
     {
         _gmx_selparser_error(scanner, "Cannot match 'group %s'", name);
         return SelectionTreeElementPointer();
     }
-    sel->setName(sel->u.cgrp.name);
+    sel->setName(foundName);
     return sel;
 }
 
@@ -966,12 +966,13 @@ _gmx_sel_init_group_by_id(int id, yyscan_t scanner)
     }
     SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_CONST));
     _gmx_selelem_set_vtype(sel, GROUP_VALUE);
-    if (!gmx_ana_indexgrps_extract(&sel->u.cgrp, grps, id))
+    std::string                 foundName;
+    if (!gmx_ana_indexgrps_extract(&sel->u.cgrp, &foundName, grps, id))
     {
         _gmx_selparser_error(scanner, "Cannot match 'group %d'", id);
         return SelectionTreeElementPointer();
     }
-    sel->setName(sel->u.cgrp.name);
+    sel->setName(foundName);
     return sel;
 }
 
@@ -1033,7 +1034,6 @@ _gmx_sel_init_selection(const char                        *name,
     if (name)
     {
         root->setName(name);
-        root->u.cgrp.name = strdup(name);
     }
     /* Update the flags */
     _gmx_selelem_update_flags(root, scanner);
@@ -1059,17 +1059,13 @@ _gmx_sel_init_selection(const char                        *name,
             && child->child->child->type == SEL_CONST
             && child->child->child->v.type == GROUP_VALUE)
         {
-            const char *grpName = child->child->child->u.cgrp.name;
-            root->setName(grpName);
-            root->u.cgrp.name = strdup(grpName);
+            root->setName(child->child->child->name());
         }
     }
     /* If there still is no name, use the selection string */
     if (root->name().empty())
     {
-        const char *selStr = _gmx_sel_lexer_pselstr(scanner);
-        root->setName(selStr);
-        root->u.cgrp.name = strdup(selStr);
+        root->setName(_gmx_sel_lexer_pselstr(scanner));
     }
 
     /* Print out some information if the parser is interactive */
@@ -1125,7 +1121,6 @@ _gmx_sel_assign_variable(const char                        *name,
     /* Create the root element */
     root.reset(new SelectionTreeElement(SEL_ROOT));
     root->setName(name);
-    root->u.cgrp.name = strdup(name);
     /* Create the subexpression element */
     root->child.reset(new SelectionTreeElement(SEL_SUBEXPR));
     root->child->setName(name);
index f294187c63bc9ba1ae53d9b119a57597fac6fa38..402af13c2f6ba33a9f605a5382ce5f1469ffaf94 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -658,12 +658,10 @@ set_poscalc_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g, bool bBase)
     if ((pc->flags & POS_COMPLWHOLE) && !bBase && pc->b.nra > g->isize)
     {
         gmx_ana_index_copy(&pc->gmax, g, true);
-        sfree(pc->gmax.name);
-        pc->gmax.name  = NULL;
     }
     else
     {
-        gmx_ana_index_set(&pc->gmax, pc->b.nra, pc->b.a, NULL, 0);
+        gmx_ana_index_set(&pc->gmax, pc->b.nra, pc->b.a, 0);
     }
 }
 
@@ -728,7 +726,7 @@ should_merge(gmx_ana_poscalc_t *pc1, gmx_ana_poscalc_t *pc2,
         return false;
     }
     /* Find the overlap between the calculations. */
-    gmx_ana_index_set(&g2, pc2->b.nra, pc2->b.a, NULL, 0);
+    gmx_ana_index_set(&g2, pc2->b.nra, pc2->b.a, 0);
     gmx_ana_index_intersection(g, g1, &g2);
     /* Do not merge if there is no overlap. */
     if (g->isize == 0)
@@ -801,8 +799,8 @@ merge_to_base(gmx_ana_poscalc_t *base, gmx_ana_poscalc_t *pc)
     int              i, bi, bj, bo;
 
     base->flags |= pc->flags & (POS_VELOCITIES | POS_FORCES);
-    gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, NULL, 0);
-    gmx_ana_index_set(&gb, base->b.nra, base->b.a, NULL, 0);
+    gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, 0);
+    gmx_ana_index_set(&gb, base->b.nra, base->b.a, 0);
     isize = gmx_ana_index_difference_size(&gp, &gb);
     if (isize > 0)
     {
@@ -857,7 +855,7 @@ merge_to_base(gmx_ana_poscalc_t *base, gmx_ana_poscalc_t *pc)
         base->b.a        = g.index;
         base->b.nalloc_a = g.isize;
         /* Refresh the gmax field */
-        gmx_ana_index_set(&base->gmax, base->b.nra, base->b.a, NULL, 0);
+        gmx_ana_index_set(&base->gmax, base->b.nra, base->b.a, 0);
     }
 }
 
@@ -912,7 +910,7 @@ setup_base(gmx_ana_poscalc_t *pc)
         return;
     }
 
-    gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, NULL, 0);
+    gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, 0);
     gmx_ana_index_clear(&g);
     gmx_ana_index_reserve(&g, pc->b.nra);
     pbase = pc;
@@ -963,7 +961,7 @@ setup_base(gmx_ana_poscalc_t *pc)
                     merge_bases(pbase, base);
                 }
             }
-            gmx_ana_index_set(&gp, pbase->b.nra, pbase->b.a, NULL, 0);
+            gmx_ana_index_set(&gp, pbase->b.nra, pbase->b.a, 0);
             gmx_ana_index_reserve(&g, pc->b.nra);
         }
         /* Proceed to the next unchecked calculation */
index 23ecee290c3ec337863d638973bd51c5a1f2b5e4..bfc0a14b74eeb79bea127758c9a8a6838c72d68f 100644 (file)
@@ -223,7 +223,6 @@ SelectionData::restoreOriginalPositions(const t_topology *top)
     {
         gmx_ana_pos_t &p = rawPositions_;
         gmx_ana_index_copy(p.g, rootElement().v.u.g, false);
-        p.g->name = NULL;
         gmx_ana_indexmap_update(&p.m, p.g, hasFlag(gmx::efSelection_DynamicMask));
         p.nr = p.m.nr;
         refreshMassesAndCharges(top);
@@ -251,8 +250,8 @@ Selection::printDebugInfo(FILE *fp, int nmaxind) const
 
     fprintf(fp, "  ");
     printInfo(fp);
-    fprintf(fp, "    ");
-    gmx_ana_index_dump(fp, p.g, -1, nmaxind);
+    fprintf(fp, "    Group ");
+    gmx_ana_index_dump(fp, p.g, nmaxind);
 
     fprintf(fp, "    Block (size=%d):", p.m.mapb.nr);
     if (!p.m.mapb.index)
index e165c53c2e6b05d602fb08591990c1b14d1fd52d..656ccb78c78b52bb3001f9e756460f2b03014737 100644 (file)
@@ -305,7 +305,8 @@ void SelectionCollection::Impl::resolveExternalGroups(
 
     if (root->type == SEL_GROUPREF)
     {
-        bool bOk = true;
+        bool        bOk = true;
+        std::string foundName;
         if (grps_ == NULL)
         {
             // TODO: Improve error messages
@@ -315,20 +316,18 @@ void SelectionCollection::Impl::resolveExternalGroups(
         else if (root->u.gref.name != NULL)
         {
             char *name = root->u.gref.name;
-            if (!gmx_ana_indexgrps_find(&root->u.cgrp, grps_, name))
+            bOk = gmx_ana_indexgrps_find(&root->u.cgrp, &foundName, grps_, name);
+            sfree(name);
+            root->u.gref.name = NULL;
+            if (!bOk)
             {
                 // TODO: Improve error messages
                 errors->append("Unknown group referenced in a selection");
-                bOk = false;
-            }
-            else
-            {
-                sfree(name);
             }
         }
         else
         {
-            if (!gmx_ana_indexgrps_extract(&root->u.cgrp, grps_,
+            if (!gmx_ana_indexgrps_extract(&root->u.cgrp, &foundName, grps_,
                                            root->u.gref.id))
             {
                 // TODO: Improve error messages
@@ -339,7 +338,7 @@ void SelectionCollection::Impl::resolveExternalGroups(
         if (bOk)
         {
             root->type = SEL_CONST;
-            root->setName(root->u.cgrp.name);
+            root->setName(foundName);
         }
     }
 
@@ -453,7 +452,7 @@ SelectionCollection::setTopology(t_topology *top, int natoms)
     }
     gmx_ana_selcollection_t *sc = &impl_->sc_;
     // Do this first, as it allocates memory, while the others don't throw.
-    gmx_ana_index_init_simple(&sc->gall, natoms, NULL);
+    gmx_ana_index_init_simple(&sc->gall, natoms);
     sc->pcc.setTopology(top);
     sc->top = top;
 }
index 85a6249a3e6011db6704b65cbb1565a57850bd53..91896fa9f2b0783d89b8ecfe57b125a95dadb433 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -379,7 +379,7 @@ TrajectoryAnalysisRunnerCommon::initFirstFrame()
                                                      "Trajectory (%d atoms) does not match topology (%d atoms)",
                                                      impl_->fr->natoms, top.topology()->atoms.nr)));
         }
-        // Check index groups if they have been initialized based on the topology.
+        // TODO: Check index groups if they have been initialized based on the topology.
         /*
            if (top)
            {