Merge branch 'release-4-5-patches'
[alexxy/gromacs.git] / src / gromacs / mdlib / mdatom.c
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * GROwing Monsters And Cloning Shrimps
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "typedefs.h"
40 #include "mdatoms.h"
41 #include "smalloc.h"
42 #include "main.h"
43 #include "qmmm.h"
44 #include "mtop_util.h"
45
46 #define ALMOST_ZERO 1e-30
47
48 t_mdatoms *init_mdatoms(FILE *fp,gmx_mtop_t *mtop,gmx_bool bFreeEnergy)
49 {
50   int    mb,a,g,nmol;
51   double tmA,tmB;
52   t_atom *atom;
53   t_mdatoms *md;
54   gmx_mtop_atomloop_all_t aloop;
55   t_ilist *ilist;
56
57   snew(md,1);
58
59   md->nenergrp = mtop->groups.grps[egcENER].nr;
60   md->bVCMgrps = FALSE;
61   tmA = 0.0;
62   tmB = 0.0;
63
64   aloop = gmx_mtop_atomloop_all_init(mtop);
65   while(gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
66     if (ggrpnr(&mtop->groups,egcVCM,a) > 0)
67       md->bVCMgrps = TRUE;
68     
69     if (bFreeEnergy && PERTURBED(*atom)) {
70       md->nPerturbed++;
71       if (atom->mB != atom->m)
72         md->nMassPerturbed++;
73       if (atom->qB != atom->q)
74         md->nChargePerturbed++;
75     }
76     
77     tmA += atom->m;
78     tmB += atom->mB;
79   }
80
81   md->tmassA = tmA;
82   md->tmassB = tmB;
83   
84   if (bFreeEnergy && fp)
85     fprintf(fp,
86             "There are %d atoms and %d charges for free energy perturbation\n",
87             md->nPerturbed,md->nChargePerturbed);
88
89   md->bOrires = gmx_mtop_ftype_count(mtop,F_ORIRES);
90
91   return md;
92 }
93
94 void atoms2md(gmx_mtop_t *mtop,t_inputrec *ir,
95               int nindex,int *index,
96               int start,int homenr,
97               t_mdatoms *md)
98 {
99   t_atoms   *atoms_mol;
100   int       i,g,ag,as,ae,molb;
101   real      mA,mB,fac;
102   t_atom    *atom;
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   }
172
173   for(i=0; (i<md->nr); i++) {
174     if (index == NULL) {
175       ag = i;
176       gmx_mtop_atomnr_to_atom(mtop,ag,&atom);
177     } else {
178       ag   = index[i];
179       molb = -1;
180       ae   = 0;
181       do {
182         molb++;
183         as = ae;
184         ae = as + molblock[molb].nmol*molblock[molb].natoms_mol;
185       } while (ag >= ae);
186       atoms_mol = &mtop->moltype[molblock[molb].type].atoms;
187       atom = &atoms_mol->atom[(ag - as) % atoms_mol->nr];
188     }
189
190     if (md->cFREEZE) {
191       md->cFREEZE[i] = ggrpnr(groups,egcFREEZE,ag);
192     }
193     if (EI_ENERGY_MINIMIZATION(ir->eI)) {
194       mA = 1.0;
195       mB = 1.0;
196     } else if (ir->eI == eiBD) {
197       /* Make the mass proportional to the friction coefficient for BD.
198        * This is necessary for the constraint algorithms.
199        */
200       if (ir->bd_fric) {
201         mA = ir->bd_fric*ir->delta_t;
202         mB = ir->bd_fric*ir->delta_t;
203       } else {
204         fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
205         mA = atom->m*fac;
206         mB = atom->mB*fac;
207       }
208     } else {
209       mA = atom->m;
210       mB = atom->mB;
211     }
212     if (md->nMassPerturbed) {
213       md->massA[i]      = mA;
214       md->massB[i]      = mB;
215     }
216     md->massT[i]        = mA;
217     if (mA == 0.0) {
218       md->invmass[i]    = 0;
219     } else if (md->cFREEZE) {
220       g = md->cFREEZE[i];
221       if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
222         /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
223          * to avoid div by zero in lincs or shake.
224          * Note that constraints can still move a partially frozen particle.
225          */
226         md->invmass[i]  = ALMOST_ZERO;
227       else
228         md->invmass[i]  = 1.0/mA;
229     } else {
230       md->invmass[i]    = 1.0/mA;
231     }
232     md->chargeA[i]      = atom->q;
233     md->typeA[i]        = atom->type;
234     if (md->nPerturbed) {
235       md->chargeB[i]    = atom->qB;
236       md->typeB[i]      = atom->typeB;
237       md->bPerturbed[i] = PERTURBED(*atom);
238     }
239     md->ptype[i]        = atom->ptype;
240     if (md->cTC)
241       md->cTC[i]        = groups->grpnr[egcTC][ag];
242     md->cENER[i]        =
243       (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
244     if (md->cACC)
245       md->cACC[i]       = groups->grpnr[egcACC][ag];
246     if (md->cVCM)
247       md->cVCM[i]       = groups->grpnr[egcVCM][ag];
248     if (md->cORF)
249       md->cORF[i]       = groups->grpnr[egcORFIT][ag];
250
251     if (md->cU1)
252       md->cU1[i]        = groups->grpnr[egcUser1][ag];
253     if (md->cU2)
254       md->cU2[i]        = groups->grpnr[egcUser2][ag];
255
256     if (ir->bQMMM) {
257       if (groups->grpnr[egcQMMM] == 0 || 
258           groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1) {
259         md->bQM[i]      = TRUE;
260       } else {
261         md->bQM[i]      = FALSE;
262       }
263     }
264   }
265
266   md->start  = start;
267   md->homenr = homenr;
268   md->lambda = 0;
269 }
270
271 void update_mdatoms(t_mdatoms *md,real lambda)
272 {
273   int    al,end;
274   real   L1=1.0-lambda;
275   
276   end=md->nr;
277
278   if (md->nMassPerturbed) {
279     for(al=0; (al<end); al++) {
280       if (md->bPerturbed[al]) {
281         md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
282         if (md->invmass[al] > 1.1*ALMOST_ZERO)
283           md->invmass[al] = 1.0/md->massT[al];
284       }
285     }
286     md->tmass = L1*md->tmassA + lambda*md->tmassB;
287   } else {
288     md->tmass = md->tmassA;
289   }
290   md->lambda = lambda;
291 }