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 || ir->bExpanded)
81 state->flags |= (1<<estLAMBDA);
82 state->flags |= (1<<estFEPSTATE);
84 state->flags |= (1<<estX);
85 if (state->lambda == NULL)
87 snew(state->lambda, efptNR);
91 snew(state->x, state->nalloc);
93 if (EI_DYNAMICS(ir->eI))
95 state->flags |= (1<<estV);
98 snew(state->v, state->nalloc);
103 state->flags |= (1<<estSDX);
104 if (state->sd_X == NULL)
106 /* sd_X is not stored in the tpx file, so we need to allocate it */
107 snew(state->sd_X, state->nalloc);
112 state->flags |= (1<<estCGP);
113 if (state->cg_p == NULL)
115 /* cg_p is not stored in the tpx file, so we need to allocate it */
116 snew(state->cg_p, state->nalloc);
119 if (EI_SD(ir->eI) || ir->eI == eiBD || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
121 state->nrng = gmx_rng_n();
123 if (EI_SD(ir->eI) || ir->eI == eiBD || ETC_ANDERSEN(ir->etc))
125 /* This will be correct later with DD */
126 state->nrng *= nnodes;
127 state->nrngi *= nnodes;
129 state->flags |= ((1<<estLD_RNG) | (1<<estLD_RNGI));
130 snew(state->ld_rng, state->nrng);
131 snew(state->ld_rngi, state->nrngi);
140 state->nmcrng = gmx_rng_n();
141 snew(state->mc_rng, state->nmcrng);
142 snew(state->mc_rngi, 1);
146 if (ir->ePBC != epbcNONE)
148 state->flags |= (1<<estBOX);
149 if (PRESERVE_SHAPE(*ir))
151 state->flags |= (1<<estBOX_REL);
153 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
155 state->flags |= (1<<estBOXV);
157 if (ir->epc != epcNO)
159 if (IR_NPT_TROTTER(ir) || (IR_NPH_TROTTER(ir)))
162 state->flags |= (1<<estNHPRES_XI);
163 state->flags |= (1<<estNHPRES_VXI);
164 state->flags |= (1<<estSVIR_PREV);
165 state->flags |= (1<<estFVIR_PREV);
166 state->flags |= (1<<estVETA);
167 state->flags |= (1<<estVOL0);
171 state->flags |= (1<<estPRES_PREV);
176 if (ir->etc == etcNOSEHOOVER)
178 state->flags |= (1<<estNH_XI);
179 state->flags |= (1<<estNH_VXI);
182 if (ir->etc == etcVRESCALE)
184 state->flags |= (1<<estTC_INT);
187 init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
188 init_ekinstate(&state->ekinstate, ir);
190 init_energyhistory(&state->enerhist);
191 init_df_history(&state->dfhist, ir->fepvals->n_lambda, ir->expandedvals->init_wl_delta);
195 void init_parallel(FILE *log, t_commrec *cr, t_inputrec *inputrec,
198 bcast_ir_mtop(cr, inputrec, mtop);
200 if (inputrec->eI == eiBD || EI_SD(inputrec->eI) || ETC_ANDERSEN(inputrec->etc))
202 /* Make sure the random seeds are different on each node */
203 inputrec->ld_seed += cr->nodeid;