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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "mdebin_bar.h"
47 #include "gromacs/fileio/enxio.h"
48 #include "gromacs/mdlib/mdebin.h"
49 #include "gromacs/mdtypes/energyhistory.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/smalloc.h"
56 /* reset the delta_h list to prepare it for new values */
57 static void mde_delta_h_reset(t_mde_delta_h *dh)
63 /* initialize the delta_h list */
64 static void mde_delta_h_init(t_mde_delta_h *dh, int nbins,
65 double dx, unsigned int ndhmax,
66 int type, int derivative, int nlambda,
72 dh->derivative = derivative;
73 dh->nlambda = nlambda;
75 snew(dh->lambda, nlambda);
76 for (i = 0; i < nlambda; i++)
79 dh->lambda[i] = lambda[i];
83 snew(dh->subblock_meta_d, dh->nlambda+1);
85 dh->ndhmax = ndhmax+2;
86 for (i = 0; i < 2; i++)
91 snew(dh->dh, dh->ndhmax);
92 snew(dh->dhf, dh->ndhmax);
94 if (nbins <= 0 || dx < GMX_REAL_EPS*10)
101 /* pre-allocate the histogram */
102 dh->nhist = 2; /* energies and derivatives histogram */
105 for (i = 0; i < dh->nhist; i++)
107 snew(dh->bin[i], dh->nbins);
110 mde_delta_h_reset(dh);
113 /* Add a value to the delta_h list */
114 static void mde_delta_h_add_dh(t_mde_delta_h *dh, double delta_h)
116 if (dh->ndh >= dh->ndhmax)
118 gmx_incons("delta_h array not big enough!");
120 dh->dh[dh->ndh] = delta_h;
124 /* construct histogram with index hi */
125 static void mde_delta_h_make_hist(t_mde_delta_h *dh, int hi, gmx_bool invert)
127 double min_dh = FLT_MAX;
128 double max_dh = -FLT_MAX;
130 double max_dh_hist; /* maximum binnable dh value */
131 double min_dh_hist; /* minimum binnable dh value */
133 double f; /* energy mult. factor */
135 /* by applying a -1 scaling factor on the energies we get the same as
136 having a negative dx, but we don't need to fix the min/max values
137 beyond inverting x0 */
140 /* first find min and max */
141 for (i = 0; i < dh->ndh; i++)
143 if (f*dh->dh[i] < min_dh)
145 min_dh = f*dh->dh[i];
147 if (f*dh->dh[i] > max_dh)
149 max_dh = f*dh->dh[i];
153 /* reset the histogram */
154 for (i = 0; i < dh->nbins; i++)
160 /* The starting point of the histogram is the lowest value found:
161 that value has the highest contribution to the free energy.
163 Get this start value in number of histogram dxs from zero,
166 dh->x0[hi] = (gmx_int64_t)floor(min_dh/dx);
168 min_dh_hist = (dh->x0[hi])*dx;
169 max_dh_hist = (dh->x0[hi] + dh->nbins + 1)*dx;
171 /* and fill the histogram*/
172 for (i = 0; i < dh->ndh; i++)
176 /* Determine the bin number. If it doesn't fit into the histogram,
177 add it to the last bin.
178 We check the max_dh_int range because converting to integers
179 might lead to overflow with unpredictable results.*/
180 if ( (f*dh->dh[i] >= min_dh_hist) && (f*dh->dh[i] <= max_dh_hist ) )
182 bin = (unsigned int)( (f*dh->dh[i] - min_dh_hist)/dx );
188 /* double-check here because of possible round-off errors*/
189 if (bin >= dh->nbins)
193 if (bin > dh->maxbin[hi])
195 dh->maxbin[hi] = bin;
201 /* make sure we include a bin with 0 if we didn't use the full
202 histogram width. This can then be used as an indication that
203 all the data was binned. */
204 if (dh->maxbin[hi] < dh->nbins-1)
211 void mde_delta_h_handle_block(t_mde_delta_h *dh, t_enxblock *blk)
213 /* first check which type we should use: histogram or raw data */
218 /* We write raw data.
219 Raw data consists of 3 subblocks: an int metadata block
220 with type and derivative index, a foreign lambda block
221 and and the data itself */
222 add_subblocks_enxblock(blk, 3);
227 dh->subblock_meta_i[0] = dh->type; /* block data type */
228 dh->subblock_meta_i[1] = dh->derivative; /* derivative direction if
229 applicable (in indices
230 starting from first coord in
231 the main delta_h_coll) */
233 blk->sub[0].type = xdr_datatype_int;
234 blk->sub[0].ival = dh->subblock_meta_i;
237 for (i = 0; i < dh->nlambda; i++)
239 dh->subblock_meta_d[i] = dh->lambda[i];
241 blk->sub[1].nr = dh->nlambda;
242 blk->sub[1].type = xdr_datatype_double;
243 blk->sub[1].dval = dh->subblock_meta_d;
246 /* check if there's actual data to be written. */
252 blk->sub[2].nr = dh->ndh;
253 /* Michael commented in 2012 that this use of explicit
254 xdr_datatype_float was good for F@H for now.
255 Apparently it's still good enough. */
256 blk->sub[2].type = xdr_datatype_float;
257 for (i = 0; i < dh->ndh; i++)
259 dh->dhf[i] = (float)dh->dh[i];
261 blk->sub[2].fval = dh->dhf;
267 blk->sub[2].type = xdr_datatype_float;
268 blk->sub[2].fval = NULL;
273 int nhist_written = 0;
277 /* TODO histogram metadata */
278 /* check if there's actual data to be written. */
281 gmx_bool prev_complete = FALSE;
282 /* Make the histogram(s) */
283 for (i = 0; i < dh->nhist; i++)
287 /* the first histogram is always normal, and the
288 second one is always reverse */
289 mde_delta_h_make_hist(dh, i, i == 1);
291 /* check whether this histogram contains all data: if the
292 last bin is 0, it does */
293 if (dh->bin[i][dh->nbins-1] == 0)
295 prev_complete = TRUE;
299 prev_complete = TRUE;
306 /* A histogram consists of 2, 3 or 4 subblocks:
307 the foreign lambda value + histogram spacing, the starting point,
308 and the histogram data (0, 1 or 2 blocks). */
309 add_subblocks_enxblock(blk, nhist_written+2);
312 /* subblock 1: the lambda value + the histogram spacing */
313 if (dh->nlambda == 1)
315 /* for backward compatibility */
316 dh->subblock_meta_d[0] = dh->lambda[0];
320 dh->subblock_meta_d[0] = -1;
321 for (i = 0; i < dh->nlambda; i++)
323 dh->subblock_meta_d[2+i] = dh->lambda[i];
326 dh->subblock_meta_d[1] = dh->dx;
327 blk->sub[0].nr = 2+ ((dh->nlambda > 1) ? dh->nlambda : 0);
328 blk->sub[0].type = xdr_datatype_double;
329 blk->sub[0].dval = dh->subblock_meta_d;
331 /* subblock 2: the starting point(s) as a long integer */
332 dh->subblock_meta_l[0] = nhist_written;
333 dh->subblock_meta_l[1] = dh->type; /*dh->derivative ? 1 : 0;*/
335 for (i = 0; i < nhist_written; i++)
337 dh->subblock_meta_l[k++] = dh->x0[i];
339 /* append the derivative data */
340 dh->subblock_meta_l[k++] = dh->derivative;
342 blk->sub[1].nr = nhist_written+3;
343 blk->sub[1].type = xdr_datatype_int64;
344 blk->sub[1].lval = dh->subblock_meta_l;
346 /* subblock 3 + 4 : the histogram data */
347 for (i = 0; i < nhist_written; i++)
349 blk->sub[i+2].nr = dh->maxbin[i]+1; /* it's +1 because size=index+1
351 blk->sub[i+2].type = xdr_datatype_int;
352 blk->sub[i+2].ival = dh->bin[i];
357 /* initialize the collection*/
358 void mde_delta_h_coll_init(t_mde_delta_h_coll *dhc, const t_inputrec *ir)
362 int ndhmax = ir->nstenergy/ir->nstcalcenergy;
363 t_lambda *fep = ir->fepvals;
365 dhc->temperature = ir->opts.ref_t[0]; /* only store system temperature */
366 dhc->start_time = 0.;
367 dhc->delta_time = ir->delta_t*ir->fepvals->nstdhdl;
368 dhc->start_time_set = FALSE;
370 /* this is the compatibility lambda value. If it is >=0, it is valid,
371 and there is either an old-style lambda or a slow growth simulation. */
372 dhc->start_lambda = ir->fepvals->init_lambda;
373 /* for continuous change of lambda values */
374 dhc->delta_lambda = ir->fepvals->delta_lambda*ir->fepvals->nstdhdl;
376 if (dhc->start_lambda < 0)
378 /* create the native lambda vectors */
379 dhc->lambda_index = fep->init_fep_state;
380 dhc->n_lambda_vec = 0;
381 for (i = 0; i < efptNR; i++)
383 if (fep->separate_dvdl[i])
388 snew(dhc->native_lambda_vec, dhc->n_lambda_vec);
389 snew(dhc->native_lambda_components, dhc->n_lambda_vec);
391 for (i = 0; i < efptNR; i++)
393 if (fep->separate_dvdl[i])
395 dhc->native_lambda_components[j] = i;
396 if (fep->init_fep_state >= 0 &&
397 fep->init_fep_state < fep->n_lambda)
399 dhc->native_lambda_vec[j] =
400 fep->all_lambda[i][fep->init_fep_state];
404 dhc->native_lambda_vec[j] = -1;
412 /* don't allocate the meta-data subblocks for lambda vectors */
413 dhc->native_lambda_vec = NULL;
414 dhc->n_lambda_vec = 0;
415 dhc->native_lambda_components = 0;
416 dhc->lambda_index = -1;
418 /* allocate metadata subblocks */
419 snew(dhc->subblock_d, 5 + dhc->n_lambda_vec);
420 snew(dhc->subblock_i, 1 + dhc->n_lambda_vec);
422 /* now decide which data to write out */
425 dhc->dh_expanded = NULL;
426 dhc->dh_energy = NULL;
429 /* total number of raw data point collections in the sample */
433 gmx_bool bExpanded = FALSE;
434 gmx_bool bEnergy = FALSE;
435 gmx_bool bPV = FALSE;
436 int n_lambda_components = 0;
438 /* first count the number of states */
441 if (fep->dhdl_derivatives == edhdlderivativesYES)
443 for (i = 0; i < efptNR; i++)
445 if (ir->fepvals->separate_dvdl[i])
452 /* add the lambdas */
453 dhc->nlambda = ir->fepvals->lambda_stop_n - ir->fepvals->lambda_start_n;
454 dhc->ndh += dhc->nlambda;
455 /* another compatibility check */
456 if (dhc->start_lambda < 0)
458 /* include one more for the specification of the state, by lambda or
460 if (ir->expandedvals->elmcmove > elmcmoveNO)
465 /* whether to print energies */
466 if (ir->fepvals->edHdLPrintEnergy != edHdLPrintEnergyNO)
473 dhc->ndh += 1; /* include pressure-volume work */
478 snew(dhc->dh, dhc->ndh);
480 /* now initialize them */
481 /* the order, for now, must match that of the dhdl.xvg file because of
482 how g_energy -odh is implemented */
486 dhc->dh_expanded = dhc->dh+n;
487 mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
488 ir->fepvals->dh_hist_spacing, ndhmax,
489 dhbtEXPANDED, 0, 0, NULL);
494 dhc->dh_energy = dhc->dh+n;
495 mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
496 ir->fepvals->dh_hist_spacing, ndhmax,
501 n_lambda_components = 0;
502 if (fep->dhdl_derivatives == edhdlderivativesYES)
504 dhc->dh_dhdl = dhc->dh + n;
505 for (i = 0; i < efptNR; i++)
507 if (ir->fepvals->separate_dvdl[i])
509 /* we give it init_lambda for compatibility */
510 mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
511 ir->fepvals->dh_hist_spacing, ndhmax,
512 dhbtDHDL, n_lambda_components, 1,
513 &(fep->init_lambda));
515 n_lambda_components++;
521 for (i = 0; i < efptNR; i++)
523 if (ir->fepvals->separate_dvdl[i])
525 n_lambda_components++; /* count the components */
530 /* add the lambdas */
531 dhc->dh_du = dhc->dh + n;
532 snew(lambda_vec, n_lambda_components);
533 for (i = ir->fepvals->lambda_start_n; i < ir->fepvals->lambda_stop_n; i++)
537 for (j = 0; j < efptNR; j++)
539 if (ir->fepvals->separate_dvdl[j])
541 lambda_vec[k++] = fep->all_lambda[j][i];
545 mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
546 ir->fepvals->dh_hist_spacing, ndhmax,
547 dhbtDH, 0, n_lambda_components, lambda_vec);
553 dhc->dh_pv = dhc->dh+n;
554 mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
555 ir->fepvals->dh_hist_spacing, ndhmax,
562 /* add a bunch of samples - note fep_state is double to allow for better data storage */
563 void mde_delta_h_coll_add_dh(t_mde_delta_h_coll *dhc,
573 if (!dhc->start_time_set)
575 dhc->start_time_set = TRUE;
576 dhc->start_time = time;
579 for (i = 0; i < dhc->ndhdl; i++)
581 mde_delta_h_add_dh(dhc->dh_dhdl+i, dhdl[i]);
583 for (i = 0; i < dhc->nlambda; i++)
585 mde_delta_h_add_dh(dhc->dh_du+i, foreign_dU[i]);
587 if (dhc->dh_pv != NULL)
589 mde_delta_h_add_dh(dhc->dh_pv, pV);
591 if (dhc->dh_energy != NULL)
593 mde_delta_h_add_dh(dhc->dh_energy, energy);
595 if (dhc->dh_expanded != NULL)
597 mde_delta_h_add_dh(dhc->dh_expanded, fep_state);
602 /* write the metadata associated with all the du blocks, and call
603 handle_block to write out all the du blocks */
604 void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
605 t_enxframe *fr, int nblock)
610 /* add one block with one subblock as the collection's own data */
612 add_blocks_enxframe(fr, nblock);
613 blk = fr->block + (nblock-1);
615 /* only allocate lambda vector component blocks if they must be written out
616 for backward compatibility */
617 if (dhc->native_lambda_components != NULL)
619 add_subblocks_enxblock(blk, 2);
623 add_subblocks_enxblock(blk, 1);
626 dhc->subblock_d[0] = dhc->temperature; /* temperature */
627 dhc->subblock_d[1] = dhc->start_time; /* time of first sample */
628 dhc->subblock_d[2] = dhc->delta_time; /* time difference between samples */
629 dhc->subblock_d[3] = dhc->start_lambda; /* old-style lambda at starttime */
630 dhc->subblock_d[4] = dhc->delta_lambda; /* lambda diff. between samples */
631 /* set the lambda vector components if they exist */
632 if (dhc->native_lambda_components != NULL)
634 for (i = 0; i < dhc->n_lambda_vec; i++)
636 dhc->subblock_d[5+i] = dhc->native_lambda_vec[i];
640 blk->sub[0].nr = 5 + dhc->n_lambda_vec;
641 blk->sub[0].type = xdr_datatype_double;
642 blk->sub[0].dval = dhc->subblock_d;
644 if (dhc->native_lambda_components != NULL)
646 dhc->subblock_i[0] = dhc->lambda_index;
647 /* set the lambda vector component IDs if they exist */
648 dhc->subblock_i[1] = dhc->n_lambda_vec;
649 for (i = 0; i < dhc->n_lambda_vec; i++)
651 dhc->subblock_i[i+2] = dhc->native_lambda_components[i];
653 blk->sub[1].nr = 2 + dhc->n_lambda_vec;
654 blk->sub[1].type = xdr_datatype_int;
655 blk->sub[1].ival = dhc->subblock_i;
658 for (i = 0; i < dhc->ndh; i++)
661 add_blocks_enxframe(fr, nblock);
662 blk = fr->block + (nblock-1);
664 mde_delta_h_handle_block(dhc->dh+i, blk);
668 /* reset the data for a new round */
669 void mde_delta_h_coll_reset(t_mde_delta_h_coll *dhc)
672 for (i = 0; i < dhc->ndh; i++)
674 if (dhc->dh[i].written)
676 /* we can now throw away the data */
677 mde_delta_h_reset(dhc->dh + i);
680 dhc->start_time_set = FALSE;
683 /* set the energyhistory variables to save state */
684 void mde_delta_h_coll_update_energyhistory(const t_mde_delta_h_coll *dhc,
685 energyhistory_t *enerhist)
687 if (enerhist->deltaHForeignLambdas == nullptr)
689 enerhist->deltaHForeignLambdas.reset(new delta_h_history_t);
690 enerhist->deltaHForeignLambdas->dh.resize(dhc->ndh);
693 delta_h_history_t * const deltaH = enerhist->deltaHForeignLambdas.get();
695 GMX_RELEASE_ASSERT(deltaH->dh.size() == static_cast<size_t>(dhc->ndh), "energy history number of delta_h histograms should match inputrec's number");
697 for (int i = 0; i < dhc->ndh; i++)
699 std::vector<real> &dh = deltaH->dh[i];
700 dh.resize(dhc->dh[i].ndh);
701 std::copy(dh.begin(), dh.end(), dhc->dh[i].dh);
703 deltaH->start_time = dhc->start_time;
704 deltaH->start_lambda = dhc->start_lambda;
709 /* restore the variables from an energyhistory */
710 void mde_delta_h_coll_restore_energyhistory(t_mde_delta_h_coll *dhc,
711 const delta_h_history_t *deltaH)
713 GMX_RELEASE_ASSERT(dhc, "Should have delta_h histograms");
714 GMX_RELEASE_ASSERT(deltaH, "Should have delta_h histograms in energy history");
715 GMX_RELEASE_ASSERT(deltaH->dh.size() == static_cast<size_t>(dhc->ndh), "energy history number of delta_h histograms should match inputrec's number");
717 for (unsigned int i = 0; i < deltaH->dh.size(); i++)
719 dhc->dh[i].ndh = deltaH->dh[i].size();
720 for (unsigned int j = 0; j < dhc->dh[i].ndh; j++)
722 dhc->dh[i].dh[j] = deltaH->dh[i][j];
725 dhc->start_time = deltaH->start_time;
726 if (deltaH->start_lambda_set)
728 dhc->start_lambda = deltaH->start_lambda;
730 if (dhc->dh[0].ndh > 0)
732 dhc->start_time_set = TRUE;
736 dhc->start_time_set = FALSE;