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