1b9362b11d667e3521485b7466446d9b3dc5b194
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / modules / freevolume.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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::analysismodules::Freevolume.
38  *
39  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40  * \ingroup module_trajectoryanalysis
41  */
42 #include "gmxpre.h"
43
44 #include "freevolume.h"
45
46 #include <string>
47
48 #include "gromacs/analysisdata/analysisdata.h"
49 #include "gromacs/analysisdata/modules/average.h"
50 #include "gromacs/analysisdata/modules/plot.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/options/basicoptions.h"
55 #include "gromacs/options/filenameoption.h"
56 #include "gromacs/options/ioptionscontainer.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/random/threefry.h"
59 #include "gromacs/random/uniformrealdistribution.h"
60 #include "gromacs/selection/nbsearch.h"
61 #include "gromacs/selection/selection.h"
62 #include "gromacs/selection/selectionoption.h"
63 #include "gromacs/topology/atomprop.h"
64 #include "gromacs/topology/mtop_util.h"
65 #include "gromacs/topology/topology.h"
66 #include "gromacs/trajectory/trajectoryframe.h"
67 #include "gromacs/trajectoryanalysis/analysissettings.h"
68 #include "gromacs/trajectoryanalysis/topologyinformation.h"
69 #include "gromacs/utility/arrayref.h"
70 #include "gromacs/utility/exceptions.h"
71 #include "gromacs/utility/pleasecite.h"
72
73 namespace gmx
74 {
75
76 namespace analysismodules
77 {
78
79 namespace
80 {
81
82 /*! \brief
83  * Class used to compute free volume in a simulations box.
84  *
85  * Inherits TrajectoryAnalysisModule and all functions from there.
86  * Does not implement any new functionality.
87  *
88  * \ingroup module_trajectoryanalysis
89  */
90 class FreeVolume : public TrajectoryAnalysisModule
91 {
92 public:
93     FreeVolume();
94     ~FreeVolume() override {}
95
96     void initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings) override;
97     void initAnalysis(const TrajectoryAnalysisSettings& settings, const TopologyInformation& top) override;
98     void analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* pdata) override;
99     void finishAnalysis(int nframes) override;
100     void writeOutput() override;
101
102 private:
103     std::string                      fnFreevol_;
104     Selection                        sel_;
105     AnalysisData                     data_;
106     AnalysisDataAverageModulePointer adata_;
107
108     int                      nmol_;
109     double                   mtot_;
110     double                   cutoff_;
111     double                   probeRadius_;
112     gmx::DefaultRandomEngine rng_;
113     int                      seed_, ninsert_;
114     AnalysisNeighborhood     nb_;
115     //! The van der Waals radius per atom
116     std::vector<double> vdw_radius_;
117
118     // Copy and assign disallowed by base.
119 };
120
121 // Constructor. Here it is important to initialize the pointer to
122 // subclasses that are elements of the main class. Here we have only
123 // one. The type of this depends on what kind of tool you need.
124 // Here we only have simple value/time kind of data.
125 FreeVolume::FreeVolume() : adata_(new AnalysisDataAverageModule())
126 {
127     // We only compute two numbers per frame
128     data_.setColumnCount(0, 2);
129     // Tell the analysis framework that this component exists
130     registerAnalysisDataset(&data_, "freevolume");
131     nmol_        = 0;
132     mtot_        = 0;
133     cutoff_      = 0;
134     probeRadius_ = 0;
135     seed_        = 0;
136     ninsert_     = 1000;
137 }
138
139
140 void FreeVolume::initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings)
141 {
142     static const char* const desc[] = {
143         "[THISMODULE] calculates the free volume in a box as",
144         "a function of time. The free volume is",
145         "plotted as a fraction of the total volume.",
146         "The program tries to insert a probe with a given radius,",
147         "into the simulations box and if the distance between the",
148         "probe and any atom is less than the sums of the",
149         "van der Waals radii of both atoms, the position is",
150         "considered to be occupied, i.e. non-free. By using a",
151         "probe radius of 0, the true free volume is computed.",
152         "By using a larger radius, e.g. 0.14 nm, roughly corresponding",
153         "to a water molecule, the free volume for a hypothetical",
154         "particle with that size will be produced.",
155         "Note however, that since atoms are treated as hard-spheres",
156         "these number are very approximate, and typically only",
157         "relative changes are meaningful, for instance by doing a",
158         "series of simulations at different temperature.[PAR]",
159         "The group specified by the selection is considered to",
160         "delineate non-free volume.",
161         "The number of insertions per unit of volume is important",
162         "to get a converged result. About 1000/nm^3 yields an overall",
163         "standard deviation that is determined by the fluctuations in",
164         "the trajectory rather than by the fluctuations due to the",
165         "random numbers.[PAR]",
166         "The results are critically dependent on the van der Waals radii;",
167         "we recommend to use the values due to Bondi (1964).[PAR]",
168         "The Fractional Free Volume (FFV) that some authors like to use",
169         "is given by 1 - 1.3*(1-Free Volume). This value is printed on",
170         "the terminal."
171     };
172
173     settings->setHelpText(desc);
174
175     // Add option for optional output file
176     options->addOption(FileNameOption("o")
177                                .filetype(eftPlot)
178                                .outputFile()
179                                .store(&fnFreevol_)
180                                .defaultBasename("freevolume")
181                                .description("Computed free volume"));
182
183     // Add option for selecting a subset of atoms
184     options->addOption(
185             SelectionOption("select").store(&sel_).defaultSelectionText("all").onlyAtoms().description(
186                     "Atoms that are considered as part of the excluded volume"));
187
188     // Add option for the probe radius and initialize it
189     options->addOption(
190             DoubleOption("radius").store(&probeRadius_).description("Radius of the probe to be inserted (nm, 0 yields the true free volume)"));
191
192     // Add option for the random number seed and initialize it to
193     // generate a value automatically
194     options->addOption(IntegerOption("seed").store(&seed_).description(
195             "Seed for random number generator (0 means generate)."));
196
197     // Add option to determine number of insertion trials per frame
198     options->addOption(IntegerOption("ninsert").store(&ninsert_).description(
199             "Number of probe insertions per cubic nm to try for each frame in the trajectory."));
200
201     // Control input settings
202     settings->setFlags(TrajectoryAnalysisSettings::efRequireTop | TrajectoryAnalysisSettings::efNoUserPBC);
203     settings->setPBC(true);
204 }
205
206
207 void FreeVolume::initAnalysis(const TrajectoryAnalysisSettings& settings, const TopologyInformation& top)
208 {
209     // Add the module that will contain the averaging and the time series
210     // for our calculation
211     data_.addModule(adata_);
212
213     // Add a module for plotting the data automatically at the end of
214     // the calculation. With this in place you only have to add data
215     // points to the data et.
216     AnalysisDataPlotModulePointer plotm_(new AnalysisDataPlotModule());
217     plotm_->setSettings(settings.plotSettings());
218     plotm_->setFileName(fnFreevol_);
219     plotm_->setTitle("Free Volume");
220     plotm_->setXAxisIsTime();
221     plotm_->setYLabel("Free Volume (%)");
222     plotm_->appendLegend("Free Volume");
223     plotm_->appendLegend("Volume");
224
225     data_.addModule(plotm_);
226
227     // Initiate variable
228     cutoff_               = 0;
229     int            nnovdw = 0;
230     AtomProperties aps;
231     auto           atoms = top.copyAtoms();
232
233     // Compute total mass
234     mtot_ = 0;
235     for (int i = 0; (i < atoms->nr); i++)
236     {
237         mtot_ += atoms->atom[i].m;
238     }
239
240     // Extracts number of molecules
241     nmol_ = gmx_mtop_num_molecules(*top.mtop());
242
243     // Loop over atoms in the selection using an iterator
244     const int           maxnovdw = 10;
245     ArrayRef<const int> atomind  = sel_.atomIndices();
246     for (ArrayRef<const int>::iterator ai = atomind.begin(); (ai < atomind.end()); ++ai)
247     {
248         // Dereference the iterator to obtain an atom number
249         int  i = *ai;
250         real value;
251
252         // Lookup the Van der Waals radius of this atom
253         int resnr = atoms->atom[i].resind;
254         if (aps.setAtomProperty(epropVDW, *(atoms->resinfo[resnr].name), *(atoms->atomname[i]), &value))
255         {
256             vdw_radius_.push_back(value);
257             if (value > cutoff_)
258             {
259                 cutoff_ = value;
260             }
261         }
262         else
263         {
264             nnovdw++;
265             if (nnovdw < maxnovdw)
266             {
267                 fprintf(stderr, "Could not determine VDW radius for %s-%s. Set to zero.\n",
268                         *(atoms->resinfo[resnr].name), *(atoms->atomname[i]));
269             }
270             vdw_radius_.push_back(0.0);
271         }
272     }
273
274     // Increase cutoff by proberadius to make sure we do not miss
275     // anything
276     cutoff_ += probeRadius_;
277
278     if (nnovdw >= maxnovdw)
279     {
280         fprintf(stderr, "Could not determine VDW radius for %d particles. These were set to zero.\n", nnovdw);
281     }
282
283     if (seed_ == 0)
284     {
285         seed_ = static_cast<int>(gmx::makeRandomSeed());
286     }
287
288     // Print parameters to output. Maybe should make dependent on
289     // verbosity flag?
290     printf("cutoff       = %g nm\n", cutoff_);
291     printf("probe_radius = %g nm\n", probeRadius_);
292     printf("seed         = %d\n", seed_);
293     printf("ninsert      = %d probes per nm^3\n", ninsert_);
294
295     // Initiate the random number generator
296     rng_.seed(seed_);
297
298     // Initiate the neighborsearching code
299     nb_.setCutoff(cutoff_);
300 }
301
302 void FreeVolume::analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* pdata)
303 {
304     AnalysisDataHandle                 dh  = pdata->dataHandle(data_);
305     const Selection&                   sel = pdata->parallelSelection(sel_);
306     gmx::UniformRealDistribution<real> dist;
307
308     GMX_RELEASE_ASSERT(nullptr != pbc, "You have no periodic boundary conditions");
309
310     // Analysis framework magic
311     dh.startFrame(frnr, fr.time);
312
313     // Compute volume and number of insertions to perform
314     real V       = det(fr.box);
315     int  Ninsert = static_cast<int>(ninsert_ * V);
316
317     // Use neighborsearching tools!
318     AnalysisNeighborhoodSearch nbsearch = nb_.initSearch(pbc, sel);
319
320     // Then loop over insertions
321     int NinsTot = 0;
322     for (int i = 0; (i < Ninsert); i++)
323     {
324         rvec rand, ins, dx;
325
326         for (int m = 0; (m < DIM); m++)
327         {
328             // Generate random number between 0 and 1
329             rand[m] = dist(rng_);
330         }
331         // Generate random 3D position within the box
332         mvmul(fr.box, rand, ins);
333
334         // Find the first reference position within the cutoff.
335         bool                           bOverlap = false;
336         AnalysisNeighborhoodPair       pair;
337         AnalysisNeighborhoodPairSearch pairSearch = nbsearch.startPairSearch(ins);
338         while (!bOverlap && pairSearch.findNextPair(&pair))
339         {
340             int jp = pair.refIndex();
341             // Compute distance vector to first atom in the neighborlist
342             pbc_dx(pbc, ins, sel.position(jp).x(), dx);
343
344             // See whether the distance is smaller than allowed
345             bOverlap = (norm(dx) < probeRadius_ + vdw_radius_[sel.position(jp).refId()]);
346         }
347
348         if (!bOverlap)
349         {
350             // We found some free volume!
351             NinsTot++;
352         }
353     }
354     // Compute total free volume for this frame
355     double frac = 0;
356     if (Ninsert > 0)
357     {
358         frac = (100.0 * NinsTot) / Ninsert;
359     }
360     // Add the free volume fraction to the data set in column 0
361     dh.setPoint(0, frac);
362     // Add the total volume to the data set in column 1
363     dh.setPoint(1, V);
364
365     // Magic
366     dh.finishFrame();
367 }
368
369
370 void FreeVolume::finishAnalysis(int /* nframes */)
371 {
372     please_cite(stdout, "Bondi1964a");
373     please_cite(stdout, "Lourenco2013a");
374 }
375
376 void FreeVolume::writeOutput()
377 {
378     // Final results come from statistics module in analysis framework
379     double FVaver  = adata_->average(0, 0);
380     double FVerror = adata_->standardDeviation(0, 0);
381     printf("Free volume %.2f +/- %.2f %%\n", FVaver, FVerror);
382
383     double Vaver  = adata_->average(0, 1);
384     double Verror = adata_->standardDeviation(0, 1);
385     printf("Total volume %.2f +/- %.2f nm^3\n", Vaver, Verror);
386
387     printf("Number of molecules %d total mass %.2f Dalton\n", nmol_, mtot_);
388     double RhoAver  = mtot_ / (Vaver * 1e-24 * AVOGADRO);
389     double RhoError = gmx::square(RhoAver / Vaver) * Verror;
390     printf("Average molar mass: %.2f Dalton\n", mtot_ / nmol_);
391
392     double VmAver  = Vaver / nmol_;
393     double VmError = Verror / nmol_;
394     printf("Density rho: %.2f +/- %.2f nm^3\n", RhoAver, RhoError);
395     printf("Molecular volume Vm assuming homogeneity: %.4f +/- %.4f nm^3\n", VmAver, VmError);
396
397     double VvdWaver  = (1 - FVaver / 100) * VmAver;
398     double VvdWerror = 0;
399     printf("Molecular van der Waals volume assuming homogeneity:  %.4f +/- %.4f nm^3\n", VvdWaver, VvdWerror);
400
401     double FFVaver  = 1 - 1.3 * ((100 - FVaver) / 100);
402     double FFVerror = (FVerror / FVaver) * FFVaver;
403     printf("Fractional free volume %.3f +/- %.3f\n", FFVaver, FFVerror);
404 }
405
406 } // namespace
407
408 const char FreeVolumeInfo::name[]             = "freevolume";
409 const char FreeVolumeInfo::shortDescription[] = "Calculate free volume";
410
411 TrajectoryAnalysisModulePointer FreeVolumeInfo::create()
412 {
413     return TrajectoryAnalysisModulePointer(new FreeVolume);
414 }
415
416 } // namespace analysismodules
417
418 } // namespace gmx