1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
42 #include "gmx_fatal.h"
48 /* The source code in this file should be thread-safe.
49 Please keep it that way. */
51 /* This number should be increased whenever the file format changes! */
52 static const int enx_version
= 3;
55 /* Stuff for reading pre 4.1 energy files */
57 bool bOldFileOpen
; /* Is this an open old file? */
58 bool bReadFirstStep
; /* Did we read the first step? */
59 int first_step
; /* First step in the energy file */
60 int step_prev
; /* Previous step */
61 int nsum_prev
; /* Previous step sum length */
62 t_energy
*ener_prev
; /* Previous energy sums */
74 void free_enxframe(t_enxframe
*fr
)
81 sfree(fr
->disre_rm3tav
);
84 for(b
=0; b
<fr
->nblock
; b
++)
91 static void edr_strings(XDR
*xdr
,bool bRead
,int file_version
,
92 int n
,gmx_enxnm_t
**nms
)
117 if(!xdr_string(xdr
,&(nm
->name
),STRLEN
))
119 gmx_file("Cannot write energy names to file; maybe you are out of quota?");
121 if (file_version
>= 2)
123 if(!xdr_string(xdr
,&(nm
->unit
),STRLEN
))
125 gmx_file("Cannot write energy names to file; maybe you are out of quota?");
130 nm
->unit
= strdup("kJ/mol");
135 static void gen_units(int n
,char ***units
)
142 (*units
)[i
] = strdup("kJ/mol");
146 void do_enxnms(ener_file_t ef
,int *nre
,gmx_enxnm_t
**nms
)
150 bool bRead
= gmx_fio_getread(ef
->fp
);
154 gmx_fio_select(ef
->fp
);
156 xdr
= gmx_fio_getxdr(ef
->fp
);
158 if (!xdr_int(xdr
,&magic
))
162 gmx_file("Cannot write energy names to file; maybe you are out of quota?");
169 /* Assume this is an old edr format */
172 ef
->eo
.bOldFileOpen
= TRUE
;
173 ef
->eo
.bReadFirstStep
= FALSE
;
174 srenew(ef
->eo
.ener_prev
,*nre
);
178 ef
->eo
.bOldFileOpen
=FALSE
;
182 gmx_fatal(FARGS
,"Energy names magic number mismatch, this is not a GROMACS edr file");
184 file_version
= enx_version
;
185 xdr_int(xdr
,&file_version
);
186 if (file_version
> enx_version
)
188 gmx_fatal(FARGS
,"reading tpx file (%s) version %d with version %d program",gmx_fio_getname(ef
->fp
),file_version
,enx_version
);
192 if (file_version
!= enx_version
)
194 fprintf(stderr
,"Note: enx file_version %d, software version %d\n",
195 file_version
,enx_version
);
198 edr_strings(xdr
,bRead
,file_version
,*nre
,nms
);
201 static bool do_eheader(ener_file_t ef
,int *file_version
,t_enxframe
*fr
,
202 int nre_test
,bool *bWrongPrecision
,bool *bOK
)
206 int block
,i
,zero
=0,dum
=0;
207 bool bRead
= gmx_fio_getread(ef
->fp
);
212 *bWrongPrecision
= FALSE
;
216 /* The original energy frame started with a real,
217 * so we have to use a real for compatibility.
218 * This is VERY DIRTY code, since do_eheader can be called
219 * with the wrong precision set and then we could read r > -1e10,
220 * while actually the intention was r < -1e10.
221 * When nre_test >= 0, do_eheader should therefore terminate
222 * before the number of i/o calls starts depending on what has been read
223 * (which is the case for for instance the block sizes for variable
224 * number of blocks, where this number is read before).
233 /* Assume we are reading an old format */
236 if (!do_int(dum
)) *bOK
= FALSE
;
241 if (!do_int(magic
)) *bOK
= FALSE
;
242 if (magic
!= -7777777)
244 gmx_fatal(FARGS
,"Energy header magic number mismatch, this is not a GROMACS edr file");
246 *file_version
= enx_version
;
247 if (!do_int (*file_version
)) *bOK
= FALSE
;
248 if (*bOK
&& *file_version
> enx_version
)
250 gmx_fatal(FARGS
,"reading tpx file (%s) version %d with version %d program",gmx_fio_getname(ef
->fp
),file_version
,enx_version
);
252 if (!do_double(fr
->t
)) *bOK
= FALSE
;
253 if (!do_gmx_large_int(fr
->step
)) *bOK
= FALSE
;
254 if (!bRead
&& fr
->nsum
== 1) {
255 /* Do not store sums of length 1,
256 * since this does not add information.
258 if (!do_int (zero
)) *bOK
= FALSE
;
260 if (!do_int (fr
->nsum
)) *bOK
= FALSE
;
262 if (*file_version
>= 3)
264 do_gmx_large_int(fr
->nsteps
);
268 fr
->nsteps
= max(1,fr
->nsum
);
271 if (!do_int (fr
->nre
)) *bOK
= FALSE
;
272 if (!do_int (fr
->ndisre
)) *bOK
= FALSE
;
273 if (!do_int (fr
->nblock
)) *bOK
= FALSE
;
275 /* Frames could have nre=0, so we can not rely only on the fr->nre check */
276 if (bRead
&& nre_test
>= 0 &&
277 ((fr
->nre
> 0 && fr
->nre
!= nre_test
) ||
278 fr
->nre
< 0 || fr
->ndisre
< 0 || fr
->nblock
< 0))
280 *bWrongPrecision
= TRUE
;
284 if (*bOK
&& bRead
&& fr
->nblock
>fr
->nr_alloc
)
286 srenew(fr
->nr
,fr
->nblock
);
287 srenew(fr
->b_alloc
,fr
->nblock
);
288 srenew(fr
->block
,fr
->nblock
);
289 for(i
=fr
->nr_alloc
; i
<fr
->nblock
; i
++) {
293 fr
->nr_alloc
= fr
->nblock
;
295 for(block
=0; block
<fr
->nblock
; block
++)
297 if (!do_int (fr
->nr
[block
]))
302 if (!do_int (fr
->e_size
)) *bOK
= FALSE
;
303 if (!do_int (fr
->d_size
)) *bOK
= FALSE
;
304 /* Do a dummy int to keep the format compatible with the old code */
305 if (!do_int (dum
)) *bOK
= FALSE
;
307 if (*bOK
&& *file_version
== 1 && nre_test
< 0)
310 if (fp
>= ener_old_nalloc
)
312 gmx_incons("Problem with reading old format energy files");
316 if (!ef
->eo
.bReadFirstStep
)
318 ef
->eo
.bReadFirstStep
= TRUE
;
319 ef
->eo
.first_step
= fr
->step
;
320 ef
->eo
.step_prev
= fr
->step
;
321 ef
->eo
.nsum_prev
= 0;
324 fr
->nsum
= fr
->step
- ef
->eo
.first_step
+ 1;
325 fr
->nsteps
= fr
->step
- ef
->eo
.step_prev
;
331 void free_enxnms(int n
,gmx_enxnm_t
*nms
)
344 void close_enx(ener_file_t ef
)
346 if(gmx_fio_close(ef
->fp
) != 0)
348 gmx_file("Cannot close energy file; it might be corrupt, or maybe you are out of quota?");
352 static bool empty_file(const char *fn
)
359 fp
= gmx_fio_fopen(fn
,"r");
360 ret
= fread(&dum
,sizeof(dum
),1,fp
);
368 ener_file_t
open_enx(const char *fn
,const char *mode
)
371 gmx_enxnm_t
*nms
=NULL
;
374 bool bWrongPrecision
,bDum
=TRUE
;
375 struct ener_file
*ef
;
380 ef
->fp
=gmx_fio_open(fn
,mode
);
381 gmx_fio_select(ef
->fp
);
382 gmx_fio_setprecision(ef
->fp
,FALSE
);
383 do_enxnms(ef
,&nre
,&nms
);
385 do_eheader(ef
,&file_version
,fr
,nre
,&bWrongPrecision
,&bDum
);
388 gmx_file("Cannot read energy file header. Corrupt file?");
391 /* Now check whether this file is in single precision */
392 if (!bWrongPrecision
&&
393 ((fr
->e_size
&& (fr
->nre
== nre
) &&
394 (nre
*4*(long int)sizeof(float) == fr
->e_size
)) ||
396 (fr
->ndisre
*(long int)sizeof(float)*2+(long int)sizeof(int)
399 fprintf(stderr
,"Opened %s as single precision energy file\n",fn
);
400 free_enxnms(nre
,nms
);
404 gmx_fio_rewind(ef
->fp
);
405 gmx_fio_select(ef
->fp
);
406 gmx_fio_setprecision(ef
->fp
,TRUE
);
407 do_enxnms(ef
,&nre
,&nms
);
408 do_eheader(ef
,&file_version
,fr
,nre
,&bWrongPrecision
,&bDum
);
411 gmx_file("Cannot write energy file header; maybe you are out of quota?");
414 if (((fr
->e_size
&& (fr
->nre
== nre
) &&
415 (nre
*4*(long int)sizeof(double) == fr
->e_size
)) ||
417 (fr
->ndisre
*(long int)sizeof(double)*2+
418 (long int)sizeof(int) ==
420 fprintf(stderr
,"Opened %s as double precision energy file\n",
424 gmx_fatal(FARGS
,"File %s is empty",fn
);
426 gmx_fatal(FARGS
,"Energy file %s not recognized, maybe different CPU?",
429 free_enxnms(nre
,nms
);
433 gmx_fio_rewind(ef
->fp
);
436 ef
->fp
= gmx_fio_open(fn
,mode
);
444 int enx_file_pointer(const ener_file_t ef
)
449 static void convert_full_sums(ener_old_t
*ener_old
,t_enxframe
*fr
)
453 double esum_all
,eav_all
;
459 for(i
=0; i
<fr
->nre
; i
++)
461 if (fr
->ener
[i
].e
!= 0) ne
++;
462 if (fr
->ener
[i
].esum
!= 0) ns
++;
464 if (ne
> 0 && ns
== 0)
466 /* We do not have all energy sums */
471 /* Convert old full simulation sums to sums between energy frames */
472 nstep_all
= fr
->step
- ener_old
->first_step
+ 1;
473 if (fr
->nsum
> 1 && fr
->nsum
== nstep_all
&& ener_old
->nsum_prev
> 0)
475 /* Set the new sum length: the frame step difference */
476 fr
->nsum
= fr
->step
- ener_old
->step_prev
;
477 for(i
=0; i
<fr
->nre
; i
++)
479 esum_all
= fr
->ener
[i
].esum
;
480 eav_all
= fr
->ener
[i
].eav
;
481 fr
->ener
[i
].esum
= esum_all
- ener_old
->ener_prev
[i
].esum
;
482 fr
->ener
[i
].eav
= eav_all
- ener_old
->ener_prev
[i
].eav
483 - dsqr(ener_old
->ener_prev
[i
].esum
/(nstep_all
- fr
->nsum
)
484 - esum_all
/nstep_all
)*
485 (nstep_all
- fr
->nsum
)*nstep_all
/(double)fr
->nsum
;
486 ener_old
->ener_prev
[i
].esum
= esum_all
;
487 ener_old
->ener_prev
[i
].eav
= eav_all
;
489 ener_old
->nsum_prev
= nstep_all
;
491 else if (fr
->nsum
> 0)
493 if (fr
->nsum
!= nstep_all
)
495 fprintf(stderr
,"\nWARNING: something is wrong with the energy sums, will not use exact averages\n");
496 ener_old
->nsum_prev
= 0;
500 ener_old
->nsum_prev
= nstep_all
;
502 /* Copy all sums to ener_prev */
503 for(i
=0; i
<fr
->nre
; i
++)
505 ener_old
->ener_prev
[i
].esum
= fr
->ener
[i
].esum
;
506 ener_old
->ener_prev
[i
].eav
= fr
->ener
[i
].eav
;
510 ener_old
->step_prev
= fr
->step
;
513 bool do_enx(ener_file_t ef
,t_enxframe
*fr
)
517 bool bRead
,bOK
,bOK1
,bSane
;
522 bRead
= gmx_fio_getread(ef
->fp
);
525 fr
->e_size
= fr
->nre
*sizeof(fr
->ener
[0].e
)*4;
526 fr
->d_size
= fr
->ndisre
*(sizeof(fr
->disre_rm3tav
[0]) +
527 sizeof(fr
->disre_rt
[0]));
529 gmx_fio_select(ef
->fp
);
531 if (!do_eheader(ef
,&file_version
,fr
,-1,NULL
,&bOK
))
535 fprintf(stderr
,"\rLast energy frame read %d time %8.3f ",
536 ef
->framenr
-1,ef
->frametime
);
540 "\nWARNING: Incomplete energy frame: nr %d time %8.3f\n",
546 gmx_file("Cannot write energy file header; maybe you are out of quota?");
552 if ((ef
->framenr
< 20 || ef
->framenr
% 10 == 0) &&
553 (ef
->framenr
< 200 || ef
->framenr
% 100 == 0) &&
554 (ef
->framenr
< 2000 || ef
->framenr
% 1000 == 0))
556 fprintf(stderr
,"\rReading energy frame %6d time %8.3f ",
560 ef
->frametime
= fr
->t
;
562 /* Check sanity of this header */
563 bSane
= (fr
->nre
> 0 || fr
->ndisre
> 0);
564 for(block
=0; block
<fr
->nblock
; block
++)
566 bSane
= bSane
|| (fr
->nr
[block
] > 0);
568 if (!((fr
->step
>= 0) && bSane
))
570 fprintf(stderr
,"\nWARNING: there may be something wrong with energy file %s\n",
571 gmx_fio_getname(ef
->fp
));
572 fprintf(stderr
,"Found: step=%s, nre=%d, ndisre=%d, nblock=%d, time=%g.\n"
573 "Trying to skip frame expect a crash though\n",
574 gmx_step_str(fr
->step
,buf
),fr
->nre
,fr
->ndisre
,fr
->nblock
,fr
->t
);
576 if (bRead
&& fr
->nre
> fr
->e_alloc
)
578 srenew(fr
->ener
,fr
->nre
);
579 for(i
=fr
->e_alloc
; (i
<fr
->nre
); i
++)
583 fr
->ener
[i
].esum
= 0;
585 fr
->e_alloc
= fr
->nre
;
588 for(i
=0; i
<fr
->nre
; i
++)
590 bOK
= bOK
&& do_real(fr
->ener
[i
].e
);
592 /* Do not store sums of length 1,
593 * since this does not add information.
595 if (file_version
== 1 ||
596 (bRead
&& fr
->nsum
> 0) || fr
->nsum
> 1)
598 tmp1
= fr
->ener
[i
].eav
;
599 bOK
= bOK
&& do_real(tmp1
);
601 fr
->ener
[i
].eav
= tmp1
;
603 /* This is to save only in single precision (unless compiled in DP) */
604 tmp2
= fr
->ener
[i
].esum
;
605 bOK
= bOK
&& do_real(tmp2
);
607 fr
->ener
[i
].esum
= tmp2
;
609 if (file_version
== 1)
611 /* Old, unused real */
613 bOK
= bOK
&& do_real(rdum
);
618 /* Here we can not check for file_version==1, since one could have
619 * continued an old format simulation with a new one with mdrun -append.
621 if (bRead
&& ef
->eo
.bOldFileOpen
)
623 /* Convert old full simulation sums to sums between energy frames */
624 convert_full_sums(&(ef
->eo
),fr
);
628 if (bRead
&& fr
->ndisre
>fr
->d_alloc
)
630 srenew(fr
->disre_rm3tav
,fr
->ndisre
);
631 srenew(fr
->disre_rt
,fr
->ndisre
);
632 fr
->d_alloc
= fr
->ndisre
;
634 ndo_real(fr
->disre_rm3tav
,fr
->ndisre
,bOK1
);
636 ndo_real(fr
->disre_rt
,fr
->ndisre
,bOK1
);
639 for(block
=0; block
<fr
->nblock
; block
++)
641 if (bRead
&& fr
->nr
[block
]>fr
->b_alloc
[block
])
643 srenew(fr
->block
[block
],fr
->nr
[block
]);
644 fr
->b_alloc
[block
] = fr
->nr
[block
];
646 ndo_real(fr
->block
[block
],fr
->nr
[block
],bOK1
);
652 if( gmx_fio_flush(ef
->fp
) != 0)
654 gmx_file("Cannot write energy file; maybe you are out of quota?");
662 fprintf(stderr
,"\nLast energy frame read %d",
664 fprintf(stderr
,"\nWARNING: Incomplete energy frame: nr %d time %8.3f\n",
669 gmx_fatal(FARGS
,"could not write energies");
677 static real
find_energy(const char *name
, int nre
, gmx_enxnm_t
*enm
,
684 if (strcmp(enm
[i
].name
,name
) == 0)
686 return fr
->ener
[i
].e
;
690 gmx_fatal(FARGS
,"Could not find energy term named '%s'",name
);
696 void get_enx_state(const char *fn
, real t
, gmx_groups_t
*groups
, t_inputrec
*ir
,
699 /* Should match the names in mdebin.c */
700 static const char *boxvel_nm
[] = {
701 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
702 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
705 static const char *pcouplmu_nm
[] = {
706 "Pcoupl-Mu-XX", "Pcoupl-Mu-YY", "Pcoupl-Mu-ZZ",
707 "Pcoupl-Mu-YX", "Pcoupl-Mu-ZX", "Pcoupl-Mu-ZY"
709 static const char *baro_nm
[] = {
714 int ind0
[] = { XX
,YY
,ZZ
,YY
,ZZ
,ZZ
};
715 int ind1
[] = { XX
,YY
,ZZ
,XX
,XX
,YY
};
716 int nre
,nfr
,i
,j
,ni
,npcoupl
;
719 gmx_enxnm_t
*enm
=NULL
;
723 in
= open_enx(fn
,"r");
724 do_enxnms(in
,&nre
,&enm
);
727 while ((nfr
==0 || fr
->t
!= t
) && do_enx(in
,fr
)) {
731 fprintf(stderr
,"\n");
733 if (nfr
== 0 || fr
->t
!= t
)
734 gmx_fatal(FARGS
,"Could not find frame with time %f in '%s'",t
,fn
);
736 npcoupl
= TRICLINIC(ir
->compress
) ? 6 : 3;
737 if (ir
->epc
== epcPARRINELLORAHMAN
) {
738 clear_mat(state
->boxv
);
739 for(i
=0; i
<npcoupl
; i
++) {
740 state
->boxv
[ind0
[i
]][ind1
[i
]] =
741 find_energy(boxvel_nm
[i
],nre
,enm
,fr
);
743 fprintf(stderr
,"\nREAD %d BOX VELOCITIES FROM %s\n\n",npcoupl
,fn
);
746 if (ir
->etc
== etcNOSEHOOVER
)
748 for(i
=0; i
<state
->ngtc
; i
++) {
749 ni
= groups
->grps
[egcTC
].nm_ind
[i
];
750 bufi
= *(groups
->grpname
[ni
]);
751 for(j
=0; (j
<state
->nhchainlength
); j
++)
753 sprintf(buf
,"Xi-%d-%s",j
,bufi
);
754 state
->nosehoover_xi
[i
] = find_energy(buf
,nre
,enm
,fr
);
755 sprintf(buf
,"vXi-%d-%s",j
,bufi
);
756 state
->nosehoover_vxi
[i
] = find_energy(buf
,nre
,enm
,fr
);
760 fprintf(stderr
,"\nREAD %d NOSE-HOOVER Xi chains FROM %s\n\n",state
->ngtc
,fn
);
762 if (IR_NPT_TROTTER(ir
))
764 for(i
=0; i
<state
->nnhpres
; i
++) {
765 bufi
= baro_nm
[0]; /* All barostat DOF's together for now */
766 for(j
=0; (j
<state
->nhchainlength
); j
++)
768 sprintf(buf
,"Xi-%d-%s",j
,bufi
);
769 state
->nhpres_xi
[i
] = find_energy(buf
,nre
,enm
,fr
);
770 sprintf(buf
,"vXi-%d-%s",j
,bufi
);
771 state
->nhpres_vxi
[i
] = find_energy(buf
,nre
,enm
,fr
);
774 fprintf(stderr
,"\nREAD %d NOSE-HOOVER BAROSTAT Xi chains FROM %s\n\n",state
->nnhpres
,fn
);
778 free_enxnms(nre
,enm
);