Merge release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / selection / evaluate.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements functions in evaluate.h.
38  *
39  * \todo
40  * One of the major bottlenecks for selection performance is that all the
41  * evaluation is carried out for atoms.
42  * There are several cases when the evaluation could be done for residues
43  * or molecules instead, including keywords that select by residue and
44  * cases where residue centers are used as reference positions.
45  * Implementing this would require a mechanism for recognizing whether
46  * something can be evaluated by residue/molecule instead by atom, and
47  * converting selections by residue/molecule into selections by atom
48  * when necessary.
49  *
50  * \author Teemu Murtola <teemu.murtola@gmail.com>
51  * \ingroup module_selection
52  */
53 #include <string.h>
54
55 #include "gromacs/legacyheaders/vec.h"
56
57 #include "gromacs/math/utilities.h"
58 #include "gromacs/selection/indexutil.h"
59 #include "gromacs/selection/poscalc.h"
60 #include "gromacs/selection/selection.h"
61 #include "gromacs/selection/selmethod.h"
62 #include "gromacs/utility/exceptions.h"
63 #include "gromacs/utility/gmxassert.h"
64 #include "gromacs/utility/smalloc.h"
65
66 #include "evaluate.h"
67 #include "mempool.h"
68 #include "selectioncollection-impl.h"
69 #include "selelem.h"
70
71 using gmx::SelectionTreeElement;
72 using gmx::SelectionTreeElementPointer;
73
74 namespace
75 {
76
77 /*! \brief
78  * Reserves memory for a selection element from the evaluation memory pool.
79  *
80  * This class implements RAII semantics for allocating memory for selection
81  * element values from a selection evaluation memory pool.
82  *
83  * \ingroup module_selection
84  */
85 class MempoolSelelemReserver
86 {
87     public:
88         //! Constructs a reserver without initial reservation.
89         MempoolSelelemReserver() {}
90         /*! \brief
91          * Constructs a reserver with initial reservation.
92          *
93          * \param[in,out] sel    Selection element for which to reserve.
94          * \param[in]     count  Number of values to reserve.
95          *
96          * \see reserve()
97          */
98         MempoolSelelemReserver(const SelectionTreeElementPointer &sel, int count)
99         {
100             reserve(sel, count);
101         }
102         //! Frees any memory allocated using this reserver.
103         ~MempoolSelelemReserver()
104         {
105             if (sel_)
106             {
107                 sel_->mempoolRelease();
108             }
109         }
110
111         /*! \brief
112          * Reserves memory for selection element values using this reserver.
113          *
114          * \param[in,out] sel    Selection element for which to reserve.
115          * \param[in]     count  Number of values to reserve.
116          *
117          * Allocates space to store \p count output values in \p sel from the
118          * memory pool associated with \p sel, or from the heap if there is no
119          * memory pool.  Type of values to allocate is automatically determined
120          * from \p sel.
121          */
122         void reserve(const SelectionTreeElementPointer &sel, int count)
123         {
124             GMX_RELEASE_ASSERT(!sel_,
125                                "Can only reserve one element with one instance");
126             sel->mempoolReserve(count);
127             sel_ = sel;
128         }
129
130     private:
131         SelectionTreeElementPointer     sel_;
132 };
133
134 /*! \brief
135  * Reserves memory for an index group from the evaluation memory pool.
136  *
137  * This class implements RAII semantics for allocating memory for an index
138  * group from a selection evaluation memory pool.
139  *
140  * \ingroup module_selection
141  */
142 class MempoolGroupReserver
143 {
144     public:
145         /*! \brief
146          * Creates a reserver associated with a given memory pool.
147          *
148          * \param    mp  Memory pool from which to reserve memory.
149          */
150         explicit MempoolGroupReserver(gmx_sel_mempool_t *mp)
151             : mp_(mp), g_(NULL)
152         {
153         }
154         //! Frees any memory allocated using this reserver.
155         ~MempoolGroupReserver()
156         {
157             if (g_ != NULL)
158             {
159                 _gmx_sel_mempool_free_group(mp_, g_);
160             }
161         }
162
163         /*! \brief
164          * Reserves memory for an index group using this reserver.
165          *
166          * \param[in,out] g      Index group to reserve.
167          * \param[in]     count  Number of atoms to reserve space for.
168          *
169          * Allocates memory from the memory pool to store \p count atoms in
170          * \p g.
171          */
172         void reserve(gmx_ana_index_t *g, int count)
173         {
174             GMX_RELEASE_ASSERT(g_ == NULL, "Can only reserve one element with one instance");
175             _gmx_sel_mempool_alloc_group(mp_, g, count);
176             g_ = g;
177         }
178
179     private:
180         gmx_sel_mempool_t      *mp_;
181         gmx_ana_index_t        *g_;
182 };
183
184 /*! \brief
185  * Assigns a temporary value for a selection element.
186  *
187  * This class implements RAII semantics for temporarily assigning the value
188  * pointer of a selection element to point to a different location.
189  *
190  * \ingroup module_selection
191  */
192 class SelelemTemporaryValueAssigner
193 {
194     public:
195         //! Constructs an assigner without an initial assignment.
196         SelelemTemporaryValueAssigner()
197             : old_ptr_(NULL), old_nalloc_(0)
198         {
199         }
200         /*! \brief
201          * Constructs an assigner with an initial assignment.
202          *
203          * \param[in,out] sel     Selection element for which to assign.
204          * \param[in]     vsource Element to which \p sel values will point to.
205          *
206          * \see assign()
207          */
208         SelelemTemporaryValueAssigner(const SelectionTreeElementPointer &sel,
209                                       const SelectionTreeElement        &vsource)
210         {
211             assign(sel, vsource);
212         }
213         //! Undoes any temporary assignment done using this assigner.
214         ~SelelemTemporaryValueAssigner()
215         {
216             if (sel_)
217             {
218                 _gmx_selvalue_setstore_alloc(&sel_->v, old_ptr_, old_nalloc_);
219             }
220         }
221
222         /*! \brief
223          * Assigns a temporary value pointer.
224          *
225          * \param[in,out] sel     Selection element for which to assign.
226          * \param[in]     vsource Element to which \p sel values will point to.
227          *
228          * Assigns the value pointer in \p sel to point to the values in
229          * \p vsource, i.e., any access/modification to values in \p sel
230          * actually accesses values in \p vsource.
231          */
232         void assign(const SelectionTreeElementPointer &sel,
233                     const SelectionTreeElement        &vsource)
234         {
235             GMX_RELEASE_ASSERT(!sel_,
236                                "Can only assign one element with one instance");
237             GMX_RELEASE_ASSERT(sel->v.type == vsource.v.type,
238                                "Mismatching selection value types");
239             old_ptr_    = sel->v.u.ptr;
240             old_nalloc_ = sel->v.nalloc;
241             _gmx_selvalue_setstore(&sel->v, vsource.v.u.ptr);
242             sel_ = sel;
243         }
244
245     private:
246         SelectionTreeElementPointer     sel_;
247         void                           *old_ptr_;
248         int                             old_nalloc_;
249 };
250
251 } // namespace
252
253 /*!
254  * \param[in] fp       File handle to receive the output.
255  * \param[in] evalfunc Function pointer to print.
256  */
257 void
258 _gmx_sel_print_evalfunc_name(FILE *fp, gmx::sel_evalfunc evalfunc)
259 {
260     if (!evalfunc)
261     {
262         fprintf(fp, "none");
263     }
264     else if (evalfunc == &_gmx_sel_evaluate_root)
265     {
266         fprintf(fp, "root");
267     }
268     else if (evalfunc == &_gmx_sel_evaluate_static)
269     {
270         fprintf(fp, "static");
271     }
272     else if (evalfunc == &_gmx_sel_evaluate_subexpr_simple)
273     {
274         fprintf(fp, "subexpr_simple");
275     }
276     else if (evalfunc == &_gmx_sel_evaluate_subexpr_staticeval)
277     {
278         fprintf(fp, "subexpr_staticeval");
279     }
280     else if (evalfunc == &_gmx_sel_evaluate_subexpr)
281     {
282         fprintf(fp, "subexpr");
283     }
284     else if (evalfunc == &_gmx_sel_evaluate_subexprref_simple)
285     {
286         fprintf(fp, "ref_simple");
287     }
288     else if (evalfunc == &_gmx_sel_evaluate_subexprref)
289     {
290         fprintf(fp, "ref");
291     }
292     else if (evalfunc == &_gmx_sel_evaluate_method)
293     {
294         fprintf(fp, "method");
295     }
296     else if (evalfunc == &_gmx_sel_evaluate_modifier)
297     {
298         fprintf(fp, "mod");
299     }
300     else if (evalfunc == &_gmx_sel_evaluate_not)
301     {
302         fprintf(fp, "not");
303     }
304     else if (evalfunc == &_gmx_sel_evaluate_and)
305     {
306         fprintf(fp, "and");
307     }
308     else if (evalfunc == &_gmx_sel_evaluate_or)
309     {
310         fprintf(fp, "or");
311     }
312     else if (evalfunc == &_gmx_sel_evaluate_arithmetic)
313     {
314         fprintf(fp, "arithmetic");
315     }
316     else
317     {
318         fprintf(fp, "%p", (void*)(evalfunc));
319     }
320 }
321
322 /*!
323  * \param[out] data Evaluation data structure to initialize.
324  * \param[in]  mp   Memory pool for intermediate evaluation values.
325  * \param[in]  gall Index group with all the atoms.
326  * \param[in]  top  Topology structure for evaluation.
327  * \param[in]  fr   New frame for evaluation.
328  * \param[in]  pbc  New PBC information for evaluation.
329  */
330 void
331 _gmx_sel_evaluate_init(gmx_sel_evaluate_t *data,
332                        gmx_sel_mempool_t *mp, gmx_ana_index_t *gall,
333                        t_topology *top, t_trxframe *fr, t_pbc *pbc)
334 {
335     data->mp   = mp;
336     data->gall = gall;
337     data->top  = top;
338     data->fr   = fr;
339     data->pbc  = pbc;
340 }
341
342 /*! \brief
343  * Recursively initializes the flags for evaluation.
344  *
345  * \param[in,out] sel Selection element to clear.
346  *
347  * The \ref SEL_INITFRAME flag is set for \ref SEL_EXPRESSION elements whose
348  * method defines the \p init_frame callback (see sel_framefunc()), and
349  * cleared for other elements.
350  *
351  * The \ref SEL_EVALFRAME flag is cleared for all elements.
352  */
353 static void
354 init_frame_eval(SelectionTreeElementPointer sel)
355 {
356     while (sel)
357     {
358         sel->flags &= ~(SEL_INITFRAME | SEL_EVALFRAME);
359         if (sel->type == SEL_EXPRESSION)
360         {
361             if (sel->u.expr.method && sel->u.expr.method->init_frame)
362             {
363                 sel->flags |= SEL_INITFRAME;
364             }
365         }
366         if (sel->child && sel->type != SEL_SUBEXPRREF)
367         {
368             init_frame_eval(sel->child);
369         }
370         sel = sel->next;
371     }
372 }
373
374 namespace gmx
375 {
376
377 SelectionEvaluator::SelectionEvaluator()
378 {
379 }
380
381 /*!
382  * \param[in,out] coll  The selection collection to evaluate.
383  * \param[in] fr  Frame for which the evaluation should be carried out.
384  * \param[in] pbc PBC data, or NULL if no PBC should be used.
385  * \returns   0 on successful evaluation, a non-zero error code on error.
386  *
387  * This functions sets the global variables for topology, frame and PBC,
388  * clears some information in the selection to initialize the evaluation
389  * for a new frame, and evaluates \p sel and all the selections pointed by
390  * the \p next pointers of \p sel.
391  *
392  * This is the only function that user code should call if they want to
393  * evaluate a selection for a new frame.
394  */
395 void
396 SelectionEvaluator::evaluate(SelectionCollection *coll,
397                              t_trxframe *fr, t_pbc *pbc)
398 {
399     gmx_ana_selcollection_t *sc = &coll->impl_->sc_;
400     gmx_sel_evaluate_t       data;
401
402     _gmx_sel_evaluate_init(&data, sc->mempool, &sc->gall, sc->top, fr, pbc);
403     init_frame_eval(sc->root);
404     SelectionTreeElementPointer sel = sc->root;
405     while (sel)
406     {
407         /* Clear the evaluation group of subexpressions */
408         if (sel->child && sel->child->type == SEL_SUBEXPR
409             && sel->child->evaluate != NULL)
410         {
411             sel->child->u.cgrp.isize = 0;
412             /* Not strictly necessary, because the value will be overwritten
413              * during first evaluation of the subexpression anyways, but we
414              * clear the group for clarity. Note that this is _not_ done during
415              * compilation because of some additional complexities involved
416              * (see compiler.cpp), so it should not be relied upon in
417              * _gmx_sel_evaluate_subexpr(). */
418             if (sel->child->v.type == GROUP_VALUE)
419             {
420                 sel->child->v.u.g->isize = 0;
421             }
422         }
423         if (sel->evaluate)
424         {
425             sel->evaluate(&data, sel, NULL);
426         }
427         sel = sel->next;
428     }
429     /* Update selection information */
430     SelectionDataList::const_iterator isel;
431     for (isel = sc->sel.begin(); isel != sc->sel.end(); ++isel)
432     {
433         internal::SelectionData &sel = **isel;
434         sel.refreshMassesAndCharges(sc->top);
435         sel.updateCoveredFractionForFrame();
436     }
437 }
438
439 /*!
440  * \param[in,out] coll  The selection collection to evaluate.
441  * \param[in]     nframes Total number of frames.
442  */
443 void
444 SelectionEvaluator::evaluateFinal(SelectionCollection *coll, int nframes)
445 {
446     gmx_ana_selcollection_t          *sc = &coll->impl_->sc_;
447
448     SelectionDataList::const_iterator isel;
449     for (isel = sc->sel.begin(); isel != sc->sel.end(); ++isel)
450     {
451         internal::SelectionData &sel = **isel;
452         sel.restoreOriginalPositions(sc->top);
453         sel.computeAverageCoveredFraction(nframes);
454     }
455 }
456
457 } // namespace gmx
458
459 /*!
460  * \param[in] data Data for the current frame.
461  * \param[in] sel  Selection element being evaluated.
462  * \param[in] g    Group for which \p sel should be evaluated.
463  * \returns   0 on success, a non-zero error code on error.
464  *
465  * Evaluates each child of \p sel in \p g.
466  */
467 void
468 _gmx_sel_evaluate_children(gmx_sel_evaluate_t                     *data,
469                            const gmx::SelectionTreeElementPointer &sel,
470                            gmx_ana_index_t                        *g)
471 {
472     SelectionTreeElementPointer child = sel->child;
473     while (child)
474     {
475         if (child->evaluate)
476         {
477             child->evaluate(data, child, g);
478         }
479         child = child->next;
480     }
481 }
482
483 void
484 _gmx_sel_evaluate_root(gmx_sel_evaluate_t                     *data,
485                        const gmx::SelectionTreeElementPointer &sel,
486                        gmx_ana_index_t                         * /* g */)
487 {
488     if (sel->u.cgrp.isize == 0 || !sel->child->evaluate)
489     {
490         return;
491     }
492
493     sel->child->evaluate(data, sel->child,
494                          sel->u.cgrp.isize < 0 ? NULL : &sel->u.cgrp);
495 }
496
497 void
498 _gmx_sel_evaluate_static(gmx_sel_evaluate_t                      * /* data */,
499                          const gmx::SelectionTreeElementPointer &sel,
500                          gmx_ana_index_t                        *g)
501 {
502     gmx_ana_index_intersection(sel->v.u.g, &sel->u.cgrp, g);
503 }
504
505
506 /*********************************************************************
507  * SUBEXPRESSION EVALUATION
508  *********************************************************************/
509
510 /*!
511  * \param[in] data Data for the current frame.
512  * \param[in] sel  Selection element being evaluated.
513  * \param[in] g    Group for which \p sel should be evaluated.
514  * \returns   0 on success, a non-zero error code on error.
515  *
516  * Evaluates the child element (there should be exactly one) in \p g.
517  * The compiler has taken care that the child actually stores the evaluated
518  * value in the value pointer of this element.
519  *
520  * This function is used as gmx::SelectionTreeElement::evaluate for
521  * \ref SEL_SUBEXPR elements that are used only once, and hence do not need
522  * full subexpression handling.
523  */
524 void
525 _gmx_sel_evaluate_subexpr_simple(gmx_sel_evaluate_t                     *data,
526                                  const gmx::SelectionTreeElementPointer &sel,
527                                  gmx_ana_index_t                        *g)
528 {
529     if (sel->child->evaluate)
530     {
531         sel->child->evaluate(data, sel->child, g);
532     }
533     sel->v.nr = sel->child->v.nr;
534 }
535
536 /*!
537  * \param[in] data Data for the current frame.
538  * \param[in] sel  Selection element being evaluated.
539  * \param[in] g    Group for which \p sel should be evaluated.
540  * \returns   0 on success, a non-zero error code on error.
541  *
542  * If this is the first call for this frame, evaluates the child element
543  * there should be exactly one in \p g.
544  * The compiler has taken care that the child actually stores the evaluated
545  * value in the value pointer of this element.
546  * Assumes that \p g is persistent for the duration of the whole evaluation.
547  *
548  * This function is used as gmx::SelectionTreeElement::evaluate for
549  * \ref SEL_SUBEXPR elements that have a static evaluation group, and hence do
550  * not need full subexpression handling.
551  */
552 void
553 _gmx_sel_evaluate_subexpr_staticeval(gmx_sel_evaluate_t                     *data,
554                                      const gmx::SelectionTreeElementPointer &sel,
555                                      gmx_ana_index_t                        *g)
556 {
557     if (sel->u.cgrp.isize == 0)
558     {
559         sel->child->evaluate(data, sel->child, g);
560         sel->v.nr = sel->child->v.nr;
561         if (!g)
562         {
563             sel->u.cgrp.isize = -1;
564         }
565         else
566         {
567             gmx_ana_index_set(&sel->u.cgrp, g->isize, g->index, 0);
568         }
569     }
570 }
571
572 /*!
573  * \param[in]  data  Data for the current frame.
574  * \param[in]  sel   Selection element being evaluated.
575  * \param[in]  g     Group for which \p sel should be evaluated.
576  * \returns    0 on success, a non-zero error code on error.
577  *
578  * Finds the part of \p g for which the subexpression
579  * has not yet been evaluated by comparing \p g to \p sel->u.cgrp.
580  * If the part is not empty, the child expression is evaluated for this
581  * part, and the results merged to the old values of the child.
582  * The value of \p sel itself is undefined after the call.
583  *
584  * \todo
585  * The call to gmx_ana_index_difference() can take quite a lot of unnecessary
586  * time if the subexpression is evaluated either several times for the same
587  * group or for completely distinct groups.
588  * However, in the majority of cases, these situations occur when
589  * _gmx_sel_evaluate_subexpr_staticeval() can be used, so this should not be a
590  * major problem.
591  */
592 void
593 _gmx_sel_evaluate_subexpr(gmx_sel_evaluate_t                     *data,
594                           const gmx::SelectionTreeElementPointer &sel,
595                           gmx_ana_index_t                        *g)
596 {
597     gmx_ana_index_t      gmiss;
598
599     MempoolGroupReserver gmissreserver(data->mp);
600     if (sel->u.cgrp.isize == 0)
601     {
602         {
603             SelelemTemporaryValueAssigner assigner(sel->child, *sel);
604             sel->child->evaluate(data, sel->child, g);
605         }
606         gmx_ana_index_copy(&sel->u.cgrp, g, false);
607         gmiss.isize      = 0;
608     }
609     else
610     {
611         gmissreserver.reserve(&gmiss, g->isize);
612         gmx_ana_index_difference(&gmiss, g, &sel->u.cgrp);
613     }
614     if (gmiss.isize > 0)
615     {
616         MempoolSelelemReserver reserver(sel->child, gmiss.isize);
617         /* Evaluate the missing values for the child */
618         sel->child->evaluate(data, sel->child, &gmiss);
619         /* Merge the missing values to the existing ones. */
620         if (sel->v.type == GROUP_VALUE)
621         {
622             gmx_ana_index_merge(sel->v.u.g, sel->child->v.u.g, sel->v.u.g);
623         }
624         else
625         {
626             int  i, j, k;
627
628             i = sel->u.cgrp.isize - 1;
629             j = gmiss.isize - 1;
630             /* TODO: This switch is kind of ugly, but it may be difficult to
631              * do this portably without C++ templates. */
632             switch (sel->v.type)
633             {
634                 case INT_VALUE:
635                     for (k = sel->u.cgrp.isize + gmiss.isize - 1; k >= 0; k--)
636                     {
637                         if (i < 0 || (j >= 0 && sel->u.cgrp.index[i] < gmiss.index[j]))
638                         {
639                             sel->v.u.i[k] = sel->child->v.u.i[j--];
640                         }
641                         else
642                         {
643                             sel->v.u.i[k] = sel->v.u.i[i--];
644                         }
645                     }
646                     break;
647
648                 case REAL_VALUE:
649                     for (k = sel->u.cgrp.isize + gmiss.isize - 1; k >= 0; k--)
650                     {
651                         if (i < 0 || (j >= 0 && sel->u.cgrp.index[i] < gmiss.index[j]))
652                         {
653                             sel->v.u.r[k] = sel->child->v.u.r[j--];
654                         }
655                         else
656                         {
657                             sel->v.u.r[k] = sel->v.u.r[i--];
658                         }
659                     }
660                     break;
661
662                 case STR_VALUE:
663                     // Note: with the currently allowed syntax, this case is never
664                     // reached.
665                     for (k = sel->u.cgrp.isize + gmiss.isize - 1; k >= 0; k--)
666                     {
667                         if (i < 0 || (j >= 0 && sel->u.cgrp.index[i] < gmiss.index[j]))
668                         {
669                             sel->v.u.s[k] = sel->child->v.u.s[j--];
670                         }
671                         else
672                         {
673                             sel->v.u.s[k] = sel->v.u.s[i--];
674                         }
675                     }
676                     break;
677
678                 case POS_VALUE:
679                     /* TODO: Implement this */
680                     GMX_THROW(gmx::NotImplementedError("position subexpressions not implemented properly"));
681
682                 case NO_VALUE:
683                 case GROUP_VALUE:
684                     GMX_THROW(gmx::InternalError("Invalid subexpression type"));
685             }
686         }
687         gmx_ana_index_merge(&sel->u.cgrp, &sel->u.cgrp, &gmiss);
688     }
689 }
690
691 /*!
692  * \param[in] data Data for the current frame.
693  * \param[in] sel Selection element being evaluated.
694  * \param[in] g   Group for which \p sel should be evaluated.
695  * \returns   0 for success.
696  *
697  * Sets the value pointers of the child and its child to point to the same
698  * memory as the value pointer of this element to avoid copying, and then
699  * evaluates evaluates the child.
700  *
701  * This function is used as gmx::SelectionTreeElement:evaluate for
702  * \ref SEL_SUBEXPRREF elements for which the \ref SEL_SUBEXPR does not have
703  * other references.
704  */
705 void
706 _gmx_sel_evaluate_subexprref_simple(gmx_sel_evaluate_t                     *data,
707                                     const gmx::SelectionTreeElementPointer &sel,
708                                     gmx_ana_index_t                        *g)
709 {
710     if (g)
711     {
712         _gmx_selvalue_setstore(&sel->child->v, sel->v.u.ptr);
713         _gmx_selvalue_setstore_alloc(&sel->child->child->v, sel->v.u.ptr,
714                                      sel->child->child->v.nalloc);
715         sel->child->evaluate(data, sel->child, g);
716     }
717     sel->v.nr = sel->child->v.nr;
718     if (sel->u.param)
719     {
720         sel->u.param->val.nr = sel->v.nr;
721         if (sel->u.param->nvalptr)
722         {
723             *sel->u.param->nvalptr = sel->u.param->val.nr;
724         }
725     }
726 }
727
728 /*!
729  * \param[in] data Data for the current frame.
730  * \param[in] sel Selection element being evaluated.
731  * \param[in] g   Group for which \p sel should be evaluated.
732  * \returns   0 on success, a non-zero error code on error.
733  *
734  * If the value type is \ref POS_VALUE, the value of the child is simply
735  * copied to set the value of \p sel (the child subexpression should
736  * already have been evaluated by its root).
737  * If the value type is something else, the child is evaluated for the
738  * group \p g, and the value of the child is then copied.
739  * There should be only one child element.
740  *
741  * This function is used as gmx::SelectionTreeElement::evaluate for
742  * \ref SEL_SUBEXPRREF elements.
743  */
744 void
745 _gmx_sel_evaluate_subexprref(gmx_sel_evaluate_t                     *data,
746                              const gmx::SelectionTreeElementPointer &sel,
747                              gmx_ana_index_t                        *g)
748 {
749     int        i, j;
750
751     if (g != NULL && sel->child->evaluate != NULL)
752     {
753         sel->child->evaluate(data, sel->child, g);
754     }
755     const SelectionTreeElementPointer &expr = sel->child;
756     switch (sel->v.type)
757     {
758         case INT_VALUE:
759             if (!g)
760             {
761                 sel->v.nr = expr->v.nr;
762                 memcpy(sel->v.u.i, expr->v.u.i, sel->v.nr*sizeof(*sel->v.u.i));
763             }
764             else
765             {
766                 sel->v.nr = g->isize;
767                 /* Extract the values corresponding to g */
768                 for (i = j = 0; i < g->isize; ++i, ++j)
769                 {
770                     while (sel->child->u.cgrp.index[j] < g->index[i])
771                     {
772                         ++j;
773                     }
774                     sel->v.u.i[i] = expr->v.u.i[j];
775                 }
776             }
777             break;
778
779         case REAL_VALUE:
780             if (!g)
781             {
782                 sel->v.nr = expr->v.nr;
783                 memcpy(sel->v.u.r, expr->v.u.r, sel->v.nr*sizeof(*sel->v.u.r));
784             }
785             else
786             {
787                 sel->v.nr = g->isize;
788                 /* Extract the values corresponding to g */
789                 for (i = j = 0; i < g->isize; ++i, ++j)
790                 {
791                     while (sel->child->u.cgrp.index[j] < g->index[i])
792                     {
793                         ++j;
794                     }
795                     sel->v.u.r[i] = expr->v.u.r[j];
796                 }
797             }
798             break;
799
800         case STR_VALUE:
801             if (!g)
802             {
803                 sel->v.nr = expr->v.nr;
804                 memcpy(sel->v.u.s, expr->v.u.s, sel->v.nr*sizeof(*sel->v.u.s));
805             }
806             else
807             {
808                 sel->v.nr = g->isize;
809                 /* Extract the values corresponding to g */
810                 for (i = j = 0; i < g->isize; ++i, ++j)
811                 {
812                     while (sel->child->u.cgrp.index[j] < g->index[i])
813                     {
814                         ++j;
815                     }
816                     sel->v.u.s[i] = expr->v.u.s[j];
817                 }
818             }
819             break;
820
821         case POS_VALUE:
822             /* Currently, there is no need to do anything fancy here,
823              * but some future extensions may need a more flexible
824              * implementation. */
825             gmx_ana_pos_copy(sel->v.u.p, expr->v.u.p, false);
826             break;
827
828         case GROUP_VALUE:
829             if (!g)
830             {
831                 gmx_ana_index_copy(sel->v.u.g, expr->v.u.g, false);
832             }
833             else
834             {
835                 gmx_ana_index_intersection(sel->v.u.g, expr->v.u.g, g);
836             }
837             break;
838
839         default: /* should not be reached */
840             GMX_THROW(gmx::InternalError("Invalid subexpression reference type"));
841     }
842     /* Store the number of values if needed */
843     if (sel->u.param)
844     {
845         sel->u.param->val.nr = sel->v.nr;
846         if (sel->u.param->nvalptr)
847         {
848             *sel->u.param->nvalptr = sel->u.param->val.nr;
849         }
850     }
851 }
852
853 /********************************************************************
854  * METHOD EXPRESSION EVALUATION
855  ********************************************************************/
856
857 /*!
858  * \param[in] data Data for the current frame.
859  * \param[in] sel Selection element being evaluated.
860  * \param[in] g   Group for which \p sel should be evaluated.
861  * \returns   0 on success, a non-zero error code on error.
862  *
863  * Evaluates each child of a \ref SEL_EXPRESSION element.
864  * The value of \p sel is not touched.
865  *
866  * This function is not used as gmx::SelectionTreeElement::evaluate,
867  * but is used internally.
868  */
869 void
870 _gmx_sel_evaluate_method_params(gmx_sel_evaluate_t                     *data,
871                                 const gmx::SelectionTreeElementPointer &sel,
872                                 gmx_ana_index_t                        *g)
873 {
874     SelectionTreeElementPointer child = sel->child;
875     while (child)
876     {
877         if (child->evaluate && !(child->flags & SEL_EVALFRAME))
878         {
879             if (child->flags & SEL_ATOMVAL)
880             {
881                 child->evaluate(data, child, g);
882             }
883             else
884             {
885                 child->flags |= SEL_EVALFRAME;
886                 child->evaluate(data, child, NULL);
887             }
888         }
889         child = child->next;
890     }
891 }
892
893 /*!
894  * \param[in] data Data for the current frame.
895  * \param[in] sel Selection element being evaluated.
896  * \param[in] g   Group for which \p sel should be evaluated.
897  * \returns   0 on success, a non-zero error code on error.
898  *
899  * Evaluates all child selections (using _gmx_sel_evaluate_method_params())
900  * to evaluate any parameter values.
901  * If this is the first time this expression is evaluated for
902  * the frame, sel_framefunc() callback is called if one is provided.
903  * If a reference position calculation has been initialized for this element,
904  * the positions are also updated, and sel_updatefunc_pos() is used to
905  * evaluate the value. Otherwise, sel_updatefunc() is used.
906  *
907  * This function is used as gmx::SelectionTreeElement::evaluate for
908  * \ref SEL_EXPRESSION elements.
909  */
910 void
911 _gmx_sel_evaluate_method(gmx_sel_evaluate_t                     *data,
912                          const gmx::SelectionTreeElementPointer &sel,
913                          gmx_ana_index_t                        *g)
914 {
915     _gmx_sel_evaluate_method_params(data, sel, g);
916     if (sel->flags & SEL_INITFRAME)
917     {
918         sel->flags &= ~SEL_INITFRAME;
919         sel->u.expr.method->init_frame(data->top, data->fr, data->pbc,
920                                        sel->u.expr.mdata);
921     }
922     if (sel->u.expr.pc)
923     {
924         gmx_ana_poscalc_update(sel->u.expr.pc, sel->u.expr.pos, g,
925                                data->fr, data->pbc);
926         sel->u.expr.method->pupdate(data->top, data->fr, data->pbc,
927                                     sel->u.expr.pos, &sel->v,
928                                     sel->u.expr.mdata);
929     }
930     else
931     {
932         sel->u.expr.method->update(data->top, data->fr, data->pbc, g,
933                                    &sel->v, sel->u.expr.mdata);
934     }
935 }
936
937 /*!
938  * \param[in] data Data for the current frame.
939  * \param[in] sel Selection element being evaluated.
940  * \param[in] g   Group for which \p sel should be evaluated.
941  * \returns   0 on success, a non-zero error code on error.
942  *
943  * Evaluates all child selections (using _gmx_sel_evaluate_method_params())
944  * to evaluate any parameter values.
945  * If this is the first time this expression is evaluated for
946  * the frame, sel_framefunc() callback is called if one is provided.
947  * The modifier is then evaluated using sel_updatefunc_pos().
948  *
949  * This function is used as gmx::SelectionTreeElement::evaluate for
950  * \ref SEL_MODIFIER elements.
951  */
952 void
953 _gmx_sel_evaluate_modifier(gmx_sel_evaluate_t                     *data,
954                            const gmx::SelectionTreeElementPointer &sel,
955                            gmx_ana_index_t                        *g)
956 {
957     _gmx_sel_evaluate_method_params(data, sel, g);
958     if (sel->flags & SEL_INITFRAME)
959     {
960         sel->flags &= ~SEL_INITFRAME;
961         sel->u.expr.method->init_frame(data->top, data->fr, data->pbc,
962                                        sel->u.expr.mdata);
963     }
964     GMX_RELEASE_ASSERT(sel->child != NULL,
965                        "Modifier element with a value must have a child");
966     if (sel->child->v.type != POS_VALUE)
967     {
968         GMX_THROW(gmx::NotImplementedError("Non-position valued modifiers not implemented"));
969     }
970     sel->u.expr.method->pupdate(data->top, data->fr, data->pbc,
971                                 sel->child->v.u.p,
972                                 &sel->v, sel->u.expr.mdata);
973 }
974
975
976 /********************************************************************
977  * BOOLEAN EXPRESSION EVALUATION
978  ********************************************************************/
979
980 /*!
981  * \param[in] data Data for the current frame.
982  * \param[in] sel Selection element being evaluated.
983  * \param[in] g   Group for which \p sel should be evaluated.
984  * \returns   0 on success, a non-zero error code on error.
985  *
986  * Evaluates the child element (there should be only one) in the group
987  * \p g, and then sets the value of \p sel to the complement of the
988  * child value.
989  *
990  * This function is used as gmx::SelectionTreeElement::evaluate for
991  * \ref SEL_BOOLEAN elements with \ref BOOL_NOT.
992  */
993 void
994 _gmx_sel_evaluate_not(gmx_sel_evaluate_t                     *data,
995                       const gmx::SelectionTreeElementPointer &sel,
996                       gmx_ana_index_t                        *g)
997 {
998     MempoolSelelemReserver reserver(sel->child, g->isize);
999     sel->child->evaluate(data, sel->child, g);
1000     gmx_ana_index_difference(sel->v.u.g, g, sel->child->v.u.g);
1001 }
1002
1003 /*!
1004  * \param[in] data Data for the current frame.
1005  * \param[in] sel Selection element being evaluated.
1006  * \param[in] g   Group for which \p sel should be evaluated.
1007  * \returns   0 on success, a non-zero error code on error.
1008  *
1009  * Short-circuiting evaluation of logical AND expressions.
1010  *
1011  * Starts by evaluating the first child element in the group \p g.
1012  * The each following child element is evaluated in the intersection
1013  * of all the previous values until all children have been evaluated
1014  * or the intersection becomes empty.
1015  * The value of \p sel is set to the intersection of all the (evaluated)
1016  * child values.
1017  *
1018  * If the first child does not have an evaluation function, it is skipped
1019  * and the evaluation is started at the second child.
1020  * This happens if the first child is a constant expression and during
1021  * compilation it was detected that the evaluation group is always a subset
1022  * of the constant group
1023  * (currently, the compiler never detects this).
1024  *
1025  * This function is used as gmx::SelectionTreeElement::evaluate for
1026  * \ref SEL_BOOLEAN elements with \ref BOOL_AND.
1027  */
1028 void
1029 _gmx_sel_evaluate_and(gmx_sel_evaluate_t                     *data,
1030                       const gmx::SelectionTreeElementPointer &sel,
1031                       gmx_ana_index_t                        *g)
1032 {
1033     SelectionTreeElementPointer child = sel->child;
1034     /* Skip the first child if it does not have an evaluation function. */
1035     if (!child->evaluate)
1036     {
1037         child = child->next;
1038     }
1039     /* Evaluate the first child */
1040     {
1041         MempoolSelelemReserver reserver(child, g->isize);
1042         child->evaluate(data, child, g);
1043         gmx_ana_index_copy(sel->v.u.g, child->v.u.g, false);
1044     }
1045     child = child->next;
1046     while (child && sel->v.u.g->isize > 0)
1047     {
1048         MempoolSelelemReserver reserver(child, sel->v.u.g->isize);
1049         child->evaluate(data, child, sel->v.u.g);
1050         gmx_ana_index_intersection(sel->v.u.g, sel->v.u.g, child->v.u.g);
1051         child = child->next;
1052     }
1053 }
1054
1055 /*!
1056  * \param[in] data Data for the current frame.
1057  * \param[in] sel Selection element being evaluated.
1058  * \param[in] g   Group for which \p sel should be evaluated.
1059  * \returns   0 on success, a non-zero error code on error.
1060  *
1061  * Short-circuiting evaluation of logical OR expressions.
1062  *
1063  * Starts by evaluating the first child element in the group \p g.
1064  * For each subsequent child, finds the part of \p g that is not
1065  * included the value of any previous child, and evaluates the child
1066  * in that group until the last child is evaluated or all of \p g
1067  * is included in some child value.
1068  * The value of \p sel is set to the union of all the (evaluated)
1069  * child values.
1070  *
1071  * If the first child does not have an evaluation function, its value is
1072  * used without evaluation.
1073  * This happens if the first child is a constant expression, the selection
1074  * has been compiled, and the evaluation group is the same for each frame.
1075  * In this case, the compiler has taken care of that the child value is a
1076  * subset of \p g, making it unnecessary to evaluate it.
1077  *
1078  * This function is used as gmx::SelectionTreeElement::evaluate for
1079  * \ref SEL_BOOLEAN elements with \ref BOOL_OR.
1080  */
1081 void
1082 _gmx_sel_evaluate_or(gmx_sel_evaluate_t                     *data,
1083                      const gmx::SelectionTreeElementPointer &sel,
1084                      gmx_ana_index_t                        *g)
1085 {
1086     gmx_ana_index_t             tmp, tmp2;
1087
1088     SelectionTreeElementPointer child = sel->child;
1089     if (child->evaluate)
1090     {
1091         MempoolSelelemReserver reserver(child, g->isize);
1092         child->evaluate(data, child, g);
1093         gmx_ana_index_partition(sel->v.u.g, &tmp, g, child->v.u.g);
1094     }
1095     else
1096     {
1097         gmx_ana_index_partition(sel->v.u.g, &tmp, g, child->v.u.g);
1098     }
1099     child = child->next;
1100     while (child && tmp.isize > 0)
1101     {
1102         {
1103             MempoolSelelemReserver reserver(child, tmp.isize);
1104             child->evaluate(data, child, &tmp);
1105             gmx_ana_index_partition(&tmp, &tmp2, &tmp, child->v.u.g);
1106         }
1107         sel->v.u.g->isize += tmp.isize;
1108         tmp.isize          = tmp2.isize;
1109         tmp.index          = tmp2.index;
1110         child              = child->next;
1111     }
1112     gmx_ana_index_sort(sel->v.u.g);
1113 }
1114
1115
1116 /********************************************************************
1117  * ARITHMETIC EVALUATION
1118  ********************************************************************/
1119
1120 /*!
1121  * \param[in] data Data for the current frame.
1122  * \param[in] sel  Selection element being evaluated.
1123  * \param[in] g    Group for which \p sel should be evaluated.
1124  * \returns   0 on success, a non-zero error code on error.
1125  */
1126 void
1127 _gmx_sel_evaluate_arithmetic(gmx_sel_evaluate_t                     *data,
1128                              const gmx::SelectionTreeElementPointer &sel,
1129                              gmx_ana_index_t                        *g)
1130 {
1131     int         n, i, i1, i2;
1132     real        lval, rval = 0., val = 0.;
1133
1134     const SelectionTreeElementPointer &left  = sel->child;
1135     const SelectionTreeElementPointer &right = left->next;
1136
1137     SelelemTemporaryValueAssigner      assigner;
1138     MempoolSelelemReserver             reserver;
1139     if (left->mempool)
1140     {
1141         assigner.assign(left, *sel);
1142         if (right)
1143         {
1144             reserver.reserve(right, g->isize);
1145         }
1146     }
1147     else if (right && right->mempool)
1148     {
1149         assigner.assign(right, *sel);
1150     }
1151     _gmx_sel_evaluate_children(data, sel, g);
1152
1153     n         = (sel->flags & SEL_SINGLEVAL) ? 1 : g->isize;
1154     sel->v.nr = n;
1155
1156     bool bArithNeg = (sel->u.arith.type == ARITH_NEG);
1157     GMX_ASSERT(right || bArithNeg,
1158                "Right operand cannot be null except for negations");
1159     for (i = i1 = i2 = 0; i < n; ++i)
1160     {
1161         lval = left->v.u.r[i1];
1162         if (!bArithNeg)
1163         {
1164             rval = right->v.u.r[i2];
1165         }
1166         switch (sel->u.arith.type)
1167         {
1168             case ARITH_PLUS:    val = lval + rval;     break;
1169             case ARITH_MINUS:   val = lval - rval;     break;
1170             case ARITH_NEG:     val = -lval;           break;
1171             case ARITH_MULT:    val = lval * rval;     break;
1172             case ARITH_DIV:     val = lval / rval;     break;
1173             case ARITH_EXP:     val = pow(lval, rval); break;
1174         }
1175         sel->v.u.r[i] = val;
1176         if (!(left->flags & SEL_SINGLEVAL))
1177         {
1178             ++i1;
1179         }
1180         if (!bArithNeg && !(right->flags & SEL_SINGLEVAL))
1181         {
1182             ++i2;
1183         }
1184     }
1185 }