Clean-up selection poscalc implementation.
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 4 Apr 2012 16:10:05 +0000 (19:10 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Tue, 10 Apr 2012 19:32:21 +0000 (22:32 +0300)
- Replace gmx_ana_poscalc_coll_t and related functions with a C++ class.
  This resolves a few exception-related TODOs, and makes the C++ code
  somewhat cleaner.
- Replace gmx_ana_poscalc_get_type_enum() with a static array, as the
  current uses no longer need to modify the returned array.
- Remove possibility to provide an external position calculation
  collection for selections; it was not used, and the code is simpler
  without.
- Updated Doxygen comments.

Related to #867; it's easier to document the code when it behaves
better.

Change-Id: Ia8d5ba4ad22272b14bdcd521793636892821b140

13 files changed:
src/gromacs/selection/compiler.cpp
src/gromacs/selection/parsetree.cpp
src/gromacs/selection/poscalc.cpp
src/gromacs/selection/poscalc.h
src/gromacs/selection/selectioncollection-impl.h
src/gromacs/selection/selectioncollection.cpp
src/gromacs/selection/selectioncollection.h
src/gromacs/selection/selmethod.h
src/gromacs/selection/sm_position.cpp
src/gromacs/selection/symrec.cpp
src/gromacs/selection/tests/selectioncollection.cpp
src/gromacs/selection/tests/selectionoption.cpp
src/gromacs/trajectoryanalysis/cmdlinerunner.cpp

index 1e55df0f4df1fc6e2ac96ef81c6b4ee9637ae998..c7b353104618a1e1d32ca7eb0b7b1cd15c553b93 100644 (file)
@@ -2451,7 +2451,7 @@ postprocess_item_subexpressions(t_selelem *sel)
  * the method also defines the \c gmx_ana_selmethod_t::update method.
  */
 static void
-init_item_comg(t_selelem *sel, gmx_ana_poscalc_coll_t *pcc,
+init_item_comg(t_selelem *sel, gmx::PositionCalculationCollection *pcc,
                e_poscalc_t type, int flags)
 {
     t_selelem *child;
@@ -2471,7 +2471,7 @@ init_item_comg(t_selelem *sel, gmx_ana_poscalc_coll_t *pcc,
             if (!sel->u.expr.pc)
             {
                 cflags |= flags;
-                gmx_ana_poscalc_create(&sel->u.expr.pc, pcc, type, cflags);
+                sel->u.expr.pc = pcc->createCalculation(type, cflags);
             }
             else
             {
@@ -2706,13 +2706,14 @@ SelectionCompiler::compile(SelectionCollection *coll)
      * compilation. */
     /* By default, use whole residues/molecules. */
     flags = POS_COMPLWHOLE;
-    gmx_ana_poscalc_type_from_enum(coll->_impl->_rpost.c_str(), &post, &flags);
+    PositionCalculationCollection::typeFromEnum(coll->_impl->_rpost.c_str(),
+                                                &post, &flags);
     item = sc->root;
     while (item)
     {
         init_root_item(item, &sc->gall);
         postprocess_item_subexpressions(item);
-        init_item_comg(item, sc->pcc, post, flags);
+        init_item_comg(item, &sc->pcc, post, flags);
         free_item_compilerdata(item);
         item = item->next;
     }
index bcfdecc9c61691bcf58de9141af3b24339027bc3..be3c1501855e70595f492e7fa1d85d052f63289b 100644 (file)
@@ -537,7 +537,7 @@ _gmx_selelem_init_method_params(t_selelem *sel, yyscan_t scanner)
     {
         gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
 
-        sel->u.expr.method->set_poscoll(sc->pcc, mdata);
+        sel->u.expr.method->set_poscoll(&sc->pcc, mdata);
     }
     /* Store the values */
     sel->u.expr.method->param = param;
@@ -574,8 +574,8 @@ _gmx_selelem_set_method(t_selelem *sel, gmx_ana_selmethod_t *method,
  * \returns       0 on success, a non-zero error code on error.
  */
 static int
-set_refpos_type(gmx_ana_poscalc_coll_t *pcc, t_selelem *sel, const char *rpost,
-                yyscan_t scanner)
+set_refpos_type(gmx::PositionCalculationCollection *pcc, t_selelem *sel,
+                const char *rpost, yyscan_t scanner)
 {
     if (!rpost)
     {
@@ -591,8 +591,8 @@ set_refpos_type(gmx_ana_poscalc_coll_t *pcc, t_selelem *sel, const char *rpost,
         try
         {
             /* By default, use whole residues/molecules. */
-            gmx_ana_poscalc_create_enum(&sel->u.expr.pc, pcc, rpost,
-                                        POS_COMPLWHOLE);
+            sel->u.expr.pc
+                = pcc->createCalculationFromEnum(rpost, POS_COMPLWHOLE);
         }
         catch (const gmx::GromacsException &ex)
         {
@@ -770,7 +770,7 @@ _gmx_sel_init_keyword(gmx_ana_selmethod_t *method, t_selexpr_value *args,
             goto on_error;
         }
     }
-    rc = set_refpos_type(sc->pcc, child, rpost, scanner);
+    rc = set_refpos_type(&sc->pcc, child, rpost, scanner);
     if (rc != 0)
     {
         goto on_error;
@@ -829,7 +829,7 @@ _gmx_sel_init_method(gmx_ana_selmethod_t *method, t_selexpr_param *params,
         _gmx_selelem_free(root);
         return NULL;
     }
-    rc = set_refpos_type(sc->pcc, root, rpost, scanner);
+    rc = set_refpos_type(&sc->pcc, root, rpost, scanner);
     if (rc != 0)
     {
         _gmx_selelem_free(root);
index 0d4aed700e2efe4174df93993f03da4aa1e11b05..027fbbe3f3434f134892cdd1137e3de33790f837 100644 (file)
  *
  * For more info, check our website at http://www.gromacs.org
  */
-/*! \internal
- * \page poscalcengine Position calculation engine
- *
- * The header file \ref poscalc.h defines an API for calculating positions
- * in an automated way. This is useful mostly in the selection engine, in
- * particular with dynamic selections, because the same COM/COG positions
- * may be needed in several contexts. The API makes it possible to
- * optimize the evaluation such that any heavy calculation is only done once,
- * and the results just copied if needed more than once.
- * The functions also provide a convenient interface for keeping the whole
- * \c gmx_ana_pos_t structure up-to-date.
- *
- * A new collection of position calculations is allocated with
- * gmx_ana_poscalc_coll_create().
- * Calculations within one collection should share the same topology, and
- * they are optimized. Calculations in different collections do not interact.
- * The topology for a collection can be set with
- * gmx_ana_poscalc_coll_set_topology().
- * This needs to be done before calling gmx_ana_poscalc_set_maxindex() for
- * any calculation in the collection, unless that calculation does not
- * require topology information.
- * All memory allocated for a collection and the calculations in it can be
- * freed with gmx_ana_poscalc_coll_free().
- *
- * A new calculation is created with gmx_ana_poscalc_create().
- * If flags need to be adjusted later, gmx_ana_poscalc_set_flags() can be
- * used.
- * After the flags are final, the largest possible index group for which the
- * positions are needed has to be set with gmx_ana_poscalc_set_maxindex().
- * gmx_ana_poscalc_coll_set_topology() should have been called before this
- * function is called.
- * After the above calls, gmx_ana_poscalc_init_pos() can be used to initialize
- * output to a \c gmx_ana_pos_t structure. Several different structures can be
- * initialized for the same calculation; the only requirement is that the
- * structure passed later to gmx_ana_poscalc_update() has been initialized
- * properly.
- * The memory allocated for a calculation can be freed with
- * gmx_ana_poscalc_free().
- *
- * The position evaluation is simple: gmx_ana_poscalc_init_frame() should be
- * called once for each frame, and gmx_ana_poscalc_update() can then be called
- * for each calculation that is needed for that frame.
- *
- * It is also possible to initialize the calculations based on a type provided
- * as a string.
- * The possible strings are returned by gmx_ana_poscalc_create_type_enum(),
- * and the string can be converted to the parameters for
- * gmx_ana_poscalc_create() using gmx_ana_poscalc_type_from_enum().
- * gmx_ana_poscalc_create_enum() is also provided for convenience.
- */
 /*! \internal \file
  * \brief
- * Implements functions in poscalc.h.
+ * Implements gmx::PositionCalculationCollection and functions in poscalc.h.
  *
  * \todo
  * There is probably some room for optimization in the calculation of
 #include <vec.h>
 
 #include "gromacs/fatalerror/exceptions.h"
+#include "gromacs/fatalerror/gmxassert.h"
 #include "gromacs/selection/centerofmass.h"
 #include "gromacs/selection/indexutil.h"
 #include "gromacs/selection/poscalc.h"
 #include "gromacs/selection/position.h"
 
+namespace gmx
+{
+
 /*! \internal \brief
- * Collection of \c gmx_ana_poscalc_t structures for the same topology.
+ * Private implementation class for PositionCalculationCollection.
  *
- * Calculations within the same structure are optimized to eliminate duplicate
- * calculations.
+ * \ingroup module_selection
  */
-struct gmx_ana_poscalc_coll_t
+class PositionCalculationCollection::Impl
 {
-    /*! \brief
-     * Topology data.
-     *
-     * Can be NULL if none of the calculations require topology data or if
-     * gmx_ana_poscalc_coll_set_topology() has not been called.
-     */
-    t_topology               *top;
-    /** Pointer to the first data structure. */
-    gmx_ana_poscalc_t        *first;
-    /** Pointer to the last data structure. */
-    gmx_ana_poscalc_t        *last;
-    /** Whether the collection has been initialized for evaluation. */
-    bool                      bInit;
+    public:
+        Impl();
+        ~Impl();
+
+        /*! \brief
+         * Inserts a position calculation structure into its collection.
+         *
+         * \param pc     Data structure to insert.
+         * \param before Data structure before which to insert
+         *   (NULL = insert at end).
+         *
+         * Inserts \p pc to its collection before \p before.
+         * If \p before is NULL, \p pc is appended to the list.
+         */
+        void insertCalculation(gmx_ana_poscalc_t *pc, gmx_ana_poscalc_t *before);
+        /*! \brief
+         * Removes a position calculation structure from its collection.
+         *
+         * \param pc    Data structure to remove.
+         *
+         * Removes \p pc from its collection.
+         */
+        void removeCalculation(gmx_ana_poscalc_t *pc);
+
+        /*! \copydoc PositionCalculationCollection::createCalculation()
+         *
+         * This method contains the actual implementation of the similarly
+         * named method in the parent class.
+         */
+        gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
+
+        /*! \brief
+         * Topology data.
+         *
+         * Can be NULL if none of the calculations require topology data or if
+         * setTopology() has not been called.
+         */
+        t_topology               *top_;
+        //! Pointer to the first data structure.
+        gmx_ana_poscalc_t        *first_;
+        //! Pointer to the last data structure.
+        gmx_ana_poscalc_t        *last_;
+        //! Whether the collection has been initialized for evaluation.
+        bool                      bInit_;
 };
 
+} // namespace gmx
+
 /*! \internal \brief
  * Data structure for position calculation.
  */
@@ -211,11 +197,10 @@ struct gmx_ana_poscalc_t
     /** Number of references to this structure. */
     int                       refcount;
     /** Collection this calculation belongs to. */
-    gmx_ana_poscalc_coll_t   *coll;
+    gmx::PositionCalculationCollection::Impl *coll;
 };
 
-//! Strings returned by gmx_ana_poscalc_create_type_enum().
-static const char *const poscalc_enum_strings[] = {
+const char * const gmx::PositionCalculationCollection::typeEnumValues[] = {
     "atom",
     "res_com",       "res_cog",
     "mol_com",       "mol_cog",
@@ -227,8 +212,6 @@ static const char *const poscalc_enum_strings[] = {
     "dyn_mol_com",   "dyn_mol_cog",
     NULL,
 };
-//! Number of elements in poscalc_enum_strings.
-#define NENUM asize(poscalc_enum_strings)
 
 /*! \brief
  * Returns the partition type for a given position type.
@@ -250,28 +233,14 @@ index_type_for_poscalc(e_poscalc_t type)
     return INDEX_UNKNOWN;
 }
 
-/*!
- * \param[in]     post  String (typically an enum command-line argument).
- *   Allowed values: 'atom', 'res_com', 'res_cog', 'mol_com', 'mol_cog',
- *   or one of the last four prepended by 'whole_', 'part_', or 'dyn_'.
- * \param[out]    type  \c e_poscalc_t corresponding to \p post.
- * \param[in,out] flags Flags corresponding to \p post.
- *   On input, the flags should contain the default flags.
- *   On exit, the flags \ref POS_MASS, \ref POS_COMPLMAX and
- *   \ref POS_COMPLWHOLE have been set according to \p post
- *   (the completion flags are left at the default values if no completion
- *   prefix is given).
- * \throws  InternalError  if post is not recognized.
- *
- * \attention
- * Checking is not complete, and other values than those listed above
- * may be accepted for \p post, but the results are undefined.
- */
-void
-gmx_ana_poscalc_type_from_enum(const char *post, e_poscalc_t *type, int *flags)
+namespace gmx
 {
-    const char *ptr;
 
+// static
+void
+PositionCalculationCollection::typeFromEnum(const char *post,
+                                            e_poscalc_t *type, int *flags)
+{
     if (post[0] == 'a')
     {
         *type   = POS_ATOM;
@@ -280,7 +249,7 @@ gmx_ana_poscalc_type_from_enum(const char *post, e_poscalc_t *type, int *flags)
     }
 
     /* Process the prefix */
-    ptr = post;
+    const char *ptr = post;
     if (post[0] == 'w')
     {
         *flags &= ~POS_COMPLMAX;
@@ -309,11 +278,11 @@ gmx_ana_poscalc_type_from_enum(const char *post, e_poscalc_t *type, int *flags)
     }
     else
     {
-        GMX_THROW(gmx::InternalError("Unknown position calculation type"));
+        GMX_THROW(InternalError("Unknown position calculation type"));
     }
     if (strlen(ptr) < 7)
     {
-        GMX_THROW(gmx::InternalError("Unknown position calculation type"));
+        GMX_THROW(InternalError("Unknown position calculation type"));
     }
     if (ptr[6] == 'm')
     {
@@ -325,109 +294,129 @@ gmx_ana_poscalc_type_from_enum(const char *post, e_poscalc_t *type, int *flags)
     }
     else
     {
-        GMX_THROW(gmx::InternalError("Unknown position calculation type"));
+        GMX_THROW(InternalError("Unknown position calculation type"));
     }
 }
 
-/*!
- * \param[in]  bAtom    If true, the "atom" value is included.
- * \returns    NULL-terminated array of strings that contains the string
- *   values acceptable for gmx_ana_poscalc_type_from_enum().
- *
- * The first string in the returned list is always NULL to allow the list to
- * be used with Gromacs command-line parsing.
+/********************************************************************
+ * PositionCalculationCollection::Impl
  */
-const char **
-gmx_ana_poscalc_create_type_enum(bool bAtom)
+
+PositionCalculationCollection::Impl::Impl()
+    : top_(NULL), first_(NULL), last_(NULL), bInit_(false)
 {
-    const char **pcenum;
-    size_t       i;
+}
 
-    if (bAtom)
+PositionCalculationCollection::Impl::~Impl()
+{
+    // Loop backwards, because there can be internal references in that are
+    // correctly handled by this direction.
+    while (last_ != NULL)
     {
-        snew(pcenum, NENUM+1);
-        for (i = 0; i < NENUM; ++i)
+        GMX_ASSERT(last_->refcount == 1,
+                   "Dangling references to position calculations");
+        gmx_ana_poscalc_free(last_);
+    }
+}
+
+void
+PositionCalculationCollection::Impl::insertCalculation(gmx_ana_poscalc_t *pc,
+                                                       gmx_ana_poscalc_t *before)
+{
+    GMX_RELEASE_ASSERT(pc->coll == this, "Inconsistent collections");
+    if (before == NULL)
+    {
+        pc->next = NULL;
+        pc->prev = last_;
+        if (last_ != NULL)
         {
-            pcenum[i+1] = poscalc_enum_strings[i];
+            last_->next = pc;
         }
+        last_ = pc;
     }
     else
     {
-        snew(pcenum, NENUM+1-1);
-        for (i = 1; i < NENUM; ++i)
+        pc->prev     = before->prev;
+        pc->next     = before;
+        if (before->prev)
         {
-            pcenum[i] = poscalc_enum_strings[i];
+            before->prev->next = pc;
         }
+        before->prev = pc;
+    }
+    if (pc->prev == NULL)
+    {
+        first_ = pc;
     }
-    pcenum[0] = NULL;
-    return pcenum;
 }
 
-/*!
- * \param[out] pccp   Allocated position calculation collection.
- * \returns    0 for success.
- */
-int
-gmx_ana_poscalc_coll_create(gmx_ana_poscalc_coll_t **pccp)
+void
+PositionCalculationCollection::Impl::removeCalculation(gmx_ana_poscalc_t *pc)
 {
-    gmx_ana_poscalc_coll_t *pcc;
-
-    snew(pcc, 1);
-    pcc->top   = NULL;
-    pcc->first = NULL;
-    pcc->last  = NULL;
-    pcc->bInit = false;
-    *pccp = pcc;
-    return 0;
+    GMX_RELEASE_ASSERT(pc->coll == this, "Inconsistent collections");
+    if (pc->prev != NULL)
+    {
+        pc->prev->next = pc->next;
+    }
+    else if (pc == first_)
+    {
+        first_ = pc->next;
+    }
+    if (pc->next != NULL)
+    {
+        pc->next->prev = pc->prev;
+    }
+    else if (pc == last_)
+    {
+        last_ = pc->prev;
+    }
+    pc->prev = pc->next = NULL;
 }
 
-/*!
- * \param[in,out] pcc   Position calculation collection data structure.
- * \param[in]     top   Topology data structure.
- *
- * This function should be called to set the topology before using
- * gmx_ana_poscalc_set_maxindex() for any calculation that requires
- * topology information.
- */
-void
-gmx_ana_poscalc_coll_set_topology(gmx_ana_poscalc_coll_t *pcc, t_topology *top)
+gmx_ana_poscalc_t *
+PositionCalculationCollection::Impl::createCalculation(e_poscalc_t type, int flags)
 {
-    pcc->top = top;
+    gmx_ana_poscalc_t *pc;
+
+    snew(pc, 1);
+    pc->type     = type;
+    pc->itype    = index_type_for_poscalc(type);
+    gmx_ana_poscalc_set_flags(pc, flags);
+    pc->refcount = 1;
+    pc->coll     = this;
+    insertCalculation(pc, NULL);
+    return pc;
 }
 
-/*!
- * \param[in] pcc   Position calculation collection to free.
- *
- * The pointer \p pcc is invalid after the call.
- * Any calculations in the collection are also freed, no matter how many
- * references to them are left.
+
+/********************************************************************
+ * PositionCalculationCollection
  */
+
+PositionCalculationCollection::PositionCalculationCollection()
+    : impl_(new Impl)
+{
+}
+
+PositionCalculationCollection::~PositionCalculationCollection()
+{
+}
+
 void
-gmx_ana_poscalc_coll_free(gmx_ana_poscalc_coll_t *pcc)
+PositionCalculationCollection::setTopology(t_topology *top)
 {
-    while (pcc->first)
-    {
-        gmx_ana_poscalc_free(pcc->first);
-    }
-    sfree(pcc);
+    impl_->top_ = top;
 }
 
-/*!
- * \param[in] fp    File handle to receive the output.
- * \param[in] pcc   Position calculation collection to print.
- *
- * The output is very technical, making this function mainly useful for
- * debugging purposes.
- */
 void
-gmx_ana_poscalc_coll_print_tree(FILE *fp, gmx_ana_poscalc_coll_t *pcc)
+PositionCalculationCollection::printTree(FILE *fp) const
 {
     gmx_ana_poscalc_t *pc;
     int                i, j;
 
     fprintf(fp, "Position calculations:\n");
     i  = 1;
-    pc = pcc->first;
+    pc = impl_->first_;
     while (pc)
     {
         fprintf(fp, "%2d ", i);
@@ -534,7 +523,7 @@ gmx_ana_poscalc_coll_print_tree(FILE *fp, gmx_ana_poscalc_coll_t *pcc)
 
             fprintf(fp, "   Base: ");
             j = 1;
-            base = pcc->first;
+            base = impl_->first_;
             while (base && base != pc->sbase)
             {
                 ++j;
@@ -556,74 +545,86 @@ gmx_ana_poscalc_coll_print_tree(FILE *fp, gmx_ana_poscalc_coll_t *pcc)
     }
 }
 
-/*! \brief
- * Inserts a position calculation structure into its collection.
- *
- * \param pc     Data structure to insert.
- * \param before Data structure before which to insert
- *   (NULL = insert at end).
- *
- * Inserts \p pc to its collection before \p before.
- * If \p before is NULL, \p pc is appended to the list.
- */
-static void
-insert_poscalc(gmx_ana_poscalc_t *pc, gmx_ana_poscalc_t *before)
+gmx_ana_poscalc_t *
+PositionCalculationCollection::createCalculation(e_poscalc_t type, int flags)
 {
-    if (before == NULL)
+    return impl_->createCalculation(type, flags);
+}
+
+gmx_ana_poscalc_t *
+PositionCalculationCollection::createCalculationFromEnum(const char *post, int flags)
+{
+    e_poscalc_t  type;
+    int cflags = flags;
+    typeFromEnum(post, &type, &cflags);
+    return impl_->createCalculation(type, cflags);
+}
+
+void PositionCalculationCollection::initEvaluation()
+{
+    if (impl_->bInit_)
     {
-        pc->next = NULL;
-        pc->prev = pc->coll->last;
-        if (pc->coll->last)
-        {
-            pc->coll->last->next = pc;
-        }
-        pc->coll->last = pc;
+        return;
     }
-    else
+    gmx_ana_poscalc_t *pc = impl_->first_;
+    while (pc)
     {
-        pc->prev     = before->prev;
-        pc->next     = before;
-        if (before->prev)
+        /* Initialize position storage for base calculations */
+        if (pc->p)
         {
-            before->prev->next = pc;
+            gmx_ana_poscalc_init_pos(pc, pc->p);
         }
-        before->prev = pc;
-    }
-    if (!pc->prev)
-    {
-        pc->coll->first = pc;
+        /* Construct the mapping of the base positions */
+        if (pc->sbase)
+        {
+            int                     bi, bj;
+
+            snew(pc->baseid, pc->b.nr);
+            for (bi = bj = 0; bi < pc->b.nr; ++bi, ++bj)
+            {
+                while (pc->sbase->b.a[pc->sbase->b.index[bj]] != pc->b.a[pc->b.index[bi]])
+                {
+                    ++bj;
+                }
+                pc->baseid[bi] = bj;
+            }
+        }
+        /* Free the block data for dynamic calculations */
+        if (pc->flags & POS_DYNAMIC)
+        {
+            if (pc->b.nalloc_index > 0)
+            {
+                sfree(pc->b.index);
+                pc->b.nalloc_index = 0;
+            }
+            if (pc->b.nalloc_a > 0)
+            {
+                sfree(pc->b.a);
+                pc->b.nalloc_a = 0;
+            }
+        }
+        pc = pc->next;
     }
+    impl_->bInit_ = true;
 }
 
-/*! \brief
- * Removes a position calculation structure from its collection.
- *
- * \param pc    Data structure to remove.
- *
- * Removes \p pc from its collection.
- */
-static void
-remove_poscalc(gmx_ana_poscalc_t *pc)
+void PositionCalculationCollection::initFrame()
 {
-    if (pc->prev)
+    if (!impl_->bInit_)
     {
-        pc->prev->next = pc->next;
-    }
-    else if (pc == pc->coll->first)
-    {
-        pc->coll->first = pc->next;
-    }
-    if (pc->next)
-    {
-        pc->next->prev = pc->prev;
+        initEvaluation();
     }
-    else if (pc == pc->coll->last)
+    /* Clear the evaluation flags */
+    gmx_ana_poscalc_t *pc = impl_->first_;
+    while (pc)
     {
-        pc->coll->last = pc->prev;
+        pc->bEval = false;
+        pc = pc->next;
     }
-    pc->prev = pc->next = NULL;
 }
 
+} // namespace gmx
+
 /*! \brief
  * Initializes position calculation using the maximum possible input index.
  *
@@ -636,7 +637,8 @@ remove_poscalc(gmx_ana_poscalc_t *pc)
 static void
 set_poscalc_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g, bool bBase)
 {
-    gmx_ana_index_make_block(&pc->b, pc->coll->top, g, pc->itype, pc->flags & POS_COMPLWHOLE);
+    t_topology *top = pc->coll->top_;
+    gmx_ana_index_make_block(&pc->b, top, g, pc->itype, pc->flags & POS_COMPLWHOLE);
     /* Set the type to POS_ATOM if the calculation in fact is such. */
     if (pc->b.nr == pc->b.nra)
     {
@@ -648,7 +650,7 @@ set_poscalc_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g, bool bBase)
     if (!(pc->flags & POS_COMPLWHOLE)
         && (!(pc->flags & POS_DYNAMIC) || (pc->flags & POS_COMPLMAX))
         && (pc->type == POS_RES || pc->type == POS_MOL)
-        && gmx_ana_index_has_complete_elems(g, pc->itype, pc->coll->top))
+        && gmx_ana_index_has_complete_elems(g, pc->itype, top))
     {
         pc->flags &= ~POS_COMPLMAX;
         pc->flags |= POS_COMPLWHOLE;
@@ -768,14 +770,14 @@ create_simple_base(gmx_ana_poscalc_t *pc)
     int                flags;
 
     flags = pc->flags & ~(POS_DYNAMIC | POS_MASKONLY);
-    gmx_ana_poscalc_create(&base, pc->coll, pc->type, flags);
+    base = pc->coll->createCalculation(pc->type, flags);
     set_poscalc_maxindex(base, &pc->gmax, true);
 
     snew(base->p, 1);
 
     pc->sbase = base;
-    remove_poscalc(base);
-    insert_poscalc(base, pc);
+    pc->coll->removeCalculation(base);
+    pc->coll->insertCalculation(base, pc);
 
     return base;
 }
@@ -877,9 +879,9 @@ merge_bases(gmx_ana_poscalc_t *tbase, gmx_ana_poscalc_t *mbase)
     gmx_ana_poscalc_t *pc;
 
     merge_to_base(tbase, mbase);
-    remove_poscalc(mbase);
+    mbase->coll->removeCalculation(mbase);
     /* Set tbase as the base for all calculations that had mbase */
-    pc = tbase->coll->first;
+    pc = tbase->coll->first_;
     while (pc)
     {
         if (pc->sbase == mbase)
@@ -915,7 +917,7 @@ setup_base(gmx_ana_poscalc_t *pc)
     gmx_ana_index_clear(&g);
     gmx_ana_index_reserve(&g, pc->b.nra);
     pbase = pc;
-    base = pc->coll->first;
+    base = pc->coll->first_;
     while (base)
     {
         /* Save the next calculation so that we can safely delete base */
@@ -979,57 +981,6 @@ setup_base(gmx_ana_poscalc_t *pc)
     }
 }
 
-/*!
- * \param[out] pcp   Position calculation data structure pointer to initialize.
- * \param[in,out] pcc   Position calculation collection.
- * \param[in]  type  Type of calculation.
- * \param[in]  flags Flags for setting calculation options
- *   (see \ref poscalc_flags "documentation of the flags").
- */
-void
-gmx_ana_poscalc_create(gmx_ana_poscalc_t **pcp, gmx_ana_poscalc_coll_t *pcc,
-                       e_poscalc_t type, int flags)
-{
-    gmx_ana_poscalc_t *pc;
-
-    snew(pc, 1);
-    pc->type     = type;
-    pc->itype    = index_type_for_poscalc(type);
-    gmx_ana_poscalc_set_flags(pc, flags);
-    pc->refcount = 1;
-    pc->coll     = pcc;
-    insert_poscalc(pc, NULL);
-    *pcp = pc;
-}
-
-/*!
- * \param[out] pcp   Position calculation data structure pointer to initialize.
- * \param[in,out] pcc   Position calculation collection.
- * \param[in]  post  One of the strings acceptable for
- *   gmx_ana_poscalc_type_from_enum().
- * \param[in]  flags Flags for setting calculation options
- *   (see \ref poscalc_flags "documentation of the flags").
- * \returns    0 on success, a non-zero error value on error.
- *
- * This is a convenience wrapper for gmx_ana_poscalc_create().
- * \p flags sets the default calculation options if not overridden by \p post;
- * see gmx_ana_poscalc_type_from_enum().
- *
- * \see gmx_ana_poscalc_create(), gmx_ana_poscalc_type_from_enum()
- */
-void
-gmx_ana_poscalc_create_enum(gmx_ana_poscalc_t **pcp, gmx_ana_poscalc_coll_t *pcc,
-                            const char *post, int flags)
-{
-    e_poscalc_t  type;
-    int          cflags;
-
-    cflags = flags;
-    *pcp = NULL;
-    gmx_ana_poscalc_type_from_enum(post, &type, &cflags);
-    gmx_ana_poscalc_create(pcp, pcc, type, cflags);
-}
-
 /*!
  * \param[in,out] pc    Position calculation data structure.
  * \param[in]     flags New flags.
@@ -1089,7 +1040,7 @@ gmx_ana_poscalc_set_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g)
 void
 gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p)
 {
-    gmx_ana_indexmap_init(&p->m, &pc->gmax, pc->coll->top, pc->itype);
+    gmx_ana_indexmap_init(&p->m, &pc->gmax, pc->coll->top_, pc->itype);
     /* Only do the static optimization when there is no completion */
     if (!(pc->flags & POS_DYNAMIC) && pc->b.nra == pc->gmax.isize)
     {
@@ -1127,7 +1078,7 @@ gmx_ana_poscalc_free(gmx_ana_poscalc_t *pc)
         return;
     }
 
-    remove_poscalc(pc);
+    pc->coll->removeCalculation(pc);
     if (pc->b.nalloc_index > 0)
     {
         sfree(pc->b.index);
@@ -1167,99 +1118,6 @@ gmx_ana_poscalc_requires_top(gmx_ana_poscalc_t *pc)
     return false;
 }
 
-/*!
- * \param[in,out] pcc Position calculation collection to initialize.
- *
- * This function does some final initialization of the data structures in the
- * collection to prepare them for evaluation.
- * After this function has been called, it is no longer possible to add new
- * calculations to the collection.
- *
- * This function is automatically called by gmx_ana_poscalc_init_frame() 
- * if not called by the user earlier.
- * Multiple calls to the function are ignored.
- */
-void
-gmx_ana_poscalc_init_eval(gmx_ana_poscalc_coll_t *pcc)
-{
-    gmx_ana_poscalc_t      *pc;
-    int                     bi, bj;
-
-    if (pcc->bInit)
-    {
-        return;
-    }
-    pc = pcc->first;
-    while (pc)
-    {
-        /* Initialize position storage for base calculations */
-        if (pc->p)
-        {
-            gmx_ana_poscalc_init_pos(pc, pc->p);
-        }
-        /* Construct the mapping of the base positions */
-        if (pc->sbase)
-        {
-            snew(pc->baseid, pc->b.nr);
-            for (bi = bj = 0; bi < pc->b.nr; ++bi, ++bj)
-            {
-                while (pc->sbase->b.a[pc->sbase->b.index[bj]] != pc->b.a[pc->b.index[bi]])
-                {
-                    ++bj;
-                }
-                pc->baseid[bi] = bj;
-            }
-        }
-        /* Free the block data for dynamic calculations */
-        if (pc->flags & POS_DYNAMIC)
-        {
-            if (pc->b.nalloc_index > 0)
-            {
-                sfree(pc->b.index);
-                pc->b.nalloc_index = 0;
-            }
-            if (pc->b.nalloc_a > 0)
-            {
-                sfree(pc->b.a);
-                pc->b.nalloc_a = 0;
-            }
-        }
-        pc = pc->next;
-    }
-    pcc->bInit = true;
-}
-
-/*!
- * \param[in,out] pcc Position calculation collection to initialize.
- *
- * Clears the evaluation flag for all calculations.
- * Should be called for each frame before calling gmx_ana_poscalc_update().
- *
- * This function is automatically called by gmx_ana_do() for each
- * frame, and should not be called by the user unless gmx_ana_do() is
- * not being used.
- *
- * This function calls gmx_ana_poscalc_init_eval() automatically if it has
- * not been called earlier.
- */
-void
-gmx_ana_poscalc_init_frame(gmx_ana_poscalc_coll_t *pcc)
-{
-    gmx_ana_poscalc_t      *pc;
-
-    if (!pcc->bInit)
-    {
-        gmx_ana_poscalc_init_eval(pcc);
-    }
-    /* Clear the evaluation flags */
-    pc = pcc->first;
-    while (pc)
-    {
-        pc->bEval = false;
-        pc = pc->next;
-    }
-}
-
 /*!
  * \param[in]     pc   Position calculation data.
  * \param[in,out] p    Output positions, initialized previously with
@@ -1391,6 +1249,8 @@ gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p,
         }
         /* Here, we assume that the topology has been properly initialized,
          * and do not check the return values of gmx_calc_comg*(). */
+        t_topology *top = pc->coll->top_;
+        bool bMass = pc->flags & POS_MASS;
         switch (pc->type)
         {
         case POS_ATOM:
@@ -1414,45 +1274,36 @@ gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p,
             }
             break;
         case POS_ALL:
-            gmx_calc_comg(pc->coll->top, fr->x, pc->b.nra, pc->b.a,
-                          pc->flags & POS_MASS, p->x[0]);
+            gmx_calc_comg(top, fr->x, pc->b.nra, pc->b.a, bMass, p->x[0]);
             if (p->v && fr->bV)
             {
-                gmx_calc_comg(pc->coll->top, fr->v, pc->b.nra, pc->b.a,
-                              pc->flags & POS_MASS, p->v[0]);
+                gmx_calc_comg(top, fr->v, pc->b.nra, pc->b.a, bMass, p->v[0]);
             }
             if (p->f && fr->bF)
             {
-                gmx_calc_comg_f(pc->coll->top, fr->f, pc->b.nra, pc->b.a,
-                                pc->flags & POS_MASS, p->f[0]);
+                gmx_calc_comg_f(top, fr->f, pc->b.nra, pc->b.a, bMass, p->f[0]);
             }
             break;
         case POS_ALL_PBC:
-            gmx_calc_comg_pbc(pc->coll->top, fr->x, pbc, pc->b.nra, pc->b.a,
-                              pc->flags & POS_MASS, p->x[0]);
+            gmx_calc_comg_pbc(top, fr->x, pbc, pc->b.nra, pc->b.a, bMass, p->x[0]);
             if (p->v && fr->bV)
             {
-                gmx_calc_comg(pc->coll->top, fr->v, pc->b.nra, pc->b.a,
-                              pc->flags & POS_MASS, p->v[0]);
+                gmx_calc_comg(top, fr->v, pc->b.nra, pc->b.a, bMass, p->v[0]);
             }
             if (p->f && fr->bF)
             {
-                gmx_calc_comg_f(pc->coll->top, fr->f, pc->b.nra, pc->b.a,
-                                pc->flags & POS_MASS, p->f[0]);
+                gmx_calc_comg_f(top, fr->f, pc->b.nra, pc->b.a, bMass, p->f[0]);
             }
             break;
         default:
-            gmx_calc_comg_blocka(pc->coll->top, fr->x, &pc->b,
-                                 pc->flags & POS_MASS, p->x);
+            gmx_calc_comg_blocka(top, fr->x, &pc->b, bMass, p->x);
             if (p->v && fr->bV)
             {
-                gmx_calc_comg_blocka(pc->coll->top, fr->v, &pc->b,
-                                     pc->flags & POS_MASS, p->v);
+                gmx_calc_comg_blocka(top, fr->v, &pc->b, bMass, p->v);
             }
             if (p->f && fr->bF)
             {
-                gmx_calc_comg_blocka(pc->coll->top, fr->f, &pc->b,
-                                     pc->flags & POS_MASS, p->f);
+                gmx_calc_comg_blocka(top, fr->f, &pc->b, bMass, p->f);
             }
             break;
         }
index db39e86f8b9b7e6a76abca34a8fc452e14423685..add791a2229e9df9ad59fd8c4a9a8126c42d520b 100644 (file)
  * For more info, check our website at http://www.gromacs.org
  */
 /*! \internal \file
- * \brief API for structured and optimized calculation of positions.
+ * \brief
+ * API for structured and optimized calculation of positions.
  *
- * The functions in this header are used internally by the analysis library
- * to calculate positions.
- * They can also be used in user code, but in most cases there should be no
- * need. Instead, one should write an analysis tool such that it gets all
- * positions through selections.
+ * This header declares an API for calculating positions in an automated way,
+ * for internal use by the selection engine.  This is useful in particular with
+ * dynamic selections, because the same COM/COG positions may be needed in
+ * several contexts.  The API makes it possible to optimize the evaluation such
+ * that any heavy calculation is only done once, and the results just copied if
+ * needed more than once.  The functions also provide a convenient interface
+ * for keeping the whole \c gmx_ana_pos_t structure up-to-date.
  *
- * The API is documented in more detail on a separate page:
- * \ref poscalcengine.
+ * The API is documented in more detail in gmx::PositionCalculationCollection.
  *
  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
  * \ingroup module_selection
@@ -48,6 +50,8 @@
 
 #include "../legacyheaders/typedefs.h"
 
+#include "../utility/common.h"
+
 /*! \name Flags for position calculation.
  * \anchor poscalc_flags
  */
@@ -118,42 +122,197 @@ typedef enum
     POS_ALL_PBC  /**< Calculate center for the whole group with PBC. */
 } e_poscalc_t;
 
-/** Collection of \c gmx_ana_poscalc_t structures for the same topology. */
-typedef struct gmx_ana_poscalc_coll_t gmx_ana_poscalc_coll_t;
 /** Data structure for position calculation. */
 typedef struct gmx_ana_poscalc_t gmx_ana_poscalc_t;
 
 struct gmx_ana_index_t;
 struct gmx_ana_pos_t;
 
-/** Converts a string to parameters for gmx_ana_poscalc_create(). */
-void
-gmx_ana_poscalc_type_from_enum(const char *post, e_poscalc_t *type, int *flags);
-/** Creates a list of strings for position enum parameter handling. */
-const char **
-gmx_ana_poscalc_create_type_enum(bool bAtom);
-
-/** Creates a new position calculation collection object. */
-int
-gmx_ana_poscalc_coll_create(gmx_ana_poscalc_coll_t **pccp);
-/** Sets the topology for a position calculation collection. */
-void
-gmx_ana_poscalc_coll_set_topology(gmx_ana_poscalc_coll_t *pcc, t_topology *top);
-/** Frees memory allocated for a position calculation collection. */
-void
-gmx_ana_poscalc_coll_free(gmx_ana_poscalc_coll_t *pcc);
-/** Prints information about calculations in a position calculation collection. */
-void
-gmx_ana_poscalc_coll_print_tree(FILE *fp, gmx_ana_poscalc_coll_t *pcc);
+namespace gmx
+{
+
+/*! \internal \brief
+ * Collection of \c gmx_ana_poscalc_t structures for the same topology.
+ *
+ * Calculations within one collection share the same topology, and they are
+ * optimized.  Calculations in different collections do not interact.
+ * The topology for a collection can be set with setTopology().
+ * This needs to be done before calling gmx_ana_poscalc_set_maxindex() for
+ * any calculation in the collection, unless that calculation does not
+ * require topology information.
+ *
+ * A new calculation is created with createCalculation().
+ * If flags need to be adjusted later, gmx_ana_poscalc_set_flags() can be
+ * used.
+ * After the flags are final, the largest possible index group for which the
+ * positions are needed has to be set with gmx_ana_poscalc_set_maxindex().
+ * setTopology() should have been called before this function is called.
+ * After the above calls, gmx_ana_poscalc_init_pos() can be used to initialize
+ * output to a \c gmx_ana_pos_t structure.  Several different structures can be
+ * initialized for the same calculation; the only requirement is that the
+ * structure passed later to gmx_ana_poscalc_update() has been initialized
+ * properly.
+ * The memory allocated for a calculation can be freed with
+ * gmx_ana_poscalc_free().
+ *
+ * The position evaluation is simple: initFrame() should be
+ * called once for each frame, and gmx_ana_poscalc_update() can then be called
+ * for each calculation that is needed for that frame.
+ *
+ * It is also possible to initialize the calculations based on a type provided
+ * as a string.
+ * The possible strings are listed in \ref typeEnumValues, and the string can
+ * be converted to the parameters for createCalculation() using typeFromEnum().
+ * createCalculationFromEnum() is also provided for convenience.
+ *
+ * \ingroup module_selection
+ */
+class PositionCalculationCollection
+{
+    public:
+        /*! \brief
+         * Array of strings acceptable for position calculation type enum.
+         *
+         * This array contains the acceptable values for typeFromEnum() and
+         * createCalculationFromEnum().
+         * The array contains a NULL pointer after the last item to indicate
+         * the end of the list.
+         */
+        static const char * const typeEnumValues[];
+
+        /*! \brief
+         * Converts a string to parameters for createCalculationFromEnum().
+         *
+         * \param[in]     post  String (typically an enum argument).
+         *     Allowed values: 'atom', 'res_com', 'res_cog', 'mol_com', 'mol_cog',
+         *     or one of the last four prepended by 'whole_', 'part_', or 'dyn_'.
+         * \param[out]    type  \c e_poscalc_t corresponding to \p post.
+         * \param[in,out] flags Flags corresponding to \p post.
+         *     On input, the flags should contain the default flags.
+         *     On exit, the flags \ref POS_MASS, \ref POS_COMPLMAX and
+         *     \ref POS_COMPLWHOLE have been set according to \p post
+         *     (the completion flags are left at the default values if no
+         *     completion prefix is given).
+         * \throws  InternalError  if post is not recognized.
+         *
+         * \attention
+         * Checking is not complete, and other values than those listed above
+         * may be accepted for \p post, but the results are undefined.
+         *
+         * \see typeEnumValues
+         */
+        static void typeFromEnum(const char *post, e_poscalc_t *type, int *flags);
+
+        /*! \brief
+         * Creates a new position calculation collection object.
+         *
+         * \throws  std::bad_alloc if out of memory.
+         */
+        PositionCalculationCollection();
+        /*! \brief
+         * Destroys a position calculation collection and its calculations.
+         *
+         * Any calculations in the collection are also freed, even if
+         * references to them are left.
+         */
+        ~PositionCalculationCollection();
+
+        /*! \brief
+         * Sets the topology used for the calculations.
+         *
+         * \param[in]     top   Topology data structure.
+         *
+         * This function should be called to set the topology before using
+         * gmx_ana_poscalc_set_maxindex() for any calculation that requires
+         * topology information.
+         *
+         * Does not throw.
+         */
+        void setTopology(t_topology *top);
+        /*! \brief
+         * Prints information about calculations.
+         *
+         * \param[in] fp    File handle to receive the output.
+         *
+         * The output is very technical, making this function mainly useful for
+         * debugging purposes.
+         *
+         * Does not throw.
+         */
+        void printTree(FILE *fp) const;
+
+        /*! \brief
+         * Creates a new position calculation.
+         *
+         * \param[in]  type  Type of calculation.
+         * \param[in]  flags Flags for setting calculation options
+         *   (see \ref poscalc_flags "documentation of the flags").
+         *
+         * Does not throw currently, but may throw std::bad_alloc in the
+         * future.
+         */
+        gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
+        /*! \brief
+         * Creates a new position calculation based on an enum value.
+         *
+         * \param[in]  post  One of the strings acceptable for
+         *      typeFromEnum().
+         * \param[in]  flags Flags for setting calculation options
+         *      (see \ref poscalc_flags "documentation of the flags").
+         * \throws     InternalError  if post is not recognized.
+         *
+         * This is a convenience wrapper for createCalculation().
+         * \p flags sets the default calculation options if not overridden by
+         * \p post; see typeFromEnum().
+         *
+         * May also throw std::bad_alloc in the future.
+         *
+         * \see createCalculation(), typeFromEnum()
+         */
+        gmx_ana_poscalc_t *createCalculationFromEnum(const char *post, int flags);
+
+        /*! \brief
+         * Initializes evaluation for a position calculation collection.
+         *
+         * This function does some final initialization of the data structures
+         * in the collection to prepare them for evaluation.
+         * After this function has been called, it is no longer possible to add
+         * new calculations to the collection.
+         *
+         * Multiple calls to the function are ignored.
+         *
+         * Does not throw currently, but may throw std::bad_alloc in the
+         * future.
+         */
+        void initEvaluation();
+        /*! \brief
+         * Initializes a position calculation collection for a new frame.
+         *
+         * Clears the evaluation flag for all calculations.
+         * Should be called for each frame before calling
+         * gmx_ana_poscalc_update().
+         *
+         * This function calls initEvaluation() automatically if it has not
+         * been called earlier.
+         *
+         * Does not throw, but may throw if initEvaluation() is changed to
+         * throw.
+         */
+        void initFrame();
+
+    private:
+        class Impl;
+
+        PrivateImplPointer<Impl> impl_;
+
+        /*! \brief
+         * Needed to access the implementation class from the C code.
+         */
+        friend struct ::gmx_ana_poscalc_t;
+};
+
+} // namespace gmx
 
-/** Creates a new position calculation. */
-void
-gmx_ana_poscalc_create(gmx_ana_poscalc_t **pcp, gmx_ana_poscalc_coll_t *pcc,
-                       e_poscalc_t type, int flags);
-/** Creates a new position calculation based on an enum value. */
-void
-gmx_ana_poscalc_create_enum(gmx_ana_poscalc_t **pcp, gmx_ana_poscalc_coll_t *pcc,
-                            const char *post, int flags);
 /** Sets the flags for position calculation. */
 void
 gmx_ana_poscalc_set_flags(gmx_ana_poscalc_t *pc, int flags);
@@ -170,12 +329,6 @@ gmx_ana_poscalc_free(gmx_ana_poscalc_t *pc);
 bool
 gmx_ana_poscalc_requires_top(gmx_ana_poscalc_t *pc);
 
-/** Initializes evaluation for a position calculation collection. */
-void
-gmx_ana_poscalc_init_eval(gmx_ana_poscalc_coll_t *pcc);
-/** Initializes a position calculation collection for a new frame. */
-void
-gmx_ana_poscalc_init_frame(gmx_ana_poscalc_coll_t *pcc);
 /** Updates a single COM/COG structure for a frame. */
 void
 gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc,
index d70c4dc978e9f58d4470bb2f337953489108d883..2f21c515aa5c530ae591e6af0391a93c533a9fd9 100644 (file)
@@ -47,9 +47,9 @@
 #include "../legacyheaders/typedefs.h"
 
 #include "../options/options.h"
-#include "../utility/flags.h"
 #include "../utility/uniqueptr.h"
 #include "indexutil.h"
+#include "poscalc.h"
 #include "selection.h" // For gmx::SelectionList
 #include "selectioncollection.h"
 
@@ -88,7 +88,7 @@ struct gmx_ana_selcollection_t
     /** Index group that contains all the atoms. */
     struct gmx_ana_index_t         gall;
     /** Position calculation collection used for selection position evaluation. */
-    struct gmx_ana_poscalc_coll_t *pcc;
+    gmx::PositionCalculationCollection  pcc;
     /** Memory pool used for selection evaluation. */
     struct gmx_sel_mempool_t      *mempool;
     /** Parser symbol table. */
@@ -147,21 +147,14 @@ class SelectionCollection::Impl
                 RequestList    *requests_;
         };
 
-        //! Possible flags for the selection collection.
-        enum Flag
-        {
-            efOwnPositionCollection = 1<<0,
-            efExternalGroupsSet     = 1<<1
-        };
-        //! Holds a collection of Flag values.
-        typedef FlagsTemplate<Flag> Flags;
-
-        //! Creates a new selection collection.
-        explicit Impl(gmx_ana_poscalc_coll_t *pcc);
+        /*! \brief
+         * Creates a new selection collection.
+         *
+         * \throws  std::bad_alloc if out of memory.
+         */
+        Impl();
         ~Impl();
 
-        //! Returns true if the given flag has been set.
-        bool hasFlag(Flag flag) const { return _flags.test(flag); }
         //! Clears the symbol table of the selection collection.
         void clearSymbolTable();
         /*! \brief
@@ -224,8 +217,8 @@ class SelectionCollection::Impl
          *  - 4: combine 2 and 3
          */
         int                     _debugLevel;
-        //! Flags for various properties of the collection.
-        Flags                   _flags;
+        //! Whether external groups have been set for the collection.
+        bool                    _bExternalGroupsSet;
         //! External index groups (can be NULL).
         gmx_ana_indexgrps_t    *_grps;
         //! List of selections requested for later parsing.
index ed720559bff14ce495aa6e7f69cd2c4e2ae0fc1d..22a30826cbdf1b2adb69b5616617871f62d59e7f 100644 (file)
@@ -46,9 +46,6 @@
 #include <string2.h>
 #include <xvgr.h>
 
-#include "poscalc.h"
-#include "selmethod.h"
-
 #include "gromacs/fatalerror/exceptions.h"
 #include "gromacs/fatalerror/gmxassert.h"
 #include "gromacs/fatalerror/messagestringcollector.h"
@@ -59,6 +56,7 @@
 
 #include "compiler.h"
 #include "mempool.h"
+#include "poscalc.h"
 #include "scanner.h"
 #include "selectioncollection-impl.h"
 #include "selectionoptionstorage.h"
@@ -85,30 +83,18 @@ int SelectionCollection::Impl::SelectionRequest::count() const
     return storage->maxValueCount();
 }
 
-SelectionCollection::Impl::Impl(gmx_ana_poscalc_coll_t *pcc)
+SelectionCollection::Impl::Impl()
     : _options("selection", "Common selection control"),
-      _debugLevel(0), _grps(NULL)
+      _debugLevel(0), _bExternalGroupsSet(false), _grps(NULL)
 {
     _sc.root      = NULL;
     _sc.nvars     = 0;
     _sc.varstrs   = NULL;
     _sc.top       = NULL;
     gmx_ana_index_clear(&_sc.gall);
-    _sc.pcc       = pcc;
     _sc.mempool   = NULL;
     _sc.symtab    = NULL;
 
-    // TODO: This is not exception-safe if any called function throws.
-    if (_sc.pcc == NULL)
-    {
-        int rc = gmx_ana_poscalc_coll_create(&_sc.pcc);
-        if (rc != 0)
-        {
-            // TODO: A more reasonable error
-            GMX_THROW(InternalError("Failed to create position calculation collection"));
-        }
-        _flags.set(Impl::efOwnPositionCollection);
-    }
     _gmx_sel_symtab_create(&_sc.symtab);
     gmx_ana_selmethod_register_defaults(_sc.symtab);
 }
@@ -128,10 +114,6 @@ SelectionCollection::Impl::~Impl()
     {
         _gmx_sel_mempool_destroy(_sc.mempool);
     }
-    if (hasFlag(efOwnPositionCollection))
-    {
-        gmx_ana_poscalc_coll_free(_sc.pcc);
-    }
     clearSymbolTable();
 }
 
@@ -250,8 +232,8 @@ void SelectionCollection::Impl::resolveExternalGroups(
  * SelectionCollection
  */
 
-SelectionCollection::SelectionCollection(gmx_ana_poscalc_coll_t *pcc)
-    : _impl(new Impl(pcc))
+SelectionCollection::SelectionCollection()
+    : _impl(new Impl)
 {
 }
 
@@ -277,17 +259,12 @@ SelectionCollection::initOptions()
     */
 
     Options &options = _impl->_options;
-    const char **postypes = gmx_ana_poscalc_create_type_enum(true);
-    if (postypes == NULL)
-    {
-        // TODO: Use an out-of-memory exception here
-        GMX_THROW(InternalError("Could not create position calculation enum"));
-    }
-    options.addOption(StringOption("selrpos").enumValue(postypes + 1)
-                          .store(&_impl->_rpost).defaultValue(postypes[1])
+    const char *const *postypes = PositionCalculationCollection::typeEnumValues;
+    options.addOption(StringOption("selrpos").enumValue(postypes)
+                          .store(&_impl->_rpost).defaultValue(postypes[0])
                           .description("Selection reference positions"));
-    options.addOption(StringOption("seltype").enumValue(postypes + 1)
-                          .store(&_impl->_spost).defaultValue(postypes[1])
+    options.addOption(StringOption("seltype").enumValue(postypes)
+                          .store(&_impl->_spost).defaultValue(postypes[0])
                           .description("Default selection output positions"));
     GMX_RELEASE_ASSERT(_impl->_debugLevel >= 0 && _impl->_debugLevel <= 4,
                        "Debug level out of range");
@@ -296,7 +273,6 @@ SelectionCollection::initOptions()
                           .defaultValue(debug_levels[_impl->_debugLevel])
                           .storeEnumIndex(&_impl->_debugLevel)
                           .description("Print out selection trees for debugging"));
-    sfree(postypes);
 
     return _impl->_options;
 }
@@ -306,10 +282,10 @@ void
 SelectionCollection::setReferencePosType(const char *type)
 {
     GMX_RELEASE_ASSERT(type != NULL, "Cannot assign NULL position type");
-    //! Check that the type is valid.
+    // Check that the type is valid, throw if it is not.
     e_poscalc_t  dummytype;
     int          dummyflags;
-    gmx_ana_poscalc_type_from_enum(type, &dummytype, &dummyflags);
+    PositionCalculationCollection::typeFromEnum(type, &dummytype, &dummyflags);
     _impl->_rpost = type;
 }
 
@@ -318,10 +294,10 @@ void
 SelectionCollection::setOutputPosType(const char *type)
 {
     GMX_RELEASE_ASSERT(type != NULL, "Cannot assign NULL position type");
-    //! Check that the type is valid.
+    // Check that the type is valid, throw if it is not.
     e_poscalc_t  dummytype;
     int          dummyflags;
-    gmx_ana_poscalc_type_from_enum(type, &dummytype, &dummyflags);
+    PositionCalculationCollection::typeFromEnum(type, &dummytype, &dummyflags);
     _impl->_spost = type;
 }
 
@@ -337,7 +313,7 @@ void
 SelectionCollection::setTopology(t_topology *top, int natoms)
 {
     gmx_ana_selcollection_t *sc = &_impl->_sc;
-    gmx_ana_poscalc_coll_set_topology(sc->pcc, top);
+    sc->pcc.setTopology(top);
     sc->top = top;
 
     /* Get the number of atoms from the topology if it is not given */
@@ -357,10 +333,10 @@ SelectionCollection::setTopology(t_topology *top, int natoms)
 void
 SelectionCollection::setIndexGroups(gmx_ana_indexgrps_t *grps)
 {
-    GMX_RELEASE_ASSERT(grps == NULL || !_impl->hasFlag(Impl::efExternalGroupsSet),
+    GMX_RELEASE_ASSERT(grps == NULL || !_impl->_bExternalGroupsSet,
                        "Can only set external groups once or clear them afterwards");
     _impl->_grps = grps;
-    _impl->_flags.set(Impl::efExternalGroupsSet);
+    _impl->_bExternalGroupsSet = true;
 
     MessageStringCollector errors;
     t_selelem *root = _impl->_sc.root;
@@ -387,7 +363,8 @@ SelectionCollection::requiresTopology() const
     {
         flags = 0;
         // Should not throw, because has been checked earlier.
-        gmx_ana_poscalc_type_from_enum(_impl->_rpost.c_str(), &type, &flags);
+        PositionCalculationCollection::typeFromEnum(_impl->_rpost.c_str(),
+                                                    &type, &flags);
         if (type != POS_ATOM)
         {
             return true;
@@ -397,7 +374,8 @@ SelectionCollection::requiresTopology() const
     {
         flags = 0;
         // Should not throw, because has been checked earlier.
-        gmx_ana_poscalc_type_from_enum(_impl->_spost.c_str(), &type, &flags);
+        PositionCalculationCollection::typeFromEnum(_impl->_spost.c_str(),
+                                                    &type, &flags);
         if (type != POS_ATOM)
         {
             return true;
@@ -501,7 +479,7 @@ SelectionCollection::parseFromStdin(int nr, bool bInteractive,
     yyscan_t scanner;
 
     _gmx_sel_init_lexer(&scanner, &_impl->_sc, bInteractive, nr,
-                        _impl->hasFlag(Impl::efExternalGroupsSet),
+                        _impl->_bExternalGroupsSet,
                         _impl->_grps);
     /* We don't set the lexer input here, which causes it to use a special
      * internal implementation for reading from stdin. */
@@ -517,7 +495,7 @@ SelectionCollection::parseFromFile(const std::string &filename,
     FILE *fp;
 
     _gmx_sel_init_lexer(&scanner, &_impl->_sc, false, -1,
-                        _impl->hasFlag(Impl::efExternalGroupsSet),
+                        _impl->_bExternalGroupsSet,
                         _impl->_grps);
     fp = ffopen(filename.c_str(), "r");
     _gmx_sel_set_lex_input_file(scanner, fp);
@@ -542,7 +520,7 @@ SelectionCollection::parseFromString(const std::string &str,
     yyscan_t scanner;
 
     _gmx_sel_init_lexer(&scanner, &_impl->_sc, false, -1,
-                        _impl->hasFlag(Impl::efExternalGroupsSet),
+                        _impl->_bExternalGroupsSet,
                         _impl->_grps);
     _gmx_sel_set_lex_input_str(scanner, str.c_str());
     _impl->runParser(scanner, -1, output);
@@ -556,7 +534,7 @@ SelectionCollection::compile()
     {
         GMX_THROW(InconsistentInputError("Selection requires topology information, but none provided"));
     }
-    if (!_impl->hasFlag(Impl::efExternalGroupsSet))
+    if (!_impl->_bExternalGroupsSet)
     {
         setIndexGroups(NULL);
     }
@@ -568,22 +546,19 @@ SelectionCollection::compile()
     SelectionCompiler compiler;
     compiler.compile(this);
 
-    if (_impl->hasFlag(Impl::efOwnPositionCollection))
+    if (_impl->_debugLevel >= 1)
     {
-        if (_impl->_debugLevel >= 1)
-        {
-            std::fprintf(stderr, "\n");
-            printTree(stderr, false);
-            std::fprintf(stderr, "\n");
-            gmx_ana_poscalc_coll_print_tree(stderr, _impl->_sc.pcc);
-            std::fprintf(stderr, "\n");
-        }
-        gmx_ana_poscalc_init_eval(_impl->_sc.pcc);
-        if (_impl->_debugLevel >= 1)
-        {
-            gmx_ana_poscalc_coll_print_tree(stderr, _impl->_sc.pcc);
-            std::fprintf(stderr, "\n");
-        }
+        std::fprintf(stderr, "\n");
+        printTree(stderr, false);
+        std::fprintf(stderr, "\n");
+        _impl->_sc.pcc.printTree(stderr);
+        std::fprintf(stderr, "\n");
+    }
+    _impl->_sc.pcc.initEvaluation();
+    if (_impl->_debugLevel >= 1)
+    {
+        _impl->_sc.pcc.printTree(stderr);
+        std::fprintf(stderr, "\n");
     }
 }
 
@@ -591,10 +566,7 @@ SelectionCollection::compile()
 void
 SelectionCollection::evaluate(t_trxframe *fr, t_pbc *pbc)
 {
-    if (_impl->hasFlag(Impl::efOwnPositionCollection))
-    {
-        gmx_ana_poscalc_init_frame(_impl->_sc.pcc);
-    }
+    _impl->_sc.pcc.initFrame();
 
     SelectionEvaluator evaluator;
     evaluator.evaluate(this, fr, pbc);
index e873af9237e8aebd5758532d1d50ee03f9808ddb..decc46e1873986033171b4bffdbdabf6d28207b4 100644 (file)
@@ -48,7 +48,6 @@
 #include "selection.h" // For gmx::SelectionList
 
 struct gmx_ana_indexgrps_t;
-struct gmx_ana_poscalc_coll_t;
 
 namespace gmx
 {
@@ -91,13 +90,9 @@ class SelectionCollection
         /*! \brief
          * Creates an empty selection collection.
          *
-         * \param[in] pcc  Position calculation collection to use for selection
-         *      evaluation.
-         *
-         * If \p pcc is NULL, an internal collection is created and managed by
-         * the object.
+         * \throws  std::bad_alloc if out of memory.
          */
-        explicit SelectionCollection(gmx_ana_poscalc_coll_t *pcc);
+        SelectionCollection();
         ~SelectionCollection();
 
         /*! \brief
@@ -110,7 +105,8 @@ class SelectionCollection
          * collection.
          *
          * \param[in]     type      Default selection reference position type
-         *   (one of the strings acceptable for gmx_ana_poscalc_type_from_enum()).
+         *     (one of the strings acceptable for
+         *     PositionCalculationCollection::typeFromEnum()).
          *
          * Should be called before calling the parser functions, unless
          * initOptions() has been called.  In the latter case, can still be
@@ -123,7 +119,8 @@ class SelectionCollection
          * collection.
          *
          * \param[in]     type      Default selection output position type
-         *   (one of the strings acceptable for gmx_ana_poscalc_type_from_enum()).
+         *     (one of the strings acceptable for
+         *     PositionCalculationCollection::typeFromEnum()).
          *
          * Should be called before calling the parser functions, unless
          * initOptions() has been called.  In the latter case, can still be
index f5290a0d4beea3c7656f760db07d140208361339..fd9b5d71e1ca448e208ee0d828e942a164f49639 100644 (file)
 #include "selparam.h"
 #include "selvalue.h"
 
+namespace gmx
+{
+class PositionCalculationCollection;
+} // namespace gmx
+
 struct gmx_ana_pos_t;
-struct gmx_ana_poscalc_coll_t;
 struct gmx_ana_selcollection_t;
 
 /*! \name Selection method flags
@@ -389,7 +393,7 @@ typedef void *(*sel_datafunc)(int npar, gmx_ana_selparam_t *param);
  * The pointer \p pcc should then be stored and used for initialization for
  * any position calculation structures.
  */
-typedef void  (*sel_posfunc)(struct gmx_ana_poscalc_coll_t *pcc, void *data);
+typedef void  (*sel_posfunc)(gmx::PositionCalculationCollection *pcc, void *data);
 /*! \brief
  * Does initialization based on topology and/or parameter values.
  *
index 18f262e724dc8555ba181ca5de1b676cfa1dca26..99af3684188f0fc876fedd4783d4f41354c6f087 100644 (file)
@@ -57,7 +57,7 @@
 typedef struct
 {
     /** Position calculation collection to use. */
-    gmx_ana_poscalc_coll_t *pcc;
+    gmx::PositionCalculationCollection *pcc;
     /** Index group for which the center should be evaluated. */
     gmx_ana_index_t    g;
     /** Position evaluation data structure. */
@@ -75,7 +75,7 @@ static void *
 init_data_pos(int npar, gmx_ana_selparam_t *param);
 /** Sets the position calculation collection for position evaluation selection methods. */
 static void
-set_poscoll_pos(gmx_ana_poscalc_coll_t *pcc, void *data);
+set_poscoll_pos(gmx::PositionCalculationCollection *pcc, void *data);
 /** Initializes position evaluation keywords. */
 static void
 init_kwpos(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
@@ -186,7 +186,7 @@ init_data_pos(int npar, gmx_ana_selparam_t *param)
  * \param[in,out] data  Should point to \c t_methoddata_pos.
  */
 static void
-set_poscoll_pos(gmx_ana_poscalc_coll_t *pcc, void *data)
+set_poscoll_pos(gmx::PositionCalculationCollection *pcc, void *data)
 {
     ((t_methoddata_pos *)data)->pcc = pcc;
 }
@@ -194,7 +194,7 @@ set_poscoll_pos(gmx_ana_poscalc_coll_t *pcc, void *data)
 /*!
  * \param[in,out] sel   Selection element to initialize.
  * \param[in]     type  One of the enum values acceptable for
- *   gmx_ana_poscalc_type_from_enum().
+ *     PositionCalculationCollection::typeFromEnum().
  *
  * Initializes the reference position type for position evaluation.
  * If called multiple times, the first setting takes effect, and later calls
@@ -224,7 +224,7 @@ _gmx_selelem_set_kwpos_type(t_selelem *sel, const char *type)
 /*!
  * \param[in,out] sel   Selection element to initialize.
  * \param[in]     flags Default completion flags
- *   (see gmx_ana_poscalc_type_from_enum()).
+ *     (see PositionCalculationCollection::typeFromEnum()).
  *
  * Initializes the flags for position evaluation.
  * If called multiple times, the first setting takes effect, and later calls
@@ -269,7 +269,7 @@ init_kwpos(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
     {
         d->flags |= POS_DYNAMIC;
     }
-    gmx_ana_poscalc_create_enum(&d->pc, d->pcc, d->type, d->flags);
+    d->pc = d->pcc->createCalculationFromEnum(d->type, d->flags);
     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
 }
 
@@ -286,8 +286,7 @@ init_cog(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
     t_methoddata_pos *d = (t_methoddata_pos *)data;
 
     d->flags = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
-    gmx_ana_poscalc_create(&d->pc, d->pcc, d->bPBC ? POS_ALL_PBC : POS_ALL,
-                           d->flags);
+    d->pc = d->pcc->createCalculation(d->bPBC ? POS_ALL_PBC : POS_ALL, d->flags);
     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
 }
 
@@ -305,8 +304,7 @@ init_com(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
 
     d->flags  = (param[0].flags & SPAR_DYNAMIC) ? POS_DYNAMIC : 0;
     d->flags |= POS_MASS;
-    gmx_ana_poscalc_create(&d->pc, d->pcc, d->bPBC ? POS_ALL_PBC : POS_ALL,
-                           d->flags);
+    d->pc = d->pcc->createCalculation(d->bPBC ? POS_ALL_PBC : POS_ALL, d->flags);
     gmx_ana_poscalc_set_maxindex(d->pc, &d->g);
 }
 
index c0bd9edd33cb0a761a9c14d00cd6c90a7e4dd0fd..4fba00d570f79b3300263d0eb0c62cde40df6d0c 100644 (file)
@@ -189,18 +189,18 @@ add_reserved_symbols(gmx_sel_symtab_t *tab)
 static void
 add_position_symbols(gmx_sel_symtab_t *tab)
 {
-    const char       **postypes;
     gmx_sel_symrec_t  *sym;
     gmx_sel_symrec_t  *last;
     int                i;
 
-    postypes = gmx_ana_poscalc_create_type_enum(true);
+    const char *const *postypes
+        = gmx::PositionCalculationCollection::typeEnumValues;
     last = tab->first;
     while (last && last->next)
     {
         last = last->next;
     }
-    for (i = 1; postypes[i] != NULL; ++i)
+    for (i = 0; postypes[i] != NULL; ++i)
     {
         snew(sym, 1);
         sym->name = strdup(postypes[i]);
@@ -216,7 +216,6 @@ add_position_symbols(gmx_sel_symtab_t *tab)
         }
         last = sym;
     }
-    sfree(postypes);
 }
 
 /*!
index 54225a928af6c62bb75a5946bc026b58024a7e1c..c9936196d3766bfca4555d1aa946c0b265f0f6a7 100644 (file)
@@ -83,7 +83,7 @@ class SelectionCollectionTest : public ::testing::Test
 
 
 SelectionCollectionTest::SelectionCollectionTest()
-    : _sc(NULL), _top(NULL), _frame(NULL)
+    : _top(NULL), _frame(NULL)
 {
     _sc.setReferencePosType("atom");
     _sc.setOutputPosType("atom");
index 86cd6d5224af044467de12efb45301c6fb4ed9f7..705080b58151eae5934f0a28b67240a2b930f0d4 100644 (file)
@@ -60,7 +60,7 @@ class SelectionOptionTest : public ::testing::Test
 };
 
 SelectionOptionTest::SelectionOptionTest()
-    : _sc(NULL), _options(NULL, NULL)
+    : _options(NULL, NULL)
 {
     _sc.setReferencePosType("atom");
     _sc.setOutputPosType("atom");
index d4fb89e40dcce10da04a95bb9e7d1c4cadf2dfd9..a39f6a4182958387b044f7beaee9f1ef950c2b21 100644 (file)
@@ -192,7 +192,7 @@ TrajectoryAnalysisCommandLineRunner::run(int argc, char *argv[])
 
     CopyRight(stderr, argv[0]);
 
-    SelectionCollection  selections(NULL);
+    SelectionCollection  selections;
     selections.setDebugLevel(_impl->_debugLevel);
 
     TrajectoryAnalysisSettings  settings;