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