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 * GROwing Monsters And Cloning Shrimps
47 #include "gmx_fatal.h"
55 #include "gmx_random.h"
61 #define NOT_FINISHED(l1,l2) \
62 printf("not finished yet: lines %d .. %d in %s\n",l1,l2,__FILE__)
64 static char *int_title(const char *title
,int nodeid
,char buf
[], int size
)
66 sprintf(buf
,"%s (%d)",title
,nodeid
);
71 void set_state_entries(t_state
*state
,const t_inputrec
*ir
,int nnodes
)
75 /* The entries in the state in the tpx file might not correspond
76 * with what is needed, so we correct this here.
79 if (ir
->efep
!= efepNO
)
80 state
->flags
|= (1<<estLAMBDA
);
81 state
->flags
|= (1<<estX
);
83 snew(state
->x
,state
->nalloc
);
84 if (EI_DYNAMICS(ir
->eI
)) {
85 state
->flags
|= (1<<estV
);
87 snew(state
->v
,state
->nalloc
);
89 if (ir
->eI
== eiSD2
) {
90 state
->flags
|= (1<<estSDX
);
91 if (state
->sd_X
== NULL
) {
92 /* sd_X is not stored in the tpx file, so we need to allocate it */
93 snew(state
->sd_X
,state
->nalloc
);
98 state
->flags
|= (1<<estCGP
);
99 if (state
->cg_p
== NULL
)
101 /* cg_p is not stored in the tpx file, so we need to allocate it */
102 snew(state
->cg_p
,state
->nalloc
);
105 if (EI_SD(ir
->eI
) || ir
->eI
== eiBD
|| ir
->etc
== etcVRESCALE
) {
106 state
->nrng
= gmx_rng_n();
108 if (EI_SD(ir
->eI
) || ir
->eI
== eiBD
) {
109 /* This will be correct later with DD */
110 state
->nrng
*= nnodes
;
111 state
->nrngi
*= nnodes
;
113 state
->flags
|= ((1<<estLD_RNG
) | (1<<estLD_RNGI
));
114 snew(state
->ld_rng
, state
->nrng
);
115 snew(state
->ld_rngi
,state
->nrngi
);
120 if (ir
->ePBC
!= epbcNONE
) {
121 state
->flags
|= (1<<estBOX
);
122 if (PRESERVE_SHAPE(*ir
)) {
123 state
->flags
|= (1<<estBOX_REL
);
125 if ((ir
->epc
== epcPARRINELLORAHMAN
) || (ir
->epc
== epcMTTK
)) {
126 state
->flags
|= (1<<estBOXV
);
128 if (ir
->epc
!= epcNO
) {
129 if (IR_NPT_TROTTER(ir
)) {
131 state
->flags
|= (1<<estNHPRES_XI
);
132 state
->flags
|= (1<<estNHPRES_VXI
);
133 state
->flags
|= (1<<estSVIR_PREV
);
134 state
->flags
|= (1<<estFVIR_PREV
);
135 state
->flags
|= (1<<estVETA
);
136 state
->flags
|= (1<<estVOL0
);
138 state
->flags
|= (1<<estPRES_PREV
);
143 if (ir
->etc
== etcNOSEHOOVER
) {
144 state
->flags
|= (1<<estNH_XI
);
145 state
->flags
|= (1<<estNH_VXI
);
148 if (ir
->etc
== etcVRESCALE
) {
149 state
->flags
|= (1<<estTC_INT
);
152 init_gtc_state(state
,state
->ngtc
,state
->nnhpres
,ir
->opts
.nhchainlength
); /* allocate the space for nose-hoover chains */
153 init_ekinstate(&state
->ekinstate
,ir
);
155 init_energyhistory(&state
->enerhist
);
159 void init_parallel(FILE *log
, t_commrec
*cr
, t_inputrec
*inputrec
,
162 bcast_ir_mtop(cr
,inputrec
,mtop
);
165 if (inputrec
->eI
== eiBD
|| EI_SD(inputrec
->eI
)) {
166 /* Make sure the random seeds are different on each node */
167 inputrec
->ld_seed
+= cr
->nodeid
;