3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
33 * Implements gmx::PositionCalculationCollection and functions in poscalc.h.
36 * There is probably some room for optimization in the calculation of
37 * positions with bases.
38 * In particular, the current implementation may do a lot of unnecessary
40 * The interface would need to be changed to make it possible to use the
41 * same output positions for several calculations.
44 * The current algorithm for setting up base calculations could be improved
45 * in cases when there are calculations that cannot use a common base but
46 * still overlap partially (e.g., with three calculations A, B, and C
47 * such that A could use both B and C as a base, but B and C cannot use the
49 * Setting up the bases in an optimal manner in every possible situation can
50 * be quite difficult unless several bases are allowed for one calculation,
51 * but better heuristics could probably be implemented.
52 * For best results, the setup should probably be postponed (at least
53 * partially) to gmx_ana_poscalc_init_eval().
55 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
56 * \ingroup module_selection
60 #include "gromacs/legacyheaders/smalloc.h"
61 #include "gromacs/legacyheaders/typedefs.h"
62 #include "gromacs/legacyheaders/pbc.h"
63 #include "gromacs/legacyheaders/vec.h"
65 #include "gromacs/selection/centerofmass.h"
66 #include "gromacs/selection/indexutil.h"
67 #include "gromacs/selection/poscalc.h"
68 #include "gromacs/selection/position.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/gmxassert.h"
76 * Private implementation class for PositionCalculationCollection.
78 * \ingroup module_selection
80 class PositionCalculationCollection::Impl
87 * Inserts a position calculation structure into its collection.
89 * \param pc Data structure to insert.
90 * \param before Data structure before which to insert
91 * (NULL = insert at end).
93 * Inserts \p pc to its collection before \p before.
94 * If \p before is NULL, \p pc is appended to the list.
96 void insertCalculation(gmx_ana_poscalc_t *pc, gmx_ana_poscalc_t *before);
98 * Removes a position calculation structure from its collection.
100 * \param pc Data structure to remove.
102 * Removes \p pc from its collection.
104 void removeCalculation(gmx_ana_poscalc_t *pc);
106 /*! \copydoc PositionCalculationCollection::createCalculation()
108 * This method contains the actual implementation of the similarly
109 * named method in the parent class.
111 gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
116 * Can be NULL if none of the calculations require topology data or if
117 * setTopology() has not been called.
120 //! Pointer to the first data structure.
121 gmx_ana_poscalc_t *first_;
122 //! Pointer to the last data structure.
123 gmx_ana_poscalc_t *last_;
124 //! Whether the collection has been initialized for evaluation.
131 * Data structure for position calculation.
133 struct gmx_ana_poscalc_t
136 * Type of calculation.
138 * This field may differ from the type requested by the user, because
139 * it is changed internally to the most effective calculation.
140 * For example, if the user requests a COM calculation for residues
141 * consisting of single atoms, it is simply set to POS_ATOM.
142 * To provide a consistent interface to the user, the field \p itype
143 * should be used when information should be given out.
147 * Flags for calculation options.
149 * See \ref poscalc_flags "documentation of the flags".
154 * Type for the created indices.
156 * This field always agrees with the type that the user requested, but
157 * may differ from \p type.
161 * Block data for the calculation.
165 * Mapping from the blocks to the blocks of \p sbase.
167 * If \p sbase is NULL, this field is also.
171 * Maximum evaluation group.
173 gmx_ana_index_t gmax;
175 /** Position storage for calculations that are used as a base. */
178 /** true if the positions have been evaluated for the current frame. */
181 * Base position data for this calculation.
183 * If not NULL, the centers required by this calculation have
184 * already been calculated in \p sbase.
185 * The structure pointed by \p sbase is always a static calculation.
187 struct gmx_ana_poscalc_t *sbase;
188 /** Next structure in the linked list of calculations. */
189 struct gmx_ana_poscalc_t *next;
190 /** Previous structure in the linked list of calculations. */
191 struct gmx_ana_poscalc_t *prev;
192 /** Number of references to this structure. */
194 /** Collection this calculation belongs to. */
195 gmx::PositionCalculationCollection::Impl *coll;
198 const char * const gmx::PositionCalculationCollection::typeEnumValues[] = {
200 "res_com", "res_cog",
201 "mol_com", "mol_cog",
202 "whole_res_com", "whole_res_cog",
203 "whole_mol_com", "whole_mol_cog",
204 "part_res_com", "part_res_cog",
205 "part_mol_com", "part_mol_cog",
206 "dyn_res_com", "dyn_res_cog",
207 "dyn_mol_com", "dyn_mol_cog",
212 * Returns the partition type for a given position type.
214 * \param [in] type \c e_poscalc_t value to convert.
215 * \returns Corresponding \c e_indet_t.
218 index_type_for_poscalc(e_poscalc_t type)
222 case POS_ATOM: return INDEX_ATOM;
223 case POS_RES: return INDEX_RES;
224 case POS_MOL: return INDEX_MOL;
225 case POS_ALL: return INDEX_ALL;
226 case POS_ALL_PBC: return INDEX_ALL;
228 return INDEX_UNKNOWN;
236 PositionCalculationCollection::typeFromEnum(const char *post,
237 e_poscalc_t *type, int *flags)
242 *flags &= ~(POS_MASS | POS_COMPLMAX | POS_COMPLWHOLE);
246 /* Process the prefix */
247 const char *ptr = post;
250 *flags &= ~POS_COMPLMAX;
251 *flags |= POS_COMPLWHOLE;
254 else if (post[0] == 'p')
256 *flags &= ~POS_COMPLWHOLE;
257 *flags |= POS_COMPLMAX;
260 else if (post[0] == 'd')
262 *flags &= ~(POS_COMPLMAX | POS_COMPLWHOLE);
270 else if (ptr[0] == 'm')
276 GMX_THROW(InternalError("Unknown position calculation type"));
280 GMX_THROW(InternalError("Unknown position calculation type"));
286 else if (ptr[6] == 'g')
292 GMX_THROW(InternalError("Unknown position calculation type"));
296 /********************************************************************
297 * PositionCalculationCollection::Impl
300 PositionCalculationCollection::Impl::Impl()
301 : top_(NULL), first_(NULL), last_(NULL), bInit_(false)
305 PositionCalculationCollection::Impl::~Impl()
307 // Loop backwards, because there can be internal references in that are
308 // correctly handled by this direction.
309 while (last_ != NULL)
311 GMX_ASSERT(last_->refcount == 1,
312 "Dangling references to position calculations");
313 gmx_ana_poscalc_free(last_);
318 PositionCalculationCollection::Impl::insertCalculation(gmx_ana_poscalc_t *pc,
319 gmx_ana_poscalc_t *before)
321 GMX_RELEASE_ASSERT(pc->coll == this, "Inconsistent collections");
334 pc->prev = before->prev;
338 before->prev->next = pc;
342 if (pc->prev == NULL)
349 PositionCalculationCollection::Impl::removeCalculation(gmx_ana_poscalc_t *pc)
351 GMX_RELEASE_ASSERT(pc->coll == this, "Inconsistent collections");
352 if (pc->prev != NULL)
354 pc->prev->next = pc->next;
356 else if (pc == first_)
360 if (pc->next != NULL)
362 pc->next->prev = pc->prev;
364 else if (pc == last_)
368 pc->prev = pc->next = NULL;
372 PositionCalculationCollection::Impl::createCalculation(e_poscalc_t type, int flags)
374 gmx_ana_poscalc_t *pc;
378 pc->itype = index_type_for_poscalc(type);
379 gmx_ana_poscalc_set_flags(pc, flags);
382 insertCalculation(pc, NULL);
387 /********************************************************************
388 * PositionCalculationCollection
391 PositionCalculationCollection::PositionCalculationCollection()
396 PositionCalculationCollection::~PositionCalculationCollection()
401 PositionCalculationCollection::setTopology(t_topology *top)
407 PositionCalculationCollection::printTree(FILE *fp) const
409 gmx_ana_poscalc_t *pc;
412 fprintf(fp, "Position calculations:\n");
417 fprintf(fp, "%2d ", i);
420 case POS_ATOM: fprintf(fp, "ATOM"); break;
421 case POS_RES: fprintf(fp, "RES"); break;
422 case POS_MOL: fprintf(fp, "MOL"); break;
423 case POS_ALL: fprintf(fp, "ALL"); break;
424 case POS_ALL_PBC: fprintf(fp, "ALL_PBC"); break;
426 if (pc->itype != index_type_for_poscalc(pc->type))
431 case INDEX_UNKNOWN: fprintf(fp, "???"); break;
432 case INDEX_ATOM: fprintf(fp, "ATOM"); break;
433 case INDEX_RES: fprintf(fp, "RES"); break;
434 case INDEX_MOL: fprintf(fp, "MOL"); break;
435 case INDEX_ALL: fprintf(fp, "ALL"); break;
439 fprintf(fp, " flg=");
440 if (pc->flags & POS_MASS)
444 if (pc->flags & POS_DYNAMIC)
448 if (pc->flags & POS_MASKONLY)
452 if (pc->flags & POS_COMPLMAX)
456 if (pc->flags & POS_COMPLWHOLE)
464 fprintf(fp, " nr=%d nra=%d", pc->b.nr, pc->b.nra);
465 fprintf(fp, " refc=%d", pc->refcount);
467 if (pc->gmax.nalloc_index > 0)
469 fprintf(fp, " Group: ");
470 if (pc->gmax.isize > 20)
472 fprintf(fp, " %d atoms", pc->gmax.isize);
476 for (j = 0; j < pc->gmax.isize; ++j)
478 fprintf(fp, " %d", pc->gmax.index[j] + 1);
483 if (pc->b.nalloc_a > 0)
485 fprintf(fp, " Atoms: ");
488 fprintf(fp, " %d atoms", pc->b.nra);
492 for (j = 0; j < pc->b.nra; ++j)
494 fprintf(fp, " %d", pc->b.a[j] + 1);
499 if (pc->b.nalloc_index > 0)
501 fprintf(fp, " Blocks:");
504 fprintf(fp, " %d pcs", pc->b.nr);
508 for (j = 0; j <= pc->b.nr; ++j)
510 fprintf(fp, " %d", pc->b.index[j]);
517 gmx_ana_poscalc_t *base;
519 fprintf(fp, " Base: ");
521 base = impl_->first_;
522 while (base && base != pc->sbase)
527 fprintf(fp, "%d", j);
528 if (pc->baseid && pc->b.nr <= 20)
531 for (j = 0; j < pc->b.nr; ++j)
533 fprintf(fp, " %d", pc->baseid[j]+1);
544 PositionCalculationCollection::createCalculation(e_poscalc_t type, int flags)
546 return impl_->createCalculation(type, flags);
550 PositionCalculationCollection::createCalculationFromEnum(const char *post, int flags)
554 typeFromEnum(post, &type, &cflags);
555 return impl_->createCalculation(type, cflags);
558 void PositionCalculationCollection::initEvaluation()
564 gmx_ana_poscalc_t *pc = impl_->first_;
567 /* Initialize position storage for base calculations */
570 gmx_ana_poscalc_init_pos(pc, pc->p);
572 /* Construct the mapping of the base positions */
577 snew(pc->baseid, pc->b.nr);
578 for (bi = bj = 0; bi < pc->b.nr; ++bi, ++bj)
580 while (pc->sbase->b.a[pc->sbase->b.index[bj]] != pc->b.a[pc->b.index[bi]])
587 /* Free the block data for dynamic calculations */
588 if (pc->flags & POS_DYNAMIC)
590 if (pc->b.nalloc_index > 0)
593 pc->b.nalloc_index = 0;
595 if (pc->b.nalloc_a > 0)
603 impl_->bInit_ = true;
606 void PositionCalculationCollection::initFrame()
612 /* Clear the evaluation flags */
613 gmx_ana_poscalc_t *pc = impl_->first_;
624 * Initializes position calculation using the maximum possible input index.
626 * \param[in,out] pc Position calculation data structure.
627 * \param[in] g Maximum index group for the calculation.
628 * \param[in] bBase Whether \p pc will be used as a base or not.
630 * \p bBase affects on how the \p pc->gmax field is initialized.
633 set_poscalc_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g, bool bBase)
635 t_topology *top = pc->coll->top_;
636 gmx_ana_index_make_block(&pc->b, top, g, pc->itype, pc->flags & POS_COMPLWHOLE);
637 /* Set the type to POS_ATOM if the calculation in fact is such. */
638 if (pc->b.nr == pc->b.nra)
641 pc->flags &= ~(POS_MASS | POS_COMPLMAX | POS_COMPLWHOLE);
643 /* Set the POS_COMPLWHOLE flag if the calculation in fact always uses
644 * complete residues and molecules. */
645 if (!(pc->flags & POS_COMPLWHOLE)
646 && (!(pc->flags & POS_DYNAMIC) || (pc->flags & POS_COMPLMAX))
647 && (pc->type == POS_RES || pc->type == POS_MOL)
648 && gmx_ana_index_has_complete_elems(g, pc->itype, top))
650 pc->flags &= ~POS_COMPLMAX;
651 pc->flags |= POS_COMPLWHOLE;
653 /* Setup the gmax field */
654 if ((pc->flags & POS_COMPLWHOLE) && !bBase && pc->b.nra > g->isize)
656 gmx_ana_index_copy(&pc->gmax, g, true);
657 sfree(pc->gmax.name);
658 pc->gmax.name = NULL;
662 gmx_ana_index_set(&pc->gmax, pc->b.nra, pc->b.a, NULL, 0);
667 * Checks whether a position calculation should use a base at all.
669 * \param[in] pc Position calculation data to check.
670 * \returns true if \p pc can use a base and gets some benefit out of it,
674 can_use_base(gmx_ana_poscalc_t *pc)
676 /* For atoms, it should be faster to do a simple copy, so don't use a
678 if (pc->type == POS_ATOM)
682 /* For dynamic selections that do not use completion, it is not possible
684 if ((pc->type == POS_RES || pc->type == POS_MOL)
685 && (pc->flags & POS_DYNAMIC) && !(pc->flags & (POS_COMPLMAX | POS_COMPLWHOLE)))
689 /* Dynamic calculations for a single position cannot use a base. */
690 if ((pc->type == POS_ALL || pc->type == POS_ALL_PBC)
691 && (pc->flags & POS_DYNAMIC))
699 * Checks whether two position calculations should use a common base.
701 * \param[in] pc1 Calculation 1 to check for.
702 * \param[in] pc2 Calculation 2 to check for.
703 * \param[in] g1 Index group structure that contains the atoms from
705 * \param[in,out] g Working space, should have enough allocated memory to
706 * contain the intersection of the atoms in \p pc1 and \p pc2.
707 * \returns true if the two calculations should be merged to use a common
708 * base, false otherwise.
711 should_merge(gmx_ana_poscalc_t *pc1, gmx_ana_poscalc_t *pc2,
712 gmx_ana_index_t *g1, gmx_ana_index_t *g)
716 /* Do not merge calculations with different mass weighting. */
717 if ((pc1->flags & POS_MASS) != (pc2->flags & POS_MASS))
721 /* Avoid messing up complete calculations. */
722 if ((pc1->flags & POS_COMPLWHOLE) != (pc2->flags & POS_COMPLWHOLE))
726 /* Find the overlap between the calculations. */
727 gmx_ana_index_set(&g2, pc2->b.nra, pc2->b.a, NULL, 0);
728 gmx_ana_index_intersection(g, g1, &g2);
729 /* Do not merge if there is no overlap. */
734 /* Full completion calculations always match if the type is correct. */
735 if ((pc1->flags & POS_COMPLWHOLE) && (pc2->flags & POS_COMPLWHOLE)
736 && pc1->type == pc2->type)
740 /* The calculations also match if the intersection consists of full
742 if (gmx_ana_index_has_full_ablocks(g, &pc1->b)
743 && gmx_ana_index_has_full_ablocks(g, &pc2->b))
751 * Creates a static base for position calculation.
753 * \param pc Data structure to copy.
754 * \returns Pointer to a newly allocated base for \p pc.
756 * Creates and returns a deep copy of \p pc, but clears the
757 * \ref POS_DYNAMIC and \ref POS_MASKONLY flags.
758 * The newly created structure is set as the base (\c gmx_ana_poscalc_t::sbase)
759 * of \p pc and inserted in the collection before \p pc.
761 static gmx_ana_poscalc_t *
762 create_simple_base(gmx_ana_poscalc_t *pc)
764 gmx_ana_poscalc_t *base;
767 flags = pc->flags & ~(POS_DYNAMIC | POS_MASKONLY);
768 base = pc->coll->createCalculation(pc->type, flags);
769 set_poscalc_maxindex(base, &pc->gmax, true);
774 pc->coll->removeCalculation(base);
775 pc->coll->insertCalculation(base, pc);
781 * Merges a calculation into another calculation such that the new calculation
782 * can be used as a base.
784 * \param[in,out] base Base calculation to merge to.
785 * \param[in,out] pc Position calculation to merge to \p base.
787 * After the call, \p base can be used as a base for \p pc (or any calculation
788 * that used it as a base).
789 * It is assumed that any overlap between \p base and \p pc is in complete
790 * blocks, i.e., that the merge is possible.
793 merge_to_base(gmx_ana_poscalc_t *base, gmx_ana_poscalc_t *pc)
795 gmx_ana_index_t gp, gb, g;
799 base->flags |= pc->flags & (POS_VELOCITIES | POS_FORCES);
800 gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, NULL, 0);
801 gmx_ana_index_set(&gb, base->b.nra, base->b.a, NULL, 0);
802 isize = gmx_ana_index_difference_size(&gp, &gb);
805 gmx_ana_index_clear(&g);
806 gmx_ana_index_reserve(&g, base->b.nra + isize);
807 /* Find the new blocks */
808 gmx_ana_index_difference(&g, &gp, &gb);
809 /* Count the blocks in g */
813 while (pc->b.a[pc->b.index[bi]] != g.index[i])
817 i += pc->b.index[bi+1] - pc->b.index[bi];
821 /* Merge the atoms into a temporary structure */
822 gmx_ana_index_merge(&g, &gb, &g);
823 /* Merge the blocks */
824 srenew(base->b.index, base->b.nr + bnr + 1);
828 bo = base->b.nr + bnr - 1;
829 base->b.index[bo+1] = i + 1;
832 if (bi < 0 || base->b.a[base->b.index[bi+1]-1] != g.index[i])
834 i -= pc->b.index[bj+1] - pc->b.index[bj];
839 if (bj >= 0 && pc->b.a[pc->b.index[bj+1]-1] == g.index[i])
843 i -= base->b.index[bi+1] - base->b.index[bi];
846 base->b.index[bo] = i + 1;
850 base->b.nalloc_index += bnr;
852 base->b.nra = g.isize;
854 base->b.nalloc_a = g.isize;
855 /* Refresh the gmax field */
856 gmx_ana_index_set(&base->gmax, base->b.nra, base->b.a, NULL, 0);
861 * Merges two bases into one.
863 * \param[in,out] tbase Base calculation to merge to.
864 * \param[in] mbase Base calculation to merge to \p tbase.
866 * After the call, \p mbase has been freed and \p tbase is used as the base
867 * for all calculations that previously had \p mbase as their base.
868 * It is assumed that any overlap between \p tbase and \p mbase is in complete
869 * blocks, i.e., that the merge is possible.
872 merge_bases(gmx_ana_poscalc_t *tbase, gmx_ana_poscalc_t *mbase)
874 gmx_ana_poscalc_t *pc;
876 merge_to_base(tbase, mbase);
877 mbase->coll->removeCalculation(mbase);
878 /* Set tbase as the base for all calculations that had mbase */
879 pc = tbase->coll->first_;
882 if (pc->sbase == mbase)
891 gmx_ana_poscalc_free(mbase);
895 * Setups the static base calculation for a position calculation.
897 * \param[in,out] pc Position calculation to setup the base for.
900 setup_base(gmx_ana_poscalc_t *pc)
902 gmx_ana_poscalc_t *base, *pbase, *next;
903 gmx_ana_index_t gp, g;
905 /* Exit immediately if pc should not have a base. */
906 if (!can_use_base(pc))
911 gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, NULL, 0);
912 gmx_ana_index_clear(&g);
913 gmx_ana_index_reserve(&g, pc->b.nra);
915 base = pc->coll->first_;
918 /* Save the next calculation so that we can safely delete base */
920 /* Skip pc, calculations that already have a base (we should match the
921 * base instead), as well as calculations that should not have a base.
922 * If the above conditions are met, check whether we should do a
925 if (base != pc && !base->sbase && can_use_base(base)
926 && should_merge(pbase, base, &gp, &g))
928 /* Check whether this is the first base found */
931 /* Create a real base if one is not present */
934 pbase = create_simple_base(base);
940 /* Make it a base for pc as well */
941 merge_to_base(pbase, pc);
945 else /* This was not the first base */
949 /* If it is not a real base, just make the new base as
950 * the base for it as well. */
951 merge_to_base(pbase, base);
957 /* If base is a real base, merge it with the new base
959 merge_bases(pbase, base);
962 gmx_ana_index_set(&gp, pbase->b.nra, pbase->b.a, NULL, 0);
963 gmx_ana_index_reserve(&g, pc->b.nra);
965 /* Proceed to the next unchecked calculation */
969 gmx_ana_index_deinit(&g);
971 /* If no base was found, create one if one is required */
972 if (!pc->sbase && (pc->flags & POS_DYNAMIC)
973 && (pc->flags & (POS_COMPLMAX | POS_COMPLWHOLE)))
975 create_simple_base(pc);
980 * \param[in,out] pc Position calculation data structure.
981 * \param[in] flags New flags.
983 * \p flags are added to the old flags.
984 * If calculation type is \ref POS_ATOM, \ref POS_MASS is automatically
986 * If both \ref POS_DYNAMIC and \ref POS_MASKONLY are provided,
987 * \ref POS_DYNAMIC is cleared.
988 * If calculation type is not \ref POS_RES or \ref POS_MOL,
989 * \ref POS_COMPLMAX and \ref POS_COMPLWHOLE are automatically cleared.
992 gmx_ana_poscalc_set_flags(gmx_ana_poscalc_t *pc, int flags)
994 if (pc->type == POS_ATOM)
998 if (flags & POS_MASKONLY)
1000 flags &= ~POS_DYNAMIC;
1002 if (pc->type != POS_RES && pc->type != POS_MOL)
1004 flags &= ~(POS_COMPLMAX | POS_COMPLWHOLE);
1010 * \param[in,out] pc Position calculation data structure.
1011 * \param[in] g Maximum index group for the calculation.
1013 * Subsequent calls to gmx_ana_poscalc_update() should use only subsets of
1014 * \p g for evaluation.
1016 * The topology should have been set for the collection of which \p pc is
1020 gmx_ana_poscalc_set_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g)
1022 set_poscalc_maxindex(pc, g, false);
1027 * \param[in] pc Position calculation data structure.
1028 * \param[out] p Output positions.
1030 * Calls to gmx_ana_poscalc_update() using \p pc should use only positions
1031 * initialized with this function.
1032 * The \c p->g pointer is initialized to point to an internal group that
1033 * contains the maximum index group set with gmx_ana_poscalc_set_maxindex().
1036 gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p)
1038 gmx_ana_indexmap_init(&p->m, &pc->gmax, pc->coll->top_, pc->itype);
1039 /* Only do the static optimization when there is no completion */
1040 if (!(pc->flags & POS_DYNAMIC) && pc->b.nra == pc->gmax.isize)
1042 gmx_ana_indexmap_set_static(&p->m, &pc->b);
1044 gmx_ana_pos_reserve(p, p->m.nr, 0);
1045 if (pc->flags & POS_VELOCITIES)
1047 gmx_ana_pos_reserve_velocities(p);
1049 if (pc->flags & POS_FORCES)
1051 gmx_ana_pos_reserve_forces(p);
1053 gmx_ana_pos_set_nr(p, p->m.nr);
1054 gmx_ana_pos_set_evalgrp(p, &pc->gmax);
1058 * \param pc Position calculation data to be freed.
1060 * The \p pc pointer is invalid after the call.
1063 gmx_ana_poscalc_free(gmx_ana_poscalc_t *pc)
1071 if (pc->refcount > 0)
1076 pc->coll->removeCalculation(pc);
1077 if (pc->b.nalloc_index > 0)
1081 if (pc->b.nalloc_a > 0)
1085 if (pc->flags & POS_COMPLWHOLE)
1087 gmx_ana_index_deinit(&pc->gmax);
1091 gmx_ana_pos_free(pc->p);
1095 gmx_ana_poscalc_free(pc->sbase);
1102 * \param[in] pc Position calculation data to query.
1103 * \returns true if \p pc requires topology for initialization and/or
1104 * evaluation, false otherwise.
1107 gmx_ana_poscalc_requires_top(gmx_ana_poscalc_t *pc)
1109 if ((pc->flags & POS_MASS) || pc->type == POS_RES || pc->type == POS_MOL)
1117 * \param[in] pc Position calculation data.
1118 * \param[in,out] p Output positions, initialized previously with
1119 * gmx_ana_poscalc_init_pos() using \p pc.
1120 * \param[in] g Index group to use for the update.
1121 * \param[in] fr Current frame.
1122 * \param[in] pbc PBC data, or NULL if no PBC should be used.
1124 * gmx_ana_poscalc_init_frame() should be called for each frame before calling
1128 gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p,
1129 gmx_ana_index_t *g, t_trxframe *fr, t_pbc *pbc)
1133 if (pc->bEval == true && !(pc->flags & POS_MASKONLY))
1139 gmx_ana_poscalc_update(pc->sbase, NULL, NULL, fr, pbc);
1149 gmx_ana_pos_set_evalgrp(p, g);
1151 /* Update the index map */
1152 if (pc->flags & POS_DYNAMIC)
1154 gmx_ana_indexmap_update(&p->m, g, false);
1157 else if (pc->flags & POS_MASKONLY)
1159 gmx_ana_indexmap_update(&p->m, g, true);
1163 if (!(pc->flags & POS_DYNAMIC))
1168 /* Evaluate the positions */
1171 /* TODO: It might be faster to evaluate the positions within this
1172 * loop instead of in the beginning. */
1173 if (pc->flags & POS_DYNAMIC)
1175 for (bi = 0; bi < p->nr; ++bi)
1177 bj = pc->baseid[p->m.refid[bi]];
1178 copy_rvec(pc->sbase->p->x[bj], p->x[bi]);
1182 for (bi = 0; bi < p->nr; ++bi)
1184 bj = pc->baseid[p->m.refid[bi]];
1185 copy_rvec(pc->sbase->p->v[bj], p->v[bi]);
1190 for (bi = 0; bi < p->nr; ++bi)
1192 bj = pc->baseid[p->m.refid[bi]];
1193 copy_rvec(pc->sbase->p->f[bj], p->f[bi]);
1199 for (bi = 0; bi < p->nr; ++bi)
1201 bj = pc->baseid[bi];
1202 copy_rvec(pc->sbase->p->x[bj], p->x[bi]);
1206 for (bi = 0; bi < p->nr; ++bi)
1208 bj = pc->baseid[bi];
1209 copy_rvec(pc->sbase->p->v[bj], p->v[bi]);
1214 for (bi = 0; bi < p->nr; ++bi)
1216 bj = pc->baseid[bi];
1217 copy_rvec(pc->sbase->p->f[bj], p->f[bi]);
1222 else /* pc->sbase is NULL */
1224 if (pc->flags & POS_DYNAMIC)
1226 pc->b.nr = p->m.mapb.nr;
1227 pc->b.index = p->m.mapb.index;
1228 pc->b.nra = g->isize;
1231 if (p->v && !fr->bV)
1233 for (i = 0; i < pc->b.nra; ++i)
1235 clear_rvec(p->v[i]);
1238 if (p->f && !fr->bF)
1240 for (i = 0; i < pc->b.nra; ++i)
1242 clear_rvec(p->f[i]);
1245 /* Here, we assume that the topology has been properly initialized,
1246 * and do not check the return values of gmx_calc_comg*(). */
1247 t_topology *top = pc->coll->top_;
1248 bool bMass = pc->flags & POS_MASS;
1252 for (i = 0; i < pc->b.nra; ++i)
1254 copy_rvec(fr->x[pc->b.a[i]], p->x[i]);
1258 for (i = 0; i < pc->b.nra; ++i)
1260 copy_rvec(fr->v[pc->b.a[i]], p->v[i]);
1265 for (i = 0; i < pc->b.nra; ++i)
1267 copy_rvec(fr->f[pc->b.a[i]], p->f[i]);
1272 gmx_calc_comg(top, fr->x, pc->b.nra, pc->b.a, bMass, p->x[0]);
1275 gmx_calc_comg(top, fr->v, pc->b.nra, pc->b.a, bMass, p->v[0]);
1279 gmx_calc_comg_f(top, fr->f, pc->b.nra, pc->b.a, bMass, p->f[0]);
1283 gmx_calc_comg_pbc(top, fr->x, pbc, pc->b.nra, pc->b.a, bMass, p->x[0]);
1286 gmx_calc_comg(top, fr->v, pc->b.nra, pc->b.a, bMass, p->v[0]);
1290 gmx_calc_comg_f(top, fr->f, pc->b.nra, pc->b.a, bMass, p->f[0]);
1294 gmx_calc_comg_blocka(top, fr->x, &pc->b, bMass, p->x);
1297 gmx_calc_comg_blocka(top, fr->v, &pc->b, bMass, p->v);
1301 gmx_calc_comg_blocka(top, fr->f, &pc->b, bMass, p->f);