Convert some mdrun utility code to C++
[alexxy/gromacs.git] / src / gromacs / mdlib / mdebin_bar.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "mdebin_bar.h"
40
41 #include <float.h>
42 #include <math.h>
43 #include <string.h>
44
45 #include "gromacs/fileio/enxio.h"
46 #include "gromacs/legacyheaders/mdebin.h"
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/smalloc.h"
50
51 /* reset the delta_h list to prepare it for new values */
52 static void mde_delta_h_reset(t_mde_delta_h *dh)
53 {
54     dh->ndh     = 0;
55     dh->written = FALSE;
56 }
57
58 /* initialize the delta_h list */
59 static void mde_delta_h_init(t_mde_delta_h *dh, int nbins,
60                              double dx, unsigned int  ndhmax,
61                              int type, int derivative, int nlambda,
62                              double *lambda)
63 {
64     int i;
65
66     dh->type       = type;
67     dh->derivative = derivative;
68     dh->lambda     = lambda;
69     dh->nlambda    = nlambda;
70
71     snew(dh->lambda, nlambda);
72     for (i = 0; i < nlambda; i++)
73     {
74         dh->lambda[i] = lambda[i];
75     }
76
77
78     snew(dh->subblock_meta_d, dh->nlambda+1);
79
80     dh->ndhmax = ndhmax+2;
81     for (i = 0; i < 2; i++)
82     {
83         dh->bin[i] = NULL;
84     }
85
86     snew(dh->dh, dh->ndhmax);
87     snew(dh->dhf, dh->ndhmax);
88
89     if (nbins <= 0 || dx < GMX_REAL_EPS*10)
90     {
91         dh->nhist = 0;
92     }
93     else
94     {
95         int i;
96         /* pre-allocate the histogram */
97         dh->nhist = 2; /* energies and derivatives histogram */
98         dh->dx    = dx;
99         dh->nbins = nbins;
100         for (i = 0; i < dh->nhist; i++)
101         {
102             snew(dh->bin[i], dh->nbins);
103         }
104     }
105     mde_delta_h_reset(dh);
106 }
107
108 /* Add a value to the delta_h list */
109 static void mde_delta_h_add_dh(t_mde_delta_h *dh, double delta_h)
110 {
111     if (dh->ndh >= dh->ndhmax)
112     {
113         gmx_incons("delta_h array not big enough!");
114     }
115     dh->dh[dh->ndh] = delta_h;
116     dh->ndh++;
117 }
118
119 /* construct histogram with index hi */
120 static void mde_delta_h_make_hist(t_mde_delta_h *dh, int hi, gmx_bool invert)
121 {
122     double       min_dh = FLT_MAX;
123     double       max_dh = -FLT_MAX;
124     unsigned int i;
125     double       max_dh_hist; /* maximum binnable dh value */
126     double       min_dh_hist; /* minimum binnable dh value */
127     double       dx = dh->dx;
128     double       f;           /* energy mult. factor */
129
130     /* by applying a -1 scaling factor on the energies we get the same as
131        having a negative dx, but we don't need to fix the min/max values
132        beyond inverting x0 */
133     f = invert ? -1 : 1;
134
135     /* first find min and max */
136     for (i = 0; i < dh->ndh; i++)
137     {
138         if (f*dh->dh[i] < min_dh)
139         {
140             min_dh = f*dh->dh[i];
141         }
142         if (f*dh->dh[i] > max_dh)
143         {
144             max_dh = f*dh->dh[i];
145         }
146     }
147
148     /* reset the histogram */
149     for (i = 0; i < dh->nbins; i++)
150     {
151         dh->bin[hi][i] = 0;
152     }
153     dh->maxbin[hi] = 0;
154
155     /* The starting point of the histogram is the lowest value found:
156        that value has the highest contribution to the free energy.
157
158        Get this start value in number of histogram dxs from zero,
159        as an integer.*/
160
161     dh->x0[hi] = (gmx_int64_t)floor(min_dh/dx);
162
163     min_dh_hist = (dh->x0[hi])*dx;
164     max_dh_hist = (dh->x0[hi] + dh->nbins + 1)*dx;
165
166     /* and fill the histogram*/
167     for (i = 0; i < dh->ndh; i++)
168     {
169         unsigned int bin;
170
171         /* Determine the bin number. If it doesn't fit into the histogram,
172            add it to the last bin.
173            We check the max_dh_int range because converting to integers
174            might lead to overflow with unpredictable results.*/
175         if ( (f*dh->dh[i] >= min_dh_hist) && (f*dh->dh[i] <= max_dh_hist ) )
176         {
177             bin = (unsigned int)( (f*dh->dh[i] - min_dh_hist)/dx );
178         }
179         else
180         {
181             bin = dh->nbins-1;
182         }
183         /* double-check here because of possible round-off errors*/
184         if (bin >= dh->nbins)
185         {
186             bin = dh->nbins-1;
187         }
188         if (bin > dh->maxbin[hi])
189         {
190             dh->maxbin[hi] = bin;
191         }
192
193         dh->bin[hi][bin]++;
194     }
195
196     /* make sure we include a bin with 0 if we didn't use the full
197        histogram width. This can then be used as an indication that
198        all the data was binned. */
199     if (dh->maxbin[hi] < dh->nbins-1)
200     {
201         dh->maxbin[hi] += 1;
202     }
203 }
204
205
206 void mde_delta_h_handle_block(t_mde_delta_h *dh, t_enxblock *blk)
207 {
208     /* first check which type we should use: histogram or raw data */
209     if (dh->nhist == 0)
210     {
211         int i;
212
213         /* We write raw data.
214            Raw data consists of 3 subblocks: an int metadata block
215            with type and derivative index, a foreign lambda block
216            and and the data itself */
217         add_subblocks_enxblock(blk, 3);
218
219         blk->id = enxDH;
220
221         /* subblock 1 */
222         dh->subblock_meta_i[0] = dh->type;       /* block data type */
223         dh->subblock_meta_i[1] = dh->derivative; /* derivative direction if
224                                                     applicable (in indices
225                                                     starting from first coord in
226                                                     the main delta_h_coll) */
227         blk->sub[0].nr         = 2;
228         blk->sub[0].type       = xdr_datatype_int;
229         blk->sub[0].ival       = dh->subblock_meta_i;
230
231         /* subblock 2 */
232         for (i = 0; i < dh->nlambda; i++)
233         {
234             dh->subblock_meta_d[i] = dh->lambda[i];
235         }
236         blk->sub[1].nr   = dh->nlambda;
237         blk->sub[1].type = xdr_datatype_double;
238         blk->sub[1].dval = dh->subblock_meta_d;
239
240         /* subblock 3 */
241         /* check if there's actual data to be written. */
242         /*if (dh->ndh > 1)*/
243         if (dh->ndh > 0)
244         {
245             unsigned int i;
246
247             blk->sub[2].nr = dh->ndh;
248 /* For F@H for now. */
249 #undef GMX_DOUBLE
250 #ifndef GMX_DOUBLE
251             blk->sub[2].type = xdr_datatype_float;
252             for (i = 0; i < dh->ndh; i++)
253             {
254                 dh->dhf[i] = (float)dh->dh[i];
255             }
256             blk->sub[2].fval = dh->dhf;
257 #else
258             blk->sub[2].type = xdr_datatype_double;
259             blk->sub[2].dval = dh->dh;
260 #endif
261             dh->written = TRUE;
262         }
263         else
264         {
265             blk->sub[2].nr = 0;
266 #ifndef GMX_DOUBLE
267             blk->sub[2].type = xdr_datatype_float;
268             blk->sub[2].fval = NULL;
269 #else
270             blk->sub[2].type = xdr_datatype_double;
271             blk->sub[2].dval = NULL;
272 #endif
273         }
274     }
275     else
276     {
277         int nhist_written = 0;
278         int i;
279         int k;
280
281         /* TODO histogram metadata */
282         /* check if there's actual data to be written. */
283         if (dh->ndh > 1)
284         {
285             gmx_bool prev_complete = FALSE;
286             /* Make the histogram(s) */
287             for (i = 0; i < dh->nhist; i++)
288             {
289                 if (!prev_complete)
290                 {
291                     /* the first histogram is always normal, and the
292                        second one is always reverse */
293                     mde_delta_h_make_hist(dh, i, i == 1);
294                     nhist_written++;
295                     /* check whether this histogram contains all data: if the
296                        last bin is 0, it does */
297                     if (dh->bin[i][dh->nbins-1] == 0)
298                     {
299                         prev_complete = TRUE;
300                     }
301                     if (!dh->derivative)
302                     {
303                         prev_complete = TRUE;
304                     }
305                 }
306             }
307             dh->written = TRUE;
308         }
309
310         /* A histogram consists of 2, 3 or 4 subblocks:
311            the foreign lambda value + histogram spacing, the starting point,
312            and the histogram data (0, 1 or 2 blocks). */
313         add_subblocks_enxblock(blk, nhist_written+2);
314         blk->id = enxDHHIST;
315
316         /* subblock 1: the lambda value + the histogram spacing */
317         if (dh->nlambda == 1)
318         {
319             /* for backward compatibility */
320             dh->subblock_meta_d[0] = dh->lambda[0];
321         }
322         else
323         {
324             dh->subblock_meta_d[0] = -1;
325             for (i = 0; i < dh->nlambda; i++)
326             {
327                 dh->subblock_meta_d[2+i] = dh->lambda[i];
328             }
329         }
330         dh->subblock_meta_d[1] = dh->dx;
331         blk->sub[0].nr         = 2+ ((dh->nlambda > 1) ? dh->nlambda : 0);
332         blk->sub[0].type       = xdr_datatype_double;
333         blk->sub[0].dval       = dh->subblock_meta_d;
334
335         /* subblock 2: the starting point(s) as a long integer */
336         dh->subblock_meta_l[0] = nhist_written;
337         dh->subblock_meta_l[1] = dh->type; /*dh->derivative ? 1 : 0;*/
338         k = 2;
339         for (i = 0; i < nhist_written; i++)
340         {
341             dh->subblock_meta_l[k++] = dh->x0[i];
342         }
343         /* append the derivative data */
344         dh->subblock_meta_l[k++] = dh->derivative;
345
346         blk->sub[1].nr   = nhist_written+3;
347         blk->sub[1].type = xdr_datatype_int64;
348         blk->sub[1].lval = dh->subblock_meta_l;
349
350         /* subblock 3 + 4 : the histogram data */
351         for (i = 0; i < nhist_written; i++)
352         {
353             blk->sub[i+2].nr   = dh->maxbin[i]+1; /* it's +1 because size=index+1
354                                                      in C */
355             blk->sub[i+2].type = xdr_datatype_int;
356             blk->sub[i+2].ival = dh->bin[i];
357         }
358     }
359 }
360
361 /* initialize the collection*/
362 void mde_delta_h_coll_init(t_mde_delta_h_coll *dhc, const t_inputrec *ir)
363 {
364     int       i, j, n;
365     double   *lambda_vec;
366     int       ndhmax = ir->nstenergy/ir->nstcalcenergy;
367     t_lambda *fep    = ir->fepvals;
368
369     dhc->temperature    = ir->opts.ref_t[0]; /* only store system temperature */
370     dhc->start_time     = 0.;
371     dhc->delta_time     = ir->delta_t*ir->fepvals->nstdhdl;
372     dhc->start_time_set = FALSE;
373
374     /* this is the compatibility lambda value. If it is >=0, it is valid,
375        and there is either an old-style lambda or a slow growth simulation. */
376     dhc->start_lambda = ir->fepvals->init_lambda;
377     /* for continuous change of lambda values */
378     dhc->delta_lambda = ir->fepvals->delta_lambda*ir->fepvals->nstdhdl;
379
380     if (dhc->start_lambda < 0)
381     {
382         /* create the native lambda vectors */
383         dhc->lambda_index = fep->init_fep_state;
384         dhc->n_lambda_vec = 0;
385         for (i = 0; i < efptNR; i++)
386         {
387             if (fep->separate_dvdl[i])
388             {
389                 dhc->n_lambda_vec++;
390             }
391         }
392         snew(dhc->native_lambda_vec, dhc->n_lambda_vec);
393         snew(dhc->native_lambda_components, dhc->n_lambda_vec);
394         j = 0;
395         for (i = 0; i < efptNR; i++)
396         {
397             if (fep->separate_dvdl[i])
398             {
399                 dhc->native_lambda_components[j] = i;
400                 if (fep->init_fep_state >= 0 &&
401                     fep->init_fep_state < fep->n_lambda)
402                 {
403                     dhc->native_lambda_vec[j] =
404                         fep->all_lambda[i][fep->init_fep_state];
405                 }
406                 else
407                 {
408                     dhc->native_lambda_vec[j] = -1;
409                 }
410                 j++;
411             }
412         }
413     }
414     else
415     {
416         /* don't allocate the meta-data subblocks for lambda vectors */
417         dhc->native_lambda_vec        = NULL;
418         dhc->n_lambda_vec             = 0;
419         dhc->native_lambda_components = 0;
420         dhc->lambda_index             = -1;
421     }
422     /* allocate metadata subblocks */
423     snew(dhc->subblock_d, 5 + dhc->n_lambda_vec);
424     snew(dhc->subblock_i, 1 + dhc->n_lambda_vec);
425
426     /* now decide which data to write out */
427     dhc->nlambda     = 0;
428     dhc->ndhdl       = 0;
429     dhc->dh_expanded = NULL;
430     dhc->dh_energy   = NULL;
431     dhc->dh_pv       = NULL;
432
433     /* total number of raw data point collections in the sample */
434     dhc->ndh = 0;
435
436     {
437         gmx_bool bExpanded           = FALSE;
438         gmx_bool bEnergy             = FALSE;
439         gmx_bool bPV                 = FALSE;
440         int      n_lambda_components = 0;
441
442         /* first count the number of states */
443
444         /* add the dhdl's */
445         if (fep->dhdl_derivatives == edhdlderivativesYES)
446         {
447             for (i = 0; i < efptNR; i++)
448             {
449                 if (ir->fepvals->separate_dvdl[i])
450                 {
451                     dhc->ndh   += 1;
452                     dhc->ndhdl += 1;
453                 }
454             }
455         }
456         /* add the lambdas */
457         dhc->nlambda = ir->fepvals->lambda_stop_n - ir->fepvals->lambda_start_n;
458         dhc->ndh    += dhc->nlambda;
459         /* another compatibility check */
460         if (dhc->start_lambda < 0)
461         {
462             /* include one more for the specification of the state, by lambda or
463                fep_state*/
464             if (ir->expandedvals->elmcmove > elmcmoveNO)
465             {
466                 dhc->ndh += 1;
467                 bExpanded = TRUE;
468             }
469             /* whether to print energies */
470             if (ir->fepvals->edHdLPrintEnergy != edHdLPrintEnergyNO)
471             {
472                 dhc->ndh += 1;
473                 bEnergy   = TRUE;
474             }
475             if (ir->epc > epcNO)
476             {
477                 dhc->ndh += 1;  /* include pressure-volume work */
478                 bPV       = TRUE;
479             }
480         }
481         /* allocate them */
482         snew(dhc->dh, dhc->ndh);
483
484         /* now initialize them */
485         /* the order, for now, must match that of the dhdl.xvg file because of
486            how g_energy -odh is implemented */
487         n = 0;
488         if (bExpanded)
489         {
490             dhc->dh_expanded = dhc->dh+n;
491             mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
492                              ir->fepvals->dh_hist_spacing, ndhmax,
493                              dhbtEXPANDED, 0, 0, NULL);
494             n++;
495         }
496         if (bEnergy)
497         {
498             dhc->dh_energy = dhc->dh+n;
499             mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
500                              ir->fepvals->dh_hist_spacing, ndhmax,
501                              dhbtEN, 0, 0, NULL);
502             n++;
503         }
504         /* add the dhdl's */
505         n_lambda_components = 0;
506         if (fep->dhdl_derivatives == edhdlderivativesYES)
507         {
508             dhc->dh_dhdl = dhc->dh + n;
509             for (i = 0; i < efptNR; i++)
510             {
511                 if (ir->fepvals->separate_dvdl[i])
512                 {
513                     /* we give it init_lambda for compatibility */
514                     mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
515                                      ir->fepvals->dh_hist_spacing, ndhmax,
516                                      dhbtDHDL, n_lambda_components, 1,
517                                      &(fep->init_lambda));
518                     n++;
519                     n_lambda_components++;
520                 }
521             }
522         }
523         else
524         {
525             for (i = 0; i < efptNR; i++)
526             {
527                 if (ir->fepvals->separate_dvdl[i])
528                 {
529                     n_lambda_components++; /* count the components */
530                 }
531             }
532
533         }
534         /* add the lambdas */
535         dhc->dh_du = dhc->dh + n;
536         snew(lambda_vec, n_lambda_components);
537         for (i = ir->fepvals->lambda_start_n; i < ir->fepvals->lambda_stop_n; i++)
538         {
539             int k = 0;
540
541             for (j = 0; j < efptNR; j++)
542             {
543                 if (ir->fepvals->separate_dvdl[j])
544                 {
545                     lambda_vec[k++] = fep->all_lambda[j][i];
546                 }
547             }
548
549             mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
550                              ir->fepvals->dh_hist_spacing, ndhmax,
551                              dhbtDH, 0, n_lambda_components, lambda_vec);
552             n++;
553         }
554         sfree(lambda_vec);
555         if (bPV)
556         {
557             dhc->dh_pv = dhc->dh+n;
558             mde_delta_h_init(dhc->dh+n, ir->fepvals->dh_hist_size,
559                              ir->fepvals->dh_hist_spacing, ndhmax,
560                              dhbtPV, 0, 0, NULL);
561             n++;
562         }
563     }
564 }
565
566 /* add a bunch of samples - note fep_state is double to allow for better data storage */
567 void mde_delta_h_coll_add_dh(t_mde_delta_h_coll *dhc,
568                              double              fep_state,
569                              double              energy,
570                              double              pV,
571                              double             *dhdl,
572                              double             *foreign_dU,
573                              double              time)
574 {
575     int i;
576
577     if (!dhc->start_time_set)
578     {
579         dhc->start_time_set = TRUE;
580         dhc->start_time     = time;
581     }
582
583     for (i = 0; i < dhc->ndhdl; i++)
584     {
585         mde_delta_h_add_dh(dhc->dh_dhdl+i, dhdl[i]);
586     }
587     for (i = 0; i < dhc->nlambda; i++)
588     {
589         mde_delta_h_add_dh(dhc->dh_du+i, foreign_dU[i]);
590     }
591     if (dhc->dh_pv != NULL)
592     {
593         mde_delta_h_add_dh(dhc->dh_pv, pV);
594     }
595     if (dhc->dh_energy != NULL)
596     {
597         mde_delta_h_add_dh(dhc->dh_energy, energy);
598     }
599     if (dhc->dh_expanded != NULL)
600     {
601         mde_delta_h_add_dh(dhc->dh_expanded, fep_state);
602     }
603
604 }
605
606 /* write the metadata associated with all the du blocks, and call
607    handle_block to write out all the du blocks */
608 void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
609                                    t_enxframe *fr, int nblock)
610 {
611     int         i;
612     t_enxblock *blk;
613
614     /* add one block with one subblock as the collection's own data */
615     nblock++;
616     add_blocks_enxframe(fr, nblock);
617     blk = fr->block + (nblock-1);
618
619     /* only allocate lambda vector component blocks if they must be written out
620        for backward compatibility */
621     if (dhc->native_lambda_components != NULL)
622     {
623         add_subblocks_enxblock(blk, 2);
624     }
625     else
626     {
627         add_subblocks_enxblock(blk, 1);
628     }
629
630     dhc->subblock_d[0] = dhc->temperature;  /* temperature */
631     dhc->subblock_d[1] = dhc->start_time;   /* time of first sample */
632     dhc->subblock_d[2] = dhc->delta_time;   /* time difference between samples */
633     dhc->subblock_d[3] = dhc->start_lambda; /* old-style lambda at starttime */
634     dhc->subblock_d[4] = dhc->delta_lambda; /* lambda diff. between samples */
635     /* set the lambda vector components if they exist */
636     if (dhc->native_lambda_components != NULL)
637     {
638         for (i = 0; i < dhc->n_lambda_vec; i++)
639         {
640             dhc->subblock_d[5+i] = dhc->native_lambda_vec[i];
641         }
642     }
643     blk->id          = enxDHCOLL;
644     blk->sub[0].nr   = 5 + dhc->n_lambda_vec;
645     blk->sub[0].type = xdr_datatype_double;
646     blk->sub[0].dval = dhc->subblock_d;
647
648     if (dhc->native_lambda_components != NULL)
649     {
650         dhc->subblock_i[0] = dhc->lambda_index;
651         /* set the lambda vector component IDs if they exist */
652         dhc->subblock_i[1] = dhc->n_lambda_vec;
653         for (i = 0; i < dhc->n_lambda_vec; i++)
654         {
655             dhc->subblock_i[i+2] = dhc->native_lambda_components[i];
656         }
657         blk->sub[1].nr   = 2 + dhc->n_lambda_vec;
658         blk->sub[1].type = xdr_datatype_int;
659         blk->sub[1].ival = dhc->subblock_i;
660     }
661
662     for (i = 0; i < dhc->ndh; i++)
663     {
664         nblock++;
665         add_blocks_enxframe(fr, nblock);
666         blk = fr->block + (nblock-1);
667
668         mde_delta_h_handle_block(dhc->dh+i, blk);
669     }
670 }
671
672 /* reset the data for a new round */
673 void mde_delta_h_coll_reset(t_mde_delta_h_coll *dhc)
674 {
675     int i;
676     for (i = 0; i < dhc->ndh; i++)
677     {
678         if (dhc->dh[i].written)
679         {
680             /* we can now throw away the data */
681             mde_delta_h_reset(dhc->dh + i);
682         }
683     }
684     dhc->start_time_set = FALSE;
685 }
686
687 /* set the energyhistory variables to save state */
688 void mde_delta_h_coll_update_energyhistory(t_mde_delta_h_coll *dhc,
689                                            energyhistory_t    *enerhist)
690 {
691     int i;
692     if (!enerhist->dht)
693     {
694         snew(enerhist->dht, 1);
695         snew(enerhist->dht->ndh, dhc->ndh);
696         snew(enerhist->dht->dh, dhc->ndh);
697         enerhist->dht->nndh = dhc->ndh;
698     }
699     else
700     {
701         if (enerhist->dht->nndh != dhc->ndh)
702         {
703             gmx_incons("energy history number of delta_h histograms != inputrec's number");
704         }
705     }
706     for (i = 0; i < dhc->ndh; i++)
707     {
708         enerhist->dht->dh[i]  = dhc->dh[i].dh;
709         enerhist->dht->ndh[i] = dhc->dh[i].ndh;
710     }
711     enerhist->dht->start_time   = dhc->start_time;
712     enerhist->dht->start_lambda = dhc->start_lambda;
713 }
714
715
716
717 /* restore the variables from an energyhistory */
718 void mde_delta_h_coll_restore_energyhistory(t_mde_delta_h_coll *dhc,
719                                             energyhistory_t    *enerhist)
720 {
721     int          i;
722     unsigned int j;
723
724     if (!dhc)
725     {
726         gmx_incons("No delta_h histograms found");
727     }
728     if (!enerhist->dht)
729     {
730         gmx_incons("No delta_h histograms found in energy history");
731     }
732     if (enerhist->dht->nndh != dhc->ndh)
733     {
734         gmx_incons("energy history number of delta_h histograms != inputrec's number");
735     }
736
737     for (i = 0; i < enerhist->dht->nndh; i++)
738     {
739         dhc->dh[i].ndh = enerhist->dht->ndh[i];
740         for (j = 0; j < dhc->dh[i].ndh; j++)
741         {
742             dhc->dh[i].dh[j] = enerhist->dht->dh[i][j];
743         }
744     }
745     dhc->start_time = enerhist->dht->start_time;
746     if (enerhist->dht->start_lambda_set)
747     {
748         dhc->start_lambda = enerhist->dht->start_lambda;
749     }
750     if (dhc->dh[0].ndh > 0)
751     {
752         dhc->start_time_set = TRUE;
753     }
754     else
755     {
756         dhc->start_time_set = FALSE;
757     }
758 }