Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / selection / poscalc.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013, 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 gmx::PositionCalculationCollection and functions in poscalc.h.
38  *
39  * \todo
40  * There is probably some room for optimization in the calculation of
41  * positions with bases.
42  * In particular, the current implementation may do a lot of unnecessary
43  * copying.
44  * The interface would need to be changed to make it possible to use the
45  * same output positions for several calculations.
46  *
47  * \todo
48  * The current algorithm for setting up base calculations could be improved
49  * in cases when there are calculations that cannot use a common base but
50  * still overlap partially (e.g., with three calculations A, B, and C
51  * such that A could use both B and C as a base, but B and C cannot use the
52  * same base).
53  * Setting up the bases in an optimal manner in every possible situation can
54  * be quite difficult unless several bases are allowed for one calculation,
55  * but better heuristics could probably be implemented.
56  * For best results, the setup should probably be postponed (at least
57  * partially) to gmx_ana_poscalc_init_eval().
58  *
59  * \author Teemu Murtola <teemu.murtola@gmail.com>
60  * \ingroup module_selection
61  */
62 #include <string.h>
63
64 #include "gromacs/legacyheaders/smalloc.h"
65 #include "gromacs/legacyheaders/typedefs.h"
66 #include "gromacs/legacyheaders/pbc.h"
67 #include "gromacs/legacyheaders/vec.h"
68
69 #include "gromacs/selection/centerofmass.h"
70 #include "gromacs/selection/indexutil.h"
71 #include "gromacs/selection/poscalc.h"
72 #include "gromacs/selection/position.h"
73 #include "gromacs/utility/exceptions.h"
74 #include "gromacs/utility/gmxassert.h"
75
76 namespace gmx
77 {
78
79 /*! \internal \brief
80  * Private implementation class for PositionCalculationCollection.
81  *
82  * \ingroup module_selection
83  */
84 class PositionCalculationCollection::Impl
85 {
86     public:
87         Impl();
88         ~Impl();
89
90         /*! \brief
91          * Inserts a position calculation structure into its collection.
92          *
93          * \param pc     Data structure to insert.
94          * \param before Data structure before which to insert
95          *   (NULL = insert at end).
96          *
97          * Inserts \p pc to its collection before \p before.
98          * If \p before is NULL, \p pc is appended to the list.
99          */
100         void insertCalculation(gmx_ana_poscalc_t *pc, gmx_ana_poscalc_t *before);
101         /*! \brief
102          * Removes a position calculation structure from its collection.
103          *
104          * \param pc    Data structure to remove.
105          *
106          * Removes \p pc from its collection.
107          */
108         void removeCalculation(gmx_ana_poscalc_t *pc);
109
110         /*! \copydoc PositionCalculationCollection::createCalculation()
111          *
112          * This method contains the actual implementation of the similarly
113          * named method in the parent class.
114          */
115         gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
116
117         /*! \brief
118          * Topology data.
119          *
120          * Can be NULL if none of the calculations require topology data or if
121          * setTopology() has not been called.
122          */
123         t_topology               *top_;
124         //! Pointer to the first data structure.
125         gmx_ana_poscalc_t        *first_;
126         //! Pointer to the last data structure.
127         gmx_ana_poscalc_t        *last_;
128         //! Whether the collection has been initialized for evaluation.
129         bool                      bInit_;
130 };
131
132 } // namespace gmx
133
134 /*! \internal \brief
135  * Data structure for position calculation.
136  */
137 struct gmx_ana_poscalc_t
138 {
139     /*! \brief
140      * Type of calculation.
141      *
142      * This field may differ from the type requested by the user, because
143      * it is changed internally to the most effective calculation.
144      * For example, if the user requests a COM calculation for residues
145      * consisting of single atoms, it is simply set to POS_ATOM.
146      * To provide a consistent interface to the user, the field \p itype
147      * should be used when information should be given out.
148      */
149     e_poscalc_t               type;
150     /*! \brief
151      * Flags for calculation options.
152      *
153      * See \ref poscalc_flags "documentation of the flags".
154      */
155     int                       flags;
156
157     /*! \brief
158      * Type for the created indices.
159      *
160      * This field always agrees with the type that the user requested, but
161      * may differ from \p type.
162      */
163     e_index_t                 itype;
164     /*! \brief
165      * Block data for the calculation.
166      */
167     t_blocka                  b;
168     /*! \brief
169      * Mapping from the blocks to the blocks of \p sbase.
170      *
171      * If \p sbase is NULL, this field is also.
172      */
173     int                      *baseid;
174     /*! \brief
175      * Maximum evaluation group.
176      */
177     gmx_ana_index_t           gmax;
178
179     /** Position storage for calculations that are used as a base. */
180     gmx_ana_pos_t            *p;
181
182     /** true if the positions have been evaluated for the current frame. */
183     bool                      bEval;
184     /*! \brief
185      * Base position data for this calculation.
186      *
187      * If not NULL, the centers required by this calculation have
188      * already been calculated in \p sbase.
189      * The structure pointed by \p sbase is always a static calculation.
190      */
191     struct gmx_ana_poscalc_t                 *sbase;
192     /** Next structure in the linked list of calculations. */
193     struct gmx_ana_poscalc_t                 *next;
194     /** Previous structure in the linked list of calculations. */
195     struct gmx_ana_poscalc_t                 *prev;
196     /** Number of references to this structure. */
197     int                                       refcount;
198     /** Collection this calculation belongs to. */
199     gmx::PositionCalculationCollection::Impl *coll;
200 };
201
202 const char * const gmx::PositionCalculationCollection::typeEnumValues[] = {
203     "atom",
204     "res_com",       "res_cog",
205     "mol_com",       "mol_cog",
206     "whole_res_com", "whole_res_cog",
207     "whole_mol_com", "whole_mol_cog",
208     "part_res_com",  "part_res_cog",
209     "part_mol_com",  "part_mol_cog",
210     "dyn_res_com",   "dyn_res_cog",
211     "dyn_mol_com",   "dyn_mol_cog",
212     NULL,
213 };
214
215 /*! \brief
216  * Returns the partition type for a given position type.
217  *
218  * \param [in] type  \c e_poscalc_t value to convert.
219  * \returns    Corresponding \c e_indet_t.
220  */
221 static e_index_t
222 index_type_for_poscalc(e_poscalc_t type)
223 {
224     switch (type)
225     {
226         case POS_ATOM:    return INDEX_ATOM;
227         case POS_RES:     return INDEX_RES;
228         case POS_MOL:     return INDEX_MOL;
229         case POS_ALL:     return INDEX_ALL;
230         case POS_ALL_PBC: return INDEX_ALL;
231     }
232     return INDEX_UNKNOWN;
233 }
234
235 namespace gmx
236 {
237
238 // static
239 void
240 PositionCalculationCollection::typeFromEnum(const char *post,
241                                             e_poscalc_t *type, int *flags)
242 {
243     if (post[0] == 'a')
244     {
245         *type   = POS_ATOM;
246         *flags &= ~(POS_MASS | POS_COMPLMAX | POS_COMPLWHOLE);
247         return;
248     }
249
250     /* Process the prefix */
251     const char *ptr = post;
252     if (post[0] == 'w')
253     {
254         *flags &= ~POS_COMPLMAX;
255         *flags |= POS_COMPLWHOLE;
256         ptr     = post + 6;
257     }
258     else if (post[0] == 'p')
259     {
260         *flags &= ~POS_COMPLWHOLE;
261         *flags |= POS_COMPLMAX;
262         ptr     = post + 5;
263     }
264     else if (post[0] == 'd')
265     {
266         *flags &= ~(POS_COMPLMAX | POS_COMPLWHOLE);
267         ptr     = post + 4;
268     }
269
270     if (ptr[0] == 'r')
271     {
272         *type = POS_RES;
273     }
274     else if (ptr[0] == 'm')
275     {
276         *type = POS_MOL;
277     }
278     else
279     {
280         GMX_THROW(InternalError("Unknown position calculation type"));
281     }
282     if (strlen(ptr) < 7)
283     {
284         GMX_THROW(InternalError("Unknown position calculation type"));
285     }
286     if (ptr[6] == 'm')
287     {
288         *flags |= POS_MASS;
289     }
290     else if (ptr[6] == 'g')
291     {
292         *flags &= ~POS_MASS;
293     }
294     else
295     {
296         GMX_THROW(InternalError("Unknown position calculation type"));
297     }
298 }
299
300 /********************************************************************
301  * PositionCalculationCollection::Impl
302  */
303
304 PositionCalculationCollection::Impl::Impl()
305     : top_(NULL), first_(NULL), last_(NULL), bInit_(false)
306 {
307 }
308
309 PositionCalculationCollection::Impl::~Impl()
310 {
311     // Loop backwards, because there can be internal references in that are
312     // correctly handled by this direction.
313     while (last_ != NULL)
314     {
315         GMX_ASSERT(last_->refcount == 1,
316                    "Dangling references to position calculations");
317         gmx_ana_poscalc_free(last_);
318     }
319 }
320
321 void
322 PositionCalculationCollection::Impl::insertCalculation(gmx_ana_poscalc_t *pc,
323                                                        gmx_ana_poscalc_t *before)
324 {
325     GMX_RELEASE_ASSERT(pc->coll == this, "Inconsistent collections");
326     if (before == NULL)
327     {
328         pc->next = NULL;
329         pc->prev = last_;
330         if (last_ != NULL)
331         {
332             last_->next = pc;
333         }
334         last_ = pc;
335     }
336     else
337     {
338         pc->prev     = before->prev;
339         pc->next     = before;
340         if (before->prev)
341         {
342             before->prev->next = pc;
343         }
344         before->prev = pc;
345     }
346     if (pc->prev == NULL)
347     {
348         first_ = pc;
349     }
350 }
351
352 void
353 PositionCalculationCollection::Impl::removeCalculation(gmx_ana_poscalc_t *pc)
354 {
355     GMX_RELEASE_ASSERT(pc->coll == this, "Inconsistent collections");
356     if (pc->prev != NULL)
357     {
358         pc->prev->next = pc->next;
359     }
360     else if (pc == first_)
361     {
362         first_ = pc->next;
363     }
364     if (pc->next != NULL)
365     {
366         pc->next->prev = pc->prev;
367     }
368     else if (pc == last_)
369     {
370         last_ = pc->prev;
371     }
372     pc->prev = pc->next = NULL;
373 }
374
375 gmx_ana_poscalc_t *
376 PositionCalculationCollection::Impl::createCalculation(e_poscalc_t type, int flags)
377 {
378     gmx_ana_poscalc_t *pc;
379
380     snew(pc, 1);
381     pc->type     = type;
382     pc->itype    = index_type_for_poscalc(type);
383     gmx_ana_poscalc_set_flags(pc, flags);
384     pc->refcount = 1;
385     pc->coll     = this;
386     insertCalculation(pc, NULL);
387     return pc;
388 }
389
390
391 /********************************************************************
392  * PositionCalculationCollection
393  */
394
395 PositionCalculationCollection::PositionCalculationCollection()
396     : impl_(new Impl)
397 {
398 }
399
400 PositionCalculationCollection::~PositionCalculationCollection()
401 {
402 }
403
404 void
405 PositionCalculationCollection::setTopology(t_topology *top)
406 {
407     impl_->top_ = top;
408 }
409
410 void
411 PositionCalculationCollection::printTree(FILE *fp) const
412 {
413     gmx_ana_poscalc_t *pc;
414     int                i, j;
415
416     fprintf(fp, "Position calculations:\n");
417     i  = 1;
418     pc = impl_->first_;
419     while (pc)
420     {
421         fprintf(fp, "%2d ", i);
422         switch (pc->type)
423         {
424             case POS_ATOM:    fprintf(fp, "ATOM");    break;
425             case POS_RES:     fprintf(fp, "RES");     break;
426             case POS_MOL:     fprintf(fp, "MOL");     break;
427             case POS_ALL:     fprintf(fp, "ALL");     break;
428             case POS_ALL_PBC: fprintf(fp, "ALL_PBC"); break;
429         }
430         if (pc->itype != index_type_for_poscalc(pc->type))
431         {
432             fprintf(fp, " (");
433             switch (pc->itype)
434             {
435                 case INDEX_UNKNOWN: fprintf(fp, "???");  break;
436                 case INDEX_ATOM:    fprintf(fp, "ATOM"); break;
437                 case INDEX_RES:     fprintf(fp, "RES");  break;
438                 case INDEX_MOL:     fprintf(fp, "MOL");  break;
439                 case INDEX_ALL:     fprintf(fp, "ALL");  break;
440             }
441             fprintf(fp, ")");
442         }
443         fprintf(fp, " flg=");
444         if (pc->flags & POS_MASS)
445         {
446             fprintf(fp, "M");
447         }
448         if (pc->flags & POS_DYNAMIC)
449         {
450             fprintf(fp, "D");
451         }
452         if (pc->flags & POS_MASKONLY)
453         {
454             fprintf(fp, "A");
455         }
456         if (pc->flags & POS_COMPLMAX)
457         {
458             fprintf(fp, "Cm");
459         }
460         if (pc->flags & POS_COMPLWHOLE)
461         {
462             fprintf(fp, "Cw");
463         }
464         if (!pc->flags)
465         {
466             fprintf(fp, "0");
467         }
468         fprintf(fp, " nr=%d nra=%d", pc->b.nr, pc->b.nra);
469         fprintf(fp, " refc=%d", pc->refcount);
470         fprintf(fp, "\n");
471         if (pc->gmax.nalloc_index > 0)
472         {
473             fprintf(fp, "   Group: ");
474             if (pc->gmax.isize > 20)
475             {
476                 fprintf(fp, " %d atoms", pc->gmax.isize);
477             }
478             else
479             {
480                 for (j = 0; j < pc->gmax.isize; ++j)
481                 {
482                     fprintf(fp, " %d", pc->gmax.index[j] + 1);
483                 }
484             }
485             fprintf(fp, "\n");
486         }
487         if (pc->b.nalloc_a > 0)
488         {
489             fprintf(fp, "   Atoms: ");
490             if (pc->b.nra > 20)
491             {
492                 fprintf(fp, " %d atoms", pc->b.nra);
493             }
494             else
495             {
496                 for (j = 0; j < pc->b.nra; ++j)
497                 {
498                     fprintf(fp, " %d", pc->b.a[j] + 1);
499                 }
500             }
501             fprintf(fp, "\n");
502         }
503         if (pc->b.nalloc_index > 0)
504         {
505             fprintf(fp, "   Blocks:");
506             if (pc->b.nr > 20)
507             {
508                 fprintf(fp, " %d pcs", pc->b.nr);
509             }
510             else
511             {
512                 for (j = 0; j <= pc->b.nr; ++j)
513                 {
514                     fprintf(fp, " %d", pc->b.index[j]);
515                 }
516             }
517             fprintf(fp, "\n");
518         }
519         if (pc->sbase)
520         {
521             gmx_ana_poscalc_t *base;
522
523             fprintf(fp, "   Base: ");
524             j    = 1;
525             base = impl_->first_;
526             while (base && base != pc->sbase)
527             {
528                 ++j;
529                 base = base->next;
530             }
531             fprintf(fp, "%d", j);
532             if (pc->baseid && pc->b.nr <= 20)
533             {
534                 fprintf(fp, " id:");
535                 for (j = 0; j < pc->b.nr; ++j)
536                 {
537                     fprintf(fp, " %d", pc->baseid[j]+1);
538                 }
539             }
540             fprintf(fp, "\n");
541         }
542         ++i;
543         pc = pc->next;
544     }
545 }
546
547 gmx_ana_poscalc_t *
548 PositionCalculationCollection::createCalculation(e_poscalc_t type, int flags)
549 {
550     return impl_->createCalculation(type, flags);
551 }
552
553 gmx_ana_poscalc_t *
554 PositionCalculationCollection::createCalculationFromEnum(const char *post, int flags)
555 {
556     e_poscalc_t  type;
557     int          cflags = flags;
558     typeFromEnum(post, &type, &cflags);
559     return impl_->createCalculation(type, cflags);
560 }
561
562 void PositionCalculationCollection::initEvaluation()
563 {
564     if (impl_->bInit_)
565     {
566         return;
567     }
568     gmx_ana_poscalc_t *pc = impl_->first_;
569     while (pc)
570     {
571         /* Initialize position storage for base calculations */
572         if (pc->p)
573         {
574             gmx_ana_poscalc_init_pos(pc, pc->p);
575         }
576         /* Construct the mapping of the base positions */
577         if (pc->sbase)
578         {
579             int                     bi, bj;
580
581             snew(pc->baseid, pc->b.nr);
582             for (bi = bj = 0; bi < pc->b.nr; ++bi, ++bj)
583             {
584                 while (pc->sbase->b.a[pc->sbase->b.index[bj]] != pc->b.a[pc->b.index[bi]])
585                 {
586                     ++bj;
587                 }
588                 pc->baseid[bi] = bj;
589             }
590         }
591         /* Free the block data for dynamic calculations */
592         if (pc->flags & POS_DYNAMIC)
593         {
594             if (pc->b.nalloc_index > 0)
595             {
596                 sfree(pc->b.index);
597                 pc->b.nalloc_index = 0;
598             }
599             if (pc->b.nalloc_a > 0)
600             {
601                 sfree(pc->b.a);
602                 pc->b.nalloc_a = 0;
603             }
604         }
605         pc = pc->next;
606     }
607     impl_->bInit_ = true;
608 }
609
610 void PositionCalculationCollection::initFrame()
611 {
612     if (!impl_->bInit_)
613     {
614         initEvaluation();
615     }
616     /* Clear the evaluation flags */
617     gmx_ana_poscalc_t *pc = impl_->first_;
618     while (pc)
619     {
620         pc->bEval = false;
621         pc        = pc->next;
622     }
623 }
624
625 } // namespace gmx
626
627 /*! \brief
628  * Initializes position calculation using the maximum possible input index.
629  *
630  * \param[in,out] pc  Position calculation data structure.
631  * \param[in]     g   Maximum index group for the calculation.
632  * \param[in]     bBase Whether \p pc will be used as a base or not.
633  *
634  * \p bBase affects on how the \p pc->gmax field is initialized.
635  */
636 static void
637 set_poscalc_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g, bool bBase)
638 {
639     t_topology *top = pc->coll->top_;
640     gmx_ana_index_make_block(&pc->b, top, g, pc->itype, pc->flags & POS_COMPLWHOLE);
641     /* Set the type to POS_ATOM if the calculation in fact is such. */
642     if (pc->b.nr == pc->b.nra)
643     {
644         pc->type   = POS_ATOM;
645         pc->flags &= ~(POS_MASS | POS_COMPLMAX | POS_COMPLWHOLE);
646     }
647     /* Set the POS_COMPLWHOLE flag if the calculation in fact always uses
648      * complete residues and molecules. */
649     if (!(pc->flags & POS_COMPLWHOLE)
650         && (!(pc->flags & POS_DYNAMIC) || (pc->flags & POS_COMPLMAX))
651         && (pc->type == POS_RES || pc->type == POS_MOL)
652         && gmx_ana_index_has_complete_elems(g, pc->itype, top))
653     {
654         pc->flags &= ~POS_COMPLMAX;
655         pc->flags |= POS_COMPLWHOLE;
656     }
657     /* Setup the gmax field */
658     if ((pc->flags & POS_COMPLWHOLE) && !bBase && pc->b.nra > g->isize)
659     {
660         gmx_ana_index_copy(&pc->gmax, g, true);
661     }
662     else
663     {
664         gmx_ana_index_set(&pc->gmax, pc->b.nra, pc->b.a, 0);
665     }
666 }
667
668 /*! \brief
669  * Checks whether a position calculation should use a base at all.
670  *
671  * \param[in] pc   Position calculation data to check.
672  * \returns   true if \p pc can use a base and gets some benefit out of it,
673  *   false otherwise.
674  */
675 static bool
676 can_use_base(gmx_ana_poscalc_t *pc)
677 {
678     /* For atoms, it should be faster to do a simple copy, so don't use a
679      * base. */
680     if (pc->type == POS_ATOM)
681     {
682         return false;
683     }
684     /* For dynamic selections that do not use completion, it is not possible
685      * to use a base. */
686     if ((pc->type == POS_RES || pc->type == POS_MOL)
687         && (pc->flags & POS_DYNAMIC) && !(pc->flags & (POS_COMPLMAX | POS_COMPLWHOLE)))
688     {
689         return false;
690     }
691     /* Dynamic calculations for a single position cannot use a base. */
692     if ((pc->type == POS_ALL || pc->type == POS_ALL_PBC)
693         && (pc->flags & POS_DYNAMIC))
694     {
695         return false;
696     }
697     return true;
698 }
699
700 /*! \brief
701  * Checks whether two position calculations should use a common base.
702  *
703  * \param[in]     pc1 Calculation 1 to check for.
704  * \param[in]     pc2 Calculation 2 to check for.
705  * \param[in]     g1  Index group structure that contains the atoms from
706  *   \p pc1.
707  * \param[in,out] g   Working space, should have enough allocated memory to
708  *   contain the intersection of the atoms in \p pc1 and \p pc2.
709  * \returns   true if the two calculations should be merged to use a common
710  *   base, false otherwise.
711  */
712 static bool
713 should_merge(gmx_ana_poscalc_t *pc1, gmx_ana_poscalc_t *pc2,
714              gmx_ana_index_t *g1, gmx_ana_index_t *g)
715 {
716     gmx_ana_index_t  g2;
717
718     /* Do not merge calculations with different mass weighting. */
719     if ((pc1->flags & POS_MASS) != (pc2->flags & POS_MASS))
720     {
721         return false;
722     }
723     /* Avoid messing up complete calculations. */
724     if ((pc1->flags & POS_COMPLWHOLE) != (pc2->flags & POS_COMPLWHOLE))
725     {
726         return false;
727     }
728     /* Find the overlap between the calculations. */
729     gmx_ana_index_set(&g2, pc2->b.nra, pc2->b.a, 0);
730     gmx_ana_index_intersection(g, g1, &g2);
731     /* Do not merge if there is no overlap. */
732     if (g->isize == 0)
733     {
734         return false;
735     }
736     /* Full completion calculations always match if the type is correct. */
737     if ((pc1->flags & POS_COMPLWHOLE) && (pc2->flags & POS_COMPLWHOLE)
738         && pc1->type == pc2->type)
739     {
740         return true;
741     }
742     /* The calculations also match if the intersection consists of full
743      * blocks. */
744     if (gmx_ana_index_has_full_ablocks(g, &pc1->b)
745         && gmx_ana_index_has_full_ablocks(g, &pc2->b))
746     {
747         return true;
748     }
749     return false;
750 }
751
752 /*! \brief
753  * Creates a static base for position calculation.
754  *
755  * \param     pc  Data structure to copy.
756  * \returns   Pointer to a newly allocated base for \p pc.
757  *
758  * Creates and returns a deep copy of \p pc, but clears the
759  * \ref POS_DYNAMIC and \ref POS_MASKONLY flags.
760  * The newly created structure is set as the base (\c gmx_ana_poscalc_t::sbase)
761  * of \p pc and inserted in the collection before \p pc.
762  */
763 static gmx_ana_poscalc_t *
764 create_simple_base(gmx_ana_poscalc_t *pc)
765 {
766     gmx_ana_poscalc_t *base;
767     int                flags;
768
769     flags = pc->flags & ~(POS_DYNAMIC | POS_MASKONLY);
770     base  = pc->coll->createCalculation(pc->type, flags);
771     set_poscalc_maxindex(base, &pc->gmax, true);
772
773     base->p = new gmx_ana_pos_t();
774
775     pc->sbase = base;
776     pc->coll->removeCalculation(base);
777     pc->coll->insertCalculation(base, pc);
778
779     return base;
780 }
781
782 /*! \brief
783  * Merges a calculation into another calculation such that the new calculation
784  * can be used as a base.
785  *
786  * \param[in,out] base Base calculation to merge to.
787  * \param[in,out] pc   Position calculation to merge to \p base.
788  *
789  * After the call, \p base can be used as a base for \p pc (or any calculation
790  * that used it as a base).
791  * It is assumed that any overlap between \p base and \p pc is in complete
792  * blocks, i.e., that the merge is possible.
793  */
794 static void
795 merge_to_base(gmx_ana_poscalc_t *base, gmx_ana_poscalc_t *pc)
796 {
797     gmx_ana_index_t  gp, gb, g;
798     int              isize, bnr;
799     int              i, bi, bj, bo;
800
801     base->flags |= pc->flags & (POS_VELOCITIES | POS_FORCES);
802     gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, 0);
803     gmx_ana_index_set(&gb, base->b.nra, base->b.a, 0);
804     isize = gmx_ana_index_difference_size(&gp, &gb);
805     if (isize > 0)
806     {
807         gmx_ana_index_clear(&g);
808         gmx_ana_index_reserve(&g, base->b.nra + isize);
809         /* Find the new blocks */
810         gmx_ana_index_difference(&g, &gp, &gb);
811         /* Count the blocks in g */
812         i = bi = bnr = 0;
813         while (i < g.isize)
814         {
815             while (pc->b.a[pc->b.index[bi]] != g.index[i])
816             {
817                 ++bi;
818             }
819             i += pc->b.index[bi+1] - pc->b.index[bi];
820             ++bnr;
821             ++bi;
822         }
823         /* Merge the atoms into a temporary structure */
824         gmx_ana_index_merge(&g, &gb, &g);
825         /* Merge the blocks */
826         srenew(base->b.index, base->b.nr + bnr + 1);
827         i                   = g.isize - 1;
828         bi                  = base->b.nr - 1;
829         bj                  = pc->b.nr - 1;
830         bo                  = base->b.nr + bnr - 1;
831         base->b.index[bo+1] = i + 1;
832         while (bo >= 0)
833         {
834             if (bi < 0 || base->b.a[base->b.index[bi+1]-1] != g.index[i])
835             {
836                 i -= pc->b.index[bj+1] - pc->b.index[bj];
837                 --bj;
838             }
839             else
840             {
841                 if (bj >= 0 && pc->b.a[pc->b.index[bj+1]-1] == g.index[i])
842                 {
843                     --bj;
844                 }
845                 i -= base->b.index[bi+1] - base->b.index[bi];
846                 --bi;
847             }
848             base->b.index[bo] = i + 1;
849             --bo;
850         }
851         base->b.nr           += bnr;
852         base->b.nalloc_index += bnr;
853         sfree(base->b.a);
854         base->b.nra      = g.isize;
855         base->b.a        = g.index;
856         base->b.nalloc_a = g.isize;
857         /* Refresh the gmax field */
858         gmx_ana_index_set(&base->gmax, base->b.nra, base->b.a, 0);
859     }
860 }
861
862 /*! \brief
863  * Merges two bases into one.
864  *
865  * \param[in,out] tbase Base calculation to merge to.
866  * \param[in]     mbase Base calculation to merge to \p tbase.
867  *
868  * After the call, \p mbase has been freed and \p tbase is used as the base
869  * for all calculations that previously had \p mbase as their base.
870  * It is assumed that any overlap between \p tbase and \p mbase is in complete
871  * blocks, i.e., that the merge is possible.
872  */
873 static void
874 merge_bases(gmx_ana_poscalc_t *tbase, gmx_ana_poscalc_t *mbase)
875 {
876     gmx_ana_poscalc_t *pc;
877
878     merge_to_base(tbase, mbase);
879     mbase->coll->removeCalculation(mbase);
880     /* Set tbase as the base for all calculations that had mbase */
881     pc = tbase->coll->first_;
882     while (pc)
883     {
884         if (pc->sbase == mbase)
885         {
886             pc->sbase = tbase;
887             tbase->refcount++;
888         }
889         pc = pc->next;
890     }
891     /* Free mbase */
892     mbase->refcount = 0;
893     gmx_ana_poscalc_free(mbase);
894 }
895
896 /*! \brief
897  * Setups the static base calculation for a position calculation.
898  *
899  * \param[in,out] pc   Position calculation to setup the base for.
900  */
901 static void
902 setup_base(gmx_ana_poscalc_t *pc)
903 {
904     gmx_ana_poscalc_t *base, *pbase, *next;
905     gmx_ana_index_t    gp, g;
906
907     /* Exit immediately if pc should not have a base. */
908     if (!can_use_base(pc))
909     {
910         return;
911     }
912
913     gmx_ana_index_set(&gp, pc->b.nra, pc->b.a, 0);
914     gmx_ana_index_clear(&g);
915     gmx_ana_index_reserve(&g, pc->b.nra);
916     pbase = pc;
917     base  = pc->coll->first_;
918     while (base)
919     {
920         /* Save the next calculation so that we can safely delete base */
921         next = base->next;
922         /* Skip pc, calculations that already have a base (we should match the
923          * base instead), as well as calculations that should not have a base.
924          * If the above conditions are met, check whether we should do a
925          * merge.
926          */
927         if (base != pc && !base->sbase && can_use_base(base)
928             && should_merge(pbase, base, &gp, &g))
929         {
930             /* Check whether this is the first base found */
931             if (pbase == pc)
932             {
933                 /* Create a real base if one is not present */
934                 if (!base->p)
935                 {
936                     pbase = create_simple_base(base);
937                 }
938                 else
939                 {
940                     pbase = base;
941                 }
942                 /* Make it a base for pc as well */
943                 merge_to_base(pbase, pc);
944                 pc->sbase = pbase;
945                 pbase->refcount++;
946             }
947             else /* This was not the first base */
948             {
949                 if (!base->p)
950                 {
951                     /* If it is not a real base, just make the new base as
952                      * the base for it as well. */
953                     merge_to_base(pbase, base);
954                     base->sbase = pbase;
955                     pbase->refcount++;
956                 }
957                 else
958                 {
959                     /* If base is a real base, merge it with the new base
960                      * and delete it. */
961                     merge_bases(pbase, base);
962                 }
963             }
964             gmx_ana_index_set(&gp, pbase->b.nra, pbase->b.a, 0);
965             gmx_ana_index_reserve(&g, pc->b.nra);
966         }
967         /* Proceed to the next unchecked calculation */
968         base = next;
969     }
970
971     gmx_ana_index_deinit(&g);
972
973     /* If no base was found, create one if one is required */
974     if (!pc->sbase && (pc->flags & POS_DYNAMIC)
975         && (pc->flags & (POS_COMPLMAX | POS_COMPLWHOLE)))
976     {
977         create_simple_base(pc);
978     }
979 }
980
981 /*!
982  * \param[in,out] pc    Position calculation data structure.
983  * \param[in]     flags New flags.
984  *
985  * \p flags are added to the old flags.
986  * If calculation type is \ref POS_ATOM, \ref POS_MASS is automatically
987  * cleared.
988  * If both \ref POS_DYNAMIC and \ref POS_MASKONLY are provided,
989  * \ref POS_DYNAMIC is cleared.
990  * If calculation type is not \ref POS_RES or \ref POS_MOL,
991  * \ref POS_COMPLMAX and \ref POS_COMPLWHOLE are automatically cleared.
992  */
993 void
994 gmx_ana_poscalc_set_flags(gmx_ana_poscalc_t *pc, int flags)
995 {
996     if (pc->type == POS_ATOM)
997     {
998         flags &= ~POS_MASS;
999     }
1000     if (flags & POS_MASKONLY)
1001     {
1002         flags &= ~POS_DYNAMIC;
1003     }
1004     if (pc->type != POS_RES && pc->type != POS_MOL)
1005     {
1006         flags &= ~(POS_COMPLMAX | POS_COMPLWHOLE);
1007     }
1008     pc->flags |= flags;
1009 }
1010
1011 /*!
1012  * \param[in,out] pc  Position calculation data structure.
1013  * \param[in]     g   Maximum index group for the calculation.
1014  *
1015  * Subsequent calls to gmx_ana_poscalc_update() should use only subsets of
1016  * \p g for evaluation.
1017  *
1018  * The topology should have been set for the collection of which \p pc is
1019  * a member.
1020  */
1021 void
1022 gmx_ana_poscalc_set_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g)
1023 {
1024     set_poscalc_maxindex(pc, g, false);
1025     setup_base(pc);
1026 }
1027
1028 /*!
1029  * \param[in]  pc  Position calculation data structure.
1030  * \param[out] p   Output positions.
1031  *
1032  * Calls to gmx_ana_poscalc_update() using \p pc should use only positions
1033  * initialized with this function.
1034  * The \c p->g pointer is initialized to point to an internal group that
1035  * contains the maximum index group set with gmx_ana_poscalc_set_maxindex().
1036  */
1037 void
1038 gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p)
1039 {
1040     gmx_ana_indexmap_init(&p->m, &pc->gmax, pc->coll->top_, pc->itype);
1041     /* Only do the static optimization when there is no completion */
1042     if (!(pc->flags & POS_DYNAMIC) && pc->b.nra == pc->gmax.isize)
1043     {
1044         gmx_ana_indexmap_set_static(&p->m, &pc->b);
1045     }
1046     gmx_ana_pos_reserve(p, p->m.mapb.nr, 0);
1047     if (pc->flags & POS_VELOCITIES)
1048     {
1049         gmx_ana_pos_reserve_velocities(p);
1050     }
1051     if (pc->flags & POS_FORCES)
1052     {
1053         gmx_ana_pos_reserve_forces(p);
1054     }
1055 }
1056
1057 /*!
1058  * \param  pc  Position calculation data to be freed.
1059  *
1060  * The \p pc pointer is invalid after the call.
1061  */
1062 void
1063 gmx_ana_poscalc_free(gmx_ana_poscalc_t *pc)
1064 {
1065     if (!pc)
1066     {
1067         return;
1068     }
1069
1070     pc->refcount--;
1071     if (pc->refcount > 0)
1072     {
1073         return;
1074     }
1075
1076     pc->coll->removeCalculation(pc);
1077     if (pc->b.nalloc_index > 0)
1078     {
1079         sfree(pc->b.index);
1080     }
1081     if (pc->b.nalloc_a > 0)
1082     {
1083         sfree(pc->b.a);
1084     }
1085     if (pc->flags & POS_COMPLWHOLE)
1086     {
1087         gmx_ana_index_deinit(&pc->gmax);
1088     }
1089     delete pc->p;
1090     if (pc->sbase)
1091     {
1092         gmx_ana_poscalc_free(pc->sbase);
1093         sfree(pc->baseid);
1094     }
1095     sfree(pc);
1096 }
1097
1098 /*!
1099  * \param[in] pc  Position calculation data to query.
1100  * \returns   true if \p pc requires topology for initialization and/or
1101  *   evaluation, false otherwise.
1102  */
1103 bool
1104 gmx_ana_poscalc_requires_top(gmx_ana_poscalc_t *pc)
1105 {
1106     if ((pc->flags & POS_MASS) || pc->type == POS_RES || pc->type == POS_MOL
1107         || ((pc->flags & POS_FORCES) && pc->type != POS_ATOM))
1108     {
1109         return true;
1110     }
1111     return false;
1112 }
1113
1114 /*!
1115  * \param[in]     pc   Position calculation data.
1116  * \param[in,out] p    Output positions, initialized previously with
1117  *   gmx_ana_poscalc_init_pos() using \p pc.
1118  * \param[in]     g    Index group to use for the update.
1119  * \param[in]     fr   Current frame.
1120  * \param[in]     pbc  PBC data, or NULL if no PBC should be used.
1121  *
1122  * gmx_ana_poscalc_init_frame() should be called for each frame before calling
1123  * this function.
1124  */
1125 void
1126 gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p,
1127                        gmx_ana_index_t *g, t_trxframe *fr, t_pbc *pbc)
1128 {
1129     int  i, bi, bj;
1130
1131     if (pc->bEval == true && !(pc->flags & POS_MASKONLY))
1132     {
1133         return;
1134     }
1135     if (pc->sbase)
1136     {
1137         gmx_ana_poscalc_update(pc->sbase, NULL, NULL, fr, pbc);
1138     }
1139     if (!p)
1140     {
1141         p = pc->p;
1142     }
1143     if (!g)
1144     {
1145         g = &pc->gmax;
1146     }
1147
1148     /* Update the index map */
1149     if (pc->flags & POS_DYNAMIC)
1150     {
1151         gmx_ana_indexmap_update(&p->m, g, false);
1152     }
1153     else if (pc->flags & POS_MASKONLY)
1154     {
1155         gmx_ana_indexmap_update(&p->m, g, true);
1156         if (pc->bEval)
1157         {
1158             return;
1159         }
1160     }
1161     if (!(pc->flags & POS_DYNAMIC))
1162     {
1163         pc->bEval = true;
1164     }
1165
1166     /* Evaluate the positions */
1167     if (pc->sbase)
1168     {
1169         /* TODO: It might be faster to evaluate the positions within this
1170          * loop instead of in the beginning. */
1171         if (pc->flags & POS_DYNAMIC)
1172         {
1173             for (bi = 0; bi < p->count(); ++bi)
1174             {
1175                 bj = pc->baseid[p->m.refid[bi]];
1176                 copy_rvec(pc->sbase->p->x[bj], p->x[bi]);
1177             }
1178             if (p->v)
1179             {
1180                 for (bi = 0; bi < p->count(); ++bi)
1181                 {
1182                     bj = pc->baseid[p->m.refid[bi]];
1183                     copy_rvec(pc->sbase->p->v[bj], p->v[bi]);
1184                 }
1185             }
1186             if (p->f)
1187             {
1188                 for (bi = 0; bi < p->count(); ++bi)
1189                 {
1190                     bj = pc->baseid[p->m.refid[bi]];
1191                     copy_rvec(pc->sbase->p->f[bj], p->f[bi]);
1192                 }
1193             }
1194         }
1195         else
1196         {
1197             for (bi = 0; bi < p->count(); ++bi)
1198             {
1199                 bj = pc->baseid[bi];
1200                 copy_rvec(pc->sbase->p->x[bj], p->x[bi]);
1201             }
1202             if (p->v)
1203             {
1204                 for (bi = 0; bi < p->count(); ++bi)
1205                 {
1206                     bj = pc->baseid[bi];
1207                     copy_rvec(pc->sbase->p->v[bj], p->v[bi]);
1208                 }
1209             }
1210             if (p->f)
1211             {
1212                 for (bi = 0; bi < p->count(); ++bi)
1213                 {
1214                     bj = pc->baseid[bi];
1215                     copy_rvec(pc->sbase->p->f[bj], p->f[bi]);
1216                 }
1217             }
1218         }
1219     }
1220     else /* pc->sbase is NULL */
1221     {
1222         if (pc->flags & POS_DYNAMIC)
1223         {
1224             pc->b.nr    = p->m.mapb.nr;
1225             pc->b.index = p->m.mapb.index;
1226             pc->b.nra   = g->isize;
1227             pc->b.a     = g->index;
1228         }
1229         if (p->v && !fr->bV)
1230         {
1231             for (i = 0; i < pc->b.nra; ++i)
1232             {
1233                 clear_rvec(p->v[i]);
1234             }
1235         }
1236         if (p->f && !fr->bF)
1237         {
1238             for (i = 0; i < pc->b.nra; ++i)
1239             {
1240                 clear_rvec(p->f[i]);
1241             }
1242         }
1243         /* Here, we assume that the topology has been properly initialized,
1244          * and do not check the return values of gmx_calc_comg*(). */
1245         t_topology *top   = pc->coll->top_;
1246         bool        bMass = pc->flags & POS_MASS;
1247         switch (pc->type)
1248         {
1249             case POS_ATOM:
1250                 for (i = 0; i < pc->b.nra; ++i)
1251                 {
1252                     copy_rvec(fr->x[pc->b.a[i]], p->x[i]);
1253                 }
1254                 if (p->v && fr->bV)
1255                 {
1256                     for (i = 0; i < pc->b.nra; ++i)
1257                     {
1258                         copy_rvec(fr->v[pc->b.a[i]], p->v[i]);
1259                     }
1260                 }
1261                 if (p->f && fr->bF)
1262                 {
1263                     for (i = 0; i < pc->b.nra; ++i)
1264                     {
1265                         copy_rvec(fr->f[pc->b.a[i]], p->f[i]);
1266                     }
1267                 }
1268                 break;
1269             case POS_ALL:
1270                 gmx_calc_comg(top, fr->x, pc->b.nra, pc->b.a, bMass, p->x[0]);
1271                 if (p->v && fr->bV)
1272                 {
1273                     gmx_calc_comg(top, fr->v, pc->b.nra, pc->b.a, bMass, p->v[0]);
1274                 }
1275                 if (p->f && fr->bF)
1276                 {
1277                     gmx_calc_comg_f(top, fr->f, pc->b.nra, pc->b.a, bMass, p->f[0]);
1278                 }
1279                 break;
1280             case POS_ALL_PBC:
1281                 gmx_calc_comg_pbc(top, fr->x, pbc, pc->b.nra, pc->b.a, bMass, p->x[0]);
1282                 if (p->v && fr->bV)
1283                 {
1284                     gmx_calc_comg(top, fr->v, pc->b.nra, pc->b.a, bMass, p->v[0]);
1285                 }
1286                 if (p->f && fr->bF)
1287                 {
1288                     gmx_calc_comg_f(top, fr->f, pc->b.nra, pc->b.a, bMass, p->f[0]);
1289                 }
1290                 break;
1291             default:
1292                 gmx_calc_comg_blocka(top, fr->x, &pc->b, bMass, p->x);
1293                 if (p->v && fr->bV)
1294                 {
1295                     gmx_calc_comg_blocka(top, fr->v, &pc->b, bMass, p->v);
1296                 }
1297                 if (p->f && fr->bF)
1298                 {
1299                     gmx_calc_comg_f_blocka(top, fr->f, &pc->b, bMass, p->f);
1300                 }
1301                 break;
1302         }
1303     }
1304 }