2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 #include "gromacs/utility/arrayref.h"
41 #include "mdebin_bar.h"
50 #include "gromacs/fileio/enxio.h"
51 #include "gromacs/mdlib/energyoutput.h"
52 #include "gromacs/mdtypes/energyhistory.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/trajectory/energyframe.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/smalloc.h"
60 /* Number of entries in subblock_d preceding the lambda components */
61 constexpr int c_subblockDNumPreEntries = 5;
62 /* Number of entries in subblock_i preceding the lambda components */
63 constexpr int c_subblockINumPreEntries = 2;
65 /* reset the delta_h list to prepare it for new values */
66 static void mde_delta_h_reset(t_mde_delta_h* dh)
72 /* initialize the delta_h list */
73 static void mde_delta_h_init(t_mde_delta_h* dh,
85 dh->derivative = derivative;
86 dh->nlambda = nlambda;
88 snew(dh->lambda, nlambda);
89 for (i = 0; i < nlambda; i++)
92 dh->lambda[i] = lambda[i];
96 snew(dh->subblock_meta_d, dh->nlambda + 1);
98 dh->ndhmax = ndhmax + 2;
99 for (i = 0; i < 2; i++)
101 dh->bin[i] = nullptr;
104 snew(dh->dh, dh->ndhmax);
105 snew(dh->dhf, dh->ndhmax);
107 if (nbins <= 0 || dx < GMX_REAL_EPS * 10)
114 /* pre-allocate the histogram */
115 dh->nhist = 2; /* energies and derivatives histogram */
118 for (i = 0; i < dh->nhist; i++)
120 snew(dh->bin[i], dh->nbins);
123 mde_delta_h_reset(dh);
126 static void done_mde_delta_h(t_mde_delta_h* dh)
129 sfree(dh->subblock_meta_d);
132 for (int i = 0; i < dh->nhist; i++)
138 /* Add a value to the delta_h list */
139 static void mde_delta_h_add_dh(t_mde_delta_h* dh, double delta_h)
141 if (dh->ndh >= dh->ndhmax)
143 gmx_incons("delta_h array not big enough!");
145 dh->dh[dh->ndh] = delta_h;
149 /* construct histogram with index hi */
150 static void mde_delta_h_make_hist(t_mde_delta_h* dh, int hi, gmx_bool invert)
152 double min_dh = FLT_MAX;
153 double max_dh = -FLT_MAX;
155 double max_dh_hist; /* maximum binnable dh value */
156 double min_dh_hist; /* minimum binnable dh value */
158 double f; /* energy mult. factor */
160 /* by applying a -1 scaling factor on the energies we get the same as
161 having a negative dx, but we don't need to fix the min/max values
162 beyond inverting x0 */
165 /* first find min and max */
166 for (i = 0; i < dh->ndh; i++)
168 if (f * dh->dh[i] < min_dh)
170 min_dh = f * dh->dh[i];
172 if (f * dh->dh[i] > max_dh)
174 max_dh = f * dh->dh[i];
178 /* reset the histogram */
179 for (i = 0; i < dh->nbins; i++)
185 /* The starting point of the histogram is the lowest value found:
186 that value has the highest contribution to the free energy.
188 Get this start value in number of histogram dxs from zero,
191 dh->x0[hi] = static_cast<int64_t>(floor(min_dh / dx));
193 min_dh_hist = (dh->x0[hi]) * dx;
194 max_dh_hist = (dh->x0[hi] + dh->nbins + 1) * dx;
196 /* and fill the histogram*/
197 for (i = 0; i < dh->ndh; i++)
201 /* Determine the bin number. If it doesn't fit into the histogram,
202 add it to the last bin.
203 We check the max_dh_int range because converting to integers
204 might lead to overflow with unpredictable results.*/
205 if ((f * dh->dh[i] >= min_dh_hist) && (f * dh->dh[i] <= max_dh_hist))
207 bin = static_cast<unsigned int>((f * dh->dh[i] - min_dh_hist) / dx);
213 /* double-check here because of possible round-off errors*/
214 if (bin >= dh->nbins)
218 if (bin > dh->maxbin[hi])
220 dh->maxbin[hi] = bin;
226 /* make sure we include a bin with 0 if we didn't use the full
227 histogram width. This can then be used as an indication that
228 all the data was binned. */
229 if (dh->maxbin[hi] < dh->nbins - 1)
236 static void mde_delta_h_handle_block(t_mde_delta_h* dh, t_enxblock* blk)
238 /* first check which type we should use: histogram or raw data */
243 /* We write raw data.
244 Raw data consists of 3 subblocks: an int metadata block
245 with type and derivative index, a foreign lambda block
246 and and the data itself */
247 add_subblocks_enxblock(blk, 3);
252 dh->subblock_meta_i[0] = dh->type; /* block data type */
253 dh->subblock_meta_i[1] = dh->derivative; /* derivative direction if
254 applicable (in indices
255 starting from first coord in
256 the main delta_h_coll) */
258 blk->sub[0].type = xdr_datatype_int;
259 blk->sub[0].ival = dh->subblock_meta_i;
262 for (i = 0; i < dh->nlambda; i++)
264 dh->subblock_meta_d[i] = dh->lambda[i];
266 blk->sub[1].nr = dh->nlambda;
267 blk->sub[1].type = xdr_datatype_double;
268 blk->sub[1].dval = dh->subblock_meta_d;
271 /* check if there's actual data to be written. */
277 blk->sub[2].nr = dh->ndh;
278 /* Michael commented in 2012 that this use of explicit
279 xdr_datatype_float was good for F@H for now.
280 Apparently it's still good enough. */
281 blk->sub[2].type = xdr_datatype_float;
282 for (i = 0; i < dh->ndh; i++)
284 dh->dhf[i] = static_cast<float>(dh->dh[i]);
286 blk->sub[2].fval = dh->dhf;
292 blk->sub[2].type = xdr_datatype_float;
293 blk->sub[2].fval = nullptr;
298 int nhist_written = 0;
302 /* TODO histogram metadata */
303 /* check if there's actual data to be written. */
306 gmx_bool prev_complete = FALSE;
307 /* Make the histogram(s) */
308 for (i = 0; i < dh->nhist; i++)
312 /* the first histogram is always normal, and the
313 second one is always reverse */
314 mde_delta_h_make_hist(dh, i, i == 1);
316 /* check whether this histogram contains all data: if the
317 last bin is 0, it does */
318 if (dh->bin[i][dh->nbins - 1] == 0)
320 prev_complete = TRUE;
324 prev_complete = TRUE;
331 /* A histogram consists of 2, 3 or 4 subblocks:
332 the foreign lambda value + histogram spacing, the starting point,
333 and the histogram data (0, 1 or 2 blocks). */
334 add_subblocks_enxblock(blk, nhist_written + 2);
337 /* subblock 1: the lambda value + the histogram spacing */
338 if (dh->nlambda == 1)
340 /* for backward compatibility */
341 dh->subblock_meta_d[0] = dh->lambda[0];
345 dh->subblock_meta_d[0] = -1;
346 for (i = 0; i < dh->nlambda; i++)
348 dh->subblock_meta_d[2 + i] = dh->lambda[i];
351 dh->subblock_meta_d[1] = dh->dx;
352 blk->sub[0].nr = 2 + ((dh->nlambda > 1) ? dh->nlambda : 0);
353 blk->sub[0].type = xdr_datatype_double;
354 blk->sub[0].dval = dh->subblock_meta_d;
356 /* subblock 2: the starting point(s) as a long integer */
357 dh->subblock_meta_l[0] = nhist_written;
358 dh->subblock_meta_l[1] = dh->type; /*dh->derivative ? 1 : 0;*/
360 for (i = 0; i < nhist_written; i++)
362 dh->subblock_meta_l[k++] = dh->x0[i];
364 /* append the derivative data */
365 dh->subblock_meta_l[k++] = dh->derivative;
367 blk->sub[1].nr = nhist_written + 3;
368 blk->sub[1].type = xdr_datatype_int64;
369 blk->sub[1].lval = dh->subblock_meta_l;
371 /* subblock 3 + 4 : the histogram data */
372 for (i = 0; i < nhist_written; i++)
374 blk->sub[i + 2].nr = dh->maxbin[i] + 1; /* it's +1 because size=index+1
376 blk->sub[i + 2].type = xdr_datatype_int;
377 blk->sub[i + 2].ival = dh->bin[i];
382 /* initialize the collection*/
383 void mde_delta_h_coll_init(t_mde_delta_h_coll* dhc, const t_inputrec& inputrec)
387 int ndhmax = inputrec.nstenergy / inputrec.nstcalcenergy;
388 t_lambda* fep = inputrec.fepvals.get();
390 dhc->temperature = inputrec.opts.ref_t[0]; /* only store system temperature */
391 dhc->start_time = 0.;
392 dhc->delta_time = inputrec.delta_t * inputrec.fepvals->nstdhdl;
393 dhc->start_time_set = FALSE;
395 /* this is the compatibility lambda value. If it is >=0, it is valid,
396 and there is either an old-style lambda or a slow growth simulation. */
397 dhc->start_lambda = inputrec.fepvals->init_lambda;
398 /* for continuous change of lambda values */
399 dhc->delta_lambda = inputrec.fepvals->delta_lambda * inputrec.fepvals->nstdhdl;
401 if (dhc->start_lambda < 0)
403 /* create the native lambda vectors */
404 dhc->lambda_index = fep->init_fep_state;
405 dhc->n_lambda_vec = 0;
406 for (auto i : keysOf(fep->separate_dvdl))
408 if (fep->separate_dvdl[i])
413 snew(dhc->native_lambda_vec, dhc->n_lambda_vec);
414 snew(dhc->native_lambda_components, dhc->n_lambda_vec);
416 for (auto i : keysOf(fep->separate_dvdl))
418 if (fep->separate_dvdl[i])
420 dhc->native_lambda_components[j] = static_cast<int>(i);
421 if (fep->init_fep_state >= 0 && fep->init_fep_state < fep->n_lambda)
423 dhc->native_lambda_vec[j] = fep->all_lambda[i][fep->init_fep_state];
427 dhc->native_lambda_vec[j] = -1;
435 /* don't allocate the meta-data subblocks for lambda vectors */
436 dhc->native_lambda_vec = nullptr;
437 dhc->n_lambda_vec = 0;
438 dhc->native_lambda_components = nullptr;
439 dhc->lambda_index = -1;
441 /* allocate metadata subblocks */
442 snew(dhc->subblock_d, c_subblockDNumPreEntries + dhc->n_lambda_vec);
443 snew(dhc->subblock_i, c_subblockINumPreEntries + dhc->n_lambda_vec);
445 /* now decide which data to write out */
448 dhc->dh_expanded = nullptr;
449 dhc->dh_energy = nullptr;
450 dhc->dh_pv = nullptr;
452 /* total number of raw data point collections in the sample */
456 gmx_bool bExpanded = FALSE;
457 gmx_bool bEnergy = FALSE;
458 gmx_bool bPV = FALSE;
459 int n_lambda_components = 0;
461 /* first count the number of states */
464 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
466 for (auto i : keysOf(fep->separate_dvdl))
468 if (inputrec.fepvals->separate_dvdl[i])
475 /* add the lambdas */
476 dhc->nlambda = inputrec.fepvals->lambda_stop_n - inputrec.fepvals->lambda_start_n;
477 dhc->ndh += dhc->nlambda;
478 /* another compatibility check */
479 if (dhc->start_lambda < 0)
481 /* include one more for the specification of the state, by lambda or
483 if (inputrec.expandedvals->elmcmove > LambdaMoveCalculation::No)
488 /* whether to print energies */
489 if (inputrec.fepvals->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
494 if (inputrec.epc > PressureCoupling::No)
496 dhc->ndh += 1; /* include pressure-volume work */
501 snew(dhc->dh, dhc->ndh);
503 /* now initialize them */
504 /* the order, for now, must match that of the dhdl.xvg file because of
505 how g_energy -odh is implemented */
509 dhc->dh_expanded = dhc->dh + n;
510 mde_delta_h_init(dhc->dh + n,
511 inputrec.fepvals->dh_hist_size,
512 inputrec.fepvals->dh_hist_spacing,
522 dhc->dh_energy = dhc->dh + n;
523 mde_delta_h_init(dhc->dh + n,
524 inputrec.fepvals->dh_hist_size,
525 inputrec.fepvals->dh_hist_spacing,
534 n_lambda_components = 0;
535 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
537 dhc->dh_dhdl = dhc->dh + n;
538 for (auto i : keysOf(fep->separate_dvdl))
540 if (inputrec.fepvals->separate_dvdl[i])
542 /* we give it init_lambda for compatibility */
543 mde_delta_h_init(dhc->dh + n,
544 inputrec.fepvals->dh_hist_size,
545 inputrec.fepvals->dh_hist_spacing,
550 &(fep->init_lambda));
552 n_lambda_components++;
558 for (auto i : keysOf(fep->separate_dvdl))
560 if (fep->separate_dvdl[i])
562 n_lambda_components++; /* count the components */
566 /* add the lambdas */
567 dhc->dh_du = dhc->dh + n;
568 snew(lambda_vec, n_lambda_components);
569 for (i = inputrec.fepvals->lambda_start_n; i < inputrec.fepvals->lambda_stop_n; i++)
573 for (auto j : keysOf(fep->separate_dvdl))
575 if (fep->separate_dvdl[j])
577 lambda_vec[k++] = fep->all_lambda[j][i];
581 mde_delta_h_init(dhc->dh + n,
582 inputrec.fepvals->dh_hist_size,
583 inputrec.fepvals->dh_hist_spacing,
594 dhc->dh_pv = dhc->dh + n;
595 mde_delta_h_init(dhc->dh + n,
596 inputrec.fepvals->dh_hist_size,
597 inputrec.fepvals->dh_hist_spacing,
608 void done_mde_delta_h_coll(t_mde_delta_h_coll* dhc)
614 sfree(dhc->native_lambda_vec);
615 sfree(dhc->native_lambda_components);
616 sfree(dhc->subblock_d);
617 sfree(dhc->subblock_i);
618 for (int i = 0; i < dhc->ndh; ++i)
620 done_mde_delta_h(&dhc->dh[i]);
626 /* add a bunch of samples - note fep_state is double to allow for better data storage */
627 void mde_delta_h_coll_add_dh(t_mde_delta_h_coll* dhc,
631 gmx::ArrayRef<double> dhdl,
637 if (!dhc->start_time_set)
639 dhc->start_time_set = TRUE;
640 dhc->start_time = time;
643 for (i = 0; i < dhc->ndhdl; i++)
645 mde_delta_h_add_dh(dhc->dh_dhdl + i, dhdl[i]);
647 for (i = 0; i < dhc->nlambda; i++)
649 mde_delta_h_add_dh(dhc->dh_du + i, foreign_dU[i]);
651 if (dhc->dh_pv != nullptr)
653 mde_delta_h_add_dh(dhc->dh_pv, pV);
655 if (dhc->dh_energy != nullptr)
657 mde_delta_h_add_dh(dhc->dh_energy, energy);
659 if (dhc->dh_expanded != nullptr)
661 mde_delta_h_add_dh(dhc->dh_expanded, fep_state);
665 /* write the metadata associated with all the du blocks, and call
666 handle_block to write out all the du blocks */
667 void mde_delta_h_coll_handle_block(t_mde_delta_h_coll* dhc, t_enxframe* fr, int nblock)
672 /* add one block with one subblock as the collection's own data */
674 add_blocks_enxframe(fr, nblock);
675 blk = fr->block + (nblock - 1);
677 /* only allocate lambda vector component blocks if they must be written out
678 for backward compatibility */
679 if (dhc->native_lambda_components != nullptr)
681 add_subblocks_enxblock(blk, 2);
685 add_subblocks_enxblock(blk, 1);
688 dhc->subblock_d[0] = dhc->temperature; /* temperature */
689 dhc->subblock_d[1] = dhc->start_time; /* time of first sample */
690 dhc->subblock_d[2] = dhc->delta_time; /* time difference between samples */
691 dhc->subblock_d[3] = dhc->start_lambda; /* old-style lambda at starttime */
692 dhc->subblock_d[4] = dhc->delta_lambda; /* lambda diff. between samples */
693 /* set the lambda vector components if they exist */
694 if (dhc->native_lambda_components != nullptr)
696 for (i = 0; i < dhc->n_lambda_vec; i++)
698 dhc->subblock_d[c_subblockDNumPreEntries + i] = dhc->native_lambda_vec[i];
702 blk->sub[0].nr = c_subblockDNumPreEntries + dhc->n_lambda_vec;
703 blk->sub[0].type = xdr_datatype_double;
704 blk->sub[0].dval = dhc->subblock_d;
706 if (dhc->native_lambda_components != nullptr)
708 dhc->subblock_i[0] = dhc->lambda_index;
709 /* set the lambda vector component IDs if they exist */
710 dhc->subblock_i[1] = dhc->n_lambda_vec;
711 for (i = 0; i < dhc->n_lambda_vec; i++)
713 dhc->subblock_i[c_subblockINumPreEntries + i] = dhc->native_lambda_components[i];
715 blk->sub[1].nr = c_subblockINumPreEntries + dhc->n_lambda_vec;
716 blk->sub[1].type = xdr_datatype_int;
717 blk->sub[1].ival = dhc->subblock_i;
720 for (i = 0; i < dhc->ndh; i++)
723 add_blocks_enxframe(fr, nblock);
724 blk = fr->block + (nblock - 1);
726 mde_delta_h_handle_block(dhc->dh + i, blk);
730 /* reset the data for a new round */
731 void mde_delta_h_coll_reset(t_mde_delta_h_coll* dhc)
734 for (i = 0; i < dhc->ndh; i++)
736 if (dhc->dh[i].written)
738 /* we can now throw away the data */
739 mde_delta_h_reset(dhc->dh + i);
742 dhc->start_time_set = FALSE;
745 /* set the energyhistory variables to save state */
746 void mde_delta_h_coll_update_energyhistory(const t_mde_delta_h_coll* dhc, energyhistory_t* enerhist)
748 if (enerhist->deltaHForeignLambdas == nullptr)
750 enerhist->deltaHForeignLambdas = std::make_unique<delta_h_history_t>();
751 enerhist->deltaHForeignLambdas->dh.resize(dhc->ndh);
754 delta_h_history_t* const deltaH = enerhist->deltaHForeignLambdas.get();
757 deltaH->dh.size() == static_cast<size_t>(dhc->ndh),
758 "energy history number of delta_h histograms should match inputrec's number");
760 for (int i = 0; i < dhc->ndh; i++)
762 std::vector<real>& dh = deltaH->dh[i];
763 dh.resize(dhc->dh[i].ndh);
764 std::copy(dhc->dh[i].dh, dhc->dh[i].dh + dhc->dh[i].ndh, dh.begin());
766 deltaH->start_time = dhc->start_time;
767 deltaH->start_lambda = dhc->start_lambda;
771 /* restore the variables from an energyhistory */
772 void mde_delta_h_coll_restore_energyhistory(t_mde_delta_h_coll* dhc, const delta_h_history_t* deltaH)
774 GMX_RELEASE_ASSERT(dhc, "Should have delta_h histograms");
775 GMX_RELEASE_ASSERT(deltaH, "Should have delta_h histograms in energy history");
777 deltaH->dh.size() == static_cast<size_t>(dhc->ndh),
778 "energy history number of delta_h histograms should match inputrec's number");
780 for (gmx::index i = 0; i < gmx::ssize(deltaH->dh); i++)
782 dhc->dh[i].ndh = deltaH->dh[i].size();
783 for (unsigned int j = 0; j < dhc->dh[i].ndh; j++)
785 dhc->dh[i].dh[j] = deltaH->dh[i][j];
788 dhc->start_time = deltaH->start_time;
789 if (deltaH->start_lambda_set)
791 dhc->start_lambda = deltaH->start_lambda;
793 dhc->start_time_set = (dhc->dh[0].ndh > 0);