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