Avoid using function calls in OpenMP directives
[alexxy/gromacs.git] / src / gromacs / mdlib / mdatom.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42
43 #include "typedefs.h"
44 #include "mdatoms.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "main.h"
47 #include "qmmm.h"
48 #include "mtop_util.h"
49 #include "gmx_omp_nthreads.h"
50
51 #define ALMOST_ZERO 1e-30
52
53 t_mdatoms *init_mdatoms(FILE *fp, gmx_mtop_t *mtop, gmx_bool bFreeEnergy)
54 {
55     int                     mb, a, g, nmol;
56     double                  tmA, tmB;
57     t_atom                 *atom;
58     t_mdatoms              *md;
59     gmx_mtop_atomloop_all_t aloop;
60     t_ilist                *ilist;
61
62     snew(md, 1);
63
64     md->nenergrp = mtop->groups.grps[egcENER].nr;
65     md->bVCMgrps = FALSE;
66     tmA          = 0.0;
67     tmB          = 0.0;
68
69     aloop = gmx_mtop_atomloop_all_init(mtop);
70     while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
71     {
72         if (ggrpnr(&mtop->groups, egcVCM, a) > 0)
73         {
74             md->bVCMgrps = TRUE;
75         }
76
77         if (bFreeEnergy && PERTURBED(*atom))
78         {
79             md->nPerturbed++;
80             if (atom->mB != atom->m)
81             {
82                 md->nMassPerturbed++;
83             }
84             if (atom->qB != atom->q)
85             {
86                 md->nChargePerturbed++;
87             }
88             if (atom->typeB != atom->type)
89             {
90                 md->nTypePerturbed++;
91             }
92         }
93
94         tmA += atom->m;
95         tmB += atom->mB;
96     }
97
98     md->tmassA = tmA;
99     md->tmassB = tmB;
100
101     if (bFreeEnergy && fp)
102     {
103         fprintf(fp,
104                 "There are %d atoms and %d charges for free energy perturbation\n",
105                 md->nPerturbed, md->nChargePerturbed);
106     }
107
108     md->bOrires = gmx_mtop_ftype_count(mtop, F_ORIRES);
109
110     return md;
111 }
112
113 void atoms2md(gmx_mtop_t *mtop, t_inputrec *ir,
114               int nindex, int *index,
115               int homenr,
116               t_mdatoms *md)
117 {
118     gmx_bool              bLJPME;
119     gmx_mtop_atomlookup_t alook;
120     int                   i;
121     t_grpopts            *opts;
122     gmx_groups_t         *groups;
123     gmx_molblock_t       *molblock;
124     int                   nthreads gmx_unused;
125
126     bLJPME = EVDW_PME(ir->vdwtype);
127
128     opts = &ir->opts;
129
130     groups = &mtop->groups;
131
132     molblock = mtop->molblock;
133
134     /* Index==NULL indicates no DD (unless we have a DD node with no
135      * atoms), so also check for homenr. This should be
136      * signaled properly with an extra parameter or nindex==-1.
137      */
138     if (index == NULL && (homenr > 0))
139     {
140         md->nr = mtop->natoms;
141     }
142     else
143     {
144         md->nr = nindex;
145     }
146
147     if (md->nr > md->nalloc)
148     {
149         md->nalloc = over_alloc_dd(md->nr);
150
151         if (md->nMassPerturbed)
152         {
153             srenew(md->massA, md->nalloc);
154             srenew(md->massB, md->nalloc);
155         }
156         srenew(md->massT, md->nalloc);
157         srenew(md->invmass, md->nalloc);
158         srenew(md->chargeA, md->nalloc);
159         if (bLJPME)
160         {
161             srenew(md->sqrt_c6A, md->nalloc);
162             srenew(md->sigmaA, md->nalloc);
163             srenew(md->sigma3A, md->nalloc);
164         }
165         if (md->nPerturbed)
166         {
167             srenew(md->chargeB, md->nalloc);
168             if (bLJPME)
169             {
170                 srenew(md->sqrt_c6B, md->nalloc);
171                 srenew(md->sigmaB, md->nalloc);
172                 srenew(md->sigma3B, md->nalloc);
173             }
174         }
175         srenew(md->typeA, md->nalloc);
176         if (md->nPerturbed)
177         {
178             srenew(md->typeB, md->nalloc);
179         }
180         srenew(md->ptype, md->nalloc);
181         if (opts->ngtc > 1)
182         {
183             srenew(md->cTC, md->nalloc);
184             /* We always copy cTC with domain decomposition */
185         }
186         srenew(md->cENER, md->nalloc);
187         if (opts->ngacc > 1)
188         {
189             srenew(md->cACC, md->nalloc);
190         }
191         if (opts->nFreeze &&
192             (opts->ngfrz > 1 ||
193              opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
194         {
195             srenew(md->cFREEZE, md->nalloc);
196         }
197         if (md->bVCMgrps)
198         {
199             srenew(md->cVCM, md->nalloc);
200         }
201         if (md->bOrires)
202         {
203             srenew(md->cORF, md->nalloc);
204         }
205         if (md->nPerturbed)
206         {
207             srenew(md->bPerturbed, md->nalloc);
208         }
209
210         /* Note that these user t_mdatoms array pointers are NULL
211          * when there is only one group present.
212          * Therefore, when adding code, the user should use something like:
213          * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
214          */
215         if (mtop->groups.grpnr[egcUser1] != NULL)
216         {
217             srenew(md->cU1, md->nalloc);
218         }
219         if (mtop->groups.grpnr[egcUser2] != NULL)
220         {
221             srenew(md->cU2, md->nalloc);
222         }
223
224         if (ir->bQMMM)
225         {
226             srenew(md->bQM, md->nalloc);
227         }
228         if (ir->bAdress)
229         {
230             srenew(md->wf, md->nalloc);
231             srenew(md->tf_table_index, md->nalloc);
232         }
233     }
234
235     alook = gmx_mtop_atomlookup_init(mtop);
236
237     nthreads = gmx_omp_nthreads_get(emntDefault);
238 #pragma omp parallel for num_threads(nthreads) schedule(static)
239     for (i = 0; i < md->nr; i++)
240     {
241         int      g, ag, molb;
242         real     mA, mB, fac;
243         real     c6, c12;
244         t_atom  *atom;
245
246         if (index == NULL)
247         {
248             ag = i;
249         }
250         else
251         {
252             ag   = index[i];
253         }
254         gmx_mtop_atomnr_to_atom(alook, ag, &atom);
255
256         if (md->cFREEZE)
257         {
258             md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
259         }
260         if (EI_ENERGY_MINIMIZATION(ir->eI))
261         {
262             /* Displacement is proportional to F, masses used for constraints */
263             mA = 1.0;
264             mB = 1.0;
265         }
266         else if (ir->eI == eiBD)
267         {
268             /* With BD the physical masses are irrelevant.
269              * To keep the code simple we use most of the normal MD code path
270              * for BD. Thus for constraining the masses should be proportional
271              * to the friction coefficient. We set the absolute value such that
272              * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
273              * Then if we set the (meaningless) velocity to v=dx/dt, we get the
274              * correct kinetic energy and temperature using the usual code path.
275              * Thus with BD v*dt will give the displacement and the reported
276              * temperature can signal bad integration (too large time step).
277              */
278             if (ir->bd_fric > 0)
279             {
280                 mA = 0.5*ir->bd_fric*ir->delta_t;
281                 mB = 0.5*ir->bd_fric*ir->delta_t;
282             }
283             else
284             {
285                 /* The friction coefficient is mass/tau_t */
286                 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
287                 mA  = 0.5*atom->m*fac;
288                 mB  = 0.5*atom->mB*fac;
289             }
290         }
291         else
292         {
293             mA = atom->m;
294             mB = atom->mB;
295         }
296         if (md->nMassPerturbed)
297         {
298             md->massA[i]  = mA;
299             md->massB[i]  = mB;
300         }
301         md->massT[i]    = mA;
302         if (mA == 0.0)
303         {
304             md->invmass[i]    = 0;
305         }
306         else if (md->cFREEZE)
307         {
308             g = md->cFREEZE[i];
309             if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
310             {
311                 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
312                  * to avoid div by zero in lincs or shake.
313                  * Note that constraints can still move a partially frozen particle.
314                  */
315                 md->invmass[i]  = ALMOST_ZERO;
316             }
317             else
318             {
319                 md->invmass[i]  = 1.0/mA;
320             }
321         }
322         else
323         {
324             md->invmass[i]    = 1.0/mA;
325         }
326         md->chargeA[i]      = atom->q;
327         md->typeA[i]        = atom->type;
328         if (bLJPME)
329         {
330             c6                = mtop->ffparams.iparams[atom->type*(mtop->ffparams.atnr+1)].lj.c6;
331             c12               = mtop->ffparams.iparams[atom->type*(mtop->ffparams.atnr+1)].lj.c12;
332             md->sqrt_c6A[i]   = sqrt(c6);
333             if (c6 == 0.0 || c12 == 0)
334             {
335                 md->sigmaA[i] = 1.0;
336             }
337             else
338             {
339                 md->sigmaA[i] = pow(c12/c6, 1.0/6.0);
340             }
341             md->sigma3A[i]    = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
342         }
343         if (md->nPerturbed)
344         {
345             md->bPerturbed[i] = PERTURBED(*atom);
346             md->chargeB[i]    = atom->qB;
347             md->typeB[i]      = atom->typeB;
348             if (bLJPME)
349             {
350                 c6                = mtop->ffparams.iparams[atom->typeB*(mtop->ffparams.atnr+1)].lj.c6;
351                 c12               = mtop->ffparams.iparams[atom->typeB*(mtop->ffparams.atnr+1)].lj.c12;
352                 md->sqrt_c6B[i]   = sqrt(c6);
353                 if (c6 == 0.0 || c12 == 0)
354                 {
355                     md->sigmaB[i] = 1.0;
356                 }
357                 else
358                 {
359                     md->sigmaB[i] = pow(c12/c6, 1.0/6.0);
360                 }
361                 md->sigma3B[i]    = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
362             }
363         }
364         md->ptype[i]    = atom->ptype;
365         if (md->cTC)
366         {
367             md->cTC[i]    = groups->grpnr[egcTC][ag];
368         }
369         md->cENER[i]    =
370             (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
371         if (md->cACC)
372         {
373             md->cACC[i]   = groups->grpnr[egcACC][ag];
374         }
375         if (md->cVCM)
376         {
377             md->cVCM[i]       = groups->grpnr[egcVCM][ag];
378         }
379         if (md->cORF)
380         {
381             md->cORF[i]       = groups->grpnr[egcORFIT][ag];
382         }
383
384         if (md->cU1)
385         {
386             md->cU1[i]        = groups->grpnr[egcUser1][ag];
387         }
388         if (md->cU2)
389         {
390             md->cU2[i]        = groups->grpnr[egcUser2][ag];
391         }
392
393         if (ir->bQMMM)
394         {
395             if (groups->grpnr[egcQMMM] == 0 ||
396                 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
397             {
398                 md->bQM[i]      = TRUE;
399             }
400             else
401             {
402                 md->bQM[i]      = FALSE;
403             }
404         }
405         /* Initialize AdResS weighting functions to adressw */
406         if (ir->bAdress)
407         {
408             md->wf[i]           = 1.0;
409             /* if no tf table groups specified, use default table */
410             md->tf_table_index[i] = DEFAULT_TF_TABLE;
411             if (ir->adress->n_tf_grps > 0)
412             {
413                 /* if tf table groups specified, tf is only applied to thoose energy groups*/
414                 md->tf_table_index[i] = NO_TF_TABLE;
415                 /* check wether atom is in one of the relevant energy groups and assign a table index */
416                 for (g = 0; g < ir->adress->n_tf_grps; g++)
417                 {
418                     if (md->cENER[i] == ir->adress->tf_table_index[g])
419                     {
420                         md->tf_table_index[i] = g;
421                     }
422                 }
423             }
424         }
425     }
426
427     gmx_mtop_atomlookup_destroy(alook);
428
429     md->homenr = homenr;
430     md->lambda = 0;
431 }
432
433 void update_mdatoms(t_mdatoms *md, real lambda)
434 {
435     int    al, end;
436     real   L1 = 1.0-lambda;
437
438     end = md->nr;
439
440     if (md->nMassPerturbed)
441     {
442         for (al = 0; (al < end); al++)
443         {
444             if (md->bPerturbed[al])
445             {
446                 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
447                 if (md->invmass[al] > 1.1*ALMOST_ZERO)
448                 {
449                     md->invmass[al] = 1.0/md->massT[al];
450                 }
451             }
452         }
453         md->tmass = L1*md->tmassA + lambda*md->tmassB;
454     }
455     else
456     {
457         md->tmass = md->tmassA;
458     }
459     md->lambda = lambda;
460 }