Merge release-4-5-patches into release-4-6
[alexxy/gromacs.git] / src / gmxlib / checkpoint.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  * This file is part of Gromacs        Copyright (c) 1991-2008
5  * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
6  *
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.
11  *
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
14  * 
15  * And Hey:
16  * Gnomes, ROck Monsters And Chili Sauce
17  */
18
19 /* The source code in this file should be thread-safe. 
20  Please keep it that way. */
21
22
23 #ifdef HAVE_CONFIG_H
24 #include <config.h>
25 #endif
26 #include "gmx_header_config.h"
27
28 #include <string.h>
29 #include <time.h>
30
31 #ifdef HAVE_SYS_TIME_H
32 #include <sys/time.h>
33 #endif
34
35
36 #ifdef GMX_NATIVE_WINDOWS
37 /* _chsize_s */
38 #include <io.h>
39 #include <sys/locking.h>
40 #endif
41
42
43 #include "filenm.h"
44 #include "names.h"
45 #include "typedefs.h"
46 #include "smalloc.h"
47 #include "gmxfio.h"
48 #include "xdrf.h"
49 #include "statutil.h"
50 #include "txtdump.h"
51 #include "vec.h"
52 #include "network.h"
53 #include "gmx_random.h"
54 #include "checkpoint.h"
55 #include "futil.h"
56 #include "string2.h"
57 #include <fcntl.h>
58
59
60 #ifdef GMX_FAHCORE
61 #include "corewrap.h"
62 #endif
63
64
65 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
66 char *
67 gmx_ctime_r(const time_t *clock,char *buf, int n);
68
69
70 #define CPT_MAGIC1 171817
71 #define CPT_MAGIC2 171819
72
73 #ifdef GMX_DOUBLE
74 #define GMX_CPT_BUILD_DP 1
75 #else
76 #define GMX_CPT_BUILD_DP 0
77 #endif
78
79 /* cpt_version should normally only be changed
80  * when the header of footer format changes.
81  * The state data format itself is backward and forward compatible.
82  * But old code can not read a new entry that is present in the file
83  * (but can read a new format when new entries are not present).
84  */
85 static const int cpt_version = 13;
86
87
88 const char *est_names[estNR]=
89 {
90     "FE-lambda",
91     "box", "box-rel", "box-v", "pres_prev",
92     "nosehoover-xi", "thermostat-integral",
93     "x", "v", "SDx", "CGp", "LD-rng", "LD-rng-i",
94     "disre_initf", "disre_rm3tav",
95     "orire_initf", "orire_Dtav",
96     "svir_prev", "nosehoover-vxi", "v_eta", "vol0", "nhpres_xi", "nhpres_vxi", "fvir_prev",
97 };
98
99 enum { eeksEKIN_N, eeksEKINH, eeksDEKINDL, eeksMVCOS, eeksEKINF, eeksEKINO, eeksEKINSCALEF, eeksEKINSCALEH, eeksVSCALE, eeksEKINTOTAL, eeksNR };
100
101 const char *eeks_names[eeksNR]=
102 {
103     "Ekin_n", "Ekinh", "dEkindlambda", "mv_cos",
104     "Ekinf", "Ekinh_old", "EkinScaleF_NHC", "EkinScaleH_NHC","Vscale_NHC","Ekin_Total"
105 };
106
107 enum { eenhENERGY_N, eenhENERGY_AVER, eenhENERGY_SUM, eenhENERGY_NSUM,
108        eenhENERGY_SUM_SIM, eenhENERGY_NSUM_SIM,
109        eenhENERGY_NSTEPS, eenhENERGY_NSTEPS_SIM, 
110        eenhENERGY_DELTA_H_NN,
111        eenhENERGY_DELTA_H_LIST, 
112        eenhENERGY_DELTA_H_STARTTIME, 
113        eenhENERGY_DELTA_H_STARTLAMBDA, 
114        eenhNR };
115
116 const char *eenh_names[eenhNR]=
117 {
118     "energy_n", "energy_aver", "energy_sum", "energy_nsum",
119     "energy_sum_sim", "energy_nsum_sim",
120     "energy_nsteps", "energy_nsteps_sim", 
121     "energy_delta_h_nn",
122     "energy_delta_h_list", 
123     "energy_delta_h_start_time", 
124     "energy_delta_h_start_lambda"
125 };
126
127
128
129 #ifdef GMX_NATIVE_WINDOWS
130 static int
131 gmx_wintruncate(const char *filename, __int64 size)
132 {
133 #ifdef GMX_FAHCORE
134     /*we do this elsewhere*/
135     return 0;
136 #else
137     FILE *fp;
138     int   rc;
139     
140     fp=fopen(filename,"rb+");
141     
142     if(fp==NULL)
143     {
144         return -1;
145     }
146     
147     return _chsize_s( fileno(fp), size);
148 #endif
149 }
150 #endif
151
152
153 enum { ecprREAL, ecprRVEC, ecprMATRIX };
154
155 static const char *st_names(int cptp,int ecpt)
156 {
157     switch (cptp)
158     {
159     case 0: return est_names [ecpt]; break;
160     case 1: return eeks_names[ecpt]; break;
161     case 2: return eenh_names[ecpt]; break;
162     }
163
164     return NULL;
165 }
166
167 static void cp_warning(FILE *fp)
168 {
169     fprintf(fp,"\nWARNING: Checkpoint file is corrupted or truncated\n\n");
170 }
171
172 static void cp_error()
173 {
174     gmx_fatal(FARGS,"Checkpoint file corrupted/truncated, or maybe you are out of disk space?");
175 }
176
177 static void do_cpt_string_err(XDR *xd,gmx_bool bRead,const char *desc,char **s,FILE *list)
178 {
179 #define CPTSTRLEN 1024
180     bool_t res=0;
181     
182     if (bRead)
183     {
184         snew(*s,CPTSTRLEN);
185     }
186     res = xdr_string(xd,s,CPTSTRLEN);
187     if (res == 0)
188     {
189         cp_error();
190     }
191     if (list)
192     {
193         fprintf(list,"%s = %s\n",desc,*s);
194         sfree(*s);
195     }
196 }
197
198 static int do_cpt_int(XDR *xd,const char *desc,int *i,FILE *list)
199 {
200     bool_t res=0;
201     
202     res = xdr_int(xd,i);
203     if (res == 0)
204     {
205         return -1;
206     }
207     if (list)
208     {
209         fprintf(list,"%s = %d\n",desc,*i);
210     }
211     return 0;
212 }
213
214 static int do_cpt_u_chars(XDR *xd,const char *desc,int n,unsigned char *i,FILE *list)
215 {
216     bool_t res=1;
217     int j;
218     if (list)
219     {
220         fprintf(list,"%s = ",desc);
221     }
222     for (j=0; j<n && res; j++)
223     {
224         res &= xdr_u_char(xd,&i[j]);
225         if (list)
226         {
227             fprintf(list,"%02x",i[j]);
228         }
229     }
230     if (list)
231     {
232         fprintf(list,"\n");
233     }
234     if (res == 0)
235     {
236         return -1;
237     }
238
239     return 0;
240 }
241
242 static void do_cpt_int_err(XDR *xd,const char *desc,int *i,FILE *list)
243 {
244     if (do_cpt_int(xd,desc,i,list) < 0)
245     {
246         cp_error();
247     }
248 }
249
250 static void do_cpt_step_err(XDR *xd,const char *desc,gmx_large_int_t *i,FILE *list)
251 {
252     bool_t res=0;
253     char   buf[STEPSTRSIZE];
254
255     res = xdr_gmx_large_int(xd,i,"reading checkpoint file");
256     if (res == 0)
257     {
258         cp_error();
259     }
260     if (list)
261     {
262         fprintf(list,"%s = %s\n",desc,gmx_step_str(*i,buf));
263     }
264 }
265
266 static void do_cpt_double_err(XDR *xd,const char *desc,double *f,FILE *list)
267 {
268     bool_t res=0;
269     
270     res = xdr_double(xd,f);
271     if (res == 0)
272     {
273         cp_error();
274     }
275     if (list)
276     {
277         fprintf(list,"%s = %f\n",desc,*f);
278     }
279 }
280
281 /* If nval >= 0, nval is used; on read this should match the passed value.
282  * If nval n<0, *nptr is used; on read the value is stored in nptr
283  */
284 static int do_cpte_reals_low(XDR *xd,int cptp,int ecpt,int sflags,
285                              int nval,int *nptr,real **v,
286                              FILE *list,int erealtype)
287 {
288     bool_t res=0;
289 #ifndef GMX_DOUBLE
290     int  dtc=xdr_datatype_float; 
291 #else
292     int  dtc=xdr_datatype_double;
293 #endif
294     real *vp,*va=NULL;
295     float  *vf;
296     double *vd;
297     int  nf,dt,i;
298     
299     if (list == NULL)
300     {
301         if (nval >= 0)
302         {
303             nf = nval;
304         }
305         else
306         {
307         if (nptr == NULL)
308         {
309             gmx_incons("*ntpr=NULL in do_cpte_reals_low");
310         }
311         nf = *nptr;
312         }
313     }
314     res = xdr_int(xd,&nf);
315     if (res == 0)
316     {
317         return -1;
318     }
319     if (list == NULL)
320     {
321         if (nval >= 0)
322         {
323             if (nf != nval)
324             {
325                 gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),nval,nf);
326             }
327         }
328         else
329         {
330             *nptr = nf;
331         }
332     }
333     dt = dtc;
334     res = xdr_int(xd,&dt);
335     if (res == 0)
336     {
337         return -1;
338     }
339     if (dt != dtc)
340     {
341         fprintf(stderr,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
342                 st_names(cptp,ecpt),xdr_datatype_names[dtc],
343                 xdr_datatype_names[dt]);
344     }
345     if (list || !(sflags & (1<<ecpt)))
346     {
347         snew(va,nf);
348         vp = va;
349     }
350     else
351     {
352         if (*v == NULL)
353         {
354             snew(*v,nf);
355         }
356         vp = *v;
357     }
358     if (dt == xdr_datatype_float)
359     {
360         if (dtc == xdr_datatype_float)
361         {
362             vf = (float *)vp;
363         }
364         else
365         {
366             snew(vf,nf);
367         }
368         res = xdr_vector(xd,(char *)vf,nf,
369                          (unsigned int)sizeof(float),(xdrproc_t)xdr_float);
370         if (res == 0)
371         {
372             return -1;
373         }
374         if (dtc != xdr_datatype_float)
375         {
376             for(i=0; i<nf; i++)
377             {
378                 vp[i] = vf[i];
379             }
380             sfree(vf);
381         }
382     }
383     else
384     {
385         if (dtc == xdr_datatype_double)
386         {
387             vd = (double *)vp;
388         }
389         else
390         {
391             snew(vd,nf);
392         }
393         res = xdr_vector(xd,(char *)vd,nf,
394                          (unsigned int)sizeof(double),(xdrproc_t)xdr_double);
395         if (res == 0)
396         {
397             return -1;
398         }
399         if (dtc != xdr_datatype_double)
400         {
401             for(i=0; i<nf; i++)
402             {
403                 vp[i] = vd[i];
404             }
405             sfree(vd);
406         }
407     }
408     
409     if (list)
410     {
411         switch (erealtype)
412         {
413         case ecprREAL:
414             pr_reals(list,0,st_names(cptp,ecpt),vp,nf);
415             break;
416         case ecprRVEC:
417             pr_rvecs(list,0,st_names(cptp,ecpt),(rvec *)vp,nf/3);
418             break;
419         default:
420             gmx_incons("Unknown checkpoint real type");
421         }
422     }
423     if (va)
424     {
425         sfree(va);
426     }
427
428     return 0;
429 }
430
431
432 /* This function stores n along with the reals for reading,
433  * but on reading it assumes that n matches the value in the checkpoint file,
434  * a fatal error is generated when this is not the case.
435  */
436 static int do_cpte_reals(XDR *xd,int cptp,int ecpt,int sflags,
437                          int n,real **v,FILE *list)
438 {
439     return do_cpte_reals_low(xd,cptp,ecpt,sflags,n,NULL,v,list,ecprREAL);
440 }
441
442 /* This function does the same as do_cpte_reals,
443  * except that on reading it ignores the passed value of *n
444  * and stored the value read from the checkpoint file in *n.
445  */
446 static int do_cpte_n_reals(XDR *xd,int cptp,int ecpt,int sflags,
447                            int *n,real **v,FILE *list)
448 {
449     return do_cpte_reals_low(xd,cptp,ecpt,sflags,-1,n,v,list,ecprREAL);
450 }
451
452 static int do_cpte_real(XDR *xd,int cptp,int ecpt,int sflags,
453                         real *r,FILE *list)
454 {
455     int n;
456
457     return do_cpte_reals_low(xd,cptp,ecpt,sflags,1,NULL,&r,list,ecprREAL);
458 }
459
460 static int do_cpte_ints(XDR *xd,int cptp,int ecpt,int sflags,
461                         int n,int **v,FILE *list)
462 {
463     bool_t res=0;
464     int  dtc=xdr_datatype_int;
465     int *vp,*va=NULL;
466     int  nf,dt,i;
467     
468     nf = n;
469     res = xdr_int(xd,&nf);
470     if (res == 0)
471     {
472         return -1;
473     }
474     if (list == NULL && v != NULL && nf != n)
475     {
476         gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
477     }
478     dt = dtc;
479     res = xdr_int(xd,&dt);
480     if (res == 0)
481     {
482         return -1;
483     }
484     if (dt != dtc)
485     {
486         gmx_fatal(FARGS,"Type mismatch for state entry %s, code type is %s, file type is %s\n",
487                   st_names(cptp,ecpt),xdr_datatype_names[dtc],
488                   xdr_datatype_names[dt]);
489     }
490     if (list || !(sflags & (1<<ecpt)) || v == NULL)
491     {
492         snew(va,nf);
493         vp = va;
494     }
495     else
496     {
497         if (*v == NULL)
498         {
499             snew(*v,nf);
500         }
501         vp = *v;
502     }
503     res = xdr_vector(xd,(char *)vp,nf,
504                      (unsigned int)sizeof(int),(xdrproc_t)xdr_int);
505     if (res == 0)
506     {
507         return -1;
508     }
509     if (list)
510     {
511         pr_ivec(list,0,st_names(cptp,ecpt),vp,nf,TRUE);
512     }
513     if (va)
514     {
515         sfree(va);
516     }
517
518     return 0;
519 }
520
521 static int do_cpte_int(XDR *xd,int cptp,int ecpt,int sflags,
522                        int *i,FILE *list)
523 {
524     return do_cpte_ints(xd,cptp,ecpt,sflags,1,&i,list);
525 }
526
527 static int do_cpte_doubles(XDR *xd,int cptp,int ecpt,int sflags,
528                            int n,double **v,FILE *list)
529 {
530     bool_t res=0;
531     int  dtc=xdr_datatype_double;
532     double *vp,*va=NULL;
533     int  nf,dt,i;
534     
535     nf = n;
536     res = xdr_int(xd,&nf);
537     if (res == 0)
538     {
539         return -1;
540     }
541     if (list == NULL && nf != n)
542     {
543         gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
544     }
545     dt = dtc;
546     res = xdr_int(xd,&dt);
547     if (res == 0)
548     {
549         return -1;
550     }
551     if (dt != dtc)
552     {
553         gmx_fatal(FARGS,"Precision mismatch for state entry %s, code precision is %s, file precision is %s\n",
554                   st_names(cptp,ecpt),xdr_datatype_names[dtc],
555                   xdr_datatype_names[dt]);
556     }
557     if (list || !(sflags & (1<<ecpt)))
558     {
559         snew(va,nf);
560         vp = va;
561     }
562     else
563     {
564         if (*v == NULL)
565         {
566             snew(*v,nf);
567         }
568         vp = *v;
569     }
570     res = xdr_vector(xd,(char *)vp,nf,
571                      (unsigned int)sizeof(double),(xdrproc_t)xdr_double);
572     if (res == 0)
573     {
574         return -1;
575     }
576     if (list)
577     {
578         pr_doubles(list,0,st_names(cptp,ecpt),vp,nf);
579     }
580     if (va)
581     {
582         sfree(va);
583     }
584
585     return 0;
586 }
587
588 static int do_cpte_double(XDR *xd,int cptp,int ecpt,int sflags,
589                           double *r,FILE *list)
590 {
591     return do_cpte_doubles(xd,cptp,ecpt,sflags,1,&r,list);
592 }
593
594
595 static int do_cpte_rvecs(XDR *xd,int cptp,int ecpt,int sflags,
596                          int n,rvec **v,FILE *list)
597 {
598     int n3;
599
600     return do_cpte_reals_low(xd,cptp,ecpt,sflags,
601                              n*DIM,NULL,(real **)v,list,ecprRVEC);
602 }
603
604 static int do_cpte_matrix(XDR *xd,int cptp,int ecpt,int sflags,
605                           matrix v,FILE *list)
606 {
607     real *vr;
608     real ret;
609
610     vr = (real *)&(v[0][0]);
611     ret = do_cpte_reals_low(xd,cptp,ecpt,sflags,
612                             DIM*DIM,NULL,&vr,NULL,ecprMATRIX);
613     
614     if (list && ret == 0)
615     {
616         pr_rvecs(list,0,st_names(cptp,ecpt),v,DIM);
617     }
618     
619     return ret;
620 }
621
622 static int do_cpte_matrices(XDR *xd,int cptp,int ecpt,int sflags,
623                             int n,matrix **v,FILE *list)
624 {
625     bool_t res=0;
626     matrix *vp,*va=NULL;
627     real *vr;
628     int  nf,i,j,k;
629     int  ret;
630
631     nf = n;
632     res = xdr_int(xd,&nf);
633     if (res == 0)
634     {
635         return -1;
636     }
637     if (list == NULL && nf != n)
638     {
639         gmx_fatal(FARGS,"Count mismatch for state entry %s, code count is %d, file count is %d\n",st_names(cptp,ecpt),n,nf);
640     }
641     if (list || !(sflags & (1<<ecpt)))
642     {
643         snew(va,nf);
644         vp = va;
645     }
646     else
647     {
648         if (*v == NULL)
649         {
650             snew(*v,nf);
651         }
652         vp = *v;
653     }
654     snew(vr,nf*DIM*DIM);
655     for(i=0; i<nf; i++)
656     {
657         for(j=0; j<DIM; j++)
658         {
659             for(k=0; k<DIM; k++)
660             {
661                 vr[(i*DIM+j)*DIM+k] = vp[i][j][k];
662             }
663         }
664     }
665     ret = do_cpte_reals_low(xd,cptp,ecpt,sflags,
666                             nf*DIM*DIM,NULL,&vr,NULL,ecprMATRIX);
667     for(i=0; i<nf; i++)
668     {
669         for(j=0; j<DIM; j++)
670         {
671             for(k=0; k<DIM; k++)
672             {
673                 vp[i][j][k] = vr[(i*DIM+j)*DIM+k];
674             }
675         }
676     }
677     sfree(vr);
678     
679     if (list && ret == 0)
680     {
681         for(i=0; i<nf; i++)
682         {
683             pr_rvecs(list,0,st_names(cptp,ecpt),vp[i],DIM);
684         }
685     }
686     if (va)
687     {
688         sfree(va);
689     }
690     
691     return ret;
692 }
693
694 static void do_cpt_header(XDR *xd,gmx_bool bRead,int *file_version,
695                           char **version,char **btime,char **buser,char **bmach,
696                           int *double_prec,
697                           char **fprog,char **ftime,
698                           int *eIntegrator,int *simulation_part,
699                           gmx_large_int_t *step,double *t,
700                           int *nnodes,int *dd_nc,int *npme,
701                           int *natoms,int *ngtc, int *nnhpres, int *nhchainlength,
702                           int *flags_state,int *flags_eks,int *flags_enh,
703                           FILE *list)
704 {
705     bool_t res=0;
706     int  magic;
707     int  idum=0;
708     int  i;
709     char *fhost;
710
711     if (bRead)
712     {
713         magic = -1;
714     }
715     else
716     {
717         magic = CPT_MAGIC1;
718     }
719     res = xdr_int(xd,&magic);
720     if (res == 0)
721     {
722         gmx_fatal(FARGS,"The checkpoint file is empty/corrupted, or maybe you are out of disk space?");
723     }
724     if (magic != CPT_MAGIC1)
725     {
726         gmx_fatal(FARGS,"Start of file magic number mismatch, checkpoint file has %d, should be %d\n"
727                   "The checkpoint file is corrupted or not a checkpoint file",
728                   magic,CPT_MAGIC1);
729     }
730     if (!bRead)
731     {
732         snew(fhost,255);
733 #ifdef HAVE_UNISTD_H
734         if (gethostname(fhost,255) != 0)
735         {
736             sprintf(fhost,"unknown");
737         }
738 #else
739         sprintf(fhost,"unknown");
740 #endif  
741     }
742     do_cpt_string_err(xd,bRead,"GROMACS version"           ,version,list);
743     do_cpt_string_err(xd,bRead,"GROMACS build time"        ,btime,list);
744     do_cpt_string_err(xd,bRead,"GROMACS build user"        ,buser,list);
745     do_cpt_string_err(xd,bRead,"GROMACS build machine"     ,bmach,list);
746     do_cpt_string_err(xd,bRead,"generating program"        ,fprog,list);
747     do_cpt_string_err(xd,bRead,"generation time"           ,ftime,list);
748     *file_version = cpt_version;
749     do_cpt_int_err(xd,"checkpoint file version",file_version,list);
750     if (*file_version > cpt_version)
751     {
752         gmx_fatal(FARGS,"Attempting to read a checkpoint file of version %d with code of version %d\n",*file_version,cpt_version);
753     }
754     if (*file_version >= 13)
755     {
756         do_cpt_int_err(xd,"GROMACS double precision",double_prec,list);
757     }
758     else
759     {
760         *double_prec = -1;
761     }
762     if (*file_version >= 12)
763     {
764         do_cpt_string_err(xd,bRead,"generating host"           ,&fhost,list);
765         if (list == NULL)
766         {
767             sfree(fhost);
768         }
769     }
770     do_cpt_int_err(xd,"#atoms"            ,natoms     ,list);
771     do_cpt_int_err(xd,"#T-coupling groups",ngtc       ,list);
772     if (*file_version >= 10) 
773     {
774         do_cpt_int_err(xd,"#Nose-Hoover T-chains",nhchainlength,list);
775     }
776     else
777     {
778         *nhchainlength = 1;
779     }
780     if (*file_version >= 11)
781     {
782         do_cpt_int_err(xd,"#Nose-Hoover T-chains for barostat ",nnhpres,list);
783     }
784     else
785     {
786         *nnhpres = 0;
787     }
788     do_cpt_int_err(xd,"integrator"        ,eIntegrator,list);
789         if (*file_version >= 3)
790         {
791                 do_cpt_int_err(xd,"simulation part #", simulation_part,list);
792         }
793         else
794         {
795                 *simulation_part = 1;
796         }
797     if (*file_version >= 5)
798     {
799         do_cpt_step_err(xd,"step"         ,step       ,list);
800     }
801     else
802     {
803         do_cpt_int_err(xd,"step"          ,&idum      ,list);
804         *step = idum;
805     }
806     do_cpt_double_err(xd,"t"              ,t          ,list);
807     do_cpt_int_err(xd,"#PP-nodes"         ,nnodes     ,list);
808     idum = 1;
809     do_cpt_int_err(xd,"dd_nc[x]",dd_nc ? &(dd_nc[0]) : &idum,list);
810     do_cpt_int_err(xd,"dd_nc[y]",dd_nc ? &(dd_nc[1]) : &idum,list);
811     do_cpt_int_err(xd,"dd_nc[z]",dd_nc ? &(dd_nc[2]) : &idum,list);
812     do_cpt_int_err(xd,"#PME-only nodes",npme,list);
813     do_cpt_int_err(xd,"state flags",flags_state,list);
814         if (*file_version >= 4)
815     {
816         do_cpt_int_err(xd,"ekin data flags",flags_eks,list);
817         do_cpt_int_err(xd,"energy history flags",flags_enh,list);
818     }
819     else
820     {
821         *flags_eks  = 0;
822         *flags_enh   = (*flags_state >> (estORIRE_DTAV+1));
823         *flags_state = (*flags_state & ~((1<<(estORIRE_DTAV+1)) |
824                                          (1<<(estORIRE_DTAV+2)) |
825                                          (1<<(estORIRE_DTAV+3))));
826     }
827 }
828
829 static int do_cpt_footer(XDR *xd,gmx_bool bRead,int file_version)
830 {
831     bool_t res=0;
832     int  magic;
833     
834     if (file_version >= 2)
835     {
836         magic = CPT_MAGIC2;
837         res = xdr_int(xd,&magic);
838         if (res == 0)
839         {
840             cp_error();
841         }
842         if (magic != CPT_MAGIC2)
843         {
844             return -1;
845         }
846     }
847
848     return 0;
849 }
850
851 static int do_cpt_state(XDR *xd,gmx_bool bRead,
852                         int fflags,t_state *state,
853                         gmx_bool bReadRNG,FILE *list)
854 {
855     int  sflags;
856     int  **rng_p,**rngi_p;
857     int  i;
858     int  ret;
859     int  nnht,nnhtp;
860
861     ret = 0;
862     
863     nnht = state->nhchainlength*state->ngtc;
864     nnhtp = state->nhchainlength*state->nnhpres;
865
866     if (bReadRNG)
867     {
868         rng_p  = (int **)&state->ld_rng;
869         rngi_p = &state->ld_rngi;
870     }
871     else
872     {
873         /* Do not read the RNG data */
874         rng_p  = NULL;
875         rngi_p = NULL;
876     }
877
878     sflags = state->flags;
879     for(i=0; (i<estNR && ret == 0); i++)
880     {
881         if (fflags & (1<<i))
882         {
883             switch (i)
884             {
885             case estLAMBDA:  ret = do_cpte_real(xd,0,i,sflags,&state->lambda,list); break;
886             case estBOX:     ret = do_cpte_matrix(xd,0,i,sflags,state->box,list); break;
887             case estBOX_REL: ret = do_cpte_matrix(xd,0,i,sflags,state->box_rel,list); break;
888             case estBOXV:    ret = do_cpte_matrix(xd,0,i,sflags,state->boxv,list); break;
889             case estPRES_PREV: ret = do_cpte_matrix(xd,0,i,sflags,state->pres_prev,list); break;
890             case estSVIR_PREV:  ret = do_cpte_matrix(xd,0,i,sflags,state->svir_prev,list); break;
891             case estFVIR_PREV:  ret = do_cpte_matrix(xd,0,i,sflags,state->fvir_prev,list); break;
892             case estNH_XI:   ret = do_cpte_doubles(xd,0,i,sflags,nnht,&state->nosehoover_xi,list); break;
893             case estNH_VXI:  ret = do_cpte_doubles(xd,0,i,sflags,nnht,&state->nosehoover_vxi,list); break;
894             case estNHPRES_XI:   ret = do_cpte_doubles(xd,0,i,sflags,nnhtp,&state->nhpres_xi,list); break;
895             case estNHPRES_VXI:  ret = do_cpte_doubles(xd,0,i,sflags,nnhtp,&state->nhpres_vxi,list); break;
896             case estTC_INT:  ret = do_cpte_doubles(xd,0,i,sflags,state->ngtc,&state->therm_integral,list); break;
897             case estVETA:    ret = do_cpte_real(xd,0,i,sflags,&state->veta,list); break;
898             case estVOL0:    ret = do_cpte_real(xd,0,i,sflags,&state->vol0,list); break;
899             case estX:       ret = do_cpte_rvecs(xd,0,i,sflags,state->natoms,&state->x,list); break;
900             case estV:       ret = do_cpte_rvecs(xd,0,i,sflags,state->natoms,&state->v,list); break;
901             case estSDX:     ret = do_cpte_rvecs(xd,0,i,sflags,state->natoms,&state->sd_X,list); break;
902             case estLD_RNG:  ret = do_cpte_ints(xd,0,i,sflags,state->nrng,rng_p,list); break;
903             case estLD_RNGI: ret = do_cpte_ints(xd,0,i,sflags,state->nrngi,rngi_p,list); break;
904             case estDISRE_INITF:  ret = do_cpte_real (xd,0,i,sflags,&state->hist.disre_initf,list); break;
905             case estDISRE_RM3TAV: ret = do_cpte_n_reals(xd,0,i,sflags,&state->hist.ndisrepairs,&state->hist.disre_rm3tav,list); break;
906             case estORIRE_INITF:  ret = do_cpte_real (xd,0,i,sflags,&state->hist.orire_initf,list); break;
907             case estORIRE_DTAV:   ret = do_cpte_n_reals(xd,0,i,sflags,&state->hist.norire_Dtav,&state->hist.orire_Dtav,list); break;
908             default:
909                 gmx_fatal(FARGS,"Unknown state entry %d\n"
910                           "You are probably reading a new checkpoint file with old code",i);
911             }
912         }
913     }
914     
915     return ret;
916 }
917
918 static int do_cpt_ekinstate(XDR *xd,gmx_bool bRead,
919                             int fflags,ekinstate_t *ekins,
920                             FILE *list)
921 {
922     int  i;
923     int  ret;
924
925     ret = 0;
926
927     for(i=0; (i<eeksNR && ret == 0); i++)
928     {
929         if (fflags & (1<<i))
930         {
931             switch (i)
932             {
933                 
934                         case eeksEKIN_N:     ret = do_cpte_int(xd,1,i,fflags,&ekins->ekin_n,list); break;
935                         case eeksEKINH :     ret = do_cpte_matrices(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinh,list); break;
936                         case eeksEKINF:      ret = do_cpte_matrices(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinf,list); break;
937                         case eeksEKINO:      ret = do_cpte_matrices(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinh_old,list); break;
938             case eeksEKINTOTAL:  ret = do_cpte_matrix(xd,1,i,fflags,ekins->ekin_total,list); break;
939             case eeksEKINSCALEF: ret = do_cpte_doubles(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinscalef_nhc,list); break;
940             case eeksVSCALE:     ret = do_cpte_doubles(xd,1,i,fflags,ekins->ekin_n,&ekins->vscale_nhc,list); break;
941             case eeksEKINSCALEH: ret = do_cpte_doubles(xd,1,i,fflags,ekins->ekin_n,&ekins->ekinscaleh_nhc,list); break;
942                         case eeksDEKINDL :   ret = do_cpte_real(xd,1,i,fflags,&ekins->dekindl,list); break;
943             case eeksMVCOS:      ret = do_cpte_real(xd,1,i,fflags,&ekins->mvcos,list); break;                   
944             default:
945                 gmx_fatal(FARGS,"Unknown ekin data state entry %d\n"
946                           "You are probably reading a new checkpoint file with old code",i);
947             }
948         }
949     }
950     
951     return ret;
952 }
953
954
955 static int do_cpt_enerhist(XDR *xd,gmx_bool bRead,
956                            int fflags,energyhistory_t *enerhist,
957                            FILE *list)
958 {
959     int  i;
960     int  j;
961     int  ret;
962
963     ret = 0;
964
965     if (bRead)
966     {
967         enerhist->nsteps     = 0;
968         enerhist->nsum       = 0;
969         enerhist->nsteps_sim = 0;
970         enerhist->nsum_sim   = 0;
971         enerhist->dht        = NULL;
972
973         if (fflags & (1<< eenhENERGY_DELTA_H_NN) )
974         {
975             snew(enerhist->dht,1);
976             enerhist->dht->ndh = NULL;
977             enerhist->dht->dh = NULL;
978             enerhist->dht->start_lambda_set=FALSE;
979         }
980     }
981
982     for(i=0; (i<eenhNR && ret == 0); i++)
983     {
984         if (fflags & (1<<i))
985         {
986             switch (i)
987             {
988                 case eenhENERGY_N:     ret = do_cpte_int(xd,2,i,fflags,&enerhist->nener,list); break;
989                 case eenhENERGY_AVER:  ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_ave,list); break;
990                 case eenhENERGY_SUM:   ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_sum,list); break;
991                 case eenhENERGY_NSUM:  do_cpt_step_err(xd,eenh_names[i],&enerhist->nsum,list); break;
992                 case eenhENERGY_SUM_SIM: ret = do_cpte_doubles(xd,2,i,fflags,enerhist->nener,&enerhist->ener_sum_sim,list); break;
993                 case eenhENERGY_NSUM_SIM:   do_cpt_step_err(xd,eenh_names[i],&enerhist->nsum_sim,list); break;
994                 case eenhENERGY_NSTEPS:     do_cpt_step_err(xd,eenh_names[i],&enerhist->nsteps,list); break;
995                 case eenhENERGY_NSTEPS_SIM: do_cpt_step_err(xd,eenh_names[i],&enerhist->nsteps_sim,list); break;
996                 case eenhENERGY_DELTA_H_NN: do_cpt_int_err(xd,eenh_names[i], &(enerhist->dht->nndh), list); 
997                     if (bRead) /* now allocate memory for it */
998                     {
999                         snew(enerhist->dht->dh, enerhist->dht->nndh);
1000                         snew(enerhist->dht->ndh, enerhist->dht->nndh);
1001                         for(j=0;j<enerhist->dht->nndh;j++)
1002                         {
1003                             enerhist->dht->ndh[j] = 0;
1004                             enerhist->dht->dh[j] = NULL;
1005                         }
1006                     }
1007                 break;
1008                 case eenhENERGY_DELTA_H_LIST: 
1009                     for(j=0;j<enerhist->dht->nndh;j++)
1010                     {
1011                         ret=do_cpte_n_reals(xd, 2, i, fflags, &enerhist->dht->ndh[j], &(enerhist->dht->dh[j]), list); 
1012                     }
1013                     break;
1014                 case eenhENERGY_DELTA_H_STARTTIME: 
1015                     ret=do_cpte_double(xd, 2, i, fflags, &(enerhist->dht->start_time), list); break;
1016                 case eenhENERGY_DELTA_H_STARTLAMBDA: 
1017                     enerhist->dht->start_lambda_set=TRUE;
1018                     ret=do_cpte_double(xd, 2, i, fflags, &(enerhist->dht->start_lambda), list); break;
1019                 default:
1020                     gmx_fatal(FARGS,"Unknown energy history entry %d\n"
1021                               "You are probably reading a new checkpoint file with old code",i);
1022             }
1023         }
1024     }
1025
1026     if ((fflags & (1<<eenhENERGY_SUM)) && !(fflags & (1<<eenhENERGY_SUM_SIM)))
1027     {
1028         /* Assume we have an old file format and copy sum to sum_sim */
1029         srenew(enerhist->ener_sum_sim,enerhist->nener);
1030         for(i=0; i<enerhist->nener; i++)
1031         {
1032             enerhist->ener_sum_sim[i] = enerhist->ener_sum[i];
1033         }
1034         fflags |= (1<<eenhENERGY_SUM_SIM);
1035     }
1036     
1037     if ( (fflags & (1<<eenhENERGY_NSUM)) &&
1038         !(fflags & (1<<eenhENERGY_NSTEPS)))
1039     {
1040         /* Assume we have an old file format and copy nsum to nsteps */
1041         enerhist->nsteps = enerhist->nsum;
1042         fflags |= (1<<eenhENERGY_NSTEPS);
1043     }
1044     if ( (fflags & (1<<eenhENERGY_NSUM_SIM)) &&
1045         !(fflags & (1<<eenhENERGY_NSTEPS_SIM)))
1046     {
1047         /* Assume we have an old file format and copy nsum to nsteps */
1048         enerhist->nsteps_sim = enerhist->nsum_sim;
1049         fflags |= (1<<eenhENERGY_NSTEPS_SIM);
1050     }
1051
1052     return ret;
1053 }
1054
1055 static int do_cpt_files(XDR *xd, gmx_bool bRead, 
1056                         gmx_file_position_t **p_outputfiles, int *nfiles, 
1057                         FILE *list, int file_version)
1058 {
1059     int    i,j;
1060     gmx_off_t  offset;
1061     gmx_off_t  mask = 0xFFFFFFFFL;
1062     int    offset_high,offset_low;
1063     char   *buf;
1064     gmx_file_position_t *outputfiles;
1065
1066     if (do_cpt_int(xd,"number of output files",nfiles,list) != 0)
1067     {
1068         return -1;
1069     }
1070
1071     if(bRead)
1072     {
1073         snew(*p_outputfiles,*nfiles);
1074     }
1075
1076     outputfiles = *p_outputfiles;
1077
1078     for(i=0;i<*nfiles;i++)
1079     {
1080         /* 64-bit XDR numbers are not portable, so it is stored as separate high/low fractions */
1081         if(bRead)
1082         {
1083             do_cpt_string_err(xd,bRead,"output filename",&buf,list);
1084             strncpy(outputfiles[i].filename,buf,CPTSTRLEN-1);
1085             if(list==NULL)
1086             {
1087                 sfree(buf);                     
1088             }
1089
1090             if (do_cpt_int(xd,"file_offset_high",&offset_high,list) != 0)
1091             {
1092                 return -1;
1093             }
1094             if (do_cpt_int(xd,"file_offset_low",&offset_low,list) != 0)
1095             {
1096                 return -1;
1097             }
1098 #if (SIZEOF_GMX_OFF_T > 4)
1099             outputfiles[i].offset = ( ((gmx_off_t) offset_high) << 32 ) | ( (gmx_off_t) offset_low & mask );
1100 #else
1101             outputfiles[i].offset = offset_low;
1102 #endif
1103         }
1104         else
1105         {
1106             buf = outputfiles[i].filename;
1107             do_cpt_string_err(xd,bRead,"output filename",&buf,list);
1108             /* writing */
1109             offset      = outputfiles[i].offset;
1110             if (offset == -1)
1111             {
1112                 offset_low  = -1;
1113                 offset_high = -1;
1114             }
1115             else
1116             {
1117 #if (SIZEOF_GMX_OFF_T > 4)
1118                 offset_low  = (int) (offset & mask);
1119                 offset_high = (int) ((offset >> 32) & mask);
1120 #else
1121                 offset_low  = offset;
1122                 offset_high = 0;
1123 #endif
1124             }
1125             if (do_cpt_int(xd,"file_offset_high",&offset_high,list) != 0)
1126             {
1127                 return -1;
1128             }
1129             if (do_cpt_int(xd,"file_offset_low",&offset_low,list) != 0)
1130             {
1131                 return -1;
1132             }
1133         }
1134         if (file_version >= 8)
1135         {
1136             if (do_cpt_int(xd,"file_checksum_size",&(outputfiles[i].chksum_size),
1137                            list) != 0)
1138             {
1139                 return -1;
1140             }
1141             if (do_cpt_u_chars(xd,"file_checksum",16,outputfiles[i].chksum,list) != 0)
1142             {
1143                 return -1;
1144             }
1145         } 
1146         else 
1147         {
1148             outputfiles[i].chksum_size = -1;
1149         }
1150     }
1151     return 0;
1152 }
1153
1154
1155 void write_checkpoint(const char *fn,gmx_bool bNumberAndKeep,
1156                       FILE *fplog,t_commrec *cr,
1157                       int eIntegrator,int simulation_part,
1158                       gmx_large_int_t step,double t,t_state *state)
1159 {
1160     t_fileio *fp;
1161     int  file_version;
1162     char *version;
1163     char *btime;
1164     char *buser;
1165     char *bmach;
1166     int  double_prec;
1167     char *fprog;
1168     char *fntemp; /* the temporary checkpoint file name */
1169     time_t now;
1170     char timebuf[STRLEN];
1171     int  nppnodes,npmenodes,flag_64bit;
1172     char buf[1024],suffix[5+STEPSTRSIZE],sbuf[STEPSTRSIZE];
1173     gmx_file_position_t *outputfiles;
1174     int  noutputfiles;
1175     char *ftime;
1176     int  flags_eks,flags_enh,i;
1177     t_fileio *ret;
1178                 
1179     if (PAR(cr))
1180     {
1181         if (DOMAINDECOMP(cr))
1182         {
1183             nppnodes  = cr->dd->nnodes;
1184             npmenodes = cr->npmenodes;
1185         }
1186         else
1187         {
1188             nppnodes  = cr->nnodes;
1189             npmenodes = 0;
1190         }
1191     }
1192     else
1193     {
1194         nppnodes  = 1;
1195         npmenodes = 0;
1196     }
1197
1198     /* make the new temporary filename */
1199     snew(fntemp, strlen(fn)+5+STEPSTRSIZE);
1200     strcpy(fntemp,fn);
1201     fntemp[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1202     sprintf(suffix,"_%s%s","step",gmx_step_str(step,sbuf));
1203     strcat(fntemp,suffix);
1204     strcat(fntemp,fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1205    
1206     time(&now);
1207     gmx_ctime_r(&now,timebuf,STRLEN);
1208
1209     if (fplog)
1210     { 
1211         fprintf(fplog,"Writing checkpoint, step %s at %s\n\n",
1212                 gmx_step_str(step,buf),timebuf);
1213     }
1214     
1215     /* Get offsets for open files */
1216     gmx_fio_get_output_file_positions(&outputfiles, &noutputfiles);
1217
1218     fp = gmx_fio_open(fntemp,"w");
1219         
1220     if (state->ekinstate.bUpToDate)
1221     {
1222         flags_eks =
1223             ((1<<eeksEKIN_N) | (1<<eeksEKINH) | (1<<eeksEKINF) | 
1224              (1<<eeksEKINO) | (1<<eeksEKINSCALEF) | (1<<eeksEKINSCALEH) | 
1225              (1<<eeksVSCALE) | (1<<eeksDEKINDL) | (1<<eeksMVCOS));
1226     }
1227     else
1228     {
1229         flags_eks = 0;
1230     }
1231
1232     flags_enh = 0;
1233     if (state->enerhist.nsum > 0 || state->enerhist.nsum_sim > 0)
1234     {
1235         flags_enh |= (1<<eenhENERGY_N);
1236         if (state->enerhist.nsum > 0)
1237         {
1238             flags_enh |= ((1<<eenhENERGY_AVER) | (1<<eenhENERGY_SUM) |
1239                           (1<<eenhENERGY_NSTEPS) | (1<<eenhENERGY_NSUM));
1240         }
1241         if (state->enerhist.nsum_sim > 0)
1242         {
1243             flags_enh |= ((1<<eenhENERGY_SUM_SIM) | (1<<eenhENERGY_NSTEPS_SIM) |
1244                           (1<<eenhENERGY_NSUM_SIM));
1245         }
1246         if (state->enerhist.dht)
1247         {
1248             flags_enh |= ( (1<< eenhENERGY_DELTA_H_NN) |
1249                            (1<< eenhENERGY_DELTA_H_LIST) | 
1250                            (1<< eenhENERGY_DELTA_H_STARTTIME) |
1251                            (1<< eenhENERGY_DELTA_H_STARTLAMBDA) );
1252         }
1253     }
1254
1255     
1256     version = gmx_strdup(VERSION);
1257     btime   = gmx_strdup(BUILD_TIME);
1258     buser   = gmx_strdup(BUILD_USER);
1259     bmach   = gmx_strdup(BUILD_MACHINE);
1260     double_prec = GMX_CPT_BUILD_DP;
1261     fprog   = gmx_strdup(Program());
1262
1263     ftime   = &(timebuf[0]);
1264     
1265     do_cpt_header(gmx_fio_getxdr(fp),FALSE,&file_version,
1266                   &version,&btime,&buser,&bmach,&double_prec,&fprog,&ftime,
1267                   &eIntegrator,&simulation_part,&step,&t,&nppnodes,
1268                   DOMAINDECOMP(cr) ? cr->dd->nc : NULL,&npmenodes,
1269                   &state->natoms,&state->ngtc,&state->nnhpres,
1270                   &state->nhchainlength, &state->flags,&flags_eks,&flags_enh,
1271                   NULL);
1272     
1273     sfree(version);
1274     sfree(btime);
1275     sfree(buser);
1276     sfree(bmach);
1277     sfree(fprog);
1278
1279     if((do_cpt_state(gmx_fio_getxdr(fp),FALSE,state->flags,state,TRUE,NULL) < 0)        ||
1280        (do_cpt_ekinstate(gmx_fio_getxdr(fp),FALSE,flags_eks,&state->ekinstate,NULL) < 0)||
1281        (do_cpt_enerhist(gmx_fio_getxdr(fp),FALSE,flags_enh,&state->enerhist,NULL) < 0)  ||
1282        (do_cpt_files(gmx_fio_getxdr(fp),FALSE,&outputfiles,&noutputfiles,NULL,
1283                      file_version) < 0))
1284     {
1285         gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1286     }
1287
1288     do_cpt_footer(gmx_fio_getxdr(fp),FALSE,file_version);
1289
1290     /* we really, REALLY, want to make sure to physically write the checkpoint, 
1291        and all the files it depends on, out to disk. Because we've
1292        opened the checkpoint with gmx_fio_open(), it's in our list
1293        of open files.  */
1294     ret=gmx_fio_all_output_fsync();
1295
1296     if (ret)
1297     {
1298         char buf[STRLEN];
1299         sprintf(buf,
1300                 "Cannot fsync '%s'; maybe you are out of disk space?",
1301                 gmx_fio_getname(ret));
1302
1303         if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV)==NULL)
1304         {
1305             gmx_file(buf);
1306         }
1307         else
1308         {
1309             gmx_warning(buf);
1310         }
1311     }
1312
1313     if( gmx_fio_close(fp) != 0)
1314     {
1315         gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1316     }
1317
1318     /* we don't move the checkpoint if the user specified they didn't want it,
1319        or if the fsyncs failed */
1320     if (!bNumberAndKeep && !ret)
1321     {
1322         if (gmx_fexist(fn))
1323         {
1324             /* Rename the previous checkpoint file */
1325             strcpy(buf,fn);
1326             buf[strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
1327             strcat(buf,"_prev");
1328             strcat(buf,fn+strlen(fn) - strlen(ftp2ext(fn2ftp(fn))) - 1);
1329 #ifndef GMX_FAHCORE
1330             /* we copy here so that if something goes wrong between now and
1331              * the rename below, there's always a state.cpt.
1332              * If renames are atomic (such as in POSIX systems),
1333              * this copying should be unneccesary.
1334              */
1335             gmx_file_copy(fn, buf, FALSE);
1336             /* We don't really care if this fails: 
1337              * there's already a new checkpoint.
1338              */
1339 #else
1340             gmx_file_rename(fn, buf);
1341 #endif
1342         }
1343         if (gmx_file_rename(fntemp, fn) != 0)
1344         {
1345             gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
1346         }
1347     }
1348
1349     sfree(outputfiles);
1350     sfree(fntemp);
1351
1352 #ifdef GMX_FAHCORE
1353     /*code for alternate checkpointing scheme.  moved from top of loop over 
1354       steps */
1355     fcRequestCheckPoint();
1356     if ( fcCheckPointParallel( cr->nodeid, NULL,0) == 0 ) {
1357         gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", step );
1358     }
1359 #endif /* end GMX_FAHCORE block */
1360 }
1361
1362 static void print_flag_mismatch(FILE *fplog,int sflags,int fflags)
1363 {
1364     int i;
1365     
1366     fprintf(fplog,"\nState entry mismatch between the simulation and the checkpoint file\n");
1367     fprintf(fplog,"Entries which are not present in the checkpoint file will not be updated\n");
1368     fprintf(fplog,"  %24s    %11s    %11s\n","","simulation","checkpoint");
1369     for(i=0; i<estNR; i++)
1370     {
1371         if ((sflags & (1<<i)) || (fflags & (1<<i)))
1372         {
1373             fprintf(fplog,"  %24s    %11s    %11s\n",
1374                     est_names[i],
1375                     (sflags & (1<<i)) ? "  present  " : "not present",
1376                     (fflags & (1<<i)) ? "  present  " : "not present");
1377         }
1378     }
1379 }
1380
1381 static void check_int(FILE *fplog,const char *type,int p,int f,gmx_bool *mm)
1382 {
1383         FILE *fp = fplog ? fplog : stderr;
1384
1385     if (p != f)
1386     {
1387                 fprintf(fp,"  %s mismatch,\n",type);
1388                 fprintf(fp,"    current program: %d\n",p);
1389                 fprintf(fp,"    checkpoint file: %d\n",f);
1390                 fprintf(fp,"\n");
1391         *mm = TRUE;
1392     }
1393 }
1394
1395 static void check_string(FILE *fplog,const char *type,const char *p,
1396                          const char *f,gmx_bool *mm)
1397 {
1398         FILE *fp = fplog ? fplog : stderr;
1399         
1400     if (strcmp(p,f) != 0)
1401     {
1402                 fprintf(fp,"  %s mismatch,\n",type);
1403                 fprintf(fp,"    current program: %s\n",p);
1404                 fprintf(fp,"    checkpoint file: %s\n",f);
1405                 fprintf(fp,"\n");
1406         *mm = TRUE;
1407     }
1408 }
1409
1410 static void check_match(FILE *fplog,
1411                         char *version,
1412                         char *btime,char *buser,char *bmach,int double_prec,
1413                         char *fprog,
1414                         t_commrec *cr,gmx_bool bPartDecomp,int npp_f,int npme_f,
1415                         ivec dd_nc,ivec dd_nc_f)
1416 {
1417     int  npp;
1418     gmx_bool mm;
1419     
1420     mm = FALSE;
1421     
1422     check_string(fplog,"Version"      ,VERSION      ,version,&mm);
1423     check_string(fplog,"Build time"   ,BUILD_TIME   ,btime  ,&mm);
1424     check_string(fplog,"Build user"   ,BUILD_USER   ,buser  ,&mm);
1425     check_string(fplog,"Build machine",BUILD_MACHINE,bmach  ,&mm);
1426     check_int   (fplog,"Double prec." ,GMX_CPT_BUILD_DP,double_prec,&mm);
1427     check_string(fplog,"Program name" ,Program()    ,fprog  ,&mm);
1428     
1429     check_int   (fplog,"#nodes"       ,cr->nnodes   ,npp_f+npme_f ,&mm);
1430     if (bPartDecomp)
1431     {
1432         dd_nc[XX] = 1;
1433         dd_nc[YY] = 1;
1434         dd_nc[ZZ] = 1;
1435     }
1436     if (cr->nnodes > 1)
1437     {
1438         check_int (fplog,"#PME-nodes"  ,cr->npmenodes,npme_f     ,&mm);
1439
1440         npp = cr->nnodes;
1441         if (cr->npmenodes >= 0)
1442         {
1443             npp -= cr->npmenodes;
1444         }
1445         if (npp == npp_f)
1446         {
1447             check_int (fplog,"#DD-cells[x]",dd_nc[XX]    ,dd_nc_f[XX],&mm);
1448             check_int (fplog,"#DD-cells[y]",dd_nc[YY]    ,dd_nc_f[YY],&mm);
1449             check_int (fplog,"#DD-cells[z]",dd_nc[ZZ]    ,dd_nc_f[ZZ],&mm);
1450         }
1451     }
1452     
1453     if (mm)
1454     {
1455                 fprintf(stderr,
1456                                 "Gromacs binary or parallel settings not identical to previous run.\n"
1457                                 "Continuation is exact, but is not guaranteed to be binary identical%s.\n\n",
1458                                 fplog ? ",\n see the log file for details" : "");
1459                 
1460         if (fplog)
1461         {
1462                         fprintf(fplog,
1463                                         "Gromacs binary or parallel settings not identical to previous run.\n"
1464                                         "Continuation is exact, but is not guaranteed to be binary identical.\n\n");
1465                 }
1466     }
1467 }
1468
1469 static void read_checkpoint(const char *fn,FILE **pfplog,
1470                             t_commrec *cr,gmx_bool bPartDecomp,ivec dd_nc,
1471                             int eIntegrator,gmx_large_int_t *step,double *t,
1472                             t_state *state,gmx_bool *bReadRNG,gmx_bool *bReadEkin,
1473                             int *simulation_part,
1474                             gmx_bool bAppendOutputFiles,gmx_bool bForceAppend)
1475 {
1476     t_fileio *fp;
1477     int  i,j,rc;
1478     int  file_version;
1479     char *version,*btime,*buser,*bmach,*fprog,*ftime;
1480     int  double_prec;
1481         char filename[STRLEN],buf[STEPSTRSIZE];
1482     int  nppnodes,eIntegrator_f,nppnodes_f,npmenodes_f;
1483     ivec dd_nc_f;
1484     int  natoms,ngtc,nnhpres,nhchainlength,fflags,flags_eks,flags_enh;
1485     int  d;
1486     int  ret;
1487     gmx_file_position_t *outputfiles;
1488     int  nfiles;
1489     t_fileio *chksum_file;
1490     FILE* fplog = *pfplog;
1491     unsigned char digest[16];
1492 #ifndef GMX_NATIVE_WINDOWS
1493     struct flock fl;  /* don't initialize here: the struct order is OS 
1494                          dependent! */
1495 #endif
1496
1497     const char *int_warn=
1498               "WARNING: The checkpoint file was generated with integrator %s,\n"
1499               "         while the simulation uses integrator %s\n\n";
1500     const char *sd_note=
1501         "NOTE: The checkpoint file was for %d nodes doing SD or BD,\n"
1502         "      while the simulation uses %d SD or BD nodes,\n"
1503         "      continuation will be exact, except for the random state\n\n";
1504     
1505 #ifndef GMX_NATIVE_WINDOWS
1506     fl.l_type=F_WRLCK;
1507     fl.l_whence=SEEK_SET;
1508     fl.l_start=0;
1509     fl.l_len=0;
1510     fl.l_pid=0;
1511 #endif
1512
1513     if (PARTDECOMP(cr))
1514     {
1515         gmx_fatal(FARGS,
1516                   "read_checkpoint not (yet) supported with particle decomposition");
1517     }
1518     
1519     fp = gmx_fio_open(fn,"r");
1520     do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1521                   &version,&btime,&buser,&bmach,&double_prec,&fprog,&ftime,
1522                   &eIntegrator_f,simulation_part,step,t,
1523                   &nppnodes_f,dd_nc_f,&npmenodes_f,
1524                   &natoms,&ngtc,&nnhpres,&nhchainlength,
1525                   &fflags,&flags_eks,&flags_enh,NULL);
1526
1527     if (bAppendOutputFiles &&
1528         file_version >= 13 && double_prec != GMX_CPT_BUILD_DP)
1529     {
1530         gmx_fatal(FARGS,"Output file appending requested, but the code and checkpoint file precision (single/double) don't match");
1531     }
1532     
1533     if (cr == NULL || MASTER(cr))
1534     {
1535         fprintf(stderr,"\nReading checkpoint file %s generated: %s\n\n",
1536                 fn,ftime);
1537     }
1538         
1539         /* This will not be written if we do appending, since fplog is still NULL then */
1540     if (fplog)
1541     {
1542         fprintf(fplog,"\n");
1543         fprintf(fplog,"Reading checkpoint file %s\n",fn);
1544         fprintf(fplog,"  file generated by:     %s\n",fprog);  
1545         fprintf(fplog,"  file generated at:     %s\n",ftime);  
1546         fprintf(fplog,"  GROMACS build time:    %s\n",btime);  
1547         fprintf(fplog,"  GROMACS build user:    %s\n",buser);  
1548         fprintf(fplog,"  GROMACS build machine: %s\n",bmach);  
1549         fprintf(fplog,"  GROMACS double prec.:  %d\n",double_prec);
1550         fprintf(fplog,"  simulation part #:     %d\n",*simulation_part);
1551         fprintf(fplog,"  step:                  %s\n",gmx_step_str(*step,buf));
1552         fprintf(fplog,"  time:                  %f\n",*t);  
1553         fprintf(fplog,"\n");
1554     }
1555     
1556     if (natoms != state->natoms)
1557     {
1558         gmx_fatal(FARGS,"Checkpoint file is for a system of %d atoms, while the current system consists of %d atoms",natoms,state->natoms);
1559     }
1560     if (ngtc != state->ngtc)
1561     {
1562         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);
1563     }
1564     if (nnhpres != state->nnhpres)
1565     {
1566         gmx_fatal(FARGS,"Checkpoint file is for a system of %d NH-pressure-coupling variables, while the current system consists of %d NH-pressure-coupling variables",nnhpres,state->nnhpres);
1567     }
1568
1569     init_gtc_state(state,state->ngtc,state->nnhpres,nhchainlength); /* need to keep this here to keep the tpr format working */
1570     /* write over whatever was read; we use the number of Nose-Hoover chains from the checkpoint */
1571     
1572     if (eIntegrator_f != eIntegrator)
1573     {
1574         if (MASTER(cr))
1575         {
1576             fprintf(stderr,int_warn,EI(eIntegrator_f),EI(eIntegrator));
1577         }
1578                 if(bAppendOutputFiles)
1579                 {
1580                         gmx_fatal(FARGS,
1581                                           "Output file appending requested, but input/checkpoint integrators do not match.\n"
1582                                           "Stopping the run to prevent you from ruining all your data...\n"
1583                                           "If you _really_ know what you are doing, try with the -noappend option.\n");
1584                 }
1585         if (fplog)
1586         {
1587             fprintf(fplog,int_warn,EI(eIntegrator_f),EI(eIntegrator));
1588         }
1589     }
1590
1591     if (!PAR(cr))
1592     {
1593         nppnodes = 1;
1594         cr->npmenodes = 0;
1595     }
1596     else if (bPartDecomp)
1597     {
1598         nppnodes = cr->nnodes;
1599         cr->npmenodes = 0;
1600     }
1601     else if (cr->nnodes == nppnodes_f + npmenodes_f)
1602     {
1603         if (cr->npmenodes < 0)
1604         {
1605             cr->npmenodes = npmenodes_f;
1606         }
1607         nppnodes = cr->nnodes - cr->npmenodes;
1608         if (nppnodes == nppnodes_f)
1609         {
1610             for(d=0; d<DIM; d++)
1611             {
1612                 if (dd_nc[d] == 0)
1613                 {
1614                     dd_nc[d] = dd_nc_f[d];
1615                 }
1616             }
1617         }
1618     }
1619     else
1620     {
1621         /* The number of PP nodes has not been set yet */
1622         nppnodes = -1;
1623     }
1624
1625     if ((EI_SD(eIntegrator) || eIntegrator == eiBD) && nppnodes > 0)
1626     {
1627         /* Correct the RNG state size for the number of PP nodes.
1628          * Such assignments should all be moved to one central function.
1629          */
1630         state->nrng  = nppnodes*gmx_rng_n();
1631         state->nrngi = nppnodes;
1632     }
1633     
1634     *bReadRNG = TRUE;
1635     if (fflags != state->flags)
1636     {
1637                 
1638         if (MASTER(cr))
1639         {
1640                         if(bAppendOutputFiles)
1641                         {
1642                                 gmx_fatal(FARGS,
1643                                                   "Output file appending requested, but input and checkpoint states are not identical.\n"
1644                                                   "Stopping the run to prevent you from ruining all your data...\n"
1645                                                   "You can try with the -noappend option, and get more info in the log file.\n");
1646                         }
1647                         
1648             if (getenv("GMX_ALLOW_CPT_MISMATCH") == NULL)
1649             {
1650                 gmx_fatal(FARGS,"You seem to have switched ensemble, integrator, T and/or P-coupling algorithm between the cpt and tpr file. The recommended way of doing this is passing the cpt file to grompp (with option -t) instead of to mdrun. If you know what you are doing, you can override this error by setting the env.var. GMX_ALLOW_CPT_MISMATCH");
1651             }
1652             else
1653             {
1654                 fprintf(stderr,
1655                         "WARNING: The checkpoint state entries do not match the simulation,\n"
1656                         "         see the log file for details\n\n");
1657             }
1658         }
1659                 
1660                 if(fplog)
1661                 {
1662                         print_flag_mismatch(fplog,state->flags,fflags);
1663                 }
1664     }
1665     else
1666     {
1667         if ((EI_SD(eIntegrator) || eIntegrator == eiBD) &&
1668             nppnodes != nppnodes_f)
1669         {
1670             *bReadRNG = FALSE;
1671             if (MASTER(cr))
1672             {
1673                 fprintf(stderr,sd_note,nppnodes_f,nppnodes);
1674             }
1675             if (fplog)
1676             {
1677                 fprintf(fplog ,sd_note,nppnodes_f,nppnodes);
1678             }
1679         }
1680         if (MASTER(cr))
1681         {
1682             check_match(fplog,version,btime,buser,bmach,double_prec,fprog,
1683                         cr,bPartDecomp,nppnodes_f,npmenodes_f,dd_nc,dd_nc_f);
1684         }
1685     }
1686     ret = do_cpt_state(gmx_fio_getxdr(fp),TRUE,fflags,state,*bReadRNG,NULL);
1687     if (ret)
1688     {
1689         cp_error();
1690     }
1691     ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1692                            flags_eks,&state->ekinstate,NULL);
1693     if (ret)
1694     {
1695         cp_error();
1696     }
1697     *bReadEkin = ((flags_eks & (1<<eeksEKINH)) || (flags_eks & (1<<eeksEKINF)) || (flags_eks & (1<<eeksEKINO)) ||
1698                   ((flags_eks & (1<<eeksEKINSCALEF)) | (flags_eks & (1<<eeksEKINSCALEH)) | (flags_eks & (1<<eeksVSCALE))));
1699     
1700     ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1701                           flags_enh,&state->enerhist,NULL);
1702     if (ret)
1703     {
1704         cp_error();
1705     }
1706
1707     if (file_version < 6)
1708     {
1709         const char *warn="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 averages stored in the energy file will be incorrect.";
1710
1711         fprintf(stderr,"\nWARNING: %s\n\n",warn);
1712         if (fplog)
1713         {
1714             fprintf(fplog,"\nWARNING: %s\n\n",warn);
1715         }
1716         state->enerhist.nsum     = *step;
1717         state->enerhist.nsum_sim = *step;
1718     }
1719
1720         ret = do_cpt_files(gmx_fio_getxdr(fp),TRUE,&outputfiles,&nfiles,NULL,file_version);
1721         if (ret)
1722         {
1723                 cp_error();
1724         }
1725                                            
1726     ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1727     if (ret)
1728     {
1729         cp_error();
1730     }
1731     if( gmx_fio_close(fp) != 0)
1732         {
1733         gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1734         }
1735     
1736     sfree(fprog);
1737     sfree(ftime);
1738     sfree(btime);
1739     sfree(buser);
1740     sfree(bmach);
1741         
1742         /* If the user wants to append to output files,
1743      * we use the file pointer positions of the output files stored
1744      * in the checkpoint file and truncate the files such that any frames
1745      * written after the checkpoint time are removed.
1746      * All files are md5sum checked such that we can be sure that
1747      * we do not truncate other (maybe imprortant) files.
1748          */
1749     if (bAppendOutputFiles)
1750     {
1751         if (fn2ftp(outputfiles[0].filename)!=efLOG)
1752         {
1753             /* make sure first file is log file so that it is OK to use it for 
1754              * locking
1755              */
1756             gmx_fatal(FARGS,"The first output file should always be the log "
1757                       "file but instead is: %s. Cannot do appending because of this condition.", outputfiles[0].filename);
1758         }
1759         for(i=0;i<nfiles;i++)
1760         {
1761             if (outputfiles[i].offset < 0)
1762             {
1763                 gmx_fatal(FARGS,"The original run wrote a file called '%s' which "
1764                     "is larger than 2 GB, but mdrun did not support large file"
1765                     " offsets. Can not append. Run mdrun with -noappend",
1766                     outputfiles[i].filename);
1767             }
1768 #ifdef GMX_FAHCORE
1769             chksum_file=gmx_fio_open(outputfiles[i].filename,"a");
1770
1771 #else
1772             chksum_file=gmx_fio_open(outputfiles[i].filename,"r+");
1773
1774             /* lock log file */                
1775             if (i==0)
1776             {
1777                 /* Note that there are systems where the lock operation
1778                  * will succeed, but a second process can also lock the file.
1779                  * We should probably try to detect this.
1780                  */
1781 #ifndef GMX_NATIVE_WINDOWS
1782                 if (fcntl(fileno(gmx_fio_getfp(chksum_file)), F_SETLK, &fl)
1783                     ==-1)
1784 #else
1785                 if (_locking(fileno(gmx_fio_getfp(chksum_file)), _LK_NBLCK, LONG_MAX)==-1)
1786 #endif
1787                 {
1788                     if (errno == ENOSYS)
1789                     {
1790                         if (!bForceAppend)
1791                         {
1792                             gmx_fatal(FARGS,"File locking is not supported on this system. Use -noappend or specify -append explicitly to append anyhow.");
1793                         }
1794                         else
1795                         {
1796                             fprintf(stderr,"\nNOTE: File locking is not supported on this system, will not lock %s\n\n",outputfiles[i].filename);
1797                             if (fplog)
1798                             {
1799                                 fprintf(fplog,"\nNOTE: File locking not supported on this system, will not lock %s\n\n",outputfiles[i].filename);
1800                             }
1801                         }
1802                     }
1803                     else if (errno == EACCES || errno == EAGAIN)
1804                     {
1805                         gmx_fatal(FARGS,"Failed to lock: %s. Already running "
1806                                   "simulation?", outputfiles[i].filename);
1807                     }
1808                     else
1809                     {
1810                         gmx_fatal(FARGS,"Failed to lock: %s. %s.",
1811                                   outputfiles[i].filename, strerror(errno));
1812                     }
1813                 }
1814             }
1815             
1816             /* compute md5 chksum */ 
1817             if (outputfiles[i].chksum_size != -1)
1818             {
1819                 if (gmx_fio_get_file_md5(chksum_file,outputfiles[i].offset,
1820                                      digest) != outputfiles[i].chksum_size)  /*at the end of the call the file position is at the end of the file*/
1821                 {
1822                     gmx_fatal(FARGS,"Can't read %d bytes of '%s' to compute checksum. The file has been replaced or its contents have been modified. Cannot do appending because of this condition.",
1823                               outputfiles[i].chksum_size, 
1824                               outputfiles[i].filename);
1825                 }
1826             } 
1827             if (i==0)  /*log file needs to be seeked in case we need to truncate (other files are truncated below)*/
1828             {
1829                 if (gmx_fio_seek(chksum_file,outputfiles[i].offset))
1830                 {
1831                         gmx_fatal(FARGS,"Seek error! Failed to truncate log-file: %s.", strerror(errno));
1832                 }
1833             }
1834 #endif
1835             
1836             if (i==0) /*open log file here - so that lock is never lifted 
1837                         after chksum is calculated */
1838             {
1839                 *pfplog = gmx_fio_getfp(chksum_file);
1840             }
1841             else
1842             {
1843                 gmx_fio_close(chksum_file);
1844             }
1845 #ifndef GMX_FAHCORE            
1846             /* compare md5 chksum */
1847             if (outputfiles[i].chksum_size != -1 &&
1848                 memcmp(digest,outputfiles[i].chksum,16)!=0) 
1849             {
1850                 if (debug)
1851                 {
1852                     fprintf(debug,"chksum for %s: ",outputfiles[i].filename);
1853                     for (j=0; j<16; j++)
1854                     {
1855                         fprintf(debug,"%02x",digest[j]);
1856                     }
1857                     fprintf(debug,"\n");
1858                 }
1859                 gmx_fatal(FARGS,"Checksum wrong for '%s'. The file has been replaced or its contents have been modified. Cannot do appending because of this condition.",
1860                           outputfiles[i].filename);
1861             }
1862 #endif        
1863
1864               
1865             if (i!=0) /*log file is already seeked to correct position */
1866             {
1867 #ifdef GMX_NATIVE_WINDOWS
1868                 rc = gmx_wintruncate(outputfiles[i].filename,outputfiles[i].offset);
1869 #else            
1870                 rc = truncate(outputfiles[i].filename,outputfiles[i].offset);
1871 #endif
1872                 if(rc!=0)
1873                 {
1874                     gmx_fatal(FARGS,"Truncation of file %s failed. Cannot do appending because of this failure.",outputfiles[i].filename);
1875                 }
1876             }
1877         }
1878     }
1879
1880     sfree(outputfiles);
1881 }
1882
1883
1884 void load_checkpoint(const char *fn,FILE **fplog,
1885                      t_commrec *cr,gmx_bool bPartDecomp,ivec dd_nc,
1886                      t_inputrec *ir,t_state *state,
1887                      gmx_bool *bReadRNG,gmx_bool *bReadEkin,
1888                      gmx_bool bAppend,gmx_bool bForceAppend)
1889 {
1890     gmx_large_int_t step;
1891     double t;
1892
1893     if (SIMMASTER(cr)) {
1894       /* Read the state from the checkpoint file */
1895       read_checkpoint(fn,fplog,
1896                       cr,bPartDecomp,dd_nc,
1897                       ir->eI,&step,&t,state,bReadRNG,bReadEkin,
1898                       &ir->simulation_part,bAppend,bForceAppend);
1899     }
1900     if (PAR(cr)) {
1901       gmx_bcast(sizeof(cr->npmenodes),&cr->npmenodes,cr);
1902       gmx_bcast(DIM*sizeof(dd_nc[0]),dd_nc,cr);
1903       gmx_bcast(sizeof(step),&step,cr);
1904       gmx_bcast(sizeof(*bReadRNG),bReadRNG,cr);
1905       gmx_bcast(sizeof(*bReadEkin),bReadEkin,cr);
1906     }
1907     ir->bContinuation    = TRUE;
1908     if (ir->nsteps >= 0)
1909     {
1910         ir->nsteps          += ir->init_step - step;
1911     }
1912     ir->init_step        = step;
1913         ir->simulation_part += 1;
1914 }
1915
1916 static void read_checkpoint_data(t_fileio *fp,int *simulation_part,
1917                                  gmx_large_int_t *step,double *t,t_state *state,
1918                                  gmx_bool bReadRNG,
1919                                  int *nfiles,gmx_file_position_t **outputfiles)
1920 {
1921     int  file_version;
1922     char *version,*btime,*buser,*bmach,*fprog,*ftime;
1923     int  double_prec;
1924     int  eIntegrator;
1925     int  nppnodes,npme;
1926     ivec dd_nc;
1927     int  flags_eks,flags_enh;
1928     int  nfiles_loc;
1929     gmx_file_position_t *files_loc=NULL;
1930     int  ret;
1931         
1932     do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
1933                   &version,&btime,&buser,&bmach,&double_prec,&fprog,&ftime,
1934                   &eIntegrator,simulation_part,step,t,&nppnodes,dd_nc,&npme,
1935                   &state->natoms,&state->ngtc,&state->nnhpres,&state->nhchainlength,
1936                   &state->flags,&flags_eks,&flags_enh,NULL);
1937     ret =
1938         do_cpt_state(gmx_fio_getxdr(fp),TRUE,state->flags,state,bReadRNG,NULL);
1939     if (ret)
1940     {
1941         cp_error();
1942     }
1943     ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
1944                            flags_eks,&state->ekinstate,NULL);
1945     if (ret)
1946     {
1947         cp_error();
1948     }
1949     ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
1950                           flags_enh,&state->enerhist,NULL);
1951     if (ret)
1952     {
1953         cp_error();
1954     }
1955
1956     ret = do_cpt_files(gmx_fio_getxdr(fp),TRUE,
1957                        outputfiles != NULL ? outputfiles : &files_loc,
1958                        outputfiles != NULL ? nfiles : &nfiles_loc,
1959                        NULL,file_version);
1960     if (files_loc != NULL)
1961     {
1962         sfree(files_loc);
1963     }
1964         
1965     if (ret)
1966     {
1967         cp_error();
1968     }
1969         
1970     ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
1971     if (ret)
1972     {
1973         cp_error();
1974     }
1975
1976     sfree(fprog);
1977     sfree(ftime);
1978     sfree(btime);
1979     sfree(buser);
1980     sfree(bmach);
1981 }
1982
1983 void 
1984 read_checkpoint_state(const char *fn,int *simulation_part,
1985                       gmx_large_int_t *step,double *t,t_state *state)
1986 {
1987     t_fileio *fp;
1988     
1989     fp = gmx_fio_open(fn,"r");
1990     read_checkpoint_data(fp,simulation_part,step,t,state,FALSE,NULL,NULL);
1991     if( gmx_fio_close(fp) != 0)
1992         {
1993         gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
1994         }
1995 }
1996
1997 void read_checkpoint_trxframe(t_fileio *fp,t_trxframe *fr)
1998 {
1999     t_state state;
2000     int simulation_part;
2001     gmx_large_int_t step;
2002     double t;
2003     
2004     init_state(&state,0,0,0,0);
2005     
2006     read_checkpoint_data(fp,&simulation_part,&step,&t,&state,FALSE,NULL,NULL);
2007     
2008     fr->natoms  = state.natoms;
2009     fr->bTitle  = FALSE;
2010     fr->bStep   = TRUE;
2011     fr->step    = gmx_large_int_to_int(step,
2012                                     "conversion of checkpoint to trajectory");
2013     fr->bTime   = TRUE;
2014     fr->time    = t;
2015     fr->bLambda = TRUE;
2016     fr->lambda  = state.lambda;
2017     fr->bAtoms  = FALSE;
2018     fr->bX      = (state.flags & (1<<estX));
2019     if (fr->bX)
2020     {
2021         fr->x     = state.x;
2022         state.x   = NULL;
2023     }
2024     fr->bV      = (state.flags & (1<<estV));
2025     if (fr->bV)
2026     {
2027         fr->v     = state.v;
2028         state.v   = NULL;
2029     }
2030     fr->bF      = FALSE;
2031     fr->bBox    = (state.flags & (1<<estBOX));
2032     if (fr->bBox)
2033     {
2034         copy_mat(state.box,fr->box);
2035     }
2036     done_state(&state);
2037 }
2038
2039 void list_checkpoint(const char *fn,FILE *out)
2040 {
2041     t_fileio *fp;
2042     int  file_version;
2043     char *version,*btime,*buser,*bmach,*fprog,*ftime;
2044     int  double_prec;
2045     int  eIntegrator,simulation_part,nppnodes,npme;
2046     gmx_large_int_t step;
2047     double t;
2048     ivec dd_nc;
2049     t_state state;
2050     int  flags_eks,flags_enh;
2051     int  indent;
2052     int  i,j;
2053     int  ret;
2054     gmx_file_position_t *outputfiles;
2055         int  nfiles;
2056         
2057     init_state(&state,-1,-1,-1,-1);
2058
2059     fp = gmx_fio_open(fn,"r");
2060     do_cpt_header(gmx_fio_getxdr(fp),TRUE,&file_version,
2061                   &version,&btime,&buser,&bmach,&double_prec,&fprog,&ftime,
2062                   &eIntegrator,&simulation_part,&step,&t,&nppnodes,dd_nc,&npme,
2063                   &state.natoms,&state.ngtc,&state.nnhpres,&state.nhchainlength,
2064                   &state.flags,&flags_eks,&flags_enh,out);
2065     ret = do_cpt_state(gmx_fio_getxdr(fp),TRUE,state.flags,&state,TRUE,out);
2066     if (ret)
2067     {
2068         cp_error();
2069     }
2070     ret = do_cpt_ekinstate(gmx_fio_getxdr(fp),TRUE,
2071                            flags_eks,&state.ekinstate,out);
2072     if (ret)
2073     {
2074         cp_error();
2075     }
2076     ret = do_cpt_enerhist(gmx_fio_getxdr(fp),TRUE,
2077                           flags_enh,&state.enerhist,out);
2078
2079     if (ret == 0)
2080     {
2081                 do_cpt_files(gmx_fio_getxdr(fp),TRUE,&outputfiles,&nfiles,out,file_version);
2082         }
2083         
2084     if (ret == 0)
2085     {
2086         ret = do_cpt_footer(gmx_fio_getxdr(fp),TRUE,file_version);
2087     }
2088         
2089     if (ret)
2090     {
2091         cp_warning(out);
2092     }
2093     if( gmx_fio_close(fp) != 0)
2094         {
2095         gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2096         }
2097     
2098     done_state(&state);
2099 }
2100
2101
2102 static gmx_bool exist_output_file(const char *fnm_cp,int nfile,const t_filenm fnm[])
2103 {
2104     int i;
2105
2106     /* Check if the output file name stored in the checkpoint file
2107      * is one of the output file names of mdrun.
2108      */
2109     i = 0;
2110     while (i < nfile &&
2111            !(is_output(&fnm[i]) && strcmp(fnm_cp,fnm[i].fns[0]) == 0))
2112     {
2113         i++;
2114     }
2115     
2116     return (i < nfile && gmx_fexist(fnm_cp));
2117 }
2118
2119 /* This routine cannot print tons of data, since it is called before the log file is opened. */
2120 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
2121                                      gmx_large_int_t *cpt_step,t_commrec *cr,
2122                                      gmx_bool bAppendReq,
2123                                      int nfile,const t_filenm fnm[],
2124                                      const char *part_suffix,gmx_bool *bAddPart)
2125 {
2126     t_fileio *fp;
2127     gmx_large_int_t step=0;
2128         double t;
2129     t_state state;
2130     int  nfiles;
2131     gmx_file_position_t *outputfiles;
2132     int  nexist,f;
2133     gmx_bool bAppend;
2134     char *fn,suf_up[STRLEN];
2135
2136     bAppend = FALSE;
2137
2138     if (SIMMASTER(cr)) {
2139         if(!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename,"r")) ))
2140         {
2141             *simulation_part = 0;
2142         }
2143         else 
2144         {
2145             init_state(&state,0,0,0,0);
2146
2147             read_checkpoint_data(fp,simulation_part,&step,&t,&state,FALSE,
2148                                  &nfiles,&outputfiles);
2149             if( gmx_fio_close(fp) != 0)
2150             {
2151                 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
2152             }
2153             done_state(&state);
2154
2155             if (bAppendReq)
2156             {
2157                 nexist = 0;
2158                 for(f=0; f<nfiles; f++)
2159                 {
2160                     if (exist_output_file(outputfiles[f].filename,nfile,fnm))
2161                     {
2162                         nexist++;
2163                     }
2164                 }
2165                 if (nexist == nfiles)
2166                 {
2167                     bAppend = bAppendReq;
2168                 }
2169                 else if (nexist > 0)
2170                 {
2171                     fprintf(stderr,
2172                             "Output file appending has been requested,\n"
2173                             "but some output files listed in the checkpoint file %s\n"
2174                             "are not present or are named differently by the current program:\n",
2175                             filename);
2176                     fprintf(stderr,"output files present:");
2177                     for(f=0; f<nfiles; f++)
2178                     {
2179                         if (exist_output_file(outputfiles[f].filename,
2180                                               nfile,fnm))
2181                         {
2182                             fprintf(stderr," %s",outputfiles[f].filename);
2183                         }
2184                     }
2185                     fprintf(stderr,"\n");
2186                     fprintf(stderr,"output files not present or named differently:");
2187                     for(f=0; f<nfiles; f++)
2188                     {
2189                         if (!exist_output_file(outputfiles[f].filename,
2190                                                nfile,fnm))
2191                         {
2192                             fprintf(stderr," %s",outputfiles[f].filename);
2193                         }
2194                     }
2195                     fprintf(stderr,"\n");
2196                     
2197                     gmx_fatal(FARGS,"File appending requested, but only %d of the %d output files are present",nexist,nfiles);
2198                 }
2199             }
2200             
2201             if (bAppend)
2202             {
2203                 if (nfiles == 0)
2204                 {
2205                     gmx_fatal(FARGS,"File appending requested, but no output file information is stored in the checkpoint file");
2206                 }
2207                 fn = outputfiles[0].filename;
2208                 if (strlen(fn) < 4 ||
2209                     gmx_strcasecmp(fn+strlen(fn)-4,ftp2ext(efLOG)) == 0)
2210                 {
2211                     gmx_fatal(FARGS,"File appending requested, but the log file is not the first file listed in the checkpoint file");
2212                 }
2213                 /* Set bAddPart to whether the suffix string '.part' is present
2214                  * in the log file name.
2215                  */
2216                 strcpy(suf_up,part_suffix);
2217                 upstring(suf_up);
2218                 *bAddPart = (strstr(fn,part_suffix) != NULL ||
2219                              strstr(fn,suf_up) != NULL);
2220             }
2221
2222             sfree(outputfiles);
2223         }
2224     }
2225     if (PAR(cr))
2226     {
2227         gmx_bcast(sizeof(*simulation_part),simulation_part,cr);
2228
2229         if (*simulation_part > 0 && bAppendReq)
2230         {
2231             gmx_bcast(sizeof(bAppend),&bAppend,cr);
2232             gmx_bcast(sizeof(*bAddPart),bAddPart,cr);
2233         }
2234     }
2235     if (NULL != cpt_step)
2236     {
2237         *cpt_step = step;
2238     }
2239
2240     return bAppend;
2241 }