ecadf3ab48b40ab7b274c07c6892bc74047fe11e
[alexxy/gromacs.git] / src / gromacs / mdlib / vcm.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) 2013,2014,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "vcm.h"
42
43 #include "gromacs/math/functions.h"
44 #include "gromacs/math/invertmatrix.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/math/vecdump.h"
47 #include "gromacs/mdlib/gmx_omp_nthreads.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/topology/topology.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/gmxomp.h"
55 #include "gromacs/utility/smalloc.h"
56
57 t_vcm::t_vcm(const SimulationGroups& groups, const t_inputrec& ir) :
58     integratorConservesMomentum(!EI_RANDOM(ir.eI))
59 {
60     mode     = (ir.nstcomm > 0) ? ir.comm_mode : ecmNO;
61     ndim     = ndof_com(&ir);
62     timeStep = ir.nstcomm * ir.delta_t;
63
64     if (mode == ecmANGULAR && ndim < 3)
65     {
66         gmx_fatal(FARGS, "Can not have angular comm removal with pbc=%s",
67                   c_pbcTypeNames[ir.pbcType].c_str());
68     }
69
70     if (mode != ecmNO)
71     {
72         nr = groups.groups[SimulationAtomGroupType::MassCenterVelocityRemoval].size();
73         /* Allocate one extra for a possible rest group */
74         size = nr + 1;
75         /* We need vcm->nr+1 elements per thread, but to avoid cache
76          * invalidation we add 2 elements to get a 152 byte separation.
77          */
78         stride = nr + 3;
79         if (mode == ecmANGULAR)
80         {
81             snew(group_i, size);
82
83             group_j.resize(size);
84             group_x.resize(size);
85             group_w.resize(size);
86         }
87
88         group_name.resize(size);
89         group_p.resize(size);
90         group_v.resize(size);
91         group_mass.resize(size);
92         group_ndf.resize(size);
93         for (int g = 0; (g < nr); g++)
94         {
95             group_ndf[g] = ir.opts.nrdf[g];
96             group_name[g] =
97                     *groups.groupNames[groups.groups[SimulationAtomGroupType::MassCenterVelocityRemoval][g]];
98         }
99
100         thread_vcm.resize(gmx_omp_nthreads_get(emntDefault) * stride);
101     }
102
103     nFreeze = ir.opts.nFreeze;
104 }
105
106 t_vcm::~t_vcm()
107 {
108     if (mode == ecmANGULAR)
109     {
110         sfree(group_i);
111     }
112 }
113
114 void reportComRemovalInfo(FILE* fp, const t_vcm& vcm)
115 {
116
117     /* Copy pointer to group names and print it. */
118     if (fp && vcm.mode != ecmNO)
119     {
120         fprintf(fp, "Center of mass motion removal mode is %s\n", ECOM(vcm.mode));
121         fprintf(fp,
122                 "We have the following groups for center of"
123                 " mass motion removal:\n");
124
125         for (int g = 0; (g < vcm.nr); g++)
126         {
127
128             fprintf(fp, "%3d:  %s\n", g, vcm.group_name[g]);
129         }
130     }
131 }
132
133
134 static void update_tensor(const rvec x, real m0, tensor I)
135 {
136     real xy, xz, yz;
137
138     /* Compute inertia tensor contribution due to this atom */
139     xy = x[XX] * x[YY] * m0;
140     xz = x[XX] * x[ZZ] * m0;
141     yz = x[YY] * x[ZZ] * m0;
142     I[XX][XX] += x[XX] * x[XX] * m0;
143     I[YY][YY] += x[YY] * x[YY] * m0;
144     I[ZZ][ZZ] += x[ZZ] * x[ZZ] * m0;
145     I[XX][YY] += xy;
146     I[YY][XX] += xy;
147     I[XX][ZZ] += xz;
148     I[ZZ][XX] += xz;
149     I[YY][ZZ] += yz;
150     I[ZZ][YY] += yz;
151 }
152
153 /* Center of mass code for groups */
154 void calc_vcm_grp(const t_mdatoms& md, const rvec x[], const rvec v[], t_vcm* vcm)
155 {
156     if (vcm->mode == ecmNO)
157     {
158         return;
159     }
160     int nthreads = gmx_omp_nthreads_get(emntDefault);
161
162     {
163 #pragma omp parallel num_threads(nthreads) default(none) shared(x, v, vcm, md)
164         {
165             int t = gmx_omp_get_thread_num();
166             for (int g = 0; g < vcm->size; g++)
167             {
168                 /* Reset linear momentum */
169                 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
170                 vcm_t->mass         = 0;
171                 clear_rvec(vcm_t->p);
172                 if (vcm->mode == ecmANGULAR)
173                 {
174                     /* Reset angular momentum */
175                     clear_rvec(vcm_t->j);
176                     clear_rvec(vcm_t->x);
177                     clear_mat(vcm_t->i);
178                 }
179             }
180
181 #pragma omp for schedule(static)
182             for (int i = 0; i < md.homenr; i++)
183             {
184                 int  g  = 0;
185                 real m0 = md.massT[i];
186                 if (md.cVCM)
187                 {
188                     g = md.cVCM[i];
189                 }
190                 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
191                 /* Calculate linear momentum */
192                 vcm_t->mass += m0;
193                 int m;
194                 for (m = 0; (m < DIM); m++)
195                 {
196                     vcm_t->p[m] += m0 * v[i][m];
197                 }
198
199                 if (vcm->mode == ecmANGULAR)
200                 {
201                     /* Calculate angular momentum */
202                     rvec j0;
203                     cprod(x[i], v[i], j0);
204
205                     for (m = 0; (m < DIM); m++)
206                     {
207                         vcm_t->j[m] += m0 * j0[m];
208                         vcm_t->x[m] += m0 * x[i][m];
209                     }
210                     /* Update inertia tensor */
211                     update_tensor(x[i], m0, vcm_t->i);
212                 }
213             }
214         }
215         for (int g = 0; g < vcm->size; g++)
216         {
217             /* Reset linear momentum */
218             vcm->group_mass[g] = 0;
219             clear_rvec(vcm->group_p[g]);
220             if (vcm->mode == ecmANGULAR)
221             {
222                 /* Reset angular momentum */
223                 clear_rvec(vcm->group_j[g]);
224                 clear_rvec(vcm->group_x[g]);
225                 clear_rvec(vcm->group_w[g]);
226                 clear_mat(vcm->group_i[g]);
227             }
228
229             for (int t = 0; t < nthreads; t++)
230             {
231                 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
232                 vcm->group_mass[g] += vcm_t->mass;
233                 rvec_inc(vcm->group_p[g], vcm_t->p);
234                 if (vcm->mode == ecmANGULAR)
235                 {
236                     rvec_inc(vcm->group_j[g], vcm_t->j);
237                     rvec_inc(vcm->group_x[g], vcm_t->x);
238                     m_add(vcm_t->i, vcm->group_i[g], vcm->group_i[g]);
239                 }
240             }
241         }
242     }
243 }
244
245 /*! \brief Remove the COM motion velocity from the velocities
246  *
247  * \note This routine should be called from within an OpenMP parallel region.
248  *
249  * \tparam numDimensions    Correct dimensions 0 to \p numDimensions-1
250  * \param[in]     mdatoms   The atom property and group information
251  * \param[in,out] v         The velocities to correct
252  * \param[in]     vcm       VCM data
253  */
254 template<int numDimensions>
255 static void doStopComMotionLinear(const t_mdatoms& mdatoms, rvec* v, const t_vcm& vcm)
256 {
257     const int             homenr   = mdatoms.homenr;
258     const unsigned short* group_id = mdatoms.cVCM;
259
260     if (mdatoms.cFREEZE != nullptr)
261     {
262         GMX_RELEASE_ASSERT(vcm.nFreeze != nullptr, "Need freeze dimension info with freeze groups");
263
264 #pragma omp for schedule(static)
265         for (int i = 0; i < homenr; i++)
266         {
267             unsigned short vcmGroup    = (group_id == nullptr ? 0 : group_id[i]);
268             unsigned short freezeGroup = mdatoms.cFREEZE[i];
269             for (int d = 0; d < numDimensions; d++)
270             {
271                 if (vcm.nFreeze[freezeGroup][d] == 0)
272                 {
273                     v[i][d] -= vcm.group_v[vcmGroup][d];
274                 }
275             }
276         }
277     }
278     else if (group_id == nullptr)
279     {
280 #pragma omp for schedule(static)
281         for (int i = 0; i < homenr; i++)
282         {
283             for (int d = 0; d < numDimensions; d++)
284             {
285                 v[i][d] -= vcm.group_v[0][d];
286             }
287         }
288     }
289     else
290     {
291 #pragma omp for schedule(static)
292         for (int i = 0; i < homenr; i++)
293         {
294             const int g = group_id[i];
295             for (int d = 0; d < numDimensions; d++)
296             {
297                 v[i][d] -= vcm.group_v[g][d];
298             }
299         }
300     }
301 }
302
303 /*! \brief Remove the COM motion velocity from the velocities, correct the coordinates assuming constant acceleration
304  *
305  * \note This routine should be called from within an OpenMP parallel region.
306  *
307  * \tparam numDimensions    Correct dimensions 0 to \p numDimensions-1
308  * \param[in]     homenr    The number of atoms to correct
309  * \param[in]     group_id  List of VCM group ids, when nullptr is passed all atoms are assumed to be in group 0
310  * \param[in,out] x         The coordinates to correct
311  * \param[in,out] v         The velocities to correct
312  * \param[in]     vcm       VCM data
313  */
314 template<int numDimensions>
315 static void doStopComMotionAccelerationCorrection(int                   homenr,
316                                                   const unsigned short* group_id,
317                                                   rvec* gmx_restrict x,
318                                                   rvec* gmx_restrict v,
319                                                   const t_vcm&       vcm)
320 {
321     const real xCorrectionFactor = 0.5 * vcm.timeStep;
322
323     if (group_id == nullptr)
324     {
325 #pragma omp for schedule(static)
326         for (int i = 0; i < homenr; i++)
327         {
328             for (int d = 0; d < numDimensions; d++)
329             {
330                 x[i][d] -= vcm.group_v[0][d] * xCorrectionFactor;
331                 v[i][d] -= vcm.group_v[0][d];
332             }
333         }
334     }
335     else
336     {
337 #pragma omp for schedule(static)
338         for (int i = 0; i < homenr; i++)
339         {
340             const int g = group_id[i];
341             for (int d = 0; d < numDimensions; d++)
342             {
343                 x[i][d] -= vcm.group_v[g][d] * xCorrectionFactor;
344                 v[i][d] -= vcm.group_v[g][d];
345             }
346         }
347     }
348 }
349
350 static void do_stopcm_grp(const t_mdatoms& mdatoms, rvec x[], rvec v[], const t_vcm& vcm)
351 {
352     if (vcm.mode == ecmNO)
353     {
354         return;
355     }
356     {
357         const int             homenr   = mdatoms.homenr;
358         const unsigned short* group_id = mdatoms.cVCM;
359
360         int gmx_unused nth = gmx_omp_nthreads_get(emntDefault);
361         // homenr could be shared, but gcc-8 & gcc-9 don't agree how to write that...
362         // https://www.gnu.org/software/gcc/gcc-9/porting_to.html -> OpenMP data sharing
363 #pragma omp parallel num_threads(nth) default(none) shared(x, v, vcm, group_id, mdatoms) \
364         firstprivate(homenr)
365         {
366             if (vcm.mode == ecmLINEAR || vcm.mode == ecmANGULAR
367                 || (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION && x == nullptr))
368             {
369                 /* Subtract linear momentum for v */
370                 switch (vcm.ndim)
371                 {
372                     case 1: doStopComMotionLinear<1>(mdatoms, v, vcm); break;
373                     case 2: doStopComMotionLinear<2>(mdatoms, v, vcm); break;
374                     case 3: doStopComMotionLinear<3>(mdatoms, v, vcm); break;
375                 }
376             }
377             else
378             {
379                 GMX_ASSERT(vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION,
380                            "When the mode is not linear or angular, it should be acceleration "
381                            "correction");
382                 /* Subtract linear momentum for v and x*/
383                 switch (vcm.ndim)
384                 {
385                     case 1:
386                         doStopComMotionAccelerationCorrection<1>(homenr, group_id, x, v, vcm);
387                         break;
388                     case 2:
389                         doStopComMotionAccelerationCorrection<2>(homenr, group_id, x, v, vcm);
390                         break;
391                     case 3:
392                         doStopComMotionAccelerationCorrection<3>(homenr, group_id, x, v, vcm);
393                         break;
394                 }
395             }
396             if (vcm.mode == ecmANGULAR)
397             {
398                 /* Subtract angular momentum */
399                 GMX_ASSERT(x, "Need x to compute angular momentum correction");
400
401                 int g = 0;
402 #pragma omp for schedule(static)
403                 for (int i = 0; i < homenr; i++)
404                 {
405                     if (group_id)
406                     {
407                         g = group_id[i];
408                     }
409                     /* Compute the correction to the velocity for each atom */
410                     rvec dv, dx;
411                     rvec_sub(x[i], vcm.group_x[g], dx);
412                     cprod(vcm.group_w[g], dx, dv);
413                     rvec_dec(v[i], dv);
414                 }
415             }
416         }
417     }
418 }
419
420 static void get_minv(tensor A, tensor B)
421 {
422     int    m, n;
423     double fac, rfac;
424     tensor tmp;
425
426     tmp[XX][XX] = A[YY][YY] + A[ZZ][ZZ];
427     tmp[YY][XX] = -A[XX][YY];
428     tmp[ZZ][XX] = -A[XX][ZZ];
429     tmp[XX][YY] = -A[XX][YY];
430     tmp[YY][YY] = A[XX][XX] + A[ZZ][ZZ];
431     tmp[ZZ][YY] = -A[YY][ZZ];
432     tmp[XX][ZZ] = -A[XX][ZZ];
433     tmp[YY][ZZ] = -A[YY][ZZ];
434     tmp[ZZ][ZZ] = A[XX][XX] + A[YY][YY];
435
436     /* This is a hack to prevent very large determinants */
437     rfac = (tmp[XX][XX] + tmp[YY][YY] + tmp[ZZ][ZZ]) / 3;
438     if (rfac == 0.0)
439     {
440         gmx_fatal(FARGS, "Can not stop center of mass: maybe 2dimensional system");
441     }
442     fac = 1.0 / rfac;
443     for (m = 0; (m < DIM); m++)
444     {
445         for (n = 0; (n < DIM); n++)
446         {
447             tmp[m][n] *= fac;
448         }
449     }
450     gmx::invertMatrix(tmp, B);
451     for (m = 0; (m < DIM); m++)
452     {
453         for (n = 0; (n < DIM); n++)
454         {
455             B[m][n] *= fac;
456         }
457     }
458 }
459
460 /* Processes VCM after reduction over ranks and prints warning with high VMC and fp != nullptr */
461 static void process_and_check_cm_grp(FILE* fp, t_vcm* vcm, real Temp_Max)
462 {
463     int    m, g;
464     real   ekcm, ekrot, tm, tm_1, Temp_cm;
465     rvec   jcm;
466     tensor Icm;
467
468     /* First analyse the total results */
469     if (vcm->mode != ecmNO)
470     {
471         for (g = 0; (g < vcm->nr); g++)
472         {
473             tm = vcm->group_mass[g];
474             if (tm != 0)
475             {
476                 tm_1 = 1.0 / tm;
477                 svmul(tm_1, vcm->group_p[g], vcm->group_v[g]);
478             }
479             /* Else it's zero anyway! */
480         }
481         if (vcm->mode == ecmANGULAR)
482         {
483             for (g = 0; (g < vcm->nr); g++)
484             {
485                 tm = vcm->group_mass[g];
486                 if (tm != 0)
487                 {
488                     tm_1 = 1.0 / tm;
489
490                     /* Compute center of mass for this group */
491                     for (m = 0; (m < DIM); m++)
492                     {
493                         vcm->group_x[g][m] *= tm_1;
494                     }
495
496                     /* Subtract the center of mass contribution to the
497                      * angular momentum
498                      */
499                     cprod(vcm->group_x[g], vcm->group_v[g], jcm);
500                     for (m = 0; (m < DIM); m++)
501                     {
502                         vcm->group_j[g][m] -= tm * jcm[m];
503                     }
504
505                     /* Subtract the center of mass contribution from the inertia
506                      * tensor (this is not as trivial as it seems, but due to
507                      * some cancellation we can still do it, even in parallel).
508                      */
509                     clear_mat(Icm);
510                     update_tensor(vcm->group_x[g], tm, Icm);
511                     m_sub(vcm->group_i[g], Icm, vcm->group_i[g]);
512
513                     /* Compute angular velocity, using matrix operation
514                      * Since J = I w
515                      * we have
516                      * w = I^-1 J
517                      */
518                     get_minv(vcm->group_i[g], Icm);
519                     mvmul(Icm, vcm->group_j[g], vcm->group_w[g]);
520                 }
521                 /* Else it's zero anyway! */
522             }
523         }
524     }
525     for (g = 0; (g < vcm->nr); g++)
526     {
527         ekcm = 0;
528         if (vcm->group_mass[g] != 0 && vcm->group_ndf[g] > 0)
529         {
530             for (m = 0; m < vcm->ndim; m++)
531             {
532                 ekcm += gmx::square(vcm->group_v[g][m]);
533             }
534             ekcm *= 0.5 * vcm->group_mass[g];
535             Temp_cm = 2 * ekcm / vcm->group_ndf[g];
536
537             if ((Temp_cm > Temp_Max) && fp)
538             {
539                 fprintf(fp, "Large VCM(group %s): %12.5f, %12.5f, %12.5f, Temp-cm: %12.5e\n",
540                         vcm->group_name[g], vcm->group_v[g][XX], vcm->group_v[g][YY],
541                         vcm->group_v[g][ZZ], Temp_cm);
542             }
543
544             if (vcm->mode == ecmANGULAR)
545             {
546                 ekrot = 0.5 * iprod(vcm->group_j[g], vcm->group_w[g]);
547                 // TODO: Change absolute energy comparison to relative
548                 if ((ekrot > 1) && fp && vcm->integratorConservesMomentum)
549                 {
550                     /* if we have an integrator that may not conserve momenta, skip */
551                     tm = vcm->group_mass[g];
552                     fprintf(fp, "Group %s with mass %12.5e, Ekrot %12.5e Det(I) = %12.5e\n",
553                             vcm->group_name[g], tm, ekrot, det(vcm->group_i[g]));
554                     fprintf(fp, "  COM: %12.5f  %12.5f  %12.5f\n", vcm->group_x[g][XX],
555                             vcm->group_x[g][YY], vcm->group_x[g][ZZ]);
556                     fprintf(fp, "  P:   %12.5f  %12.5f  %12.5f\n", vcm->group_p[g][XX],
557                             vcm->group_p[g][YY], vcm->group_p[g][ZZ]);
558                     fprintf(fp, "  V:   %12.5f  %12.5f  %12.5f\n", vcm->group_v[g][XX],
559                             vcm->group_v[g][YY], vcm->group_v[g][ZZ]);
560                     fprintf(fp, "  J:   %12.5f  %12.5f  %12.5f\n", vcm->group_j[g][XX],
561                             vcm->group_j[g][YY], vcm->group_j[g][ZZ]);
562                     fprintf(fp, "  w:   %12.5f  %12.5f  %12.5f\n", vcm->group_w[g][XX],
563                             vcm->group_w[g][YY], vcm->group_w[g][ZZ]);
564                     pr_rvecs(fp, 0, "Inertia tensor", vcm->group_i[g], DIM);
565                 }
566             }
567         }
568     }
569 }
570
571 void process_and_stopcm_grp(FILE* fplog, t_vcm* vcm, const t_mdatoms& mdatoms, rvec x[], rvec v[])
572 {
573     if (vcm->mode != ecmNO)
574     {
575         // TODO: Replace fixed temperature of 1 by a system value
576         process_and_check_cm_grp(fplog, vcm, 1);
577
578         do_stopcm_grp(mdatoms, x, v, *vcm);
579     }
580 }