Minor fix for checkpointing.
[gromacs.git] / src / gmxlib / checkpoint.c
blob60ba8264345f2fcf0733c40a6c4f2fb51db16e98
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
15 * And Hey:
16 * Gnomes, ROck Monsters And Chili Sauce
19 #ifdef HAVE_CONFIG_H
20 #include <config.h>
21 #endif
23 #include <string.h>
24 #include <time.h>
26 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
27 /* _chsize_s */
28 #include <io.h>
29 #include <sys/locking.h>
30 #endif
33 #include "filenm.h"
34 #include "names.h"
35 #include "typedefs.h"
36 #include "smalloc.h"
37 #include "gmxfio.h"
38 #include "xdrf.h"
39 #include "statutil.h"
40 #include "txtdump.h"
41 #include "vec.h"
42 #include "network.h"
43 #include "gmx_random.h"
44 #include "checkpoint.h"
45 #include "futil.h"
46 #include <fcntl.h>
50 #ifdef GMX_FAHCORE
51 #include "corewrap.h"
52 #endif
54 #define CPT_MAGIC1 171817
55 #define CPT_MAGIC2 171819
57 /* The source code in this file should be thread-safe.
58 Please keep it that way. */
60 /* cpt_version should normally only be changed
61 * when the header of footer format changes.
62 * The state data format itself is backward and forward compatible.
63 * But old code can not read a new entry that is present in the file
64 * (but can read a new format when new entries are not present).
66 static const int cpt_version = 8;
68 enum { ecpdtINT, ecpdtFLOAT, ecpdtDOUBLE, ecpdtNR };
70 const char *ecpdt_names[ecpdtNR] = { "int", "float", "double" };
72 const char *est_names[estNR]=
74 "FE-lambda",
75 "box", "box-rel", "box-v", "pres_prev",
76 "nosehoover-xi", "thermostat-integral",
77 "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
78 "disre_initf", "disre_rm3tav",
79 "orire_initf", "orire_Dtav",
82 enum { eeksEKINH_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksNR };
84 const char *eeks_names[eeksNR]=
86 "Ekinh_n", "Ekinh", "dEkindlambda", "mv_cos"
89 enum { eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
90 eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
91 eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM, eenhNR };
93 const char *eenh_names[eenhNR]=
95 "energy_n", "energy_aver", "energy_sum", "energy_nsum",
96 "energy_sum_sim", "energy_nsum_sim",
97 "energy_nsteps", "energy_nsteps_sim"
102 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
103 static int
104 gmx_wintruncate(const char *filename, __int64 size)
106 FILE *fp;
107 int rc;
109 fp=fopen(filename,"r+");
111 if(fp==NULL)
113 return -1;
116 return _chsize_s( fileno(fp), size);
118 #endif
122 enum { ecprREAL, ecprRVEC, ecprMATRIX };
124 static const char *st_names(int cptp,int ecpt)
126 switch (cptp)
128 case 0: return est_names [ecpt]; break;
129 case 1: return eeks_names[ecpt]; break;
130 case 2: return eenh_names[ecpt]; break;
133 return NULL;
136 static void cp_warning(FILE *fp)
138 fprintf(fp,"\nWARNING: Checkpoint file is corrupted or truncated\n\n");
141 static void cp_error()
143 gmx_fatal(FARGS,"Checkpoint file corrupted/truncated, or maybe you are out of quota?");
146 static void do_cpt_string_err(XDR *xd,bool bRead,const char *desc,char **s,FILE *list)
148 #define CPTSTRLEN 1024
149 bool_t res=0;
151 if (bRead)
153 snew(*s,CPTSTRLEN);
155 res = xdr_string(xd,s,CPTSTRLEN);
156 if (res == 0)
158 cp_error();
160 if (list)
162 fprintf(list,"%s = %s\n",desc,*s);
163 sfree(*s);
167 static int do_cpt_int(XDR *xd,const char *desc,int *i,FILE *list)
169 bool_t res=0;
171 res = xdr_int(xd,i);
172 if (res == 0)
174 return -1;
176 if (list)
178 fprintf(list,"%s = %d\n",desc,*i);
180 return 0;
183 static int do_cpt_u_chars(XDR *xd,const char *desc,int n,unsigned char *i,FILE *list)
185 bool_t res=1;
186 int j;
187 if (list)
189 fprintf(list,"%s = ",desc);
191 for (j=0; j<n && res; j++)
193 res &= xdr_u_char(xd,&i[j]);
194 if (list)
196 fprintf(list,"%02x",i[j]);
199 if (list)
201 fprintf(list,"\n");
203 if (res == 0)
205 return -1;
208 return 0;
211 static void do_cpt_int_err(XDR *xd,const char *desc,int *i,FILE *list)
213 if (do_cpt_int(xd,desc,i,list) < 0)
215 cp_error();
219 static void do_cpt_step_err(XDR *xd,const char *desc,gmx_large_int_t *i,FILE *list)
221 bool_t res=0;
222 char buf[22];
224 res = xdr_gmx_large_int(xd,i,"reading checkpoint file");
225 if (res == 0)
227 cp_error();
229 if (list)
231 fprintf(list,"%s = %s\n",desc,gmx_step_str(*i,buf));
235 static void do_cpt_double_err(XDR *xd,const char *desc,double *f,FILE *list)
237 bool_t res=0;
239 res = xdr_double(xd,f);
240 if (res == 0)
242 cp_error();
244 if (list)
246 fprintf(list,"%s = %f\n",desc,*f);
251 static int do_cpte_reals_low(XDR *xd,int cptp,int ecpt,int sflags,
252 int n,real **v,
253 FILE *list,int erealtype)
255 bool_t res=0;
256 #ifndef GMX_DOUBLE
257 int dtc=ecpdtFLOAT;
258 #else
259 int dtc=ecpdtDOUBLE;
260 #endif
261 real *vp,*va=NULL;
262 float *vf;
263 double *vd;
264 int nf,dt,i;
266 nf = n;
267 res = xdr_int(xd,&nf);
268 if (res == 0)
270 return -1;
272 if (list == NULL && nf != n)
274 gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
276 dt = dtc;
277 res = xdr_int(xd,&dt);
278 if (res == 0)
280 return -1;
282 if (dt != dtc)
284 fprintf(stderr,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
285 st_names(cptp,ecpt),ecpdt_names[dtc],ecpdt_names[dt]);
287 if (list || !(sflags & (1<<ecpt)))
289 snew(va,nf);
290 vp = va;
292 else
294 if (*v == NULL)
296 snew(*v,nf);
298 vp = *v;
300 if (dt == ecpdtFLOAT)
302 if (dtc == ecpdtFLOAT)
304 vf = (float *)vp;
306 else
308 snew(vf,nf);
310 res = xdr_vector(xd,(char *)vf,nf,
311 (unsigned int)sizeof(float),(xdrproc_t)xdr_float);
312 if (res == 0)
314 return -1;
316 if (dtc != ecpdtFLOAT)
318 for(i=0; i<nf; i++)
320 vp[i] = vf[i];
322 sfree(vf);
325 else
327 if (dtc == ecpdtDOUBLE)
329 vd = (double *)vp;
331 else
333 snew(vd,nf);
335 res = xdr_vector(xd,(char *)vd,nf,
336 (unsigned int)sizeof(double),(xdrproc_t)xdr_double);
337 if (res == 0)
339 return -1;
341 if (dtc != ecpdtDOUBLE)
343 for(i=0; i<nf; i++)
345 vp[i] = vd[i];
347 sfree(vd);
351 if (list)
353 switch (erealtype)
355 case ecprREAL:
356 pr_reals(list,0,st_names(cptp,ecpt),vp,nf);
357 break;
358 case ecprRVEC:
359 pr_rvecs(list,0,st_names(cptp,ecpt),(rvec *)vp,nf/3);
360 break;
361 default:
362 gmx_incons("Unknown checkpoint real type");
365 if (va)
367 sfree(va);
370 return 0;
375 static int do_cpte_reals(XDR *xd,int cptp,int ecpt,int sflags,
376 int n,real **v,FILE *list)
378 return do_cpte_reals_low(xd,cptp,ecpt,sflags,n,v,list,ecprREAL);
381 static int do_cpte_real(XDR *xd,int cptp,int ecpt,int sflags,
382 real *r,FILE *list)
384 return do_cpte_reals_low(xd,cptp,ecpt,sflags,1,&r,list,ecprREAL);
387 static int do_cpte_ints(XDR *xd,int cptp,int ecpt,int sflags,
388 int n,int **v,FILE *list)
390 bool_t res=0;
391 int dtc=ecpdtINT;
392 int *vp,*va=NULL;
393 int nf,dt,i;
395 nf = n;
396 res = xdr_int(xd,&nf);
397 if (res == 0)
399 return -1;
401 if (list == NULL && v != NULL && nf != n)
403 gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
405 dt = dtc;
406 res = xdr_int(xd,&dt);
407 if (res == 0)
409 return -1;
411 if (dt != dtc)
413 gmx_fatal(FARGS,"Type mismatch for state entry %s, code type is %s, file type is %s\n",
414 st_names(cptp,ecpt),ecpdt_names[dtc],ecpdt_names[dt]);
416 if (list || !(sflags & (1<<ecpt)) || v == NULL)
418 snew(va,nf);
419 vp = va;
421 else
423 if (*v == NULL)
425 snew(*v,nf);
427 vp = *v;
429 res = xdr_vector(xd,(char *)vp,nf,
430 (unsigned int)sizeof(int),(xdrproc_t)xdr_int);
431 if (res == 0)
433 return -1;
435 if (list)
437 pr_ivec(list,0,st_names(cptp,ecpt),vp,nf,TRUE);
439 if (va)
441 sfree(va);
444 return 0;
447 static int do_cpte_int(XDR *xd,int cptp,int ecpt,int sflags,
448 int *i,FILE *list)
450 return do_cpte_ints(xd,cptp,ecpt,sflags,1,&i,list);
453 static int do_cpte_doubles(XDR *xd,int cptp,int ecpt,int sflags,
454 int n,double **v,FILE *list)
456 bool_t res=0;
457 int dtc=ecpdtDOUBLE;
458 double *vp,*va=NULL;
459 int nf,dt,i;
461 nf = n;
462 res = xdr_int(xd,&nf);
463 if (res == 0)
465 return -1;
467 if (list == NULL && nf != n)
469 gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
471 dt = dtc;
472 res = xdr_int(xd,&dt);
473 if (res == 0)
475 return -1;
477 if (dt != dtc)
479 gmx_fatal(FARGS,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
480 st_names(cptp,ecpt),ecpdt_names[dtc],ecpdt_names[dt]);
482 if (list || !(sflags & (1<<ecpt)))
484 snew(va,nf);
485 vp = va;
487 else
489 if (*v == NULL)
491 snew(*v,nf);
493 vp = *v;
495 res = xdr_vector(xd,(char *)vp,nf,
496 (unsigned int)sizeof(double),(xdrproc_t)xdr_double);
497 if (res == 0)
499 return -1;
501 if (list)
503 pr_doubles(list,0,st_names(cptp,ecpt),vp,nf);
505 if (va)
507 sfree(va);
510 return 0;
514 static int do_cpte_rvecs(XDR *xd,int cptp,int ecpt,int sflags,
515 int n,rvec **v,FILE *list)
517 return do_cpte_reals_low(xd,cptp,ecpt,sflags,
518 n*DIM,(real **)v,list,ecprRVEC);
521 static int do_cpte_matrix(XDR *xd,int cptp,int ecpt,int sflags,
522 matrix v,FILE *list)
524 real *vr;
525 real ret;
527 vr = (real *)&(v[0][0]);
528 ret = do_cpte_reals_low(xd,cptp,ecpt,sflags,DIM*DIM,&vr,NULL,ecprMATRIX);
530 if (list && ret == 0)
532 pr_rvecs(list,0,st_names(cptp,ecpt),v,DIM);
535 return ret;
538 static int do_cpte_matrices(XDR *xd,int cptp,int ecpt,int sflags,
539 int n,matrix **v,FILE *list)
541 bool_t res=0;
542 matrix *vp,*va=NULL;
543 real *vr;
544 int nf,i,j,k;
545 int ret;
547 nf = n;
548 res = xdr_int(xd,&nf);
549 if (res == 0)
551 return -1;
553 if (list == NULL && nf != n)
555 gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
557 if (list || !(sflags & (1<<ecpt)))
559 snew(va,nf);
560 vp = va;
562 else
564 if (*v == NULL)
566 snew(*v,nf);
568 vp = *v;
570 snew(vr,nf*DIM*DIM);
571 for(i=0; i<nf; i++)
573 for(j=0; j<DIM; j++)
575 for(k=0; k<DIM; k++)
577 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
581 ret = do_cpte_reals_low(xd,cptp,ecpt,sflags,
582 nf*DIM*DIM,&vr,NULL,ecprMATRIX);
583 for(i=0; i<nf; i++)
585 for(j=0; j<DIM; j++)
587 for(k=0; k<DIM; k++)
589 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
593 sfree(vr);
595 if (list && ret == 0)
597 for(i=0; i<nf; i++)
599 pr_rvecs(list,0,st_names(cptp,ecpt),vp[i],DIM);
602 if (va)
604 sfree(va);
607 return ret;
610 static void do_cpt_header(XDR *xd,bool bRead,int *file_version,
611 char **version,char **btime,char **buser,char **bmach,
612 char **fprog,char **ftime,
613 int *eIntegrator,int *simulation_part,
614 gmx_large_int_t *step,double *t,
615 int *nnodes,int *dd_nc,int *npme,
616 int *natoms,int *ngtc,
617 int *flags_state,int *flags_eks,int *flags_enh,
618 FILE *list)
620 bool_t res=0;
621 int magic;
622 int idum=0;
623 int i;
625 if (bRead)
627 magic = -1;
629 else
631 magic = CPT_MAGIC1;
633 res = xdr_int(xd,&magic);
634 if (res == 0)
636 gmx_fatal(FARGS,"The checkpoint file is empty/corrupted, or maybe you are out of quota?");
638 if (magic != CPT_MAGIC1)
640 gmx_fatal(FARGS,"Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
641 "The checkpoint file is corrupted or not a checkpoint file",
642 magic,CPT_MAGIC1);
644 do_cpt_string_err(xd,bRead,"GROMACS version" ,version,list);
645 do_cpt_string_err(xd,bRead,"GROMACS build time" ,btime,list);
646 do_cpt_string_err(xd,bRead,"GROMACS build user" ,buser,list);
647 do_cpt_string_err(xd,bRead,"GROMACS build machine" ,bmach,list);
648 do_cpt_string_err(xd,bRead,"generating program" ,fprog,list);
649 do_cpt_string_err(xd,bRead,"generation time" ,ftime,list);
650 *file_version = cpt_version;
651 do_cpt_int_err(xd,"checkpoint file version",file_version,list);
652 if (*file_version > cpt_version)
654 gmx_fatal(FARGS,"Attempting to read a checkpoint file of version %d with code of version %d\n",*file_version,cpt_version);
656 do_cpt_int_err(xd,"#atoms" ,natoms ,list);
657 do_cpt_int_err(xd,"#T-coupling groups",ngtc ,list);
658 do_cpt_int_err(xd,"integrator" ,eIntegrator,list);
659 if (*file_version >= 3)
661 do_cpt_int_err(xd,"simulation part #", simulation_part,list);
663 else
665 *simulation_part = 1;
667 if (*file_version >= 5)
669 do_cpt_step_err(xd,"step" ,step ,list);
671 else
673 do_cpt_int_err(xd,"step" ,&idum ,list);
674 *step = idum;
676 do_cpt_double_err(xd,"t" ,t ,list);
677 do_cpt_int_err(xd,"#PP-nodes" ,nnodes ,list);
678 idum = 1;
679 do_cpt_int_err(xd,"dd_nc[x]",dd_nc ? &(dd_nc[0]) : &idum,list);
680 do_cpt_int_err(xd,"dd_nc[y]",dd_nc ? &(dd_nc[1]) : &idum,list);
681 do_cpt_int_err(xd,"dd_nc[z]",dd_nc ? &(dd_nc[2]) : &idum,list);
682 do_cpt_int_err(xd,"#PME-only nodes",npme,list);
683 do_cpt_int_err(xd,"state flags",flags_state,list);
684 if (*file_version >= 4)
686 do_cpt_int_err(xd,"ekin data flags",flags_eks,list);
687 do_cpt_int_err(xd,"energy history flags",flags_enh,list);
689 else
691 *flags_eks = 0;
692 *flags_enh = (*flags_state >> (estORIRE_DTAV+1));
693 *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
694 (1<<(estORIRE_DTAV+2)) |
695 (1<<(estORIRE_DTAV+3))));
699 static int do_cpt_footer(XDR *xd,bool bRead,int file_version)
701 bool_t res=0;
702 int magic;
704 if (file_version >= 2)
706 magic = CPT_MAGIC2;
707 res = xdr_int(xd,&magic);
708 if (res == 0)
710 cp_error();
712 if (magic != CPT_MAGIC2)
714 return -1;
718 return 0;
721 static int do_cpt_state(XDR *xd,bool bRead,
722 int fflags,t_state *state,
723 bool bReadRNG,FILE *list)
725 int sflags;
726 int **rng_p,**rngi_p;
727 int i;
728 int ret;
730 ret = 0;
732 if (bReadRNG)
734 rng_p = (int **)&state->ld_rng;
735 rngi_p = &state->ld_rngi;
737 else
739 /* Do not read the RNG data */
740 rng_p = NULL;
741 rngi_p = NULL;
744 sflags = state->flags;
745 for(i=0; (i<estNR && ret == 0); i++)
747 if (fflags & (1<<i))
749 switch (i)
751 case estLAMBDA: ret = do_cpte_real (xd,0,i,sflags,&state->lambda,list); break;
752 case estBOX: ret = do_cpte_matrix(xd,0,i,sflags,state->box,list); break;
753 case estBOX_REL: ret = do_cpte_matrix(xd,0,i,sflags,state->box_rel,list); break;
754 case estBOXV: ret = do_cpte_matrix(xd,0,i,sflags,state->boxv,list); break;
755 case estPRES_PREV: ret = do_cpte_matrix(xd,0,i,sflags,state->pres_prev,list); break;
756 case estNH_XI: ret = do_cpte_reals (xd,0,i,sflags,state->ngtc,&state->nosehoover_xi,list); break;
757 case estTC_INT: ret = do_cpte_doubles(xd,0,i,sflags,state->ngtc,&state->therm_integral,list); break;
758 case estX: ret = do_cpte_rvecs (xd,0,i,sflags,state->natoms,&state->x,list); break;
759 case estV: ret = do_cpte_rvecs (xd,0,i,sflags,state->natoms,&state->v,list); break;
760 case estSDX: ret = do_cpte_rvecs (xd,0,i,sflags,state->natoms,&state->sd_X,list); break;
761 case estLD_RNG: ret = do_cpte_ints (xd,0,i,sflags,state->nrng,rng_p,list); break;
762 case estLD_RNGI: ret = do_cpte_ints (xd,0,i,sflags,state->nrngi,rngi_p,list); break;
763 case estDISRE_INITF: ret = do_cpte_real (xd,0,i,sflags,&state->hist.disre_initf,list); break;
764 case estDISRE_RM3TAV: ret = do_cpte_reals(xd,0,i,sflags,state->hist.ndisrepairs,&state->hist.disre_rm3tav,list); break;
765 case estORIRE_INITF: ret = do_cpte_real (xd,0,i,sflags,&state->hist.orire_initf,list); break;
766 case estORIRE_DTAV: ret = do_cpte_reals(xd,0,i,sflags,state->hist.norire_Dtav,&state->hist.orire_Dtav,list); break;
767 default:
768 gmx_fatal(FARGS,"Unknown state entry %d\n"
769 "You are probably reading a new checkpoint file with old code",i);
774 return ret;
777 static int do_cpt_ekinstate(XDR *xd,bool bRead,
778 int fflags,ekinstate_t *ekins,
779 FILE *list)
781 int i;
782 int ret;
784 ret = 0;
786 for(i=0; (i<eeksNR && ret == 0); i++)
788 if (fflags & (1<<i))
790 switch (i)
793 case eeksEKINH_N: ret = do_cpte_int(xd,1,i,fflags,&ekins->ekinh_n,list); break;
794 case eeksEKINH: ret = do_cpte_matrices(xd,1,i,fflags,ekins->ekinh_n,&ekins->ekinh,list); break;
795 case eeksDEKINDL: ret = do_cpte_real(xd,1,i,fflags,&ekins->dekindl,list); break;
796 case eeksMVCOS: ret = do_cpte_real(xd,1,i,fflags,&ekins->mvcos,list); break;
797 default:
798 gmx_fatal(FARGS,"Unknown ekin data state entry %d\n"
799 "You are probably reading a new checkpoint file with old code",i);
804 return ret;
807 static int do_cpt_enerhist(XDR *xd,bool bRead,
808 int fflags,energyhistory_t *enerhist,
809 FILE *list)
811 int i;
812 int ret;
814 ret = 0;
816 if (bRead)
818 enerhist->nsteps = 0;
819 enerhist->nsum = 0;
820 enerhist->nsteps_sim = 0;
821 enerhist->nsum_sim = 0;
824 for(i=0; (i<eenhNR && ret == 0); i++)
826 if (fflags & (1<<i))
828 switch (i)
830 case eenhENERGY_N: ret = do_cpte_int(xd,2,i,fflags,&enerhist->nener,list); break;
831 case eenhENERGY_AVER: ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_ave,list); break;
832 case eenhENERGY_SUM: ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_sum,list); break;
833 case eenhENERGY_NSUM: do_cpt_step_err(xd,eenh_names[i],&enerhist->nsum,list); break;
834 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_sum_sim,list); break;
835 case eenhENERGY_NSUM_SIM: do_cpt_step_err(xd,eenh_names[i],&enerhist->nsum_sim,list); break;
836 case eenhENERGY_NSTEPS: do_cpt_step_err(xd,eenh_names[i],&enerhist->nsteps,list); break;
837 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd,eenh_names[i],&enerhist->nsteps_sim,list); break;
838 default:
839 gmx_fatal(FARGS,"Unknown energy history entry %d\n"
840 "You are probably reading a new checkpoint file with old code",i);
845 if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
847 /* Assume we have an old file format and copy sum to sum_sim */
848 srenew(enerhist->ener_sum_sim,enerhist->nener);
849 for(i=0; i<enerhist->nener; i++)
851 enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
853 fflags |= (1<<eenhENERGY_SUM_SIM);
856 if ( (fflags & (1<<eenhENERGY_NSUM)) &&
857 !(fflags & (1<<eenhENERGY_NSTEPS)))
859 /* Assume we have an old file format and copy nsum to nsteps */
860 enerhist->nsteps = enerhist->nsum;
861 fflags |= (1<<eenhENERGY_NSTEPS);
863 if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
864 !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
866 /* Assume we have an old file format and copy nsum to nsteps */
867 enerhist->nsteps_sim = enerhist->nsum_sim;
868 fflags |= (1<<eenhENERGY_NSTEPS_SIM);
871 return ret;
874 static int do_cpt_files(XDR *xd, bool bRead,
875 gmx_file_position_t **p_outputfiles, int *nfiles,
876 FILE *list, int file_version)
878 int i,j;
879 off_t offset;
880 off_t mask = 0xFFFFFFFFL;
881 int offset_high,offset_low;
882 char *buf;
883 gmx_file_position_t *outputfiles;
885 if (do_cpt_int(xd,"number of output files",nfiles,list) != 0)
887 return -1;
890 if(bRead)
892 snew(*p_outputfiles,*nfiles);
895 outputfiles = *p_outputfiles;
897 for(i=0;i<*nfiles;i++)
899 /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
900 if(bRead)
902 do_cpt_string_err(xd,bRead,"output filename",&buf,list);
903 strncpy(outputfiles[i].filename,buf,CPTSTRLEN-1);
904 if(list==NULL)
906 sfree(buf);
909 if (do_cpt_int(xd,"file_offset_high",&offset_high,list) != 0)
911 return -1;
913 if (do_cpt_int(xd,"file_offset_low",&offset_low,list) != 0)
915 return -1;
917 #if (SIZEOF_OFF_T > 4)
918 outputfiles[i].offset = ( ((off_t) offset_high) << 32 ) | ( (off_t) offset_low & mask );
919 #else
920 outputfiles[i].offset = offset_low;
921 #endif
923 else
925 buf = outputfiles[i].filename;
926 do_cpt_string_err(xd,bRead,"output filename",&buf,list);
927 /* writing */
928 offset = outputfiles[i].offset;
929 if (offset == -1)
931 offset_low = -1;
932 offset_high = -1;
934 else
936 #if (SIZEOF_OFF_T > 4)
937 offset_low = (int) (offset & mask);
938 offset_high = (int) ((offset >> 32) & mask);
939 #else
940 offset_low = offset;
941 offset_high = 0;
942 #endif
944 if (do_cpt_int(xd,"file_offset_high",&offset_high,list) != 0)
946 return -1;
948 if (do_cpt_int(xd,"file_offset_low",&offset_low,list) != 0)
950 return -1;
953 if (file_version >= 8)
955 if (do_cpt_int(xd,"file_checksum_size",&(outputfiles[i].chksum_size),
956 list) != 0)
958 return -1;
960 if (do_cpt_u_chars(xd,"file_checksum",16,outputfiles[i].chksum,list) != 0)
962 return -1;
965 else
967 outputfiles[i].chksum_size = -1;
970 return 0;
974 void write_checkpoint(const char *fn,FILE *fplog,t_commrec *cr,
975 int eIntegrator,int simulation_part,
976 gmx_large_int_t step,double t,t_state *state)
978 int fp;
979 int file_version;
980 char *version;
981 char *btime;
982 char *buser;
983 char *bmach;
984 char *fprog;
985 char *ftime;
986 time_t now;
987 int nppnodes,npmenodes,flag_64bit;
988 char buf[1024];
989 gmx_file_position_t *outputfiles;
990 int noutputfiles;
991 int flags_eks,flags_enh,i;
993 if (PAR(cr))
995 if (DOMAINDECOMP(cr))
997 nppnodes = cr->dd->nnodes;
998 npmenodes = cr->npmenodes;
1000 else
1002 nppnodes = cr->nnodes;
1003 npmenodes = 0;
1006 else
1008 nppnodes = 1;
1009 npmenodes = 0;
1012 if (gmx_fexist(fn))
1014 /* Rename the previous checkpoint file */
1015 strcpy(buf,fn);
1016 buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1017 strcat(buf,"_prev");
1018 strcat(buf,fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1019 (void)remove(buf); /* Unix will overwrite buf if it exists, but for windows we need to remove first */
1020 if(rename(fn,buf) != 0)
1022 gmx_file("Cannot rename checkpoint file; maybe you are out of quota?");
1026 now = time(NULL);
1027 ftime = strdup(ctime(&now));
1028 ftime[strlen(ftime)-1] = '\0';
1030 /* No need to pollute stderr every time we write a checkpoint file */
1031 /* fprintf(stderr,"\nWriting checkpoint, step %d at %s\n",step,ftime); */
1032 if (fplog)
1034 fprintf(fplog,"Writing checkpoint, step %s at %s\n\n",
1035 gmx_step_str(step,buf),ftime);
1038 /* Get offsets for open files */
1039 gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1041 fp = gmx_fio_open(fn,"w");
1043 if (state->ekinstate.bUpToDate)
1045 flags_eks =
1046 ((1<<eeksEKINH_N) | (1<<eeksEKINH) | (1<<eeksDEKINDL) |
1047 (1<<eeksMVCOS));
1049 else
1051 flags_eks = 0;
1054 flags_enh = 0;
1055 if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1057 flags_enh |= (1<<eenhENERGY_N);
1058 if (state->enerhist.nsum > 0)
1060 flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1061 (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1063 if (state->enerhist.nsum_sim > 0)
1065 flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1066 (1<<eenhENERGY_NSUM_SIM));
1071 version = strdup(VERSION);
1072 btime = strdup(BUILD_TIME);
1073 buser = strdup(BUILD_USER);
1074 bmach = strdup(BUILD_MACHINE);
1075 fprog = strdup(Program());
1077 do_cpt_header(gmx_fio_getxdr(fp),FALSE,&file_version,
1078 &version,&btime,&buser,&bmach,&fprog,&ftime,
1079 &eIntegrator,&simulation_part,&step,&t,&nppnodes,
1080 DOMAINDECOMP(cr) ? cr->dd->nc : NULL,&npmenodes,
1081 &state->natoms,&state->ngtc,
1082 &state->flags,&flags_eks,&flags_enh,NULL);
1084 sfree(version);
1085 sfree(btime);
1086 sfree(buser);
1087 sfree(bmach);
1088 sfree(fprog);
1090 if( (do_cpt_state(gmx_fio_getxdr(fp),FALSE,state->flags,state,TRUE,NULL) < 0) ||
1091 (do_cpt_ekinstate(gmx_fio_getxdr(fp),FALSE,flags_eks,&state->ekinstate,NULL) < 0) ||
1092 (do_cpt_enerhist(gmx_fio_getxdr(fp),FALSE,flags_enh,&state->enerhist,NULL) < 0) ||
1093 (do_cpt_files(gmx_fio_getxdr(fp),FALSE,&outputfiles,&noutputfiles,NULL,file_version) < 0))
1095 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1098 do_cpt_footer(gmx_fio_getxdr(fp),FALSE,file_version);
1100 if( gmx_fio_close(fp) != 0)
1102 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1105 sfree(ftime);
1106 sfree(outputfiles);
1107 #ifdef GMX_FAHCORE
1108 /*code for alternate checkpointing scheme. moved from top of loop over steps */
1109 fcRequestCheckPoint();
1110 if ( fcCheckPointParallel( cr->nodeid, NULL,0) == 0 ) {
1111 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", step );
1113 #endif /* end FAHCORE block */
1117 static void print_flag_mismatch(FILE *fplog,int sflags,int fflags)
1119 int i;
1121 fprintf(fplog,"\nState entry mismatch between the simulation and the checkpoint file\n");
1122 fprintf(fplog,"Entries which are not present in the checkpoint file will not be updated\n");
1123 fprintf(fplog," %24s %11s %11s\n","","simulation","checkpoint");
1124 for(i=0; i<estNR; i++)
1126 if ((sflags & (1<<i)) || (fflags & (1<<i)))
1128 fprintf(fplog," %24s %11s %11s\n",
1129 est_names[i],
1130 (sflags & (1<<i)) ? " present " : "not present",
1131 (fflags & (1<<i)) ? " present " : "not present");
1136 static void check_int(FILE *fplog,const char *type,int p,int f,bool *mm)
1138 FILE *fp = fplog ? fplog : stderr;
1140 if (p != f)
1142 fprintf(fp," %s mismatch,\n",type);
1143 fprintf(fp," current program: %d\n",p);
1144 fprintf(fp," checkpoint file: %d\n",f);
1145 fprintf(fp,"\n");
1146 *mm = TRUE;
1150 static void check_string(FILE *fplog,const char *type,const char *p,
1151 const char *f,bool *mm)
1153 FILE *fp = fplog ? fplog : stderr;
1155 if (strcmp(p,f) != 0)
1157 fprintf(fp," %s mismatch,\n",type);
1158 fprintf(fp," current program: %s\n",p);
1159 fprintf(fp," checkpoint file: %s\n",f);
1160 fprintf(fp,"\n");
1161 *mm = TRUE;
1165 static void check_match(FILE *fplog,
1166 char *version,
1167 char *btime,char *buser,char *bmach,char *fprog,
1168 t_commrec *cr,bool bPartDecomp,int npp_f,int npme_f,
1169 ivec dd_nc,ivec dd_nc_f)
1171 int npp;
1172 bool mm;
1174 mm = FALSE;
1176 check_string(fplog,"Version" ,VERSION ,version,&mm);
1177 check_string(fplog,"Build time" ,BUILD_TIME ,btime ,&mm);
1178 check_string(fplog,"Build user" ,BUILD_USER ,buser ,&mm);
1179 check_string(fplog,"Build machine",BUILD_MACHINE,bmach ,&mm);
1180 check_string(fplog,"Program name" ,Program() ,fprog ,&mm);
1182 npp = cr->nnodes - cr->npmenodes;
1183 check_int (fplog,"#nodes" ,cr->nnodes ,npp_f+npme_f ,&mm);
1184 if (bPartDecomp)
1186 dd_nc[XX] = 1;
1187 dd_nc[YY] = 1;
1188 dd_nc[ZZ] = 1;
1190 if (npp > 1)
1192 check_int (fplog,"#PME-nodes" ,cr->npmenodes,npme_f ,&mm);
1193 if (npp == npp_f)
1195 check_int (fplog,"#DD-cells[x]",dd_nc[XX] ,dd_nc_f[XX],&mm);
1196 check_int (fplog,"#DD-cells[y]",dd_nc[YY] ,dd_nc_f[YY],&mm);
1197 check_int (fplog,"#DD-cells[z]",dd_nc[ZZ] ,dd_nc_f[ZZ],&mm);
1201 if (mm)
1203 fprintf(stderr,
1204 "Gromacs binary or parallel settings not identical to previous run.\n"
1205 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1206 fplog ? ",\n see the log file for details" : "");
1208 if (fplog)
1210 fprintf(fplog,
1211 "Gromacs binary or parallel settings not identical to previous run.\n"
1212 "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1217 static void read_checkpoint(const char *fn,FILE **pfplog,
1218 t_commrec *cr,bool bPartDecomp,ivec dd_nc,
1219 int eIntegrator,gmx_large_int_t *step,double *t,
1220 t_state *state,bool *bReadRNG,bool *bReadEkin,
1221 int *simulation_part,bool bAppendOutputFiles)
1223 int fp,i,j,rc;
1224 int file_version;
1225 char *version,*btime,*buser,*bmach,*fprog,*ftime;
1226 char filename[STRLEN],buf[22];
1227 int nppnodes,eIntegrator_f,nppnodes_f,npmenodes_f;
1228 ivec dd_nc_f;
1229 int natoms,ngtc,fflags,flags_eks,flags_enh;
1230 int d;
1231 int ret;
1232 gmx_file_position_t *outputfiles;
1233 int nfiles;
1234 int chksum_file;
1235 FILE* fplog = *pfplog;
1236 unsigned char digest[16];
1237 #if !((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1238 struct flock fl = { F_WRLCK, SEEK_SET, 0, 0, 0 };
1239 #endif
1241 const char *int_warn=
1242 "WARNING: The checkpoint file was generator with integrator %s,\n"
1243 " while the simulation uses integrator %s\n\n";
1244 const char *sd_note=
1245 "NOTE: The checkpoint file was for %d nodes doing SD or BD,\n"
1246 " while the simulation uses %d SD or BD nodes,\n"
1247 " continuation will be exact, except for the random state\n\n";
1249 if (PARTDECOMP(cr))
1251 gmx_fatal(FARGS,
1252 "read_checkpoint not (yet) supported with particle decomposition");
1255 fp = gmx_fio_open(fn,"r");
1256 do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1257 &version,&btime,&buser,&bmach,&fprog,&ftime,
1258 &eIntegrator_f,simulation_part,step,t,
1259 &nppnodes_f,dd_nc_f,&npmenodes_f,
1260 &natoms,&ngtc,
1261 &fflags,&flags_eks,&flags_enh,NULL);
1263 if (cr == NULL || MASTER(cr))
1265 fprintf(stderr,"\nReading checkpoint file %s generated: %s\n\n",
1266 fn,ftime);
1269 /* This will not be written if we do appending, since fplog is still NULL then */
1270 if (fplog)
1272 fprintf(fplog,"\n");
1273 fprintf(fplog,"Reading checkpoint file %s\n",fn);
1274 fprintf(fplog," file generated by: %s\n",fprog);
1275 fprintf(fplog," file generated at: %s\n",ftime);
1276 fprintf(fplog," GROMACS build time: %s\n",btime);
1277 fprintf(fplog," GROMACS build user: %s\n",buser);
1278 fprintf(fplog," GROMACS build machine: %s\n",bmach);
1279 fprintf(fplog," simulation part #: %d\n",*simulation_part);
1280 fprintf(fplog," step: %s\n",gmx_step_str(*step,buf));
1281 fprintf(fplog," time: %f\n",*t);
1282 fprintf(fplog,"\n");
1285 if (natoms != state->natoms)
1287 gmx_fatal(FARGS,"Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms",natoms,state->natoms);
1289 if (ngtc != state->ngtc)
1291 gmx_fatal(FARGS,"Checkpoint file is for a system of %d T-coupling groups, while the current system consists of %d T-coupling groups",ngtc,state->ngtc);
1293 if (eIntegrator_f != eIntegrator)
1295 if (MASTER(cr))
1297 fprintf(stderr,int_warn,EI(eIntegrator_f),EI(eIntegrator));
1299 if(bAppendOutputFiles)
1301 gmx_fatal(FARGS,
1302 "Output file appending requested, but input/checkpoint integrators do not match.\n"
1303 "Stopping the run to prevent you from ruining all your data...\n"
1304 "If you _really_ know what you are doing, try without the -append option.\n");
1306 if (fplog)
1308 fprintf(fplog,int_warn,EI(eIntegrator_f),EI(eIntegrator));
1312 if (!PAR(cr))
1314 nppnodes = 1;
1316 else if (bPartDecomp)
1318 nppnodes = cr->nnodes;
1320 else if (cr->nnodes == nppnodes_f + npmenodes_f)
1322 if (cr->npmenodes < 0)
1324 cr->npmenodes = npmenodes_f;
1326 nppnodes = cr->nnodes - cr->npmenodes;
1327 if (nppnodes == nppnodes_f)
1329 for(d=0; d<DIM; d++)
1331 if (dd_nc[d] == 0)
1333 dd_nc[d] = dd_nc_f[d];
1338 else
1340 /* The number of PP nodes has not been set yet */
1341 nppnodes = -1;
1344 if ((EI_SD(eIntegrator) || eIntegrator == eiBD) && nppnodes > 0)
1346 /* Correct the RNG state size for the number of PP nodes.
1347 * Such assignments should all be moved to one central function.
1349 state->nrng = nppnodes*gmx_rng_n();
1350 state->nrngi = nppnodes;
1353 *bReadRNG = TRUE;
1354 if (fflags != state->flags)
1357 if (MASTER(cr))
1359 if(bAppendOutputFiles)
1361 gmx_fatal(FARGS,
1362 "Output file appending requested, but input and checkpoint states are not identical.\n"
1363 "Stopping the run to prevent you from ruining all your data...\n"
1364 "You can try without the -append option, and get more info in the log file.\n");
1367 fprintf(stderr,
1368 "WARNING: The checkpoint state entries do not match the simulation,\n"
1369 " see the log file for details\n\n");
1372 if(fplog)
1374 print_flag_mismatch(fplog,state->flags,fflags);
1377 else
1379 if ((EI_SD(eIntegrator) || eIntegrator == eiBD) &&
1380 nppnodes != nppnodes_f)
1382 *bReadRNG = FALSE;
1383 if (MASTER(cr))
1385 fprintf(stderr,sd_note,nppnodes_f,nppnodes);
1387 if (fplog)
1389 fprintf(fplog ,sd_note,nppnodes_f,nppnodes);
1392 if (MASTER(cr))
1394 check_match(fplog,version,btime,buser,bmach,fprog,
1395 cr,bPartDecomp,nppnodes_f,npmenodes_f,dd_nc,dd_nc_f);
1398 ret = do_cpt_state(gmx_fio_getxdr(fp),TRUE,fflags,state,*bReadRNG,NULL);
1399 if (ret)
1401 cp_error();
1403 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1404 flags_eks,&state->ekinstate,NULL);
1405 if (ret)
1407 cp_error();
1409 *bReadEkin = (flags_eks & (1<<eeksEKINH));
1411 ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1412 flags_enh,&state->enerhist,NULL);
1413 if (ret)
1415 cp_error();
1418 if (file_version < 6)
1420 fprintf(stderr,"\nWARNING: Reading checkpoint file in old format, assuming that the run that generated this file started at step 0, if this is not the case the energy averages will be incorrect.\n\n");
1421 if (fplog)
1423 fprintf(fplog,"\nWARNING: Reading checkpoint file in old format, assuming that the run that generated this file started at step 0, if this is not the case the energy averages will be incorrect.\n\n");
1425 state->enerhist.nsum = *step;
1426 state->enerhist.nsum_sim = *step;
1429 ret = do_cpt_files(gmx_fio_getxdr(fp),TRUE,&outputfiles,&nfiles,NULL,file_version);
1430 if (ret)
1432 cp_error();
1435 ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1436 if (ret)
1438 cp_error();
1440 if( gmx_fio_close(fp) != 0)
1442 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1445 sfree(fprog);
1446 sfree(ftime);
1447 sfree(btime);
1448 sfree(buser);
1449 sfree(bmach);
1451 /* If the user wants to append to output files (danger), we use the file pointer
1452 * positions of the output files stored in the checkpoint file and truncate the
1453 * files such that any frames written after the checkpoint time are removed.
1455 * You will get REALLY fun problems if you use the -append option by provide
1456 * mdrun with other input files (in-frame truncation in the wrong places). Suit yourself!
1458 if (bAppendOutputFiles)
1460 if (fn2ftp(outputfiles[0].filename)!=efLOG)
1462 /* make sure first file is log file so that it is OK to use it for
1463 * locking
1465 gmx_fatal(FARGS,"The first output file should always be the log "
1466 "file but instead is: %s", outputfiles[0].filename);
1468 for(i=0;i<nfiles;i++)
1470 if (outputfiles[i].filename,outputfiles[i].offset < 0)
1472 gmx_fatal(FARGS,"The original run wrote a file called '%s' which "
1473 "is larger than 2 GB, but mdrun did not support large file"
1474 " offsets. Can not append. Run mdrun without -append",
1475 outputfiles[i].filename);
1477 #ifdef FAHCORE
1478 chksum_file=gmx_fio_open(outputfiles[i].filename,"a");
1480 #else
1481 chksum_file=gmx_fio_open(outputfiles[i].filename,"r+");
1483 /* lock log file */
1484 if (i==0)
1486 #if !((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1487 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl)
1488 ==-1)
1489 #else
1490 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX)==-1)
1491 #endif
1493 gmx_fatal(FARGS,"Failed to lock: %s. Already running "
1494 "simulation?", outputfiles[i].filename);
1498 /* compute md5 chksum */
1499 if (outputfiles[i].chksum_size != -1)
1501 if (gmx_fio_get_file_md5(chksum_file,outputfiles[i].offset,
1502 digest) != outputfiles[i].chksum_size)
1504 gmx_fatal(FARGS,"Can't read %d bytes of '%s' to compute"
1505 " checksum.", outputfiles[i].chksum_size,
1506 outputfiles[i].filename);
1509 else if (i==0) /*log file need to be seeked even when not reading md5*/
1511 gmx_fio_seek(chksum_file,outputfiles[i].offset);
1513 #endif
1515 if (i==0) /*open log file here - so that lock is never lifted
1516 after chksum is calculated */
1518 *pfplog = gmx_fio_getfp(chksum_file);
1520 else
1522 gmx_fio_close(chksum_file);
1524 #ifndef FAHCORE
1525 /* compare md5 chksum */
1526 if (outputfiles[i].chksum_size != -1 &&
1527 memcmp(digest,outputfiles[i].chksum,16)!=0)
1529 if (debug)
1531 fprintf(debug,"chksum for %s: ",outputfiles[i].filename);
1532 for (j=0; j<16; j++)
1534 fprintf(debug,"%02x",digest[j]);
1536 fprintf(debug,"\n");
1538 gmx_fatal(FARGS,"Checksum wrong for '%s'.",
1539 outputfiles[i].filename);
1541 #endif
1544 if (i!=0) /*log file is already seeked to correct position */
1546 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
1547 rc = gmx_wintruncate(outputfiles[i].filename,outputfiles[i].offset);
1548 #else
1549 rc = truncate(outputfiles[i].filename,outputfiles[i].offset);
1550 #endif
1551 if(rc!=0)
1553 gmx_fatal(FARGS,"Truncation of file %s failed.",outputfiles[i].filename);
1559 sfree(outputfiles);
1563 void load_checkpoint(const char *fn,FILE **fplog,
1564 t_commrec *cr,bool bPartDecomp,ivec dd_nc,
1565 t_inputrec *ir,t_state *state,
1566 bool *bReadRNG,bool *bReadEkin,bool bAppend)
1568 gmx_large_int_t step;
1569 double t;
1571 if (SIMMASTER(cr)) {
1572 /* Read the state from the checkpoint file */
1573 read_checkpoint(fn,fplog,
1574 cr,bPartDecomp,dd_nc,
1575 ir->eI,&step,&t,state,bReadRNG,bReadEkin,
1576 &ir->simulation_part,bAppend);
1578 if (PAR(cr)) {
1579 gmx_bcast(sizeof(cr->npmenodes),&cr->npmenodes,cr);
1580 gmx_bcast(DIM*sizeof(dd_nc[0]),dd_nc,cr);
1581 gmx_bcast(sizeof(step),&step,cr);
1582 gmx_bcast(sizeof(*bReadRNG),bReadRNG,cr);
1583 gmx_bcast(sizeof(*bReadEkin),bReadEkin,cr);
1585 ir->bContinuation = TRUE;
1586 ir->nsteps += ir->init_step - step;
1587 ir->init_step = step;
1588 ir->simulation_part += 1;
1591 static void low_read_checkpoint_state(int fp,int *simulation_part,
1592 gmx_large_int_t *step,double *t,t_state *state,
1593 bool bReadRNG)
1595 int file_version;
1596 char *version,*btime,*buser,*bmach,*fprog,*ftime;
1597 int eIntegrator;
1598 int nppnodes,npme;
1599 ivec dd_nc;
1600 int flags_eks,flags_enh;
1601 int ret;
1602 gmx_file_position_t *outputfiles;
1603 int nfiles;
1605 do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1606 &version,&btime,&buser,&bmach,&fprog,&ftime,
1607 &eIntegrator,simulation_part,step,t,&nppnodes,dd_nc,&npme,
1608 &state->natoms,&state->ngtc,
1609 &state->flags,&flags_eks,&flags_enh,NULL);
1610 ret =
1611 do_cpt_state(gmx_fio_getxdr(fp),TRUE,state->flags,state,bReadRNG,NULL);
1612 if (ret)
1614 cp_error();
1616 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1617 flags_eks,&state->ekinstate,NULL);
1618 if (ret)
1620 cp_error();
1622 ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1623 flags_enh,&state->enerhist,NULL);
1624 if (ret)
1626 cp_error();
1629 ret = do_cpt_files(gmx_fio_getxdr(fp),TRUE,&outputfiles,&nfiles,NULL,file_version);
1630 sfree(outputfiles);
1632 if (ret)
1634 cp_error();
1637 ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1638 if (ret)
1640 cp_error();
1643 sfree(fprog);
1644 sfree(ftime);
1645 sfree(btime);
1646 sfree(buser);
1647 sfree(bmach);
1650 void
1651 read_checkpoint_state(const char *fn,int *simulation_part,
1652 gmx_large_int_t *step,double *t,t_state *state)
1654 int fp;
1656 fp = gmx_fio_open(fn,"r");
1657 low_read_checkpoint_state(fp,simulation_part,step,t,state,TRUE);
1658 if( gmx_fio_close(fp) != 0)
1660 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1664 void read_checkpoint_trxframe(int fp,t_trxframe *fr)
1666 t_state state;
1667 int simulation_part;
1668 gmx_large_int_t step;
1669 double t;
1671 init_state(&state,0,0);
1673 low_read_checkpoint_state(fp,&simulation_part,&step,&t,&state,FALSE);
1675 fr->natoms = state.natoms;
1676 fr->bTitle = FALSE;
1677 fr->bStep = TRUE;
1678 fr->step = gmx_large_int_to_int(step,
1679 "conversion of checkpoint to trajectory");
1680 fr->bTime = TRUE;
1681 fr->time = t;
1682 fr->bLambda = TRUE;
1683 fr->lambda = state.lambda;
1684 fr->bAtoms = FALSE;
1685 fr->bX = (state.flags & (1<<estX));
1686 if (fr->bX)
1688 fr->x = state.x;
1689 state.x = NULL;
1691 fr->bV = (state.flags & (1<<estV));
1692 if (fr->bV)
1694 fr->v = state.v;
1695 state.v = NULL;
1697 fr->bF = FALSE;
1698 fr->bBox = (state.flags & (1<<estBOX));
1699 if (fr->bBox)
1701 copy_mat(state.box,fr->box);
1703 done_state(&state);
1706 void list_checkpoint(const char *fn,FILE *out)
1708 int fp;
1709 int file_version;
1710 char *version,*btime,*buser,*bmach,*fprog,*ftime;
1711 int eIntegrator,simulation_part,nppnodes,npme;
1712 gmx_large_int_t step;
1713 double t;
1714 ivec dd_nc;
1715 t_state state;
1716 int flags_eks,flags_enh;
1717 int indent;
1718 int i,j;
1719 int ret;
1720 gmx_file_position_t *outputfiles;
1721 int nfiles;
1723 init_state(&state,-1,-1);
1725 fp = gmx_fio_open(fn,"r");
1726 do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1727 &version,&btime,&buser,&bmach,&fprog,&ftime,
1728 &eIntegrator,&simulation_part,&step,&t,&nppnodes,dd_nc,&npme,
1729 &state.natoms,&state.ngtc,
1730 &state.flags,&flags_eks,&flags_enh,out);
1731 ret = do_cpt_state(gmx_fio_getxdr(fp),TRUE,state.flags,&state,TRUE,out);
1732 if (ret)
1734 cp_error();
1736 ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1737 flags_eks,&state.ekinstate,out);
1738 if (ret)
1740 cp_error();
1742 ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1743 flags_enh,&state.enerhist,out);
1745 if (ret == 0)
1747 do_cpt_files(gmx_fio_getxdr(fp),TRUE,&outputfiles,&nfiles,out,file_version);
1750 if (ret == 0)
1752 ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1755 if (ret)
1757 cp_warning(out);
1759 if( gmx_fio_close(fp) != 0)
1761 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1764 done_state(&state);
1768 /* This routine cannot print tons of data, since it is called before the log file is opened. */
1769 void read_checkpoint_simulation_part(const char *filename, int *simulation_part,
1770 gmx_large_int_t *cpt_step,t_commrec *cr)
1772 int fp;
1773 int file_version;
1774 char *version,*btime,*buser,*bmach,*fprog,*ftime;
1775 int eIntegrator_f,nppnodes_f,npmenodes_f;
1776 ivec dd_nc_f;
1777 gmx_large_int_t step=0;
1778 double t;
1779 int natoms,ngtc,fflags,flags_eks,flags_enh;
1781 if (SIMMASTER(cr)) {
1782 if(!gmx_fexist(filename) || ( (fp = gmx_fio_open(filename,"r")) < 0 ))
1784 *simulation_part = 0;
1786 else
1788 do_cpt_header(gmx_fio_getxdr(fp),
1789 TRUE,&file_version,
1790 &version,&btime,&buser,&bmach,&fprog,&ftime,
1791 &eIntegrator_f,simulation_part,
1792 &step,&t,&nppnodes_f,dd_nc_f,&npmenodes_f,
1793 &natoms,&ngtc,
1794 &fflags,&flags_eks,&flags_enh,NULL);
1795 if( gmx_fio_close(fp) != 0)
1797 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of quota?");
1800 sfree(version);
1801 sfree(btime);
1802 sfree(buser);
1803 sfree(bmach);
1804 sfree(fprog);
1805 sfree(ftime);
1808 if (PAR(cr)) {
1809 gmx_bcast(sizeof(*simulation_part),simulation_part,cr);
1811 if (NULL != cpt_step)
1813 *cpt_step = step;