Improved version for builds from modified trees.
[gromacs/qmmm-gamess-us.git] / src / mdlib / init.c
blob860dc0e1bfecdd7fdcfa3772aa00bd8525201973
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GROwing Monsters And Cloning Shrimps
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <stdio.h>
40 #include "typedefs.h"
41 #include "tpxio.h"
42 #include "smalloc.h"
43 #include "vec.h"
44 #include "main.h"
45 #include "mvdata.h"
46 #include "gmx_fatal.h"
47 #include "symtab.h"
48 #include "txtdump.h"
49 #include "mdatoms.h"
50 #include "mdrun.h"
51 #include "statutil.h"
52 #include "names.h"
53 #include "calcgrid.h"
54 #include "gmx_random.h"
55 #include "update.h"
56 #include "mdebin.h"
58 #define BUFSIZE 256
60 #define NOT_FINISHED(l1,l2) \
61 printf("not finished yet: lines %d .. %d in %s\n",l1,l2,__FILE__)
63 static char *int_title(const char *title,int nodeid,char buf[], int size)
65 sprintf(buf,"%s (%d)",title,nodeid);
67 return buf;
70 static void set_state_entries(t_state *state,t_inputrec *ir,int nnodes)
72 int nnhpres;
74 /* The entries in the state in the tpx file might not correspond
75 * with what is needed, so we correct this here.
77 state->flags = 0;
78 if (ir->efep != efepNO)
79 state->flags |= (1<<estLAMBDA);
80 state->flags |= (1<<estX);
81 if (state->x == NULL)
82 snew(state->x,state->nalloc);
83 if (EI_DYNAMICS(ir->eI)) {
84 state->flags |= (1<<estV);
85 if (state->v == NULL)
86 snew(state->v,state->nalloc);
88 if (ir->eI == eiSD2) {
89 state->flags |= (1<<estSDX);
90 if (state->sd_X == NULL) {
91 /* sd_X is not stored in the tpx file, so we need to allocate it */
92 snew(state->sd_X,state->nalloc);
95 if (ir->eI == eiCG) {
96 state->flags |= (1<<estCGP);
98 if (EI_SD(ir->eI) || ir->eI == eiBD || ir->etc == etcVRESCALE) {
99 state->nrng = gmx_rng_n();
100 state->nrngi = 1;
101 if (EI_SD(ir->eI) || ir->eI == eiBD) {
102 /* This will be correct later with DD */
103 state->nrng *= nnodes;
104 state->nrngi *= nnodes;
106 state->flags |= ((1<<estLD_RNG) | (1<<estLD_RNGI));
107 snew(state->ld_rng, state->nrng);
108 snew(state->ld_rngi,state->nrngi);
109 } else {
110 state->nrng = 0;
112 state->nnhpres = 0;
113 if (ir->ePBC != epbcNONE) {
114 state->flags |= (1<<estBOX);
115 if (PRESERVE_SHAPE(*ir)) {
116 state->flags |= (1<<estBOX_REL);
118 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK)) {
119 state->flags |= (1<<estBOXV);
121 if (ir->epc != epcNO) {
122 if (IR_NPT_TROTTER(ir)) {
123 state->nnhpres = 1;
124 state->flags |= (1<<estNHPRES_XI);
125 state->flags |= (1<<estNHPRES_VXI);
126 state->flags |= (1<<estSVIR_PREV);
127 state->flags |= (1<<estFVIR_PREV);
128 state->flags |= (1<<estVETA);
129 state->flags |= (1<<estVOL0);
130 } else {
131 state->flags |= (1<<estPRES_PREV);
136 if (ir->etc == etcNOSEHOOVER) {
137 state->flags |= (1<<estNH_XI);
138 state->flags |= (1<<estNH_VXI);
141 if (ir->etc == etcVRESCALE) {
142 state->flags |= (1<<estTC_INT);
145 init_gtc_state(state,state->ngtc,state->nnhpres,ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
146 init_ekinstate(&state->ekinstate,ir);
148 init_energyhistory(&state->enerhist);
151 void init_single(FILE *fplog,t_inputrec *inputrec,
152 const char *tpxfile,gmx_mtop_t *mtop,
153 t_state *state)
155 read_tpx_state(tpxfile,inputrec,state,NULL,mtop);
156 set_state_entries(state,inputrec,1);
158 if (fplog)
159 pr_inputrec(fplog,0,"Input Parameters",inputrec,FALSE);
162 void init_parallel(FILE *log,const char *tpxfile,t_commrec *cr,
163 t_inputrec *inputrec,gmx_mtop_t *mtop,
164 t_state *state,
165 int list)
167 char buf[256];
169 if (MASTER(cr)) {
170 init_inputrec(inputrec);
171 read_tpx_state(tpxfile,inputrec,state,NULL,mtop);
172 /* When we will be doing domain decomposition with separate PME nodes
173 * the rng entries will be too large, we correct for this later.
175 set_state_entries(state,inputrec,cr->nnodes);
177 bcast_ir_mtop(cr,inputrec,mtop);
179 if (inputrec->eI == eiBD || EI_SD(inputrec->eI)) {
180 /* Make sure the random seeds are different on each node */
181 inputrec->ld_seed += cr->nodeid;
184 /* Printing */
185 if (list!=0 && log!=NULL)
187 if (list&LIST_INPUTREC)
188 pr_inputrec(log,0,"parameters of the run",inputrec,FALSE);
189 if (list&LIST_X)
190 pr_rvecs(log,0,"box",state->box,DIM);
191 if (list&LIST_X)
192 pr_rvecs(log,0,"box_rel",state->box_rel,DIM);
193 if (list&LIST_V)
194 pr_rvecs(log,0,"boxv",state->boxv,DIM);
195 if (list&LIST_X)
196 pr_rvecs(log,0,int_title("x",0,buf,255),state->x,state->natoms);
197 if (list&LIST_V)
198 pr_rvecs(log,0,int_title("v",0,buf,255),state->v,state->natoms);
199 if (list&LIST_TOP)
200 pr_mtop(log,0,int_title("topology",cr->nodeid,buf,255),mtop,TRUE);
201 fflush(log);