2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * Implements functions in evaluate.h.
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
50 * \author Teemu Murtola <teemu.murtola@gmail.com>
51 * \ingroup module_selection
57 #include "gromacs/math/utilities.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/selection/indexutil.h"
60 #include "gromacs/selection/poscalc.h"
61 #include "gromacs/selection/selection.h"
62 #include "gromacs/selection/selmethod.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/smalloc.h"
69 #include "selectioncollection-impl.h"
72 using gmx::SelectionTreeElement;
73 using gmx::SelectionTreeElementPointer;
79 * Reserves memory for a selection element from the evaluation memory pool.
81 * This class implements RAII semantics for allocating memory for selection
82 * element values from a selection evaluation memory pool.
84 * \ingroup module_selection
86 class MempoolSelelemReserver
89 //! Constructs a reserver without initial reservation.
90 MempoolSelelemReserver() {}
92 * Constructs a reserver with initial reservation.
94 * \param[in,out] sel Selection element for which to reserve.
95 * \param[in] count Number of values to reserve.
99 MempoolSelelemReserver(const SelectionTreeElementPointer &sel, int count)
103 //! Frees any memory allocated using this reserver.
104 ~MempoolSelelemReserver()
108 sel_->mempoolRelease();
113 * Reserves memory for selection element values using this reserver.
115 * \param[in,out] sel Selection element for which to reserve.
116 * \param[in] count Number of values to reserve.
118 * Allocates space to store \p count output values in \p sel from the
119 * memory pool associated with \p sel, or from the heap if there is no
120 * memory pool. Type of values to allocate is automatically determined
123 void reserve(const SelectionTreeElementPointer &sel, int count)
125 GMX_RELEASE_ASSERT(!sel_,
126 "Can only reserve one element with one instance");
127 sel->mempoolReserve(count);
132 SelectionTreeElementPointer sel_;
136 * Reserves memory for an index group from the evaluation memory pool.
138 * This class implements RAII semantics for allocating memory for an index
139 * group from a selection evaluation memory pool.
141 * \ingroup module_selection
143 class MempoolGroupReserver
147 * Creates a reserver associated with a given memory pool.
149 * \param mp Memory pool from which to reserve memory.
151 explicit MempoolGroupReserver(gmx_sel_mempool_t *mp)
155 //! Frees any memory allocated using this reserver.
156 ~MempoolGroupReserver()
160 _gmx_sel_mempool_free_group(mp_, g_);
165 * Reserves memory for an index group using this reserver.
167 * \param[in,out] g Index group to reserve.
168 * \param[in] count Number of atoms to reserve space for.
170 * Allocates memory from the memory pool to store \p count atoms in
173 void reserve(gmx_ana_index_t *g, int count)
175 GMX_RELEASE_ASSERT(g_ == NULL, "Can only reserve one element with one instance");
176 _gmx_sel_mempool_alloc_group(mp_, g, count);
181 gmx_sel_mempool_t *mp_;
186 * Assigns a temporary value for a selection element.
188 * This class implements RAII semantics for temporarily assigning the value
189 * pointer of a selection element to point to a different location.
191 * \ingroup module_selection
193 class SelelemTemporaryValueAssigner
196 //! Constructs an assigner without an initial assignment.
197 SelelemTemporaryValueAssigner()
198 : old_ptr_(NULL), old_nalloc_(0)
202 * Constructs an assigner with an initial assignment.
204 * \param[in,out] sel Selection element for which to assign.
205 * \param[in] vsource Element to which \p sel values will point to.
209 SelelemTemporaryValueAssigner(const SelectionTreeElementPointer &sel,
210 const SelectionTreeElement &vsource)
212 assign(sel, vsource);
214 //! Undoes any temporary assignment done using this assigner.
215 ~SelelemTemporaryValueAssigner()
219 _gmx_selvalue_setstore_alloc(&sel_->v, old_ptr_, old_nalloc_);
224 * Assigns a temporary value pointer.
226 * \param[in,out] sel Selection element for which to assign.
227 * \param[in] vsource Element to which \p sel values will point to.
229 * Assigns the value pointer in \p sel to point to the values in
230 * \p vsource, i.e., any access/modification to values in \p sel
231 * actually accesses values in \p vsource.
233 void assign(const SelectionTreeElementPointer &sel,
234 const SelectionTreeElement &vsource)
236 GMX_RELEASE_ASSERT(!sel_,
237 "Can only assign one element with one instance");
238 GMX_RELEASE_ASSERT(sel->v.type == vsource.v.type,
239 "Mismatching selection value types");
240 old_ptr_ = sel->v.u.ptr;
241 old_nalloc_ = sel->v.nalloc;
242 _gmx_selvalue_setstore(&sel->v, vsource.v.u.ptr);
247 SelectionTreeElementPointer sel_;
255 * \param[in] fp File handle to receive the output.
256 * \param[in] evalfunc Function pointer to print.
259 _gmx_sel_print_evalfunc_name(FILE *fp, gmx::sel_evalfunc evalfunc)
265 else if (evalfunc == &_gmx_sel_evaluate_root)
269 else if (evalfunc == &_gmx_sel_evaluate_static)
271 fprintf(fp, "static");
273 else if (evalfunc == &_gmx_sel_evaluate_subexpr_simple)
275 fprintf(fp, "subexpr_simple");
277 else if (evalfunc == &_gmx_sel_evaluate_subexpr_staticeval)
279 fprintf(fp, "subexpr_staticeval");
281 else if (evalfunc == &_gmx_sel_evaluate_subexpr)
283 fprintf(fp, "subexpr");
285 else if (evalfunc == &_gmx_sel_evaluate_subexprref_simple)
287 fprintf(fp, "ref_simple");
289 else if (evalfunc == &_gmx_sel_evaluate_subexprref)
293 else if (evalfunc == &_gmx_sel_evaluate_method)
295 fprintf(fp, "method");
297 else if (evalfunc == &_gmx_sel_evaluate_modifier)
301 else if (evalfunc == &_gmx_sel_evaluate_not)
305 else if (evalfunc == &_gmx_sel_evaluate_and)
309 else if (evalfunc == &_gmx_sel_evaluate_or)
313 else if (evalfunc == &_gmx_sel_evaluate_arithmetic)
315 fprintf(fp, "arithmetic");
319 fprintf(fp, "%p", (void*)(evalfunc));
324 * \param[out] data Evaluation data structure to initialize.
325 * \param[in] mp Memory pool for intermediate evaluation values.
326 * \param[in] gall Index group with all the atoms.
327 * \param[in] top Topology structure for evaluation.
328 * \param[in] fr New frame for evaluation.
329 * \param[in] pbc New PBC information for evaluation.
332 _gmx_sel_evaluate_init(gmx_sel_evaluate_t *data,
333 gmx_sel_mempool_t *mp, gmx_ana_index_t *gall,
334 t_topology *top, t_trxframe *fr, t_pbc *pbc)
344 * Recursively initializes the flags for evaluation.
346 * \param[in,out] sel Selection element to clear.
348 * The \ref SEL_INITFRAME flag is set for \ref SEL_EXPRESSION elements whose
349 * method defines the \p init_frame callback (see sel_framefunc()), and
350 * cleared for other elements.
352 * The \ref SEL_EVALFRAME flag is cleared for all elements.
355 init_frame_eval(SelectionTreeElementPointer sel)
359 sel->flags &= ~(SEL_INITFRAME | SEL_EVALFRAME);
360 if (sel->type == SEL_EXPRESSION)
362 if (sel->u.expr.method && sel->u.expr.method->init_frame)
364 sel->flags |= SEL_INITFRAME;
367 if (sel->child && sel->type != SEL_SUBEXPRREF)
369 init_frame_eval(sel->child);
378 SelectionEvaluator::SelectionEvaluator()
383 * \param[in,out] coll The selection collection to evaluate.
384 * \param[in] fr Frame for which the evaluation should be carried out.
385 * \param[in] pbc PBC data, or NULL if no PBC should be used.
386 * \returns 0 on successful evaluation, a non-zero error code on error.
388 * This functions sets the global variables for topology, frame and PBC,
389 * clears some information in the selection to initialize the evaluation
390 * for a new frame, and evaluates \p sel and all the selections pointed by
391 * the \p next pointers of \p sel.
393 * This is the only function that user code should call if they want to
394 * evaluate a selection for a new frame.
397 SelectionEvaluator::evaluate(SelectionCollection *coll,
398 t_trxframe *fr, t_pbc *pbc)
400 gmx_ana_selcollection_t *sc = &coll->impl_->sc_;
401 gmx_sel_evaluate_t data;
403 _gmx_sel_evaluate_init(&data, sc->mempool, &sc->gall, sc->top, fr, pbc);
404 init_frame_eval(sc->root);
405 SelectionTreeElementPointer sel = sc->root;
408 /* Clear the evaluation group of subexpressions */
409 if (sel->child && sel->child->type == SEL_SUBEXPR
410 && sel->child->evaluate != NULL)
412 sel->child->u.cgrp.isize = 0;
413 /* Not strictly necessary, because the value will be overwritten
414 * during first evaluation of the subexpression anyways, but we
415 * clear the group for clarity. Note that this is _not_ done during
416 * compilation because of some additional complexities involved
417 * (see compiler.cpp), so it should not be relied upon in
418 * _gmx_sel_evaluate_subexpr(). */
419 if (sel->child->v.type == GROUP_VALUE)
421 sel->child->v.u.g->isize = 0;
426 sel->evaluate(&data, sel, NULL);
430 /* Update selection information */
431 SelectionDataList::const_iterator isel;
432 for (isel = sc->sel.begin(); isel != sc->sel.end(); ++isel)
434 internal::SelectionData &sel = **isel;
435 sel.refreshMassesAndCharges(sc->top);
436 sel.updateCoveredFractionForFrame();
441 * \param[in,out] coll The selection collection to evaluate.
442 * \param[in] nframes Total number of frames.
445 SelectionEvaluator::evaluateFinal(SelectionCollection *coll, int nframes)
447 gmx_ana_selcollection_t *sc = &coll->impl_->sc_;
449 SelectionDataList::const_iterator isel;
450 for (isel = sc->sel.begin(); isel != sc->sel.end(); ++isel)
452 internal::SelectionData &sel = **isel;
453 sel.restoreOriginalPositions(sc->top);
454 sel.computeAverageCoveredFraction(nframes);
461 * \param[in] data Data for the current frame.
462 * \param[in] sel Selection element being evaluated.
463 * \param[in] g Group for which \p sel should be evaluated.
464 * \returns 0 on success, a non-zero error code on error.
466 * Evaluates each child of \p sel in \p g.
469 _gmx_sel_evaluate_children(gmx_sel_evaluate_t *data,
470 const gmx::SelectionTreeElementPointer &sel,
473 SelectionTreeElementPointer child = sel->child;
478 child->evaluate(data, child, g);
485 _gmx_sel_evaluate_root(gmx_sel_evaluate_t *data,
486 const gmx::SelectionTreeElementPointer &sel,
487 gmx_ana_index_t * /* g */)
489 if (sel->u.cgrp.isize == 0 || !sel->child->evaluate)
494 sel->child->evaluate(data, sel->child,
495 sel->u.cgrp.isize < 0 ? NULL : &sel->u.cgrp);
499 _gmx_sel_evaluate_static(gmx_sel_evaluate_t * /* data */,
500 const gmx::SelectionTreeElementPointer &sel,
503 if (sel->flags & SEL_UNSORTED)
505 // This only works if g contains all the atoms, but that is currently
506 // the only supported case.
507 gmx_ana_index_copy(sel->v.u.g, &sel->u.cgrp, false);
511 gmx_ana_index_intersection(sel->v.u.g, &sel->u.cgrp, g);
516 /*********************************************************************
517 * SUBEXPRESSION EVALUATION
518 *********************************************************************/
521 * \param[in] data Data for the current frame.
522 * \param[in] sel Selection element being evaluated.
523 * \param[in] g Group for which \p sel should be evaluated.
524 * \returns 0 on success, a non-zero error code on error.
526 * Evaluates the child element (there should be exactly one) in \p g.
527 * The compiler has taken care that the child actually stores the evaluated
528 * value in the value pointer of this element.
530 * This function is used as gmx::SelectionTreeElement::evaluate for
531 * \ref SEL_SUBEXPR elements that are used only once, and hence do not need
532 * full subexpression handling.
535 _gmx_sel_evaluate_subexpr_simple(gmx_sel_evaluate_t *data,
536 const gmx::SelectionTreeElementPointer &sel,
539 if (sel->child->evaluate)
541 sel->child->evaluate(data, sel->child, g);
543 sel->v.nr = sel->child->v.nr;
547 * \param[in] data Data for the current frame.
548 * \param[in] sel Selection element being evaluated.
549 * \param[in] g Group for which \p sel should be evaluated.
550 * \returns 0 on success, a non-zero error code on error.
552 * If this is the first call for this frame, evaluates the child element
553 * there should be exactly one in \p g.
554 * The compiler has taken care that the child actually stores the evaluated
555 * value in the value pointer of this element.
556 * Assumes that \p g is persistent for the duration of the whole evaluation.
558 * This function is used as gmx::SelectionTreeElement::evaluate for
559 * \ref SEL_SUBEXPR elements that have a static evaluation group, and hence do
560 * not need full subexpression handling.
563 _gmx_sel_evaluate_subexpr_staticeval(gmx_sel_evaluate_t *data,
564 const gmx::SelectionTreeElementPointer &sel,
567 if (sel->u.cgrp.isize == 0)
569 sel->child->evaluate(data, sel->child, g);
570 sel->v.nr = sel->child->v.nr;
573 sel->u.cgrp.isize = -1;
577 gmx_ana_index_set(&sel->u.cgrp, g->isize, g->index, 0);
583 * \param[in] data Data for the current frame.
584 * \param[in] sel Selection element being evaluated.
585 * \param[in] g Group for which \p sel should be evaluated.
586 * \returns 0 on success, a non-zero error code on error.
588 * Finds the part of \p g for which the subexpression
589 * has not yet been evaluated by comparing \p g to \p sel->u.cgrp.
590 * If the part is not empty, the child expression is evaluated for this
591 * part, and the results merged to the old values of the child.
592 * The value of \p sel itself is undefined after the call.
595 * The call to gmx_ana_index_difference() can take quite a lot of unnecessary
596 * time if the subexpression is evaluated either several times for the same
597 * group or for completely distinct groups.
598 * However, in the majority of cases, these situations occur when
599 * _gmx_sel_evaluate_subexpr_staticeval() can be used, so this should not be a
603 _gmx_sel_evaluate_subexpr(gmx_sel_evaluate_t *data,
604 const gmx::SelectionTreeElementPointer &sel,
607 gmx_ana_index_t gmiss;
609 MempoolGroupReserver gmissreserver(data->mp);
610 if (sel->u.cgrp.isize == 0)
613 SelelemTemporaryValueAssigner assigner(sel->child, *sel);
614 sel->child->evaluate(data, sel->child, g);
616 gmx_ana_index_copy(&sel->u.cgrp, g, false);
621 gmissreserver.reserve(&gmiss, g->isize);
622 gmx_ana_index_difference(&gmiss, g, &sel->u.cgrp);
626 MempoolSelelemReserver reserver(sel->child, gmiss.isize);
627 /* Evaluate the missing values for the child */
628 sel->child->evaluate(data, sel->child, &gmiss);
629 /* Merge the missing values to the existing ones. */
630 if (sel->v.type == GROUP_VALUE)
632 gmx_ana_index_merge(sel->v.u.g, sel->child->v.u.g, sel->v.u.g);
638 i = sel->u.cgrp.isize - 1;
640 /* TODO: This switch is kind of ugly, but it may be difficult to
641 * do this portably without C++ templates. */
645 for (k = sel->u.cgrp.isize + gmiss.isize - 1; k >= 0; k--)
647 if (i < 0 || (j >= 0 && sel->u.cgrp.index[i] < gmiss.index[j]))
649 sel->v.u.i[k] = sel->child->v.u.i[j--];
653 sel->v.u.i[k] = sel->v.u.i[i--];
659 for (k = sel->u.cgrp.isize + gmiss.isize - 1; k >= 0; k--)
661 if (i < 0 || (j >= 0 && sel->u.cgrp.index[i] < gmiss.index[j]))
663 sel->v.u.r[k] = sel->child->v.u.r[j--];
667 sel->v.u.r[k] = sel->v.u.r[i--];
673 // Note: with the currently allowed syntax, this case is never
675 for (k = sel->u.cgrp.isize + gmiss.isize - 1; k >= 0; k--)
677 if (i < 0 || (j >= 0 && sel->u.cgrp.index[i] < gmiss.index[j]))
679 sel->v.u.s[k] = sel->child->v.u.s[j--];
683 sel->v.u.s[k] = sel->v.u.s[i--];
689 /* TODO: Implement this */
690 GMX_THROW(gmx::NotImplementedError("position subexpressions not implemented properly"));
694 GMX_THROW(gmx::InternalError("Invalid subexpression type"));
697 gmx_ana_index_merge(&sel->u.cgrp, &sel->u.cgrp, &gmiss);
702 * \param[in] data Data for the current frame.
703 * \param[in] sel Selection element being evaluated.
704 * \param[in] g Group for which \p sel should be evaluated.
705 * \returns 0 for success.
707 * Sets the value pointers of the child and its child to point to the same
708 * memory as the value pointer of this element to avoid copying, and then
709 * evaluates evaluates the child.
711 * This function is used as gmx::SelectionTreeElement:evaluate for
712 * \ref SEL_SUBEXPRREF elements for which the \ref SEL_SUBEXPR does not have
716 _gmx_sel_evaluate_subexprref_simple(gmx_sel_evaluate_t *data,
717 const gmx::SelectionTreeElementPointer &sel,
722 _gmx_selvalue_setstore(&sel->child->v, sel->v.u.ptr);
723 _gmx_selvalue_setstore_alloc(&sel->child->child->v, sel->v.u.ptr,
724 sel->child->child->v.nalloc);
725 sel->child->evaluate(data, sel->child, g);
727 sel->v.nr = sel->child->v.nr;
730 sel->u.param->val.nr = sel->v.nr;
731 if (sel->u.param->nvalptr)
733 *sel->u.param->nvalptr = sel->u.param->val.nr;
739 * \param[in] data Data for the current frame.
740 * \param[in] sel Selection element being evaluated.
741 * \param[in] g Group for which \p sel should be evaluated.
742 * \returns 0 on success, a non-zero error code on error.
744 * If the value type is \ref POS_VALUE, the value of the child is simply
745 * copied to set the value of \p sel (the child subexpression should
746 * already have been evaluated by its root).
747 * If the value type is something else, the child is evaluated for the
748 * group \p g, and the value of the child is then copied.
749 * There should be only one child element.
751 * This function is used as gmx::SelectionTreeElement::evaluate for
752 * \ref SEL_SUBEXPRREF elements.
755 _gmx_sel_evaluate_subexprref(gmx_sel_evaluate_t *data,
756 const gmx::SelectionTreeElementPointer &sel,
761 if (g != NULL && sel->child->evaluate != NULL)
763 sel->child->evaluate(data, sel->child, g);
765 const SelectionTreeElementPointer &expr = sel->child;
771 sel->v.nr = expr->v.nr;
772 memcpy(sel->v.u.i, expr->v.u.i, sel->v.nr*sizeof(*sel->v.u.i));
776 sel->v.nr = g->isize;
777 /* Extract the values corresponding to g */
778 for (i = j = 0; i < g->isize; ++i, ++j)
780 while (sel->child->u.cgrp.index[j] < g->index[i])
784 sel->v.u.i[i] = expr->v.u.i[j];
792 sel->v.nr = expr->v.nr;
793 memcpy(sel->v.u.r, expr->v.u.r, sel->v.nr*sizeof(*sel->v.u.r));
797 sel->v.nr = g->isize;
798 /* Extract the values corresponding to g */
799 for (i = j = 0; i < g->isize; ++i, ++j)
801 while (sel->child->u.cgrp.index[j] < g->index[i])
805 sel->v.u.r[i] = expr->v.u.r[j];
813 sel->v.nr = expr->v.nr;
814 memcpy(sel->v.u.s, expr->v.u.s, sel->v.nr*sizeof(*sel->v.u.s));
818 sel->v.nr = g->isize;
819 /* Extract the values corresponding to g */
820 for (i = j = 0; i < g->isize; ++i, ++j)
822 while (sel->child->u.cgrp.index[j] < g->index[i])
826 sel->v.u.s[i] = expr->v.u.s[j];
832 /* Currently, there is no need to do anything fancy here,
833 * but some future extensions may need a more flexible
835 gmx_ana_pos_copy(sel->v.u.p, expr->v.u.p, false);
841 gmx_ana_index_copy(sel->v.u.g, expr->v.u.g, false);
845 gmx_ana_index_intersection(sel->v.u.g, expr->v.u.g, g);
849 default: /* should not be reached */
850 GMX_THROW(gmx::InternalError("Invalid subexpression reference type"));
852 /* Store the number of values if needed */
855 sel->u.param->val.nr = sel->v.nr;
856 if (sel->u.param->nvalptr)
858 *sel->u.param->nvalptr = sel->u.param->val.nr;
863 /********************************************************************
864 * METHOD EXPRESSION EVALUATION
865 ********************************************************************/
868 * \param[in] data Data for the current frame.
869 * \param[in] sel Selection element being evaluated.
870 * \param[in] g Group for which \p sel should be evaluated.
871 * \returns 0 on success, a non-zero error code on error.
873 * Evaluates each child of a \ref SEL_EXPRESSION element.
874 * The value of \p sel is not touched.
876 * This function is not used as gmx::SelectionTreeElement::evaluate,
877 * but is used internally.
880 _gmx_sel_evaluate_method_params(gmx_sel_evaluate_t *data,
881 const gmx::SelectionTreeElementPointer &sel,
884 SelectionTreeElementPointer child = sel->child;
887 if (child->evaluate && !(child->flags & SEL_EVALFRAME))
889 if (child->flags & SEL_ATOMVAL)
891 child->evaluate(data, child, g);
895 child->flags |= SEL_EVALFRAME;
896 child->evaluate(data, child, NULL);
904 * \param[in] data Data for the current frame.
905 * \param[in] sel Selection element being evaluated.
906 * \param[in] g Group for which \p sel should be evaluated.
907 * \returns 0 on success, a non-zero error code on error.
909 * Evaluates all child selections (using _gmx_sel_evaluate_method_params())
910 * to evaluate any parameter values.
911 * If this is the first time this expression is evaluated for
912 * the frame, sel_framefunc() callback is called if one is provided.
913 * If a reference position calculation has been initialized for this element,
914 * the positions are also updated, and sel_updatefunc_pos() is used to
915 * evaluate the value. Otherwise, sel_updatefunc() is used.
917 * This function is used as gmx::SelectionTreeElement::evaluate for
918 * \ref SEL_EXPRESSION elements.
921 _gmx_sel_evaluate_method(gmx_sel_evaluate_t *data,
922 const gmx::SelectionTreeElementPointer &sel,
925 _gmx_sel_evaluate_method_params(data, sel, g);
926 if (sel->flags & SEL_INITFRAME)
928 sel->flags &= ~SEL_INITFRAME;
929 sel->u.expr.method->init_frame(data->top, data->fr, data->pbc,
934 gmx_ana_poscalc_update(sel->u.expr.pc, sel->u.expr.pos, g,
935 data->fr, data->pbc);
936 sel->u.expr.method->pupdate(data->top, data->fr, data->pbc,
937 sel->u.expr.pos, &sel->v,
942 sel->u.expr.method->update(data->top, data->fr, data->pbc, g,
943 &sel->v, sel->u.expr.mdata);
948 * \param[in] data Data for the current frame.
949 * \param[in] sel Selection element being evaluated.
950 * \param[in] g Group for which \p sel should be evaluated.
951 * \returns 0 on success, a non-zero error code on error.
953 * Evaluates all child selections (using _gmx_sel_evaluate_method_params())
954 * to evaluate any parameter values.
955 * If this is the first time this expression is evaluated for
956 * the frame, sel_framefunc() callback is called if one is provided.
957 * The modifier is then evaluated using sel_updatefunc_pos().
959 * This function is used as gmx::SelectionTreeElement::evaluate for
960 * \ref SEL_MODIFIER elements.
963 _gmx_sel_evaluate_modifier(gmx_sel_evaluate_t *data,
964 const gmx::SelectionTreeElementPointer &sel,
967 _gmx_sel_evaluate_method_params(data, sel, g);
968 if (sel->flags & SEL_INITFRAME)
970 sel->flags &= ~SEL_INITFRAME;
971 sel->u.expr.method->init_frame(data->top, data->fr, data->pbc,
974 GMX_RELEASE_ASSERT(sel->child != NULL,
975 "Modifier element with a value must have a child");
976 if (sel->child->v.type != POS_VALUE)
978 GMX_THROW(gmx::NotImplementedError("Non-position valued modifiers not implemented"));
980 sel->u.expr.method->pupdate(data->top, data->fr, data->pbc,
982 &sel->v, sel->u.expr.mdata);
986 /********************************************************************
987 * BOOLEAN EXPRESSION EVALUATION
988 ********************************************************************/
991 * \param[in] data Data for the current frame.
992 * \param[in] sel Selection element being evaluated.
993 * \param[in] g Group for which \p sel should be evaluated.
994 * \returns 0 on success, a non-zero error code on error.
996 * Evaluates the child element (there should be only one) in the group
997 * \p g, and then sets the value of \p sel to the complement of the
1000 * This function is used as gmx::SelectionTreeElement::evaluate for
1001 * \ref SEL_BOOLEAN elements with \ref BOOL_NOT.
1004 _gmx_sel_evaluate_not(gmx_sel_evaluate_t *data,
1005 const gmx::SelectionTreeElementPointer &sel,
1008 MempoolSelelemReserver reserver(sel->child, g->isize);
1009 sel->child->evaluate(data, sel->child, g);
1010 gmx_ana_index_difference(sel->v.u.g, g, sel->child->v.u.g);
1014 * \param[in] data Data for the current frame.
1015 * \param[in] sel Selection element being evaluated.
1016 * \param[in] g Group for which \p sel should be evaluated.
1017 * \returns 0 on success, a non-zero error code on error.
1019 * Short-circuiting evaluation of logical AND expressions.
1021 * Starts by evaluating the first child element in the group \p g.
1022 * The each following child element is evaluated in the intersection
1023 * of all the previous values until all children have been evaluated
1024 * or the intersection becomes empty.
1025 * The value of \p sel is set to the intersection of all the (evaluated)
1028 * If the first child does not have an evaluation function, it is skipped
1029 * and the evaluation is started at the second child.
1030 * This happens if the first child is a constant expression and during
1031 * compilation it was detected that the evaluation group is always a subset
1032 * of the constant group
1033 * (currently, the compiler never detects this).
1035 * This function is used as gmx::SelectionTreeElement::evaluate for
1036 * \ref SEL_BOOLEAN elements with \ref BOOL_AND.
1039 _gmx_sel_evaluate_and(gmx_sel_evaluate_t *data,
1040 const gmx::SelectionTreeElementPointer &sel,
1043 SelectionTreeElementPointer child = sel->child;
1044 /* Skip the first child if it does not have an evaluation function. */
1045 if (!child->evaluate)
1047 child = child->next;
1049 /* Evaluate the first child */
1051 MempoolSelelemReserver reserver(child, g->isize);
1052 child->evaluate(data, child, g);
1053 gmx_ana_index_copy(sel->v.u.g, child->v.u.g, false);
1055 child = child->next;
1056 while (child && sel->v.u.g->isize > 0)
1058 MempoolSelelemReserver reserver(child, sel->v.u.g->isize);
1059 child->evaluate(data, child, sel->v.u.g);
1060 gmx_ana_index_intersection(sel->v.u.g, sel->v.u.g, child->v.u.g);
1061 child = child->next;
1066 * \param[in] data Data for the current frame.
1067 * \param[in] sel Selection element being evaluated.
1068 * \param[in] g Group for which \p sel should be evaluated.
1069 * \returns 0 on success, a non-zero error code on error.
1071 * Short-circuiting evaluation of logical OR expressions.
1073 * Starts by evaluating the first child element in the group \p g.
1074 * For each subsequent child, finds the part of \p g that is not
1075 * included the value of any previous child, and evaluates the child
1076 * in that group until the last child is evaluated or all of \p g
1077 * is included in some child value.
1078 * The value of \p sel is set to the union of all the (evaluated)
1081 * If the first child does not have an evaluation function, its value is
1082 * used without evaluation.
1083 * This happens if the first child is a constant expression, the selection
1084 * has been compiled, and the evaluation group is the same for each frame.
1085 * In this case, the compiler has taken care of that the child value is a
1086 * subset of \p g, making it unnecessary to evaluate it.
1088 * This function is used as gmx::SelectionTreeElement::evaluate for
1089 * \ref SEL_BOOLEAN elements with \ref BOOL_OR.
1092 _gmx_sel_evaluate_or(gmx_sel_evaluate_t *data,
1093 const gmx::SelectionTreeElementPointer &sel,
1096 gmx_ana_index_t tmp, tmp2;
1098 SelectionTreeElementPointer child = sel->child;
1099 if (child->evaluate)
1101 MempoolSelelemReserver reserver(child, g->isize);
1102 child->evaluate(data, child, g);
1103 gmx_ana_index_partition(sel->v.u.g, &tmp, g, child->v.u.g);
1107 gmx_ana_index_partition(sel->v.u.g, &tmp, g, child->v.u.g);
1109 child = child->next;
1110 while (child && tmp.isize > 0)
1113 MempoolSelelemReserver reserver(child, tmp.isize);
1114 child->evaluate(data, child, &tmp);
1115 gmx_ana_index_partition(&tmp, &tmp2, &tmp, child->v.u.g);
1117 sel->v.u.g->isize += tmp.isize;
1118 tmp.isize = tmp2.isize;
1119 tmp.index = tmp2.index;
1120 child = child->next;
1122 gmx_ana_index_sort(sel->v.u.g);
1126 /********************************************************************
1127 * ARITHMETIC EVALUATION
1128 ********************************************************************/
1131 * \param[in] data Data for the current frame.
1132 * \param[in] sel Selection element being evaluated.
1133 * \param[in] g Group for which \p sel should be evaluated.
1134 * \returns 0 on success, a non-zero error code on error.
1137 _gmx_sel_evaluate_arithmetic(gmx_sel_evaluate_t *data,
1138 const gmx::SelectionTreeElementPointer &sel,
1142 real lval, rval = 0., val = 0.;
1144 const SelectionTreeElementPointer &left = sel->child;
1145 const SelectionTreeElementPointer &right = left->next;
1147 SelelemTemporaryValueAssigner assigner;
1148 MempoolSelelemReserver reserver;
1151 assigner.assign(left, *sel);
1154 reserver.reserve(right, g->isize);
1157 else if (right && right->mempool)
1159 assigner.assign(right, *sel);
1161 _gmx_sel_evaluate_children(data, sel, g);
1163 n = (sel->flags & SEL_SINGLEVAL) ? 1 : g->isize;
1166 bool bArithNeg = (sel->u.arith.type == ARITH_NEG);
1167 GMX_ASSERT(right || bArithNeg,
1168 "Right operand cannot be null except for negations");
1169 for (i = i1 = i2 = 0; i < n; ++i)
1171 lval = left->v.u.r[i1];
1174 rval = right->v.u.r[i2];
1176 switch (sel->u.arith.type)
1178 case ARITH_PLUS: val = lval + rval; break;
1179 case ARITH_MINUS: val = lval - rval; break;
1180 case ARITH_NEG: val = -lval; break;
1181 case ARITH_MULT: val = lval * rval; break;
1182 case ARITH_DIV: val = lval / rval; break;
1183 case ARITH_EXP: val = pow(lval, rval); break;
1185 sel->v.u.r[i] = val;
1186 if (!(left->flags & SEL_SINGLEVAL))
1190 if (!bArithNeg && !(right->flags & SEL_SINGLEVAL))