2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 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.
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.
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.
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.
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.
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.
37 * Implements gmx::analysismodules::Freevolume.
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40 * \ingroup module_trajectoryanalysis
44 #include "freevolume.h"
48 #include "gromacs/analysisdata/analysisdata.h"
49 #include "gromacs/analysisdata/modules/average.h"
50 #include "gromacs/analysisdata/modules/plot.h"
51 #include "gromacs/fileio/trx.h"
52 #include "gromacs/legacyheaders/copyrite.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/options/basicoptions.h"
55 #include "gromacs/options/filenameoption.h"
56 #include "gromacs/options/options.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/random/random.h"
59 #include "gromacs/selection/nbsearch.h"
60 #include "gromacs/selection/selection.h"
61 #include "gromacs/selection/selectionoption.h"
62 #include "gromacs/topology/atomprop.h"
63 #include "gromacs/topology/topology.h"
64 #include "gromacs/trajectoryanalysis/analysissettings.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/exceptions.h"
71 namespace analysismodules
78 * Class used to compute free volume in a simulations box.
80 * Inherits TrajectoryAnalysisModule and all functions from there.
81 * Does not implement any new functionality.
83 * \ingroup module_trajectoryanalysis
85 class FreeVolume : public TrajectoryAnalysisModule
92 virtual ~FreeVolume();
94 //! Set the options and setting
95 virtual void initOptions(Options *options,
96 TrajectoryAnalysisSettings *settings);
98 //! First routine called by the analysis framework
99 virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
100 const TopologyInformation &top);
102 //! Call for each frame of the trajectory
103 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
104 TrajectoryAnalysisModuleData *pdata);
106 //! Last routine called by the analysis framework
107 virtual void finishAnalysis(int nframes);
109 //! Routine to write output, that is additional over the built-in
110 virtual void writeOutput();
113 std::string fnFreevol_;
116 AnalysisDataAverageModulePointer adata_;
124 AnalysisNeighborhood nb_;
125 //! The van der Waals radius per atom
126 std::vector<double> vdw_radius_;
128 // Copy and assign disallowed by base.
131 // Constructor. Here it is important to initialize the pointer to
132 // subclasses that are elements of the main class. Here we have only
133 // one. The type of this depends on what kind of tool you need.
134 // Here we only have simple value/time kind of data.
135 FreeVolume::FreeVolume()
136 : TrajectoryAnalysisModule(FreeVolumeInfo::name, FreeVolumeInfo::shortDescription),
137 adata_(new AnalysisDataAverageModule())
139 // We only compute two numbers per frame
140 data_.setColumnCount(0, 2);
141 // Tell the analysis framework that this component exists
142 registerAnalysisDataset(&data_, "freevolume");
153 FreeVolume::~FreeVolume()
155 // Destroy C structures where there is no automatic memory release
156 // C++ takes care of memory in classes (hopefully)
159 gmx_rng_destroy(rng_);
165 FreeVolume::initOptions(Options *options,
166 TrajectoryAnalysisSettings *settings)
168 static const char *const desc[] = {
169 "[THISMODULE] calculates the free volume in a box as",
170 "a function of time. The free volume is",
171 "plotted as a fraction of the total volume.",
172 "The program tries to insert a probe with a given radius,",
173 "into the simulations box and if the distance between the",
174 "probe and any atom is less than the sums of the",
175 "van der Waals radii of both atoms, the position is",
176 "considered to be occupied, i.e. non-free. By using a",
177 "probe radius of 0, the true free volume is computed.",
178 "By using a larger radius, e.g. 0.14 nm, roughly corresponding",
179 "to a water molecule, the free volume for a hypothetical",
180 "particle with that size will be produced.",
181 "Note however, that since atoms are treated as hard-spheres",
182 "these number are very approximate, and typically only",
183 "relative changes are meaningful, for instance by doing a",
184 "series of simulations at different temperature.[PAR]",
185 "The group specified by the selection is considered to",
186 "delineate non-free volume.",
187 "The number of insertions per unit of volume is important",
188 "to get a converged result. About 1000/nm^3 yields an overall",
189 "standard deviation that is determined by the fluctuations in",
190 "the trajectory rather than by the fluctuations due to the",
191 "random numbers.[PAR]",
192 "The results are critically dependent on the van der Waals radii;",
193 "we recommend to use the values due to Bondi (1964).[PAR]",
194 "The Fractional Free Volume (FFV) that some authors like to use",
195 "is given by 1 - 1.3*(1-Free Volume). This value is printed on",
199 // Add the descriptive text (program help text) to the options
200 options->setDescription(desc);
202 // Add option for optional output file
203 options->addOption(FileNameOption("o").filetype(eftPlot).outputFile()
204 .store(&fnFreevol_).defaultBasename("freevolume")
205 .description("Computed free volume"));
207 // Add option for selecting a subset of atoms
208 options->addOption(SelectionOption("select")
209 .store(&sel_).defaultSelectionText("all")
212 // Add option for the probe radius and initialize it
213 options->addOption(DoubleOption("radius").store(&probeRadius_)
214 .description("Radius of the probe to be inserted (nm, 0 yields the true free volume)"));
216 // Add option for the random number seed and initialize it to
217 // generate a value automatically
218 options->addOption(IntegerOption("seed").store(&seed_)
219 .description("Seed for random number generator."));
221 // Add option to determine number of insertion trials per frame
222 options->addOption(IntegerOption("ninsert").store(&ninsert_)
223 .description("Number of probe insertions per cubic nm to try for each frame in the trajectory."));
225 // Control input settings
226 settings->setFlags(TrajectoryAnalysisSettings::efRequireTop |
227 TrajectoryAnalysisSettings::efNoUserPBC);
228 settings->setPBC(true);
233 FreeVolume::initAnalysis(const TrajectoryAnalysisSettings &settings,
234 const TopologyInformation &top)
236 // Add the module that will contain the averaging and the time series
237 // for our calculation
238 data_.addModule(adata_);
240 // Add a module for plotting the data automatically at the end of
241 // the calculation. With this in place you only have to add data
242 // points to the data et.
243 AnalysisDataPlotModulePointer plotm_(new AnalysisDataPlotModule());
244 plotm_->setSettings(settings.plotSettings());
245 plotm_->setFileName(fnFreevol_);
246 plotm_->setTitle("Free Volume");
247 plotm_->setXAxisIsTime();
248 plotm_->setYLabel("Free Volume (%)");
249 plotm_->appendLegend("Free Volume");
250 plotm_->appendLegend("Volume");
252 data_.addModule(plotm_);
257 gmx_atomprop_t aps = gmx_atomprop_init();
258 t_atoms *atoms = &(top.topology()->atoms);
260 // Compute total mass
262 for (int i = 0; (i < atoms->nr); i++)
264 mtot_ += atoms->atom[i].m;
267 // Extracts number of molecules
268 nmol_ = top.topology()->mols.nr;
270 // Loop over atoms in the selection using an iterator
271 const int maxnovdw = 10;
272 ConstArrayRef<int> atomind = sel_.atomIndices();
273 for (ConstArrayRef<int>::iterator ai = atomind.begin(); (ai < atomind.end()); ++ai)
275 // Dereference the iterator to obtain an atom number
279 // Lookup the Van der Waals radius of this atom
280 int resnr = atoms->atom[i].resind;
281 if (TRUE == gmx_atomprop_query(aps, epropVDW,
282 *(atoms->resinfo[resnr].name),
283 *(atoms->atomname[i]),
286 vdw_radius_.push_back(value);
295 if (nnovdw < maxnovdw)
297 fprintf(stderr, "Could not determine VDW radius for %s-%s. Set to zero.\n",
298 *(atoms->resinfo[resnr].name),
299 *(atoms->atomname[i]));
301 vdw_radius_.push_back(0.0);
304 gmx_atomprop_destroy(aps);
306 // Increase cutoff by proberadius to make sure we do not miss
308 cutoff_ += probeRadius_;
310 if (nnovdw >= maxnovdw)
312 fprintf(stderr, "Could not determine VDW radius for %d particles. These were set to zero.\n", nnovdw);
315 // Print parameters to output. Maybe should make dependent on
317 printf("cutoff = %g nm\n", cutoff_);
318 printf("probe_radius = %g nm\n", probeRadius_);
319 printf("seed = %d\n", seed_);
320 printf("ninsert = %d probes per nm^3\n", ninsert_);
322 // Initiate the random number generator
323 rng_ = gmx_rng_init(seed_);
325 // Initiate the neighborsearching code
326 nb_.setCutoff(cutoff_);
330 FreeVolume::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
331 TrajectoryAnalysisModuleData *pdata)
333 AnalysisDataHandle dh = pdata->dataHandle(data_);
334 const Selection &sel = pdata->parallelSelection(sel_);
336 GMX_RELEASE_ASSERT(NULL != pbc, "You have no periodic boundary conditions");
338 // Analysis framework magic
339 dh.startFrame(frnr, fr.time);
341 // Compute volume and number of insertions to perform
342 real V = det(fr.box);
343 int Ninsert = static_cast<int>(ninsert_*V);
345 // Use neighborsearching tools!
346 AnalysisNeighborhoodSearch nbsearch = nb_.initSearch(pbc, sel);
348 // Then loop over insertions
350 for (int i = 0; (i < Ninsert); i++)
354 for (int m = 0; (m < DIM); m++)
356 // Generate random number between 0 and 1
357 rand[m] = gmx_rng_uniform_real(rng_);
359 // Generate random 3D position within the box
360 mvmul(fr.box, rand, ins);
362 // Find the first reference position within the cutoff.
363 bool bOverlap = false;
364 AnalysisNeighborhoodPair pair;
365 AnalysisNeighborhoodPairSearch pairSearch = nbsearch.startPairSearch(ins);
366 while (!bOverlap && pairSearch.findNextPair(&pair))
368 int jp = pair.refIndex();
369 // Compute distance vector to first atom in the neighborlist
370 pbc_dx(pbc, ins, sel.position(jp).x(), dx);
372 // See whether the distance is smaller than allowed
373 bOverlap = (norm(dx) <
374 probeRadius_+vdw_radius_[sel.position(jp).refId()]);
380 // We found some free volume!
384 // Compute total free volume for this frame
388 frac = (100.0*NinsTot)/Ninsert;
390 // Add the free volume fraction to the data set in column 0
391 dh.setPoint(0, frac);
392 // Add the total volume to the data set in column 1
401 FreeVolume::finishAnalysis(int /* nframes */)
403 please_cite(stdout, "Bondi1964a");
404 please_cite(stdout, "Lourenco2013a");
408 FreeVolume::writeOutput()
410 // Final results come from statistics module in analysis framework
411 double FVaver = adata_->average(0, 0);
412 double FVerror = adata_->standardDeviation(0, 0);
413 printf("Free volume %.2f +/- %.2f %%\n", FVaver, FVerror);
415 double Vaver = adata_->average(0, 1);
416 double Verror = adata_->standardDeviation(0, 1);
417 printf("Total volume %.2f +/- %.2f nm^3\n", Vaver, Verror);
419 printf("Number of molecules %d total mass %.2f Dalton\n", nmol_, mtot_);
420 double RhoAver = mtot_ / (Vaver * 1e-24 * AVOGADRO);
421 double RhoError = sqr(RhoAver / Vaver)*Verror;
422 printf("Average molar mass: %.2f Dalton\n", mtot_/nmol_);
424 double VmAver = Vaver/nmol_;
425 double VmError = Verror/nmol_;
426 printf("Density rho: %.2f +/- %.2f nm^3\n", RhoAver, RhoError);
427 printf("Molecular volume Vm assuming homogeneity: %.4f +/- %.4f nm^3\n",
430 double VvdWaver = (1-FVaver/100)*VmAver;
431 double VvdWerror = 0;
432 printf("Molecular van der Waals volume assuming homogeneity: %.4f +/- %.4f nm^3\n",
433 VvdWaver, VvdWerror);
435 double FFVaver = 1-1.3*((100-FVaver)/100);
436 double FFVerror = (FVerror/FVaver)*FFVaver;
437 printf("Fractional free volume %.3f +/- %.3f\n", FFVaver, FFVerror);
442 const char FreeVolumeInfo::name[] = "freevolume";
443 const char FreeVolumeInfo::shortDescription[] =
444 "Calculate free volume";
446 TrajectoryAnalysisModulePointer FreeVolumeInfo::create()
448 return TrajectoryAnalysisModulePointer(new FreeVolume);
451 } // namespace analysismodules