Code beautification with uncrustify
[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     {
69         if (ggrpnr(&mtop->groups, egcVCM, a) > 0)
70         {
71             md->bVCMgrps = TRUE;
72         }
73
74         if (bFreeEnergy && PERTURBED(*atom))
75         {
76             md->nPerturbed++;
77             if (atom->mB != atom->m)
78             {
79                 md->nMassPerturbed++;
80             }
81             if (atom->qB != atom->q)
82             {
83                 md->nChargePerturbed++;
84             }
85         }
86
87         tmA += atom->m;
88         tmB += atom->mB;
89     }
90
91     md->tmassA = tmA;
92     md->tmassB = tmB;
93
94     if (bFreeEnergy && fp)
95     {
96         fprintf(fp,
97                 "There are %d atoms and %d charges for free energy perturbation\n",
98                 md->nPerturbed, md->nChargePerturbed);
99     }
100
101     md->bOrires = gmx_mtop_ftype_count(mtop, F_ORIRES);
102
103     return md;
104 }
105
106 void atoms2md(gmx_mtop_t *mtop, t_inputrec *ir,
107               int nindex, int *index,
108               int start, int homenr,
109               t_mdatoms *md)
110 {
111     gmx_mtop_atomlookup_t alook;
112     int                   i;
113     t_grpopts            *opts;
114     gmx_groups_t         *groups;
115     gmx_molblock_t       *molblock;
116
117     opts = &ir->opts;
118
119     groups = &mtop->groups;
120
121     molblock = mtop->molblock;
122
123     /* Index==NULL indicates particle decomposition,
124      * unless we have an empty DD node, so also check for homenr and start.
125      * This should be signaled properly with an extra parameter or nindex==-1.
126      */
127     if (index == NULL && (homenr > 0 || start > 0))
128     {
129         md->nr = mtop->natoms;
130     }
131     else
132     {
133         md->nr = nindex;
134     }
135
136     if (md->nr > md->nalloc)
137     {
138         md->nalloc = over_alloc_dd(md->nr);
139
140         if (md->nMassPerturbed)
141         {
142             srenew(md->massA, md->nalloc);
143             srenew(md->massB, md->nalloc);
144         }
145         srenew(md->massT, md->nalloc);
146         srenew(md->invmass, md->nalloc);
147         srenew(md->chargeA, md->nalloc);
148         if (md->nPerturbed)
149         {
150             srenew(md->chargeB, md->nalloc);
151         }
152         srenew(md->typeA, md->nalloc);
153         if (md->nPerturbed)
154         {
155             srenew(md->typeB, md->nalloc);
156         }
157         srenew(md->ptype, md->nalloc);
158         if (opts->ngtc > 1)
159         {
160             srenew(md->cTC, md->nalloc);
161             /* We always copy cTC with domain decomposition */
162         }
163         srenew(md->cENER, md->nalloc);
164         if (opts->ngacc > 1)
165         {
166             srenew(md->cACC, md->nalloc);
167         }
168         if (opts->nFreeze &&
169             (opts->ngfrz > 1 ||
170              opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
171         {
172             srenew(md->cFREEZE, md->nalloc);
173         }
174         if (md->bVCMgrps)
175         {
176             srenew(md->cVCM, md->nalloc);
177         }
178         if (md->bOrires)
179         {
180             srenew(md->cORF, md->nalloc);
181         }
182         if (md->nPerturbed)
183         {
184             srenew(md->bPerturbed, md->nalloc);
185         }
186
187         /* Note that these user t_mdatoms array pointers are NULL
188          * when there is only one group present.
189          * Therefore, when adding code, the user should use something like:
190          * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
191          */
192         if (mtop->groups.grpnr[egcUser1] != NULL)
193         {
194             srenew(md->cU1, md->nalloc);
195         }
196         if (mtop->groups.grpnr[egcUser2] != NULL)
197         {
198             srenew(md->cU2, md->nalloc);
199         }
200
201         if (ir->bQMMM)
202         {
203             srenew(md->bQM, md->nalloc);
204         }
205         if (ir->bAdress)
206         {
207             srenew(md->wf, md->nalloc);
208             srenew(md->tf_table_index, md->nalloc);
209         }
210     }
211
212     alook = gmx_mtop_atomlookup_init(mtop);
213
214 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntDefault)) schedule(static)
215     for (i = 0; i < md->nr; i++)
216     {
217         int      g, ag, molb;
218         real     mA, mB, fac;
219         t_atom  *atom;
220
221         if (index == NULL)
222         {
223             ag = i;
224         }
225         else
226         {
227             ag   = index[i];
228         }
229         gmx_mtop_atomnr_to_atom(alook, ag, &atom);
230
231         if (md->cFREEZE)
232         {
233             md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
234         }
235         if (EI_ENERGY_MINIMIZATION(ir->eI))
236         {
237             /* Displacement is proportional to F, masses used for constraints */
238             mA = 1.0;
239             mB = 1.0;
240         }
241         else if (ir->eI == eiBD)
242         {
243             /* With BD the physical masses are irrelevant.
244              * To keep the code simple we use most of the normal MD code path
245              * for BD. Thus for constraining the masses should be proportional
246              * to the friction coefficient. We set the absolute value such that
247              * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
248              * Then if we set the (meaningless) velocity to v=dx/dt, we get the
249              * correct kinetic energy and temperature using the usual code path.
250              * Thus with BD v*dt will give the displacement and the reported
251              * temperature can signal bad integration (too large time step).
252              */
253             if (ir->bd_fric > 0)
254             {
255                 mA = 0.5*ir->bd_fric*ir->delta_t;
256                 mB = 0.5*ir->bd_fric*ir->delta_t;
257             }
258             else
259             {
260                 /* The friction coefficient is mass/tau_t */
261                 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
262                 mA  = 0.5*atom->m*fac;
263                 mB  = 0.5*atom->mB*fac;
264             }
265         }
266         else
267         {
268             mA = atom->m;
269             mB = atom->mB;
270         }
271         if (md->nMassPerturbed)
272         {
273             md->massA[i]  = mA;
274             md->massB[i]  = mB;
275         }
276         md->massT[i]    = mA;
277         if (mA == 0.0)
278         {
279             md->invmass[i]    = 0;
280         }
281         else if (md->cFREEZE)
282         {
283             g = md->cFREEZE[i];
284             if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
285             {
286                 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
287                  * to avoid div by zero in lincs or shake.
288                  * Note that constraints can still move a partially frozen particle.
289                  */
290                 md->invmass[i]  = ALMOST_ZERO;
291             }
292             else
293             {
294                 md->invmass[i]  = 1.0/mA;
295             }
296         }
297         else
298         {
299             md->invmass[i]    = 1.0/mA;
300         }
301         md->chargeA[i]  = atom->q;
302         md->typeA[i]    = atom->type;
303         if (md->nPerturbed)
304         {
305             md->chargeB[i]    = atom->qB;
306             md->typeB[i]      = atom->typeB;
307             md->bPerturbed[i] = PERTURBED(*atom);
308         }
309         md->ptype[i]    = atom->ptype;
310         if (md->cTC)
311         {
312             md->cTC[i]    = groups->grpnr[egcTC][ag];
313         }
314         md->cENER[i]    =
315             (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
316         if (md->cACC)
317         {
318             md->cACC[i]   = groups->grpnr[egcACC][ag];
319         }
320         if (md->cVCM)
321         {
322             md->cVCM[i]       = groups->grpnr[egcVCM][ag];
323         }
324         if (md->cORF)
325         {
326             md->cORF[i]       = groups->grpnr[egcORFIT][ag];
327         }
328
329         if (md->cU1)
330         {
331             md->cU1[i]        = groups->grpnr[egcUser1][ag];
332         }
333         if (md->cU2)
334         {
335             md->cU2[i]        = groups->grpnr[egcUser2][ag];
336         }
337
338         if (ir->bQMMM)
339         {
340             if (groups->grpnr[egcQMMM] == 0 ||
341                 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
342             {
343                 md->bQM[i]      = TRUE;
344             }
345             else
346             {
347                 md->bQM[i]      = FALSE;
348             }
349         }
350         /* Initialize AdResS weighting functions to adressw */
351         if (ir->bAdress)
352         {
353             md->wf[i]           = 1.0;
354             /* if no tf table groups specified, use default table */
355             md->tf_table_index[i] = DEFAULT_TF_TABLE;
356             if (ir->adress->n_tf_grps > 0)
357             {
358                 /* if tf table groups specified, tf is only applied to thoose energy groups*/
359                 md->tf_table_index[i] = NO_TF_TABLE;
360                 /* check wether atom is in one of the relevant energy groups and assign a table index */
361                 for (g = 0; g < ir->adress->n_tf_grps; g++)
362                 {
363                     if (md->cENER[i] == ir->adress->tf_table_index[g])
364                     {
365                         md->tf_table_index[i] = g;
366                     }
367                 }
368             }
369         }
370     }
371
372     gmx_mtop_atomlookup_destroy(alook);
373
374     md->start  = start;
375     md->homenr = homenr;
376     md->lambda = 0;
377 }
378
379 void update_mdatoms(t_mdatoms *md, real lambda)
380 {
381     int    al, end;
382     real   L1 = 1.0-lambda;
383
384     end = md->nr;
385
386     if (md->nMassPerturbed)
387     {
388         for (al = 0; (al < end); al++)
389         {
390             if (md->bPerturbed[al])
391             {
392                 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
393                 if (md->invmass[al] > 1.1*ALMOST_ZERO)
394                 {
395                     md->invmass[al] = 1.0/md->massT[al];
396                 }
397             }
398         }
399         md->tmass = L1*md->tmassA + lambda*md->tmassB;
400     }
401     else
402     {
403         md->tmass = md->tmassA;
404     }
405     md->lambda = lambda;
406 }