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