Fixes for clang-tidy-9
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / modules / angle.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2011-2018, The GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief
38  * Implements gmx::analysismodules::Angle.
39  *
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  * \ingroup module_trajectoryanalysis
42  */
43 #include "gmxpre.h"
44
45 #include "angle.h"
46
47 #include <algorithm>
48 #include <memory>
49 #include <string>
50 #include <vector>
51
52 #include "gromacs/analysisdata/analysisdata.h"
53 #include "gromacs/analysisdata/modules/average.h"
54 #include "gromacs/analysisdata/modules/histogram.h"
55 #include "gromacs/analysisdata/modules/plot.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/options/basicoptions.h"
59 #include "gromacs/options/filenameoption.h"
60 #include "gromacs/options/ioptionscontainer.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/selection/selection.h"
63 #include "gromacs/selection/selectionoption.h"
64 #include "gromacs/trajectory/trajectoryframe.h"
65 #include "gromacs/trajectoryanalysis/analysissettings.h"
66 #include "gromacs/utility/arrayref.h"
67 #include "gromacs/utility/enumerationhelpers.h"
68 #include "gromacs/utility/exceptions.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/stringutil.h"
71
72 namespace gmx
73 {
74
75 namespace analysismodules
76 {
77
78 namespace
79 {
80
81 /********************************************************************
82  * Helper classes
83  */
84
85 /*! \brief
86  * Helper to encapsulate logic for looping over input selections.
87  *
88  * This class provides two-dimensional iteration:
89  *  - Over _angle groups_, corresponding to an input selection.  If the input
90  *    selection list contains a single selection, that selection gets used
91  *    for all angle groups.
92  *  - Within a group, over _values_, each consisting of a fixed number of
93  *    selection positions.  If there is only a single value within a selection,
94  *    that value is returned over and over again.
95  * This transparently provides the semantics of using a single selection/vector
96  * to compute angles against multiple selections/vectors as described in the
97  * tool help text.
98  *
99  * This class isn't perferctly self-contained and requires the caller to know
100  * some of the internals to use it properly, but it serves its purpose for this
101  * single analysis tool by simplifying the loops.
102  * Some methods have also been tailored to allow the caller to use it a bit
103  * more easily.
104  *
105  * \ingroup module_trajectoryanalysis
106  */
107 class AnglePositionIterator
108 {
109 public:
110     /*! \brief
111      * Creates an iterator to loop over input selection positions.
112      *
113      * \param[in] selections       List of selections.
114      * \param[in] posCountPerValue Number of selection positions that
115      *     constitute a single value for the iteration.
116      *
117      * If \p selections is empty, and/or \p posCountPerValue is zero, the
118      * iterator can still be advanced and hasValue()/hasSingleValue()
119      * called, but values cannot be accessed.
120      */
121     AnglePositionIterator(const SelectionList& selections, int posCountPerValue) :
122         selections_(selections),
123         posCountPerValue_(posCountPerValue),
124         currentSelection_(0),
125         nextPosition_(0)
126     {
127     }
128
129     //! Advances the iterator to the next group of angles.
130     void nextGroup()
131     {
132         if (selections_.size() > 1)
133         {
134             ++currentSelection_;
135         }
136         nextPosition_ = 0;
137     }
138     //! Advances the iterator to the next angle in the current group.
139     void nextValue()
140     {
141         if (!hasSingleValue())
142         {
143             nextPosition_ += posCountPerValue_;
144         }
145     }
146
147     /*! \brief
148      * Returns whether this iterator represents any values.
149      *
150      * If the return value is `false`, only nextGroup(), nextValue() and
151      * hasSingleValue() are allowed to be called.
152      */
153     bool hasValue() const { return !selections_.empty(); }
154     /*! \brief
155      * Returns whether the current selection only contains a single value.
156      *
157      * Returns `false` if hasValue() returns false, which allows cutting
158      * some corners in consistency checks.
159      */
160     bool hasSingleValue() const
161     {
162         return hasValue() && currentSelection().posCount() == posCountPerValue_;
163     }
164     //! Returns whether the current selection is dynamic.
165     bool isDynamic() const { return currentSelection().isDynamic(); }
166     /*! \brief
167      * Returns whether positions in the current value are either all
168      * selected or all unselected.
169      */
170     bool allValuesConsistentlySelected() const
171     {
172         if (posCountPerValue_ <= 1)
173         {
174             return true;
175         }
176         const bool bSelected = currentPosition(0).selected();
177         for (int i = 1; i < posCountPerValue_; ++i)
178         {
179             if (currentPosition(i).selected() != bSelected)
180             {
181                 return false;
182             }
183         }
184         return true;
185     }
186     /*! \brief
187      * Returns whether positions in the current value are selected.
188      *
189      * Only works reliably if allValuesConsistentlySelected() returns
190      * `true`.
191      */
192     bool currentValuesSelected() const
193     {
194         return selections_.empty() || currentPosition(0).selected();
195     }
196
197     //! Returns the currently active selection.
198     const Selection& currentSelection() const
199     {
200         GMX_ASSERT(currentSelection_ < ssize(selections_), "Accessing an invalid selection");
201         return selections_[currentSelection_];
202     }
203     //! Returns the `i`th position for the current value.
204     SelectionPosition currentPosition(int i) const
205     {
206         return currentSelection().position(nextPosition_ + i);
207     }
208     /*! \brief
209      * Extracts all coordinates corresponding to the current value.
210      *
211      * \param[out] x  Array to which the positions are extracted.
212      *
213      * \p x should contain at minimum the number of positions per value
214      * passed to the constructor.
215      */
216     void getCurrentPositions(rvec x[]) const
217     {
218         GMX_ASSERT(posCountPerValue_ > 0, "Accessing positions for an invalid angle type");
219         GMX_ASSERT(nextPosition_ + posCountPerValue_ <= currentSelection().posCount(),
220                    "Accessing an invalid position");
221         for (int i = 0; i < posCountPerValue_; ++i)
222         {
223             copy_rvec(currentPosition(i).x(), x[i]);
224         }
225     }
226
227 private:
228     const SelectionList& selections_;
229     const int            posCountPerValue_;
230     int                  currentSelection_;
231     int                  nextPosition_;
232
233     GMX_DISALLOW_COPY_AND_ASSIGN(AnglePositionIterator);
234 };
235
236 /********************************************************************
237  * Actual analysis module
238  */
239
240 //! How to interpret the selections in -group1.
241 enum class Group1Type : int
242 {
243     Angle,
244     Dihedral,
245     Vector,
246     Plane,
247     Count
248 };
249 //! How to interpret the selections in -group2.
250 enum class Group2Type : int
251 {
252     None,
253     Vector,
254     Plane,
255     TimeZero,
256     Z,
257     SphereNormal,
258     Count
259 };
260 //! String values corresponding to Group1Type.
261 const EnumerationArray<Group1Type, const char*> c_group1TypeEnumNames = { { "angle", "dihedral",
262                                                                             "vector", "plane" } };
263 //! String values corresponding to Group2Type.
264 const EnumerationArray<Group2Type, const char*> c_group2TypeEnumNames = {
265     { "none", "vector", "plane", "t0", "z", "sphnorm" }
266 };
267
268 class Angle : public TrajectoryAnalysisModule
269 {
270 public:
271     Angle();
272
273     void initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings) override;
274     void optionsFinished(TrajectoryAnalysisSettings* settings) override;
275     void initAnalysis(const TrajectoryAnalysisSettings& settings, const TopologyInformation& top) override;
276
277     void analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* pdata) override;
278
279     void finishAnalysis(int nframes) override;
280     void writeOutput() override;
281
282 private:
283     void initFromSelections(const SelectionList& sel1, const SelectionList& sel2);
284     void checkSelections(const SelectionList& sel1, const SelectionList& sel2) const;
285
286     SelectionList        sel1_;
287     SelectionList        sel2_;
288     SelectionOptionInfo* sel1info_;
289     SelectionOptionInfo* sel2info_;
290     std::string          fnAverage_;
291     std::string          fnAll_;
292     std::string          fnHistogram_;
293
294     Group1Type g1type_;
295     Group2Type g2type_;
296     double     binWidth_;
297
298     AnalysisData                             angles_;
299     AnalysisDataFrameAverageModulePointer    averageModule_;
300     AnalysisDataSimpleHistogramModulePointer histogramModule_;
301
302     std::vector<int>               angleCount_;
303     int                            natoms1_;
304     int                            natoms2_;
305     std::vector<std::vector<RVec>> vt0_;
306
307     // Copy and assign disallowed by base.
308 };
309
310 Angle::Angle() :
311     sel1info_(nullptr),
312     sel2info_(nullptr),
313     g1type_(Group1Type::Angle),
314     g2type_(Group2Type::None),
315     binWidth_(1.0),
316     natoms1_(0),
317     natoms2_(0)
318 {
319     averageModule_ = std::make_unique<AnalysisDataFrameAverageModule>();
320     angles_.addModule(averageModule_);
321     histogramModule_ = std::make_unique<AnalysisDataSimpleHistogramModule>();
322     angles_.addModule(histogramModule_);
323
324     registerAnalysisDataset(&angles_, "angle");
325     registerBasicDataset(averageModule_.get(), "average");
326     registerBasicDataset(&histogramModule_->averager(), "histogram");
327 }
328
329
330 void Angle::initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings)
331 {
332     static const char* const desc[] = {
333         "[THISMODULE] computes different types of angles between vectors.",
334         "It supports both vectors defined by two positions and normals of",
335         "planes defined by three positions.",
336         "The z axis or the local normal of a sphere can also be used as",
337         "one of the vectors.",
338         "There are also convenience options 'angle' and 'dihedral' for",
339         "calculating bond angles and dihedrals defined by three/four",
340         "positions.[PAR]",
341         "The type of the angle is specified with [TT]-g1[tt] and [TT]-g2[tt].",
342         "If [TT]-g1[tt] is [TT]angle[tt] or [TT]dihedral[tt], [TT]-g2[tt]",
343         "should not be specified.",
344         "In this case, [TT]-group1[tt] should specify one or more selections,",
345         "and each should contain triplets or quartets of positions that define",
346         "the angles to be calculated.[PAR]",
347         "If [TT]-g1[tt] is [TT]vector[tt] or [TT]plane[tt], [TT]-group1[tt]",
348         "should specify selections that contain either pairs ([TT]vector[tt])",
349         "or triplets ([TT]plane[tt]) of positions. For vectors, the positions",
350         "set the endpoints of the vector, and for planes, the three positions",
351         "are used to calculate the normal of the plane. In both cases,",
352         "[TT]-g2[tt] specifies the other vector to use (see below).[PAR]",
353         "With [TT]-g2 vector[tt] or [TT]-g2 plane[tt], [TT]-group2[tt] should",
354         "specify another set of vectors. [TT]-group1[tt] and [TT]-group2[tt]",
355         "should specify the same number of selections. It is also allowed to",
356         "only have a single selection for one of the options, in which case",
357         "the same selection is used with each selection in the other group.",
358         "Similarly, for each selection in [TT]-group1[tt], the corresponding",
359         "selection in [TT]-group2[tt] should specify the same number of",
360         "vectors or a single vector. In the latter case, the angle is",
361         "calculated between that single vector and each vector from the other",
362         "selection.[PAR]",
363         "With [TT]-g2 sphnorm[tt], each selection in [TT]-group2[tt] should",
364         "specify a single position that is the center of the sphere.",
365         "The second vector is calculated as the vector from the center to the",
366         "midpoint of the positions specified by [TT]-group1[tt].[PAR]",
367         "With [TT]-g2 z[tt], [TT]-group2[tt] is not necessary, and angles",
368         "between the first vectors and the positive Z axis are calculated.[PAR]",
369         "With [TT]-g2 t0[tt], [TT]-group2[tt] is not necessary, and angles",
370         "are calculated from the vectors as they are in the first frame.[PAR]",
371         "There are three options for output:",
372         "[TT]-oav[tt] writes an xvg file with the time and the average angle",
373         "for each frame.",
374         "[TT]-oall[tt] writes all the individual angles.",
375         "[TT]-oh[tt] writes a histogram of the angles. The bin width can be",
376         "set with [TT]-binw[tt].",
377         "For [TT]-oav[tt] and [TT]-oh[tt], separate average/histogram is",
378         "computed for each selection in [TT]-group1[tt]."
379     };
380
381     settings->setHelpText(desc);
382
383     options->addOption(FileNameOption("oav")
384                                .filetype(eftPlot)
385                                .outputFile()
386                                .store(&fnAverage_)
387                                .defaultBasename("angaver")
388                                .description("Average angles as a function of time"));
389     options->addOption(FileNameOption("oall")
390                                .filetype(eftPlot)
391                                .outputFile()
392                                .store(&fnAll_)
393                                .defaultBasename("angles")
394                                .description("All angles as a function of time"));
395     options->addOption(FileNameOption("oh")
396                                .filetype(eftPlot)
397                                .outputFile()
398                                .store(&fnHistogram_)
399                                .defaultBasename("anghist")
400                                .description("Histogram of the angles"));
401
402     options->addOption(EnumOption<Group1Type>("g1")
403                                .enumValue(c_group1TypeEnumNames)
404                                .store(&g1type_)
405                                .description("Type of analysis/first vector group"));
406     options->addOption(
407             EnumOption<Group2Type>("g2").enumValue(c_group2TypeEnumNames).store(&g2type_).description("Type of second vector group"));
408     options->addOption(
409             DoubleOption("binw").store(&binWidth_).description("Binwidth for -oh in degrees"));
410
411     sel1info_ = options->addOption(
412             SelectionOption("group1").required().dynamicMask().storeVector(&sel1_).multiValue().description(
413                     "First analysis/vector selection"));
414     sel2info_ = options->addOption(
415             SelectionOption("group2").dynamicMask().storeVector(&sel2_).multiValue().description(
416                     "Second analysis/vector selection"));
417 }
418
419
420 void Angle::optionsFinished(TrajectoryAnalysisSettings* /* settings */)
421 {
422     const bool bSingle = (g1type_ == Group1Type::Angle || g1type_ == Group1Type::Dihedral);
423
424     if (bSingle && g2type_ != Group2Type::None)
425     {
426         GMX_THROW(
427                 InconsistentInputError("Cannot use a second group (-g2) with "
428                                        "-g1 angle or dihedral"));
429     }
430     if (bSingle && sel2info_->isSet())
431     {
432         GMX_THROW(
433                 InconsistentInputError("Cannot provide a second selection "
434                                        "(-group2) with -g1 angle or dihedral"));
435     }
436     if (!bSingle && g2type_ == Group2Type::None)
437     {
438         GMX_THROW(
439                 InconsistentInputError("Should specify a second group (-g2) "
440                                        "if the first group is not an angle or a dihedral"));
441     }
442
443     // Set up the number of positions per angle.
444     switch (g1type_)
445     {
446         case Group1Type::Angle: natoms1_ = 3; break;
447         case Group1Type::Dihedral: natoms1_ = 4; break;
448         case Group1Type::Vector: natoms1_ = 2; break;
449         case Group1Type::Plane: natoms1_ = 3; break;
450         default: GMX_THROW(InternalError("invalid -g1 value"));
451     }
452     switch (g2type_)
453     {
454         case Group2Type::None: natoms2_ = 0; break;
455         case Group2Type::Vector: natoms2_ = 2; break;
456         case Group2Type::Plane: natoms2_ = 3; break;
457         case Group2Type::TimeZero: // Intended to fall through
458         case Group2Type::Z: natoms2_ = 0; break;
459         case Group2Type::SphereNormal: natoms2_ = 1; break;
460         default: GMX_THROW(InternalError("invalid -g2 value"));
461     }
462     if (natoms2_ == 0 && sel2info_->isSet())
463     {
464         GMX_THROW(InconsistentInputError(
465                 "Cannot provide a second selection (-group2) with -g2 t0 or z"));
466     }
467     // TODO: If bSingle is not set, the second selection option should be
468     // required.
469 }
470
471
472 void Angle::initFromSelections(const SelectionList& sel1, const SelectionList& sel2)
473 {
474     const int  angleGroups         = std::max(sel1.size(), sel2.size());
475     const bool bHasSecondSelection = natoms2_ > 0;
476
477     if (bHasSecondSelection && sel1.size() != sel2.size() && std::min(sel1.size(), sel2.size()) != 1)
478     {
479         GMX_THROW(InconsistentInputError(
480                 "-group1 and -group2 should specify the same number of selections"));
481     }
482
483     AnglePositionIterator iter1(sel1, natoms1_);
484     AnglePositionIterator iter2(sel2, natoms2_);
485     for (int g = 0; g < angleGroups; ++g, iter1.nextGroup(), iter2.nextGroup())
486     {
487         const int posCount1 = iter1.currentSelection().posCount();
488         if (natoms1_ > 1 && posCount1 % natoms1_ != 0)
489         {
490             GMX_THROW(InconsistentInputError(formatString(
491                     "Number of positions in selection %d in the first group not divisible by %d",
492                     static_cast<int>(g + 1), natoms1_)));
493         }
494         const int angleCount1 = posCount1 / natoms1_;
495         int       angleCount  = angleCount1;
496
497         if (bHasSecondSelection)
498         {
499             const int posCount2 = iter2.currentSelection().posCount();
500             if (natoms2_ > 1 && posCount2 % natoms2_ != 0)
501             {
502                 GMX_THROW(InconsistentInputError(
503                         formatString("Number of positions in selection %d in the second group not "
504                                      "divisible by %d",
505                                      static_cast<int>(g + 1), natoms2_)));
506             }
507             if (g2type_ == Group2Type::SphereNormal && posCount2 != 1)
508             {
509                 GMX_THROW(InconsistentInputError(
510                         "The second group should contain a single position with -g2 sphnorm"));
511             }
512
513             const int angleCount2 = posCount2 / natoms2_;
514             angleCount            = std::max(angleCount1, angleCount2);
515             if (angleCount1 != angleCount2 && std::min(angleCount1, angleCount2) != 1)
516             {
517                 GMX_THROW(InconsistentInputError(
518                         "Number of vectors defined by the two groups are not the same"));
519             }
520         }
521         angleCount_.push_back(angleCount);
522     }
523 }
524
525
526 void Angle::checkSelections(const SelectionList& sel1, const SelectionList& sel2) const
527 {
528     AnglePositionIterator iter1(sel1, natoms1_);
529     AnglePositionIterator iter2(sel2, natoms2_);
530     for (size_t g = 0; g < angleCount_.size(); ++g, iter1.nextGroup(), iter2.nextGroup())
531     {
532         if (iter1.isDynamic() || (iter2.hasValue() && iter2.isDynamic()))
533         {
534             for (int n = 0; n < angleCount_[g]; ++n, iter1.nextValue(), iter2.nextValue())
535             {
536                 bool bOk = true;
537                 if (!iter1.allValuesConsistentlySelected())
538                 {
539                     bOk = false;
540                 }
541                 if (!iter2.allValuesConsistentlySelected())
542                 {
543                     bOk = false;
544                 }
545                 if (angleCount_[g] > 1)
546                 {
547                     if (iter1.hasSingleValue() && !iter1.currentValuesSelected())
548                     {
549                         bOk = false;
550                     }
551                     if (iter2.hasSingleValue() && !iter2.currentValuesSelected())
552                     {
553                         bOk = false;
554                     }
555                 }
556                 if (iter2.hasValue()
557                     && (angleCount_[g] == 1 || (!iter1.hasSingleValue() && !iter2.hasSingleValue()))
558                     && iter1.currentValuesSelected() != iter2.currentValuesSelected())
559                 {
560                     bOk = false;
561                 }
562                 if (!bOk)
563                 {
564                     std::string message = formatString(
565                             "Dynamic selection %d does not select "
566                             "a consistent set of angles over the frames",
567                             static_cast<int>(g + 1));
568                     GMX_THROW(InconsistentInputError(message));
569                 }
570             }
571         }
572     }
573 }
574
575
576 void Angle::initAnalysis(const TrajectoryAnalysisSettings& settings, const TopologyInformation& /* top */)
577 {
578     initFromSelections(sel1_, sel2_);
579
580     // checkSelections() ensures that both selection lists have the same size.
581     angles_.setDataSetCount(angleCount_.size());
582     for (size_t i = 0; i < angleCount_.size(); ++i)
583     {
584         angles_.setColumnCount(i, angleCount_[i]);
585     }
586     double histogramMin = (g1type_ == Group1Type::Dihedral ? -180.0 : 0);
587     histogramModule_->init(histogramFromRange(histogramMin, 180.0).binWidth(binWidth_).includeAll());
588
589     if (g2type_ == Group2Type::TimeZero)
590     {
591         vt0_.resize(sel1_.size());
592         for (size_t g = 0; g < sel1_.size(); ++g)
593         {
594             vt0_[g].resize(sel1_[g].posCount() / natoms1_);
595         }
596     }
597
598     if (!fnAverage_.empty())
599     {
600         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
601         plotm->setFileName(fnAverage_);
602         plotm->setTitle("Average angle");
603         plotm->setXAxisIsTime();
604         plotm->setYLabel("Angle (degrees)");
605         // TODO: Consider adding information about the second selection,
606         // and/or a subtitle describing what kind of angle this is.
607         for (size_t g = 0; g < sel1_.size(); ++g)
608         {
609             plotm->appendLegend(sel1_[g].name());
610         }
611         averageModule_->addModule(plotm);
612     }
613
614     if (!fnAll_.empty())
615     {
616         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
617         plotm->setFileName(fnAll_);
618         plotm->setTitle("Angle");
619         plotm->setXAxisIsTime();
620         plotm->setYLabel("Angle (degrees)");
621         // TODO: Add legends? (there can be a massive amount of columns)
622         angles_.addModule(plotm);
623     }
624
625     if (!fnHistogram_.empty())
626     {
627         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
628         plotm->setFileName(fnHistogram_);
629         plotm->setTitle("Angle histogram");
630         plotm->setXLabel("Angle (degrees)");
631         plotm->setYLabel("Probability");
632         // TODO: Consider adding information about the second selection,
633         // and/or a subtitle describing what kind of angle this is.
634         for (size_t g = 0; g < sel1_.size(); ++g)
635         {
636             plotm->appendLegend(sel1_[g].name());
637         }
638         histogramModule_->averager().addModule(plotm);
639     }
640 }
641
642
643 //! Helper method to calculate a vector from two or three positions.
644 void calc_vec(int natoms, rvec x[], t_pbc* pbc, rvec xout, rvec cout)
645 {
646     switch (natoms)
647     {
648         case 2:
649             if (pbc)
650             {
651                 pbc_dx(pbc, x[1], x[0], xout);
652             }
653             else
654             {
655                 rvec_sub(x[1], x[0], xout);
656             }
657             svmul(0.5, xout, cout);
658             rvec_add(x[0], cout, cout);
659             break;
660         case 3:
661         {
662             rvec v1, v2;
663             if (pbc)
664             {
665                 pbc_dx(pbc, x[1], x[0], v1);
666                 pbc_dx(pbc, x[2], x[0], v2);
667             }
668             else
669             {
670                 rvec_sub(x[1], x[0], v1);
671                 rvec_sub(x[2], x[0], v2);
672             }
673             cprod(v1, v2, xout);
674             rvec_add(x[0], x[1], cout);
675             rvec_add(cout, x[2], cout);
676             svmul(1.0 / 3.0, cout, cout);
677             break;
678         }
679         default: GMX_RELEASE_ASSERT(false, "Incorrectly initialized number of atoms");
680     }
681 }
682
683
684 void Angle::analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* pdata)
685 {
686     AnalysisDataHandle   dh   = pdata->dataHandle(angles_);
687     const SelectionList& sel1 = TrajectoryAnalysisModuleData::parallelSelections(sel1_);
688     const SelectionList& sel2 = TrajectoryAnalysisModuleData::parallelSelections(sel2_);
689
690     checkSelections(sel1, sel2);
691
692     dh.startFrame(frnr, fr.time);
693
694     AnglePositionIterator iter1(sel1, natoms1_);
695     AnglePositionIterator iter2(sel2, natoms2_);
696     for (size_t g = 0; g < angleCount_.size(); ++g, iter1.nextGroup(), iter2.nextGroup())
697     {
698         rvec v1, v2;
699         rvec c1, c2;
700
701         // v2 & c2 are conditionally set in the switch statement below, and conditionally
702         // used in a different switch statement later. Apparently the clang static analyzer
703         // thinks there are cases where they can be used uninitialzed (which I cannot find),
704         // but to avoid trouble if we ever change just one of the switch statements it
705         // makes sense to clear them outside the first switch.
706
707         clear_rvec(v2);
708         clear_rvec(c2);
709
710         switch (g2type_)
711         {
712             case Group2Type::Z: v2[ZZ] = 1.0; break;
713             case Group2Type::SphereNormal: copy_rvec(sel2_[g].position(0).x(), c2); break;
714             default:
715                 // do nothing
716                 break;
717         }
718
719         dh.selectDataSet(g);
720         for (int n = 0; n < angleCount_[g]; ++n, iter1.nextValue(), iter2.nextValue())
721         {
722             rvec x[4];
723             // x[] will be assigned below based on the number of atoms used to initialize iter1,
724             // which in turn should correspond perfectly to g1type_ (which determines how many we
725             // read), but unsurprisingly the static analyzer chokes a bit on that.
726             clear_rvecs(4, x);
727
728             real angle;
729             // checkSelections() ensures that this reflects all the involved
730             // positions.
731             const bool bPresent = iter1.currentValuesSelected() && iter2.currentValuesSelected();
732             iter1.getCurrentPositions(x);
733             switch (g1type_)
734             {
735                 case Group1Type::Angle:
736                     if (pbc)
737                     {
738                         pbc_dx(pbc, x[0], x[1], v1);
739                         pbc_dx(pbc, x[2], x[1], v2);
740                     }
741                     else
742                     {
743                         rvec_sub(x[0], x[1], v1);
744                         rvec_sub(x[2], x[1], v2);
745                     }
746                     angle = gmx_angle(v1, v2);
747                     break;
748                 case Group1Type::Dihedral:
749                 {
750                     rvec dx[3];
751                     if (pbc)
752                     {
753                         pbc_dx(pbc, x[0], x[1], dx[0]);
754                         pbc_dx(pbc, x[2], x[1], dx[1]);
755                         pbc_dx(pbc, x[2], x[3], dx[2]);
756                     }
757                     else
758                     {
759                         rvec_sub(x[0], x[1], dx[0]);
760                         rvec_sub(x[2], x[1], dx[1]);
761                         rvec_sub(x[2], x[3], dx[2]);
762                     }
763                     cprod(dx[0], dx[1], v1);
764                     cprod(dx[1], dx[2], v2);
765                     angle    = gmx_angle(v1, v2);
766                     real ipr = iprod(dx[0], v2);
767                     if (ipr < 0)
768                     {
769                         angle = -angle;
770                     }
771                     break;
772                 }
773                 case Group1Type::Vector:
774                 case Group1Type::Plane:
775                     calc_vec(natoms1_, x, pbc, v1, c1);
776                     switch (g2type_)
777                     {
778                         case Group2Type::Vector:
779                         case Group2Type::Plane:
780                             iter2.getCurrentPositions(x);
781                             calc_vec(natoms2_, x, pbc, v2, c2);
782                             break;
783                         case Group2Type::TimeZero:
784                             // FIXME: This is not parallelizable.
785                             if (frnr == 0)
786                             {
787                                 copy_rvec(v1, vt0_[g][n]);
788                             }
789                             copy_rvec(vt0_[g][n], v2);
790                             break;
791                         case Group2Type::Z: c1[XX] = c1[YY] = 0.0; break;
792                         case Group2Type::SphereNormal:
793                             if (pbc)
794                             {
795                                 pbc_dx(pbc, c1, c2, v2);
796                             }
797                             else
798                             {
799                                 rvec_sub(c1, c2, v2);
800                             }
801                             break;
802                         default: GMX_THROW(InternalError("invalid -g2 value"));
803                     }
804                     angle = gmx_angle(v1, v2);
805                     break;
806                 default: GMX_THROW(InternalError("invalid -g1 value"));
807             }
808             dh.setPoint(n, angle * RAD2DEG, bPresent);
809         }
810     }
811     dh.finishFrame();
812 }
813
814
815 void Angle::finishAnalysis(int /*nframes*/)
816 {
817     AbstractAverageHistogram& averageHistogram = histogramModule_->averager();
818     averageHistogram.normalizeProbability();
819     averageHistogram.done();
820 }
821
822
823 void Angle::writeOutput() {}
824
825 } // namespace
826
827 const char AngleInfo::name[]             = "gangle";
828 const char AngleInfo::shortDescription[] = "Calculate angles";
829
830 TrajectoryAnalysisModulePointer AngleInfo::create()
831 {
832     return TrajectoryAnalysisModulePointer(new Angle);
833 }
834
835 } // namespace analysismodules
836
837 } // namespace gmx