6acf73ab1d2a344614b13634fa0c0ed342b42014
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / modules / distance.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010-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::Distance.
39  *
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  * \ingroup module_trajectoryanalysis
42  */
43 #include "gmxpre.h"
44
45 #include "distance.h"
46
47 #include <memory>
48 #include <string>
49
50 #include "gromacs/analysisdata/analysisdata.h"
51 #include "gromacs/analysisdata/modules/average.h"
52 #include "gromacs/analysisdata/modules/histogram.h"
53 #include "gromacs/analysisdata/modules/plot.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/options/basicoptions.h"
56 #include "gromacs/options/filenameoption.h"
57 #include "gromacs/options/ioptionscontainer.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/selection/selection.h"
60 #include "gromacs/selection/selectionoption.h"
61 #include "gromacs/trajectory/trajectoryframe.h"
62 #include "gromacs/trajectoryanalysis/analysissettings.h"
63 #include "gromacs/utility/arrayref.h"
64 #include "gromacs/utility/exceptions.h"
65 #include "gromacs/utility/stringutil.h"
66
67 namespace gmx
68 {
69
70 namespace analysismodules
71 {
72
73 namespace
74 {
75
76 class Distance : public TrajectoryAnalysisModule
77 {
78 public:
79     Distance();
80
81     void initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings) override;
82     void initAnalysis(const TrajectoryAnalysisSettings& settings, const TopologyInformation& top) override;
83
84     void analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* pdata) override;
85
86     void finishAnalysis(int nframes) override;
87     void writeOutput() override;
88
89 private:
90     SelectionList sel_;
91     std::string   fnAverage_;
92     std::string   fnAll_;
93     std::string   fnXYZ_;
94     std::string   fnHistogram_;
95     std::string   fnAllStats_;
96     double        meanLength_;
97     double        lengthDev_;
98     double        binWidth_;
99
100     AnalysisData                             distances_;
101     AnalysisData                             xyz_;
102     AnalysisDataAverageModulePointer         summaryStatsModule_;
103     AnalysisDataAverageModulePointer         allStatsModule_;
104     AnalysisDataFrameAverageModulePointer    averageModule_;
105     AnalysisDataSimpleHistogramModulePointer histogramModule_;
106
107     // Copy and assign disallowed by base.
108 };
109
110 Distance::Distance() :
111     meanLength_(0.1),
112     lengthDev_(1.0),
113     binWidth_(0.001),
114     summaryStatsModule_(std::make_unique<AnalysisDataAverageModule>()),
115     allStatsModule_(std::make_unique<AnalysisDataAverageModule>()),
116     averageModule_(std::make_unique<AnalysisDataFrameAverageModule>()),
117     histogramModule_(std::make_unique<AnalysisDataSimpleHistogramModule>())
118 {
119     summaryStatsModule_->setAverageDataSets(true);
120     distances_.addModule(summaryStatsModule_);
121     distances_.addModule(allStatsModule_);
122     distances_.addModule(averageModule_);
123     distances_.addModule(histogramModule_);
124
125     registerAnalysisDataset(&distances_, "dist");
126     registerAnalysisDataset(&xyz_, "xyz");
127     registerBasicDataset(summaryStatsModule_.get(), "stats");
128     registerBasicDataset(allStatsModule_.get(), "allstats");
129     registerBasicDataset(averageModule_.get(), "average");
130     registerBasicDataset(&histogramModule_->averager(), "histogram");
131 }
132
133
134 void Distance::initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings)
135 {
136     static const char* const desc[] = {
137         "[THISMODULE] calculates distances between pairs of positions",
138         "as a function of time. Each selection specifies an independent set",
139         "of distances to calculate. Each selection should consist of pairs",
140         "of positions, and the distances are computed between positions 1-2,",
141         "3-4, etc.[PAR]",
142         "[TT]-oav[tt] writes the average distance as a function of time for",
143         "each selection.",
144         "[TT]-oall[tt] writes all the individual distances.",
145         "[TT]-oxyz[tt] does the same, but the x, y, and z components of the",
146         "distance are written instead of the norm.",
147         "[TT]-oh[tt] writes a histogram of the distances for each selection.",
148         "The location of the histogram is set with [TT]-len[tt] and",
149         "[TT]-tol[tt]. Bin width is set with [TT]-binw[tt].",
150         "[TT]-oallstat[tt] writes out the average and standard deviation for",
151         "each individual distance, calculated over the frames.[PAR]",
152         "Note that [THISMODULE] calculates distances between fixed pairs",
153         "(1-2, 3-4, etc.) within a single selection.  To calculate distances",
154         "between two selections, including minimum, maximum, and pairwise",
155         "distances, use [gmx-pairdist]."
156     };
157
158     settings->setHelpText(desc);
159
160     options->addOption(FileNameOption("oav")
161                                .filetype(eftPlot)
162                                .outputFile()
163                                .store(&fnAverage_)
164                                .defaultBasename("distave")
165                                .description("Average distances as function of time"));
166     options->addOption(FileNameOption("oall")
167                                .filetype(eftPlot)
168                                .outputFile()
169                                .store(&fnAll_)
170                                .defaultBasename("dist")
171                                .description("All distances as function of time"));
172     options->addOption(FileNameOption("oxyz")
173                                .filetype(eftPlot)
174                                .outputFile()
175                                .store(&fnXYZ_)
176                                .defaultBasename("distxyz")
177                                .description("Distance components as function of time"));
178     options->addOption(FileNameOption("oh")
179                                .filetype(eftPlot)
180                                .outputFile()
181                                .store(&fnHistogram_)
182                                .defaultBasename("disthist")
183                                .description("Histogram of the distances"));
184     options->addOption(FileNameOption("oallstat")
185                                .filetype(eftPlot)
186                                .outputFile()
187                                .store(&fnAllStats_)
188                                .defaultBasename("diststat")
189                                .description("Statistics for individual distances"));
190     options->addOption(
191             SelectionOption("select").storeVector(&sel_).required().dynamicMask().multiValue().description(
192                     "Position pairs to calculate distances for"));
193     // TODO: Extend the histogramming implementation to allow automatic
194     // extension of the histograms to cover the data, removing the need for
195     // the first two options.
196     options->addOption(
197             DoubleOption("len").store(&meanLength_).description("Mean distance for histogramming"));
198     options->addOption(
199             DoubleOption("tol").store(&lengthDev_).description("Width of full distribution as fraction of [TT]-len[tt]"));
200     options->addOption(
201             DoubleOption("binw").store(&binWidth_).description("Bin width for histogramming"));
202 }
203
204
205 /*! \brief
206  * Checks that selections conform to the expectations of the tool.
207  */
208 void checkSelections(const SelectionList& sel)
209 {
210     for (size_t g = 0; g < sel.size(); ++g)
211     {
212         if (sel[g].posCount() % 2 != 0)
213         {
214             std::string message = formatString(
215                     "Selection '%s' does not evaluate into an even number of positions "
216                     "(there are %d positions)",
217                     sel[g].name(), sel[g].posCount());
218             GMX_THROW(InconsistentInputError(message));
219         }
220         if (sel[g].isDynamic())
221         {
222             for (int i = 0; i < sel[g].posCount(); i += 2)
223             {
224                 if (sel[g].position(i).selected() != sel[g].position(i + 1).selected())
225                 {
226                     std::string message = formatString(
227                             "Dynamic selection %d does not select "
228                             "a consistent set of pairs over the frames",
229                             static_cast<int>(g + 1));
230                     GMX_THROW(InconsistentInputError(message));
231                 }
232             }
233         }
234     }
235 }
236
237
238 void Distance::initAnalysis(const TrajectoryAnalysisSettings& settings, const TopologyInformation& /*top*/)
239 {
240     checkSelections(sel_);
241
242     distances_.setDataSetCount(sel_.size());
243     xyz_.setDataSetCount(sel_.size());
244     for (size_t i = 0; i < sel_.size(); ++i)
245     {
246         const int distCount = sel_[i].posCount() / 2;
247         distances_.setColumnCount(i, distCount);
248         xyz_.setColumnCount(i, distCount * 3);
249     }
250     const double histogramMin = (1.0 - lengthDev_) * meanLength_;
251     const double histogramMax = (1.0 + lengthDev_) * meanLength_;
252     histogramModule_->init(histogramFromRange(histogramMin, histogramMax).binWidth(binWidth_).includeAll());
253
254     if (!fnAverage_.empty())
255     {
256         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
257         plotm->setFileName(fnAverage_);
258         plotm->setTitle("Average distance");
259         plotm->setXAxisIsTime();
260         plotm->setYLabel("Distance (nm)");
261         for (size_t g = 0; g < sel_.size(); ++g)
262         {
263             plotm->appendLegend(sel_[g].name());
264         }
265         averageModule_->addModule(plotm);
266     }
267
268     if (!fnAll_.empty())
269     {
270         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
271         plotm->setFileName(fnAll_);
272         plotm->setTitle("Distance");
273         plotm->setXAxisIsTime();
274         plotm->setYLabel("Distance (nm)");
275         // TODO: Add legends? (there can be a massive amount of columns)
276         distances_.addModule(plotm);
277     }
278
279     if (!fnXYZ_.empty())
280     {
281         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
282         plotm->setFileName(fnXYZ_);
283         plotm->setTitle("Distance");
284         plotm->setXAxisIsTime();
285         plotm->setYLabel("Distance (nm)");
286         // TODO: Add legends? (there can be a massive amount of columns)
287         xyz_.addModule(plotm);
288     }
289
290     if (!fnHistogram_.empty())
291     {
292         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
293         plotm->setFileName(fnHistogram_);
294         plotm->setTitle("Distance histogram");
295         plotm->setXLabel("Distance (nm)");
296         plotm->setYLabel("Probability");
297         for (size_t g = 0; g < sel_.size(); ++g)
298         {
299             plotm->appendLegend(sel_[g].name());
300         }
301         histogramModule_->averager().addModule(plotm);
302     }
303
304     if (!fnAllStats_.empty())
305     {
306         AnalysisDataPlotModulePointer plotm(new AnalysisDataPlotModule(settings.plotSettings()));
307         plotm->setFileName(fnAllStats_);
308         plotm->setErrorsAsSeparateColumn(true);
309         plotm->setTitle("Statistics for individual distances");
310         plotm->setXLabel("Distance index");
311         plotm->setYLabel("Average/standard deviation (nm)");
312         for (size_t g = 0; g < sel_.size(); ++g)
313         {
314             plotm->appendLegend(std::string(sel_[g].name()) + " avg");
315             plotm->appendLegend(std::string(sel_[g].name()) + " std.dev.");
316         }
317         // TODO: Consider whether this output format is the best possible.
318         allStatsModule_->addModule(plotm);
319     }
320 }
321
322
323 void Distance::analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* pdata)
324 {
325     AnalysisDataHandle   distHandle = pdata->dataHandle(distances_);
326     AnalysisDataHandle   xyzHandle  = pdata->dataHandle(xyz_);
327     const SelectionList& sel        = TrajectoryAnalysisModuleData::parallelSelections(sel_);
328
329     checkSelections(sel);
330
331     distHandle.startFrame(frnr, fr.time);
332     xyzHandle.startFrame(frnr, fr.time);
333     for (size_t g = 0; g < sel.size(); ++g)
334     {
335         distHandle.selectDataSet(g);
336         xyzHandle.selectDataSet(g);
337         for (int i = 0, n = 0; i < sel[g].posCount(); i += 2, ++n)
338         {
339             const SelectionPosition& p1 = sel[g].position(i);
340             const SelectionPosition& p2 = sel[g].position(i + 1);
341             rvec                     dx;
342             if (pbc != nullptr)
343             {
344                 pbc_dx(pbc, p2.x(), p1.x(), dx);
345             }
346             else
347             {
348                 rvec_sub(p2.x(), p1.x(), dx);
349             }
350             real dist     = norm(dx);
351             bool bPresent = p1.selected() && p2.selected();
352             distHandle.setPoint(n, dist, bPresent);
353             xyzHandle.setPoints(n * 3, 3, dx, bPresent);
354         }
355     }
356     distHandle.finishFrame();
357     xyzHandle.finishFrame();
358 }
359
360
361 void Distance::finishAnalysis(int /*nframes*/)
362 {
363     AbstractAverageHistogram& averageHistogram = histogramModule_->averager();
364     averageHistogram.normalizeProbability();
365     averageHistogram.done();
366 }
367
368
369 void Distance::writeOutput()
370 {
371     SelectionList::const_iterator sel;
372     int                           index;
373     for (sel = sel_.begin(), index = 0; sel != sel_.end(); ++sel, ++index)
374     {
375         printf("%s:\n", sel->name());
376         printf("  Number of samples:  %d\n", summaryStatsModule_->sampleCount(index, 0));
377         printf("  Average distance:   %-8.5f nm\n", summaryStatsModule_->average(index, 0));
378         printf("  Standard deviation: %-8.5f nm\n", summaryStatsModule_->standardDeviation(index, 0));
379     }
380 }
381
382 } // namespace
383
384 const char DistanceInfo::name[]             = "distance";
385 const char DistanceInfo::shortDescription[] = "Calculate distances between pairs of positions";
386
387 TrajectoryAnalysisModulePointer DistanceInfo::create()
388 {
389     return TrajectoryAnalysisModulePointer(new Distance);
390 }
391
392 } // namespace analysismodules
393
394 } // namespace gmx