a2d998d8ae6655d8098625442026f0f913da38b1
[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/gmxlib/network.h"
44 #include "gromacs/math/functions.h"
45 #include "gromacs/math/invertmatrix.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/math/vecdump.h"
48 #include "gromacs/mdlib/gmx_omp_nthreads.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/gmxomp.h"
56 #include "gromacs/utility/smalloc.h"
57
58 t_vcm::t_vcm(const SimulationGroups& groups, const t_inputrec& ir) :
59     integratorConservesMomentum(!EI_RANDOM(ir.eI))
60 {
61     mode     = (ir.nstcomm > 0) ? ir.comm_mode : ecmNO;
62     ndim     = ndof_com(&ir);
63     timeStep = ir.nstcomm * ir.delta_t;
64
65     if (mode == ecmANGULAR && ndim < 3)
66     {
67         gmx_fatal(FARGS, "Can not have angular comm removal with pbc=%s", epbc_names[ir.ePBC]);
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     int nthreads = gmx_omp_nthreads_get(emntDefault);
157     if (vcm->mode != ecmNO)
158     {
159 #pragma omp parallel num_threads(nthreads)
160         {
161             int t = gmx_omp_get_thread_num();
162             for (int g = 0; g < vcm->size; g++)
163             {
164                 /* Reset linear momentum */
165                 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
166                 vcm_t->mass         = 0;
167                 clear_rvec(vcm_t->p);
168                 if (vcm->mode == ecmANGULAR)
169                 {
170                     /* Reset angular momentum */
171                     clear_rvec(vcm_t->j);
172                     clear_rvec(vcm_t->x);
173                     clear_mat(vcm_t->i);
174                 }
175             }
176
177 #pragma omp for schedule(static)
178             for (int i = 0; i < md.homenr; i++)
179             {
180                 int  g  = 0;
181                 real m0 = md.massT[i];
182                 if (md.cVCM)
183                 {
184                     g = md.cVCM[i];
185                 }
186                 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
187                 /* Calculate linear momentum */
188                 vcm_t->mass += m0;
189                 int m;
190                 for (m = 0; (m < DIM); m++)
191                 {
192                     vcm_t->p[m] += m0 * v[i][m];
193                 }
194
195                 if (vcm->mode == ecmANGULAR)
196                 {
197                     /* Calculate angular momentum */
198                     rvec j0;
199                     cprod(x[i], v[i], j0);
200
201                     for (m = 0; (m < DIM); m++)
202                     {
203                         vcm_t->j[m] += m0 * j0[m];
204                         vcm_t->x[m] += m0 * x[i][m];
205                     }
206                     /* Update inertia tensor */
207                     update_tensor(x[i], m0, vcm_t->i);
208                 }
209             }
210         }
211         for (int g = 0; g < vcm->size; g++)
212         {
213             /* Reset linear momentum */
214             vcm->group_mass[g] = 0;
215             clear_rvec(vcm->group_p[g]);
216             if (vcm->mode == ecmANGULAR)
217             {
218                 /* Reset angular momentum */
219                 clear_rvec(vcm->group_j[g]);
220                 clear_rvec(vcm->group_x[g]);
221                 clear_rvec(vcm->group_w[g]);
222                 clear_mat(vcm->group_i[g]);
223             }
224
225             for (int t = 0; t < nthreads; t++)
226             {
227                 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
228                 vcm->group_mass[g] += vcm_t->mass;
229                 rvec_inc(vcm->group_p[g], vcm_t->p);
230                 if (vcm->mode == ecmANGULAR)
231                 {
232                     rvec_inc(vcm->group_j[g], vcm_t->j);
233                     rvec_inc(vcm->group_x[g], vcm_t->x);
234                     m_add(vcm_t->i, vcm->group_i[g], vcm->group_i[g]);
235                 }
236             }
237         }
238     }
239 }
240
241 /*! \brief Remove the COM motion velocity from the velocities
242  *
243  * \note This routine should be called from within an OpenMP parallel region.
244  *
245  * \tparam numDimensions    Correct dimensions 0 to \p numDimensions-1
246  * \param[in]     mdatoms   The atom property and group information
247  * \param[in,out] v         The velocities to correct
248  * \param[in]     vcm       VCM data
249  */
250 template<int numDimensions>
251 static void doStopComMotionLinear(const t_mdatoms& mdatoms, rvec* v, const t_vcm& vcm)
252 {
253     const int             homenr   = mdatoms.homenr;
254     const unsigned short* group_id = mdatoms.cVCM;
255
256     if (mdatoms.cFREEZE != nullptr)
257     {
258         GMX_RELEASE_ASSERT(vcm.nFreeze != nullptr, "Need freeze dimension info with freeze groups");
259
260 #pragma omp for schedule(static)
261         for (int i = 0; i < homenr; i++)
262         {
263             unsigned short vcmGroup    = (group_id == nullptr ? 0 : group_id[i]);
264             unsigned short freezeGroup = mdatoms.cFREEZE[i];
265             for (int d = 0; d < numDimensions; d++)
266             {
267                 if (vcm.nFreeze[freezeGroup][d] == 0)
268                 {
269                     v[i][d] -= vcm.group_v[vcmGroup][d];
270                 }
271             }
272         }
273     }
274     else if (group_id == nullptr)
275     {
276 #pragma omp for schedule(static)
277         for (int i = 0; i < homenr; i++)
278         {
279             for (int d = 0; d < numDimensions; d++)
280             {
281                 v[i][d] -= vcm.group_v[0][d];
282             }
283         }
284     }
285     else
286     {
287 #pragma omp for schedule(static)
288         for (int i = 0; i < homenr; i++)
289         {
290             const int g = group_id[i];
291             for (int d = 0; d < numDimensions; d++)
292             {
293                 v[i][d] -= vcm.group_v[g][d];
294             }
295         }
296     }
297 }
298
299 /*! \brief Remove the COM motion velocity from the velocities, correct the coordinates assuming constant acceleration
300  *
301  * \note This routine should be called from within an OpenMP parallel region.
302  *
303  * \tparam numDimensions    Correct dimensions 0 to \p numDimensions-1
304  * \param[in]     homenr    The number of atoms to correct
305  * \param[in]     group_id  List of VCM group ids, when nullptr is passed all atoms are assumed to be in group 0
306  * \param[in,out] x         The coordinates to correct
307  * \param[in,out] v         The velocities to correct
308  * \param[in]     vcm       VCM data
309  */
310 template<int numDimensions>
311 static void doStopComMotionAccelerationCorrection(int                   homenr,
312                                                   const unsigned short* group_id,
313                                                   rvec* gmx_restrict x,
314                                                   rvec* gmx_restrict v,
315                                                   const t_vcm&       vcm)
316 {
317     const real xCorrectionFactor = 0.5 * vcm.timeStep;
318
319     if (group_id == nullptr)
320     {
321 #pragma omp for schedule(static)
322         for (int i = 0; i < homenr; i++)
323         {
324             for (int d = 0; d < numDimensions; d++)
325             {
326                 x[i][d] -= vcm.group_v[0][d] * xCorrectionFactor;
327                 v[i][d] -= vcm.group_v[0][d];
328             }
329         }
330     }
331     else
332     {
333 #pragma omp for schedule(static)
334         for (int i = 0; i < homenr; i++)
335         {
336             const int g = group_id[i];
337             for (int d = 0; d < numDimensions; d++)
338             {
339                 x[i][d] -= vcm.group_v[g][d] * xCorrectionFactor;
340                 v[i][d] -= vcm.group_v[g][d];
341             }
342         }
343     }
344 }
345
346 static void do_stopcm_grp(const t_mdatoms& mdatoms, rvec x[], rvec v[], const t_vcm& vcm)
347 {
348     if (vcm.mode != ecmNO)
349     {
350         const int             homenr   = mdatoms.homenr;
351         const unsigned short* group_id = mdatoms.cVCM;
352
353         int gmx_unused nth = gmx_omp_nthreads_get(emntDefault);
354 #pragma omp parallel num_threads(nth)
355         {
356             if (vcm.mode == ecmLINEAR || vcm.mode == ecmANGULAR
357                 || (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION && x == nullptr))
358             {
359                 /* Subtract linear momentum for v */
360                 switch (vcm.ndim)
361                 {
362                     case 1: doStopComMotionLinear<1>(mdatoms, v, vcm); break;
363                     case 2: doStopComMotionLinear<2>(mdatoms, v, vcm); break;
364                     case 3: doStopComMotionLinear<3>(mdatoms, v, vcm); break;
365                 }
366             }
367             else
368             {
369                 GMX_ASSERT(vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION,
370                            "When the mode is not linear or angular, it should be acceleration "
371                            "correction");
372                 /* Subtract linear momentum for v and x*/
373                 switch (vcm.ndim)
374                 {
375                     case 1:
376                         doStopComMotionAccelerationCorrection<1>(homenr, group_id, x, v, vcm);
377                         break;
378                     case 2:
379                         doStopComMotionAccelerationCorrection<2>(homenr, group_id, x, v, vcm);
380                         break;
381                     case 3:
382                         doStopComMotionAccelerationCorrection<3>(homenr, group_id, x, v, vcm);
383                         break;
384                 }
385             }
386             if (vcm.mode == ecmANGULAR)
387             {
388                 /* Subtract angular momentum */
389                 GMX_ASSERT(x, "Need x to compute angular momentum correction");
390
391                 int g = 0;
392 #pragma omp for schedule(static)
393                 for (int i = 0; i < homenr; i++)
394                 {
395                     if (group_id)
396                     {
397                         g = group_id[i];
398                     }
399                     /* Compute the correction to the velocity for each atom */
400                     rvec dv, dx;
401                     rvec_sub(x[i], vcm.group_x[g], dx);
402                     cprod(vcm.group_w[g], dx, dv);
403                     rvec_dec(v[i], dv);
404                 }
405             }
406         }
407     }
408 }
409
410 static void get_minv(tensor A, tensor B)
411 {
412     int    m, n;
413     double fac, rfac;
414     tensor tmp;
415
416     tmp[XX][XX] = A[YY][YY] + A[ZZ][ZZ];
417     tmp[YY][XX] = -A[XX][YY];
418     tmp[ZZ][XX] = -A[XX][ZZ];
419     tmp[XX][YY] = -A[XX][YY];
420     tmp[YY][YY] = A[XX][XX] + A[ZZ][ZZ];
421     tmp[ZZ][YY] = -A[YY][ZZ];
422     tmp[XX][ZZ] = -A[XX][ZZ];
423     tmp[YY][ZZ] = -A[YY][ZZ];
424     tmp[ZZ][ZZ] = A[XX][XX] + A[YY][YY];
425
426     /* This is a hack to prevent very large determinants */
427     rfac = (tmp[XX][XX] + tmp[YY][YY] + tmp[ZZ][ZZ]) / 3;
428     if (rfac == 0.0)
429     {
430         gmx_fatal(FARGS, "Can not stop center of mass: maybe 2dimensional system");
431     }
432     fac = 1.0 / rfac;
433     for (m = 0; (m < DIM); m++)
434     {
435         for (n = 0; (n < DIM); n++)
436         {
437             tmp[m][n] *= fac;
438         }
439     }
440     gmx::invertMatrix(tmp, B);
441     for (m = 0; (m < DIM); m++)
442     {
443         for (n = 0; (n < DIM); n++)
444         {
445             B[m][n] *= fac;
446         }
447     }
448 }
449
450 /* Processes VCM after reduction over ranks and prints warning with high VMC and fp != nullptr */
451 static void process_and_check_cm_grp(FILE* fp, t_vcm* vcm, real Temp_Max)
452 {
453     int    m, g;
454     real   ekcm, ekrot, tm, tm_1, Temp_cm;
455     rvec   jcm;
456     tensor Icm;
457
458     /* First analyse the total results */
459     if (vcm->mode != ecmNO)
460     {
461         for (g = 0; (g < vcm->nr); g++)
462         {
463             tm = vcm->group_mass[g];
464             if (tm != 0)
465             {
466                 tm_1 = 1.0 / tm;
467                 svmul(tm_1, vcm->group_p[g], vcm->group_v[g]);
468             }
469             /* Else it's zero anyway! */
470         }
471         if (vcm->mode == ecmANGULAR)
472         {
473             for (g = 0; (g < vcm->nr); g++)
474             {
475                 tm = vcm->group_mass[g];
476                 if (tm != 0)
477                 {
478                     tm_1 = 1.0 / tm;
479
480                     /* Compute center of mass for this group */
481                     for (m = 0; (m < DIM); m++)
482                     {
483                         vcm->group_x[g][m] *= tm_1;
484                     }
485
486                     /* Subtract the center of mass contribution to the
487                      * angular momentum
488                      */
489                     cprod(vcm->group_x[g], vcm->group_v[g], jcm);
490                     for (m = 0; (m < DIM); m++)
491                     {
492                         vcm->group_j[g][m] -= tm * jcm[m];
493                     }
494
495                     /* Subtract the center of mass contribution from the inertia
496                      * tensor (this is not as trivial as it seems, but due to
497                      * some cancellation we can still do it, even in parallel).
498                      */
499                     clear_mat(Icm);
500                     update_tensor(vcm->group_x[g], tm, Icm);
501                     m_sub(vcm->group_i[g], Icm, vcm->group_i[g]);
502
503                     /* Compute angular velocity, using matrix operation
504                      * Since J = I w
505                      * we have
506                      * w = I^-1 J
507                      */
508                     get_minv(vcm->group_i[g], Icm);
509                     mvmul(Icm, vcm->group_j[g], vcm->group_w[g]);
510                 }
511                 /* Else it's zero anyway! */
512             }
513         }
514     }
515     for (g = 0; (g < vcm->nr); g++)
516     {
517         ekcm = 0;
518         if (vcm->group_mass[g] != 0 && vcm->group_ndf[g] > 0)
519         {
520             for (m = 0; m < vcm->ndim; m++)
521             {
522                 ekcm += gmx::square(vcm->group_v[g][m]);
523             }
524             ekcm *= 0.5 * vcm->group_mass[g];
525             Temp_cm = 2 * ekcm / vcm->group_ndf[g];
526
527             if ((Temp_cm > Temp_Max) && fp)
528             {
529                 fprintf(fp, "Large VCM(group %s): %12.5f, %12.5f, %12.5f, Temp-cm: %12.5e\n",
530                         vcm->group_name[g], vcm->group_v[g][XX], vcm->group_v[g][YY],
531                         vcm->group_v[g][ZZ], Temp_cm);
532             }
533
534             if (vcm->mode == ecmANGULAR)
535             {
536                 ekrot = 0.5 * iprod(vcm->group_j[g], vcm->group_w[g]);
537                 // TODO: Change absolute energy comparison to relative
538                 if ((ekrot > 1) && fp && vcm->integratorConservesMomentum)
539                 {
540                     /* if we have an integrator that may not conserve momenta, skip */
541                     tm = vcm->group_mass[g];
542                     fprintf(fp, "Group %s with mass %12.5e, Ekrot %12.5e Det(I) = %12.5e\n",
543                             vcm->group_name[g], tm, ekrot, det(vcm->group_i[g]));
544                     fprintf(fp, "  COM: %12.5f  %12.5f  %12.5f\n", vcm->group_x[g][XX],
545                             vcm->group_x[g][YY], vcm->group_x[g][ZZ]);
546                     fprintf(fp, "  P:   %12.5f  %12.5f  %12.5f\n", vcm->group_p[g][XX],
547                             vcm->group_p[g][YY], vcm->group_p[g][ZZ]);
548                     fprintf(fp, "  V:   %12.5f  %12.5f  %12.5f\n", vcm->group_v[g][XX],
549                             vcm->group_v[g][YY], vcm->group_v[g][ZZ]);
550                     fprintf(fp, "  J:   %12.5f  %12.5f  %12.5f\n", vcm->group_j[g][XX],
551                             vcm->group_j[g][YY], vcm->group_j[g][ZZ]);
552                     fprintf(fp, "  w:   %12.5f  %12.5f  %12.5f\n", vcm->group_w[g][XX],
553                             vcm->group_w[g][YY], vcm->group_w[g][ZZ]);
554                     pr_rvecs(fp, 0, "Inertia tensor", vcm->group_i[g], DIM);
555                 }
556             }
557         }
558     }
559 }
560
561 void process_and_stopcm_grp(FILE* fplog, t_vcm* vcm, const t_mdatoms& mdatoms, rvec x[], rvec v[])
562 {
563     if (vcm->mode != ecmNO)
564     {
565         // TODO: Replace fixed temperature of 1 by a system value
566         process_and_check_cm_grp(fplog, vcm, 1);
567
568         do_stopcm_grp(mdatoms, x, v, *vcm);
569     }
570 }