2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
5 * Copyright (c) 2018,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.
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.
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.
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.
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.
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.
39 * Tool for extracting AWH (Accelerated Weight Histogram) data from energy files.
41 * \author Viveca Lindahl
55 #include "gromacs/commandline/pargs.h"
56 #include "gromacs/fileio/enxio.h"
57 #include "gromacs/fileio/gmxfio.h"
58 #include "gromacs/fileio/oenv.h"
59 #include "gromacs/fileio/tpxio.h"
60 #include "gromacs/fileio/trxio.h"
61 #include "gromacs/fileio/xvgr.h"
62 #include "gromacs/gmxana/gmx_ana.h"
63 #include "gromacs/math/units.h"
64 #include "gromacs/mdtypes/awh_params.h"
65 #include "gromacs/mdtypes/inputrec.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/topology/topology.h"
68 #include "gromacs/trajectory/energyframe.h"
69 #include "gromacs/utility/arraysize.h"
70 #include "gromacs/utility/basedefinitions.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/futil.h"
73 #include "gromacs/utility/gmxassert.h"
74 #include "gromacs/utility/smalloc.h"
75 #include "gromacs/utility/stringutil.h"
77 using gmx::AwhBiasParams
;
83 //! Enum for choosing AWH graph output
84 enum class AwhGraphSelection
86 Pmf
, //!< Only the PMF
87 All
//!< All possible AWH graphs
90 //! Energy unit options
99 /*! \brief All meta-data that is shared for one output file type for one bias */
103 /*! \brief Constructor, Set the output base file name and title.
105 * Result is a valid object, but will produce empty output files.
107 * \param[in] filename The name for output files, frame time, and possibly bias number, will be added per file/frame.
108 * \param[in] baseTitle The base title of the plot, the bias number might be added.
109 * \param[in] numBias The total number of AWH biases in the system.
110 * \param[in] biasIndex The index of this bias.
112 OutputFile(const std::string
& filename
, const std::string
& baseTitle
, int numBias
, int biasIndex
);
114 /*! \brief Initializes the output file setup for the AWH output.
116 * \param[in] subBlockStart Index of the first sub-block to write in the energy frame.
117 * \param[in] numSubBlocks The total number of sub-blocks in the framw.
118 * \param[in] awhBiasParams The AWH bias parameters.
119 * \param[in] graphSelection Selects which graphs to plot.
120 * \param[in] energyUnit Requested energy unit in output.
121 * \param[in] kTValue kB*T in kJ/mol.
123 void initializeAwhOutputFile(int subBlockStart
,
125 const AwhBiasParams
* awhBiasParams
,
126 AwhGraphSelection graphSelection
,
127 EnergyUnit energyUnit
,
130 /*! \brief Initializes the output file setup for the fricion output.
132 * \param[in] subBlockStart Index of the first sub-block to write in the energy frame.
133 * \param[in] numSubBlocks The total number of sub-blocks in the framw.
134 * \param[in] awhBiasParams The AWH bias parameters.
135 * \param[in] energyUnit Requested energy unit in output.
136 * \param[in] kTValue kB*T in kJ/mol.
138 void initializeFrictionOutputFile(int subBlockStart
,
140 const AwhBiasParams
* awhBiasParams
,
141 EnergyUnit energyUnit
,
144 /*! \brief Opens a single output file for a bias, prints title and legends.
146 * \param[in] time The time for of frame to be written.
147 * \param[in] oenv The output environment.
148 * \returns the file pointer.
150 FILE* openBiasOutputFile(double time
, const gmx_output_env_t
* oenv
) const;
152 /*! \brief Writes data selected from \p block to file.
154 * \param[in] block The energy block with the data to write.
155 * \param[in] subBlockStart The sub-block start index.
156 * \param[in] fp The file pointer.
158 void writeData(const t_enxblock
& block
, int subBlockStart
, FILE* fp
) const;
161 std::string baseFilename_
; /**< Base of the output file name. */
162 std::string title_
; /**< Title for the graph. */
163 int numDim_
; /**< Number of dimensions. */
164 int firstGraphSubBlock_
; /**< Index in the energy sub-blocks for the first graph. */
165 int numGraph_
; /**< Number of actual graphs. */
166 bool useKTForEnergy_
; /**< Whether to use kT instead of kJ/mol for energy. */
167 std::vector
<real
> scaleFactor_
; /**< Scale factors from energy file data to output for each of the numGraph quantities. */
169 std::vector
<std::string
> legend_
; /**< Legends for the output */
170 std::string xLabel_
; /**< Label for the x-axis. */
171 std::string yLabel_
; /**< Label for the y-axis. */
174 /*! \brief All meta-data that is shared for all output files for one bias */
175 struct BiasOutputSetup
178 BiasOutputSetup(int subblockStart
,
179 std::unique_ptr
<OutputFile
> awhOutputFile
,
180 std::unique_ptr
<OutputFile
> frictionOutputFile
) :
181 subblockStart_(subblockStart
),
182 awhOutputFile_(std::move(awhOutputFile
)),
183 frictionOutputFile_(std::move(frictionOutputFile
))
187 //! Return the AWH output file data.
188 const OutputFile
& awhOutputFile() const
190 GMX_RELEASE_ASSERT(awhOutputFile_
!= nullptr,
191 "awhOutputFile() called without initialized AWH output file");
193 return *awhOutputFile_
;
196 //! Return the a pointer to the friction output file data, can return nullptr
197 const OutputFile
* frictionOutputFile() const { return frictionOutputFile_
.get(); }
199 //! Return the starting subblock.
200 int subblockStart() const { return subblockStart_
; }
203 const int subblockStart_
; /**< The start index of the subblocks to read. */
204 std::unique_ptr
<OutputFile
> awhOutputFile_
; /**< The standard AWH output file data. */
205 std::unique_ptr
<OutputFile
> frictionOutputFile_
; /**< The friction/metric tensor output file data */
208 /*! \brief All options and meta-data needed for the AWH output */
213 AwhReader(const AwhParams
* awhParams
,
215 const t_filenm
* filenames
,
216 AwhGraphSelection awhGraphSelection
,
217 EnergyUnit energyUnit
,
219 const t_enxblock
* block
);
221 //! Extract and print AWH data for an AWH block of one energy frame
222 void processAwhFrame(const t_enxblock
& block
, double time
, const gmx_output_env_t
* oenv
) const;
225 std::vector
<BiasOutputSetup
> biasOutputSetups_
; /**< The output setups, one for each AWH bias. */
227 const real kT_
; /**< kB*T in kJ/mol. */
233 //! Enum for selecting output file type.
234 enum class OutputFileType
236 Awh
, //!< AWH, all data except friction tensor.
237 Friction
//!< The full friction tensor.
240 /*! \brief The maximum number of observables per subblock, not including the full friction tensor.
242 * This number, as well as the legends in make_legend() should match
243 * the output that mdrun writes. It would be better to define these
244 * values in a single location.
246 constexpr int maxAwhGraphs
= 6;
248 /*! \brief Constructs a legend for a standard awh output file */
249 std::vector
<std::string
> makeLegend(const AwhBiasParams
* awhBiasParams
,
250 OutputFileType outputFileType
,
253 const std::array
<std::string
, maxAwhGraphs
> legendBase
= {
254 { "PMF", "Coord bias", "Coord distr", "Ref value distr", "Target ref value distr",
258 std::vector
<std::string
> legend
;
259 /* Give legends to dimensions higher than the first */
260 for (int d
= 1; d
< awhBiasParams
->ndim
; d
++)
262 legend
.push_back(gmx::formatString("dim%d", d
));
265 switch (outputFileType
)
267 case OutputFileType::Awh
:
269 /* Add as many legends as possible from the "base" legend list */
270 size_t legendBaseIndex
= 0;
271 while (legend
.size() < numLegend
&& legendBaseIndex
< legendBase
.size())
273 legend
.push_back(legendBase
[legendBaseIndex
]);
278 case OutputFileType::Friction
:
279 for (int i0
= 0; i0
< awhBiasParams
->ndim
; i0
++)
281 for (int i1
= 0; i1
<= i0
; i1
++)
283 legend
.push_back(gmx::formatString("%d,%d", i0
, i1
));
289 GMX_RELEASE_ASSERT(legend
.size() == numLegend
,
290 "The number of legends requested for printing and the number generated "
291 "should be the same");
298 OutputFile::OutputFile(const std::string
& filename
, const std::string
& baseTitle
, int numBias
, int biasIndex
) :
300 firstGraphSubBlock_(0),
302 useKTForEnergy_(false)
305 baseFilename_
= filename
.substr(0, filename
.find('.'));
309 baseFilename_
+= gmx::formatString("%d", biasIndex
+ 1);
310 title_
+= gmx::formatString(" %d", biasIndex
+ 1);
314 void OutputFile::initializeAwhOutputFile(int subblockStart
,
316 const AwhBiasParams
* awhBiasParams
,
317 AwhGraphSelection graphSelection
,
318 EnergyUnit energyUnit
,
321 /* The first subblock with actual graph y-values is index 1 + ndim */
322 numDim_
= awhBiasParams
->ndim
;
323 firstGraphSubBlock_
= subblockStart
+ 1 + numDim_
;
324 if (graphSelection
== AwhGraphSelection::All
)
326 /* There are one metadata and ndim coordinate blocks */
327 numGraph_
= std::min(numSubBlocks
- 1 - numDim_
, maxAwhGraphs
);
333 useKTForEnergy_
= (energyUnit
== EnergyUnit::KT
);
334 scaleFactor_
.resize(numGraph_
, 1);
335 if (!useKTForEnergy_
)
337 /* The first two graphs are in units of energy, multiply by kT */
338 std::fill(scaleFactor_
.begin(), scaleFactor_
.begin() + std::min(2, numGraph_
), kTValue
);
340 int numLegend
= numDim_
- 1 + numGraph_
;
341 legend_
= makeLegend(awhBiasParams
, OutputFileType::Awh
, numLegend
);
342 /* We could have both length and angle coordinates in a single bias */
343 xLabel_
= "(nm, deg or lambda state)";
344 yLabel_
= useKTForEnergy_
? "(k\\sB\\NT)" : "(kJ/mol)";
345 if (graphSelection
== AwhGraphSelection::All
)
347 yLabel_
+= gmx::formatString(", (nm\\S-%d\\N or rad\\S-%d\\N), (-)", awhBiasParams
->ndim
,
348 awhBiasParams
->ndim
);
352 /*! \brief Initializes the output file setup for the fricion output (note that the filename is not set here). */
353 void OutputFile::initializeFrictionOutputFile(int subBlockStart
,
355 const AwhBiasParams
* awhBiasParams
,
356 EnergyUnit energyUnit
,
359 /* The first subblock with actual graph y-values is index 1 + ndim */
360 numDim_
= awhBiasParams
->ndim
;
361 int numTensorElements
= (numDim_
* (numDim_
+ 1)) / 2;
363 /* The friction tensor elements are always the last subblocks */
364 if (numSubBlocks
< 1 + numDim_
+ maxAwhGraphs
+ numTensorElements
)
367 "You requested friction tensor output, but the AWH data in the energy file does "
368 "not contain the friction tensor");
370 GMX_ASSERT(numSubBlocks
== 1 + numDim_
+ maxAwhGraphs
+ numTensorElements
,
371 "The number of sub-blocks per bias should be 1 + ndim + maxAwhGraphs + (ndim*(ndim "
374 firstGraphSubBlock_
= subBlockStart
+ numSubBlocks
- numTensorElements
;
375 numGraph_
= numTensorElements
;
376 useKTForEnergy_
= (energyUnit
== EnergyUnit::KT
);
377 scaleFactor_
.resize(numGraph_
, useKTForEnergy_
? 1 : kTValue
);
378 int numLegend
= numDim_
- 1 + numGraph_
;
379 legend_
= makeLegend(awhBiasParams
, OutputFileType::Friction
, numLegend
);
380 xLabel_
= "(nm, deg or lambda state)";
383 yLabel_
= "friction/k\\sB\\NT (nm\\S-2\\Nps, rad\\S-2\\Nps or ps)";
388 "friction (kJ mol\\S-1\\Nnm\\S-2\\Nps, kJ mol\\S-1\\Nrad\\S-2\\Nps or kJ "
393 AwhReader::AwhReader(const AwhParams
* awhParams
,
395 const t_filenm
* filenames
,
396 AwhGraphSelection awhGraphSelection
,
397 EnergyUnit energyUnit
,
399 const t_enxblock
* block
) :
402 GMX_RELEASE_ASSERT(block
!= nullptr, "NULL pointer passed to initAwhReader");
404 bool outputFriction
= opt2bSet("-fric", numFileOptions
, filenames
);
406 /* The first subblock of each AWH block has metadata about
407 * the number of subblocks belonging to that AWH block.
408 * (block->nsub gives only the total number of subblocks and not how
409 * they are distributed between the AWH biases).
412 /* Keep track of the first subblock of this AWH */
413 int subblockStart
= 0;
414 for (int k
= 0; k
< awhParams
->numBias
; k
++)
416 AwhBiasParams
* awhBiasParams
= &awhParams
->awhBiasParams
[k
];
418 int numSubBlocks
= static_cast<int>(block
->sub
[subblockStart
].fval
[0]);
420 std::unique_ptr
<OutputFile
> awhOutputFile(new OutputFile(
421 opt2fn("-o", numFileOptions
, filenames
), "AWH", awhParams
->numBias
, k
));
423 awhOutputFile
->initializeAwhOutputFile(subblockStart
, numSubBlocks
, awhBiasParams
,
424 awhGraphSelection
, energyUnit
, kT
);
426 std::unique_ptr
<OutputFile
> frictionOutputFile
;
429 frictionOutputFile
= std::make_unique
<OutputFile
>(
430 opt2fn("-fric", numFileOptions
, filenames
), "Friction tensor", awhParams
->numBias
, k
);
432 frictionOutputFile
->initializeFrictionOutputFile(subblockStart
, numSubBlocks
,
433 awhBiasParams
, energyUnit
, kT
);
436 biasOutputSetups_
.emplace_back(BiasOutputSetup(subblockStart
, std::move(awhOutputFile
),
437 std::move(frictionOutputFile
)));
439 subblockStart
+= numSubBlocks
;
443 FILE* OutputFile::openBiasOutputFile(double time
, const gmx_output_env_t
* oenv
) const
445 std::string filename
= baseFilename_
+ gmx::formatString("_t%g.xvg", time
);
447 FILE* fp
= xvgropen(filename
.c_str(), title_
.c_str(), xLabel_
, yLabel_
, oenv
);
448 xvgrLegend(fp
, legend_
, oenv
);
453 /*! \brief Prints data selected by \p outputFile from \p block to \p fp */
454 void OutputFile::writeData(const t_enxblock
& block
, int subBlockStart
, FILE* fp
) const
456 int numPoints
= block
.sub
[subBlockStart
+ 1].nr
;
457 for (int j
= 0; j
< numPoints
; j
++)
459 /* Print the coordinates for numDim dimensions */
460 for (int d
= 0; d
< numDim_
; d
++)
462 fprintf(fp
, " %8.4f", block
.sub
[subBlockStart
+ 1 + d
].fval
[j
]);
465 /* Print numGraph observables */
466 for (int i
= 0; i
< numGraph_
; i
++)
468 fprintf(fp
, " %g", block
.sub
[firstGraphSubBlock_
+ i
].fval
[j
] * scaleFactor_
[i
]);
475 void AwhReader::processAwhFrame(const t_enxblock
& block
, double time
, const gmx_output_env_t
* oenv
) const
477 /* We look for AWH data every energy frame and count the no of AWH frames found. We only extract every 'skip' AWH frame. */
479 for (const auto& setup
: biasOutputSetups_
)
481 const int subStart
= setup
.subblockStart();
483 /* Each frame and AWH instance extracted generates one xvg file. */
485 const OutputFile
& awhOutputFile
= setup
.awhOutputFile();
487 FILE* fpAwh
= awhOutputFile
.openBiasOutputFile(time
, oenv
);
489 /* Now do the actual printing. Metadata in first subblock is treated separately. */
490 fprintf(fpAwh
, "# AWH metadata: target error = %.2f kT = %.2f kJ/mol\n",
491 block
.sub
[subStart
].fval
[1], block
.sub
[subStart
].fval
[1] * kT_
);
493 fprintf(fpAwh
, "# AWH metadata: log sample weight = %4.2f\n", block
.sub
[subStart
].fval
[2]);
495 awhOutputFile
.writeData(block
, subStart
, fpAwh
);
500 const OutputFile
* frictionOutputFile
= setup
.frictionOutputFile();
501 if (frictionOutputFile
!= nullptr)
503 FILE* fpFriction
= frictionOutputFile
->openBiasOutputFile(time
, oenv
);
505 frictionOutputFile
->writeData(block
, subStart
, fpFriction
);
507 gmx_ffclose(fpFriction
);
512 /*! \brief The main function for the AWH tool */
513 int gmx_awh(int argc
, char* argv
[])
515 const char* desc
[] = { "[THISMODULE] extracts AWH data from an energy file.",
516 "One or two files are written per AWH bias per time frame.",
517 "The bias index, if more than one, is appended to the file, as well as",
518 "the time of the frame. By default only the PMF is printed.",
519 "With [TT]-more[tt] the bias, target and coordinate distributions",
521 "With [TT]-more[tt] the bias, target and coordinate distributions",
522 "are also printed, as well as the metric sqrt(det(friction_tensor))",
523 "normalized such that the average is 1.",
524 "Option [TT]-fric[tt] prints all components of the friction tensor",
525 "to an additional set of files." };
526 static gmx_bool moreGraphs
= FALSE
;
528 static gmx_bool kTUnit
= FALSE
;
530 { "-skip", FALSE
, etINT
, { &skip
}, "Skip number of frames between data points" },
531 { "-more", FALSE
, etBOOL
, { &moreGraphs
}, "Print more output" },
536 "Print free energy output in units of kT instead of kJ/mol" }
541 gmx_enxnm_t
* enm
= nullptr;
544 gmx_output_env_t
* oenv
;
546 t_filenm fnm
[] = { { efEDR
, "-f", nullptr, ffREAD
},
547 { efTPR
, "-s", nullptr, ffREAD
},
548 { efXVG
, "-o", "awh", ffWRITE
},
549 { efXVG
, "-fric", "friction", ffOPTWR
} };
550 const int nfile
= asize(fnm
);
551 if (!parse_common_args(&argc
, argv
, PCA_CAN_VIEW
| PCA_CAN_BEGIN
| PCA_CAN_END
, nfile
, fnm
,
552 asize(pa
), pa
, asize(desc
), desc
, 0, nullptr, &oenv
))
558 fp
= open_enx(ftp2fn(efEDR
, nfile
, fnm
), "r");
559 do_enxnms(fp
, &nre
, &enm
);
561 /* We just need the AWH parameters from inputrec. These are used to initialize
562 the AWH reader when we have a frame to read later on. */
566 read_tpx(ftp2fn(efTPR
, nfile
, fnm
), &ir
, box
, &natoms
, nullptr, nullptr, &mtop
);
570 gmx_fatal(FARGS
, "No AWH data in %s\n", opt2fn("-f", nfile
, fnm
));
573 std::unique_ptr
<AwhReader
> awhReader
;
575 /* Initiate counters */
577 int awhFrameCounter
= 0;
581 haveFrame
= do_enx(fp
, frame
);
583 bool useFrame
= false;
585 t_enxblock
* block
= nullptr;
589 timeCheck
= check_times(frame
->t
);
593 block
= find_block_id_enxframe(frame
, enxAWH
, nullptr);
595 if (block
!= nullptr)
597 if ((skip
== 0) || (awhFrameCounter
% skip
== 0))
608 /* We read a valid frame with an AWH block, so we can use it */
610 /* Currently we have the number of subblocks per AWH stored
611 * in the first subblock (we cannot get this directly from the tpr),
612 * so we have to read an AWH block before making the legends.
614 if (awhReader
== nullptr)
616 AwhGraphSelection awhGraphSelection
=
617 (moreGraphs
? AwhGraphSelection::All
: AwhGraphSelection::Pmf
);
618 EnergyUnit energyUnit
= (kTUnit
? EnergyUnit::KT
: EnergyUnit::KJPerMol
);
619 awhReader
= std::make_unique
<AwhReader
>(ir
.awhParams
, nfile
, fnm
, awhGraphSelection
,
620 energyUnit
, BOLTZ
* ir
.opts
.ref_t
[0], block
);
623 awhReader
->processAwhFrame(*block
, frame
->t
, oenv
);
625 } while (haveFrame
&& (timeCheck
<= 0));
627 fprintf(stderr
, "\n");