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/legacyheaders/copyrite.h"
50 #include "gromacs/analysisdata/analysisdata.h"
51 #include "gromacs/analysisdata/modules/average.h"
52 #include "gromacs/analysisdata/modules/plot.h"
53 #include "gromacs/fileio/trx.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/options/basicoptions.h"
56 #include "gromacs/options/filenameoption.h"
57 #include "gromacs/options/options.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/random/random.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/topology.h"
65 #include "gromacs/trajectoryanalysis/analysissettings.h"
66 #include "gromacs/utility/arrayref.h"
67 #include "gromacs/utility/exceptions.h"
72 namespace analysismodules
79 * Class used to compute free volume in a simulations box.
81 * Inherits TrajectoryAnalysisModule and all functions from there.
82 * Does not implement any new functionality.
84 * \ingroup module_trajectoryanalysis
86 class FreeVolume : public TrajectoryAnalysisModule
93 virtual ~FreeVolume();
95 //! Set the options and setting
96 virtual void initOptions(Options *options,
97 TrajectoryAnalysisSettings *settings);
99 //! First routine called by the analysis framework
100 virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
101 const TopologyInformation &top);
103 //! Call for each frame of the trajectory
104 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
105 TrajectoryAnalysisModuleData *pdata);
107 //! Last routine called by the analysis framework
108 virtual void finishAnalysis(int nframes);
110 //! Routine to write output, that is additional over the built-in
111 virtual void writeOutput();
114 std::string fnFreevol_;
117 AnalysisDataAverageModulePointer adata_;
125 AnalysisNeighborhood nb_;
126 //! The van der Waals radius per atom
127 std::vector<double> vdw_radius_;
129 // Copy and assign disallowed by base.
132 // Constructor. Here it is important to initialize the pointer to
133 // subclasses that are elements of the main class. Here we have only
134 // one. The type of this depends on what kind of tool you need.
135 // Here we only have simple value/time kind of data.
136 FreeVolume::FreeVolume()
137 : TrajectoryAnalysisModule(FreeVolumeInfo::name, FreeVolumeInfo::shortDescription),
138 adata_(new AnalysisDataAverageModule())
140 // We only compute two numbers per frame
141 data_.setColumnCount(0, 2);
142 // Tell the analysis framework that this component exists
143 registerAnalysisDataset(&data_, "freevolume");
154 FreeVolume::~FreeVolume()
156 // Destroy C structures where there is no automatic memory release
157 // C++ takes care of memory in classes (hopefully)
160 gmx_rng_destroy(rng_);
166 FreeVolume::initOptions(Options *options,
167 TrajectoryAnalysisSettings *settings)
169 static const char *const desc[] = {
170 "[THISMODULE] calculates the free volume in a box as",
171 "a function of time. The free volume is",
172 "plotted as a fraction of the total volume.",
173 "The program tries to insert a probe with a given radius,",
174 "into the simulations box and if the distance between the",
175 "probe and any atom is less than the sums of the",
176 "van der Waals radii of both atoms, the position is",
177 "considered to be occupied, i.e. non-free. By using a",
178 "probe radius of 0, the true free volume is computed.",
179 "By using a larger radius, e.g. 0.14 nm, roughly corresponding",
180 "to a water molecule, the free volume for a hypothetical",
181 "particle with that size will be produced.",
182 "Note however, that since atoms are treated as hard-spheres",
183 "these number are very approximate, and typically only",
184 "relative changes are meaningful, for instance by doing a",
185 "series of simulations at different temperature.[PAR]",
186 "The group specified by the selection is considered to",
187 "delineate non-free volume.",
188 "The number of insertions per unit of volume is important",
189 "to get a converged result. About 1000/nm^3 yields an overall",
190 "standard deviation that is determined by the fluctuations in",
191 "the trajectory rather than by the fluctuations due to the",
192 "random numbers.[PAR]",
193 "The results are critically dependent on the van der Waals radii;",
194 "we recommend to use the values due to Bondi (1964).[PAR]",
195 "The Fractional Free Volume (FFV) that some authors like to use",
196 "is given by 1 - 1.3*(1-Free Volume). This value is printed on",
200 // Add the descriptive text (program help text) to the options
201 options->setDescription(desc);
203 // Add option for optional output file
204 options->addOption(FileNameOption("o").filetype(eftPlot).outputFile()
205 .store(&fnFreevol_).defaultBasename("freevolume")
206 .description("Computed free volume"));
208 // Add option for selecting a subset of atoms
209 options->addOption(SelectionOption("select")
210 .store(&sel_).defaultSelectionText("all")
213 // Add option for the probe radius and initialize it
214 options->addOption(DoubleOption("radius").store(&probeRadius_)
215 .description("Radius of the probe to be inserted (nm, 0 yields the true free volume)"));
217 // Add option for the random number seed and initialize it to
218 // generate a value automatically
219 options->addOption(IntegerOption("seed").store(&seed_)
220 .description("Seed for random number generator."));
222 // Add option to determine number of insertion trials per frame
223 options->addOption(IntegerOption("ninsert").store(&ninsert_)
224 .description("Number of probe insertions per cubic nm to try for each frame in the trajectory."));
226 // Control input settings
227 settings->setFlags(TrajectoryAnalysisSettings::efRequireTop |
228 TrajectoryAnalysisSettings::efNoUserPBC);
229 settings->setPBC(true);
234 FreeVolume::initAnalysis(const TrajectoryAnalysisSettings &settings,
235 const TopologyInformation &top)
237 // Add the module that will contain the averaging and the time series
238 // for our calculation
239 data_.addModule(adata_);
241 // Add a module for plotting the data automatically at the end of
242 // the calculation. With this in place you only have to add data
243 // points to the data et.
244 AnalysisDataPlotModulePointer plotm_(new AnalysisDataPlotModule());
245 plotm_->setSettings(settings.plotSettings());
246 plotm_->setFileName(fnFreevol_);
247 plotm_->setTitle("Free Volume");
248 plotm_->setXAxisIsTime();
249 plotm_->setYLabel("Free Volume (%)");
250 plotm_->appendLegend("Free Volume");
251 plotm_->appendLegend("Volume");
253 data_.addModule(plotm_);
258 gmx_atomprop_t aps = gmx_atomprop_init();
259 t_atoms *atoms = &(top.topology()->atoms);
261 // Compute total mass
263 for (int i = 0; (i < atoms->nr); i++)
265 mtot_ += atoms->atom[i].m;
268 // Extracts number of molecules
269 nmol_ = top.topology()->mols.nr;
271 // Loop over atoms in the selection using an iterator
272 const int maxnovdw = 10;
273 ConstArrayRef<int> atomind = sel_.atomIndices();
274 for (ConstArrayRef<int>::iterator ai = atomind.begin(); (ai < atomind.end()); ++ai)
276 // Dereference the iterator to obtain an atom number
280 // Lookup the Van der Waals radius of this atom
281 int resnr = atoms->atom[i].resind;
282 if (TRUE == gmx_atomprop_query(aps, epropVDW,
283 *(atoms->resinfo[resnr].name),
284 *(atoms->atomname[i]),
287 vdw_radius_.push_back(value);
296 if (nnovdw < maxnovdw)
298 fprintf(stderr, "Could not determine VDW radius for %s-%s. Set to zero.\n",
299 *(atoms->resinfo[resnr].name),
300 *(atoms->atomname[i]));
302 vdw_radius_.push_back(0.0);
305 gmx_atomprop_destroy(aps);
307 // Increase cutoff by proberadius to make sure we do not miss
309 cutoff_ += probeRadius_;
311 if (nnovdw >= maxnovdw)
313 fprintf(stderr, "Could not determine VDW radius for %d particles. These were set to zero.\n", nnovdw);
316 // Print parameters to output. Maybe should make dependent on
318 printf("cutoff = %g nm\n", cutoff_);
319 printf("probe_radius = %g nm\n", probeRadius_);
320 printf("seed = %d\n", seed_);
321 printf("ninsert = %d probes per nm^3\n", ninsert_);
323 // Initiate the random number generator
324 rng_ = gmx_rng_init(seed_);
326 // Initiate the neighborsearching code
327 nb_.setCutoff(cutoff_);
331 FreeVolume::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
332 TrajectoryAnalysisModuleData *pdata)
334 AnalysisDataHandle dh = pdata->dataHandle(data_);
335 const Selection &sel = pdata->parallelSelection(sel_);
337 GMX_RELEASE_ASSERT(NULL != pbc, "You have no periodic boundary conditions");
339 // Analysis framework magic
340 dh.startFrame(frnr, fr.time);
342 // Compute volume and number of insertions to perform
343 real V = det(fr.box);
344 int Ninsert = static_cast<int>(ninsert_*V);
346 // Use neighborsearching tools!
347 AnalysisNeighborhoodSearch nbsearch = nb_.initSearch(pbc, sel);
349 // Then loop over insertions
351 for (int i = 0; (i < Ninsert); i++)
355 for (int m = 0; (m < DIM); m++)
357 // Generate random number between 0 and 1
358 rand[m] = gmx_rng_uniform_real(rng_);
360 // Generate random 3D position within the box
361 mvmul(fr.box, rand, ins);
363 // Find the first reference position within the cutoff.
364 bool bOverlap = false;
365 AnalysisNeighborhoodPair pair;
366 AnalysisNeighborhoodPairSearch pairSearch = nbsearch.startPairSearch(ins);
367 while (!bOverlap && pairSearch.findNextPair(&pair))
369 int jp = pair.refIndex();
370 // Compute distance vector to first atom in the neighborlist
371 pbc_dx(pbc, ins, sel.position(jp).x(), dx);
373 // See whether the distance is smaller than allowed
374 bOverlap = (norm(dx) <
375 probeRadius_+vdw_radius_[sel.position(jp).refId()]);
381 // We found some free volume!
385 // Compute total free volume for this frame
389 frac = (100.0*NinsTot)/Ninsert;
391 // Add the free volume fraction to the data set in column 0
392 dh.setPoint(0, frac);
393 // Add the total volume to the data set in column 1
402 FreeVolume::finishAnalysis(int /* nframes */)
404 please_cite(stdout, "Bondi1964a");
405 please_cite(stdout, "Lourenco2013a");
409 FreeVolume::writeOutput()
411 // Final results come from statistics module in analysis framework
412 double FVaver = adata_->average(0, 0);
413 double FVerror = adata_->standardDeviation(0, 0);
414 printf("Free volume %.2f +/- %.2f %%\n", FVaver, FVerror);
416 double Vaver = adata_->average(0, 1);
417 double Verror = adata_->standardDeviation(0, 1);
418 printf("Total volume %.2f +/- %.2f nm^3\n", Vaver, Verror);
420 printf("Number of molecules %d total mass %.2f Dalton\n", nmol_, mtot_);
421 double RhoAver = mtot_ / (Vaver * 1e-24 * AVOGADRO);
422 double RhoError = sqr(RhoAver / Vaver)*Verror;
423 printf("Average molar mass: %.2f Dalton\n", mtot_/nmol_);
425 double VmAver = Vaver/nmol_;
426 double VmError = Verror/nmol_;
427 printf("Density rho: %.2f +/- %.2f nm^3\n", RhoAver, RhoError);
428 printf("Molecular volume Vm assuming homogeneity: %.4f +/- %.4f nm^3\n",
431 double VvdWaver = (1-FVaver/100)*VmAver;
432 double VvdWerror = 0;
433 printf("Molecular van der Waals volume assuming homogeneity: %.4f +/- %.4f nm^3\n",
434 VvdWaver, VvdWerror);
436 double FFVaver = 1-1.3*((100-FVaver)/100);
437 double FFVerror = (FVerror/FVaver)*FFVaver;
438 printf("Fractional free volume %.3f +/- %.3f\n", FFVaver, FFVerror);
443 const char FreeVolumeInfo::name[] = "freevolume";
444 const char FreeVolumeInfo::shortDescription[] =
445 "Calculate free volume";
447 TrajectoryAnalysisModulePointer FreeVolumeInfo::create()
449 return TrajectoryAnalysisModulePointer(new FreeVolume);
452 } // namespace analysismodules