Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / mdatom.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.2.0
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.
15
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.
20  * 
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.
27  * 
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.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * GROwing Monsters And Cloning Shrimps
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include "typedefs.h"
41 #include "mdatoms.h"
42 #include "smalloc.h"
43 #include "main.h"
44 #include "qmmm.h"
45 #include "mtop_util.h"
46 #include "gmx_omp_nthreads.h"
47
48 #define ALMOST_ZERO 1e-30
49
50 t_mdatoms *init_mdatoms(FILE *fp,gmx_mtop_t *mtop,gmx_bool bFreeEnergy)
51 {
52   int    mb,a,g,nmol;
53   double tmA,tmB;
54   t_atom *atom;
55   t_mdatoms *md;
56   gmx_mtop_atomloop_all_t aloop;
57   t_ilist *ilist;
58
59   snew(md,1);
60
61   md->nenergrp = mtop->groups.grps[egcENER].nr;
62   md->bVCMgrps = FALSE;
63   tmA = 0.0;
64   tmB = 0.0;
65
66   aloop = gmx_mtop_atomloop_all_init(mtop);
67   while(gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
68     if (ggrpnr(&mtop->groups,egcVCM,a) > 0)
69       md->bVCMgrps = TRUE;
70     
71     if (bFreeEnergy && PERTURBED(*atom)) {
72       md->nPerturbed++;
73       if (atom->mB != atom->m)
74         md->nMassPerturbed++;
75       if (atom->qB != atom->q)
76         md->nChargePerturbed++;
77     }
78     
79     tmA += atom->m;
80     tmB += atom->mB;
81   }
82
83   md->tmassA = tmA;
84   md->tmassB = tmB;
85   
86   if (bFreeEnergy && fp)
87     fprintf(fp,
88             "There are %d atoms and %d charges for free energy perturbation\n",
89             md->nPerturbed,md->nChargePerturbed);
90
91   md->bOrires = gmx_mtop_ftype_count(mtop,F_ORIRES);
92
93   return md;
94 }
95
96 void atoms2md(gmx_mtop_t *mtop,t_inputrec *ir,
97               int nindex,int *index,
98               int start,int homenr,
99               t_mdatoms *md)
100 {
101   gmx_mtop_atomlookup_t alook;
102   int       i;
103   t_grpopts *opts;
104   gmx_groups_t *groups;
105   gmx_molblock_t *molblock;
106
107   opts = &ir->opts;
108
109   groups = &mtop->groups;
110
111   molblock = mtop->molblock;
112
113   /* Index==NULL indicates particle decomposition,
114    * unless we have an empty DD node, so also check for homenr and start.
115    * This should be signaled properly with an extra parameter or nindex==-1.
116    */
117   if (index == NULL && (homenr > 0 || start > 0)) {
118     md->nr = mtop->natoms;
119   } else {
120     md->nr = nindex;
121   }
122
123   if (md->nr > md->nalloc) {
124     md->nalloc = over_alloc_dd(md->nr);
125
126     if (md->nMassPerturbed) {
127       srenew(md->massA,md->nalloc);
128       srenew(md->massB,md->nalloc);
129     }
130     srenew(md->massT,md->nalloc);
131     srenew(md->invmass,md->nalloc);
132     srenew(md->chargeA,md->nalloc);
133     if (md->nPerturbed) {
134       srenew(md->chargeB,md->nalloc);
135     }
136     srenew(md->typeA,md->nalloc);
137     if (md->nPerturbed) {
138       srenew(md->typeB,md->nalloc);
139     }
140     srenew(md->ptype,md->nalloc);
141     if (opts->ngtc > 1) {
142       srenew(md->cTC,md->nalloc);
143       /* We always copy cTC with domain decomposition */
144     }
145     srenew(md->cENER,md->nalloc);
146     if (opts->ngacc > 1)
147       srenew(md->cACC,md->nalloc);
148     if (opts->nFreeze &&
149         (opts->ngfrz > 1 ||
150          opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
151       srenew(md->cFREEZE,md->nalloc);
152     if (md->bVCMgrps)
153       srenew(md->cVCM,md->nalloc);
154     if (md->bOrires)
155       srenew(md->cORF,md->nalloc);
156     if (md->nPerturbed)
157       srenew(md->bPerturbed,md->nalloc);
158     
159     /* Note that these user t_mdatoms array pointers are NULL
160      * when there is only one group present.
161      * Therefore, when adding code, the user should use something like:
162      * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
163      */
164     if (mtop->groups.grpnr[egcUser1] != NULL)
165       srenew(md->cU1,md->nalloc);
166     if (mtop->groups.grpnr[egcUser2] != NULL)
167       srenew(md->cU2,md->nalloc);
168     
169     if (ir->bQMMM)
170       srenew(md->bQM,md->nalloc);
171     if (ir->bAdress) {
172       srenew(md->wf,md->nalloc);
173       srenew(md->tf_table_index,md->nalloc);
174     }
175   }
176
177   alook = gmx_mtop_atomlookup_init(mtop);
178
179 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntDefault)) schedule(static)
180   for(i=0; i<md->nr; i++) {
181     int     g,ag,molb;
182     real    mA,mB,fac;
183     t_atom  *atom;
184
185     if (index == NULL) {
186       ag = i;
187     } else {
188       ag   = index[i];
189     }
190     gmx_mtop_atomnr_to_atom(alook,ag,&atom);
191
192     if (md->cFREEZE) {
193       md->cFREEZE[i] = ggrpnr(groups,egcFREEZE,ag);
194     }
195         if (EI_ENERGY_MINIMIZATION(ir->eI))
196         {
197             /* Displacement is proportional to F, masses used for constraints */
198             mA = 1.0;
199             mB = 1.0;
200         }
201         else if (ir->eI == eiBD)
202         {
203             /* With BD the physical masses are irrelevant.
204              * To keep the code simple we use most of the normal MD code path
205              * for BD. Thus for constraining the masses should be proportional
206              * to the friction coefficient. We set the absolute value such that
207              * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
208              * Then if we set the (meaningless) velocity to v=dx/dt, we get the
209              * correct kinetic energy and temperature using the usual code path.
210              * Thus with BD v*dt will give the displacement and the reported
211              * temperature can signal bad integration (too large time step).
212              */
213             if (ir->bd_fric > 0)
214             {
215                 mA = 0.5*ir->bd_fric*ir->delta_t;
216                 mB = 0.5*ir->bd_fric*ir->delta_t;
217             }
218             else
219             {
220                 /* The friction coefficient is mass/tau_t */
221                 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
222                 mA = 0.5*atom->m*fac;
223                 mB = 0.5*atom->mB*fac;
224             }
225         }
226         else
227         {
228             mA = atom->m;
229             mB = atom->mB;
230         }
231     if (md->nMassPerturbed) {
232       md->massA[i]      = mA;
233       md->massB[i]      = mB;
234     }
235     md->massT[i]        = mA;
236     if (mA == 0.0) {
237       md->invmass[i]    = 0;
238     } else if (md->cFREEZE) {
239       g = md->cFREEZE[i];
240       if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
241         /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
242          * to avoid div by zero in lincs or shake.
243          * Note that constraints can still move a partially frozen particle.
244          */
245         md->invmass[i]  = ALMOST_ZERO;
246       else
247         md->invmass[i]  = 1.0/mA;
248     } else {
249       md->invmass[i]    = 1.0/mA;
250     }
251     md->chargeA[i]      = atom->q;
252     md->typeA[i]        = atom->type;
253     if (md->nPerturbed) {
254       md->chargeB[i]    = atom->qB;
255       md->typeB[i]      = atom->typeB;
256       md->bPerturbed[i] = PERTURBED(*atom);
257     }
258     md->ptype[i]        = atom->ptype;
259     if (md->cTC)
260       md->cTC[i]        = groups->grpnr[egcTC][ag];
261     md->cENER[i]        =
262       (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
263     if (md->cACC)
264       md->cACC[i]       = groups->grpnr[egcACC][ag];
265     if (md->cVCM)
266       md->cVCM[i]       = groups->grpnr[egcVCM][ag];
267     if (md->cORF)
268       md->cORF[i]       = groups->grpnr[egcORFIT][ag];
269
270     if (md->cU1)
271       md->cU1[i]        = groups->grpnr[egcUser1][ag];
272     if (md->cU2)
273       md->cU2[i]        = groups->grpnr[egcUser2][ag];
274
275     if (ir->bQMMM) {
276       if (groups->grpnr[egcQMMM] == 0 || 
277           groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1) {
278         md->bQM[i]      = TRUE;
279       } else {
280         md->bQM[i]      = FALSE;
281       }
282     }
283     /* Initialize AdResS weighting functions to adressw */
284     if (ir->bAdress){
285        md->wf[i]           = 1.0;
286         /* if no tf table groups specified, use default table */
287        md->tf_table_index[i] = DEFAULT_TF_TABLE;
288        if (ir->adress->n_tf_grps > 0){
289             /* if tf table groups specified, tf is only applied to thoose energy groups*/
290             md->tf_table_index[i] = NO_TF_TABLE;
291             /* check wether atom is in one of the relevant energy groups and assign a table index */
292             for (g=0; g<ir->adress->n_tf_grps; g++){
293                 if (md->cENER[i] == ir->adress->tf_table_index[g]){
294                    md->tf_table_index[i] = g;
295                 }
296             }
297         }
298     }
299   }
300
301   gmx_mtop_atomlookup_destroy(alook);
302
303   md->start  = start;
304   md->homenr = homenr;
305   md->lambda = 0;
306 }
307
308 void update_mdatoms(t_mdatoms *md,real lambda)
309 {
310   int    al,end;
311   real   L1=1.0-lambda;
312   
313   end=md->nr;
314
315   if (md->nMassPerturbed) {
316     for(al=0; (al<end); al++) {
317       if (md->bPerturbed[al]) {
318         md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
319         if (md->invmass[al] > 1.1*ALMOST_ZERO)
320           md->invmass[al] = 1.0/md->massT[al];
321       }
322     }
323     md->tmass = L1*md->tmassA + lambda*md->tmassB;
324   } else {
325     md->tmass = md->tmassA;
326   }
327   md->lambda = lambda;
328 }