Refactor md_enums
[alexxy/gromacs.git] / src / gromacs / listed_forces / disre.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,2021, 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 "disre.h"
42
43 #include "config.h"
44
45 #include <cmath>
46 #include <cstdlib>
47 #include <cstring>
48
49 #include <algorithm>
50
51 #include "gromacs/gmxlib/network.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdrunutility/multisim.h"
55 #include "gromacs/mdtypes/commrec.h"
56 #include "gromacs/mdtypes/fcdata.h"
57 #include "gromacs/mdtypes/inputrec.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/mdtypes/state.h"
60 #include "gromacs/pbcutil/ishift.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/topology/topology.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/futil.h"
66 #include "gromacs/utility/gmxassert.h"
67 #include "gromacs/utility/pleasecite.h"
68 #include "gromacs/utility/smalloc.h"
69
70 void init_disres(FILE*                 fplog,
71                  const gmx_mtop_t*     mtop,
72                  t_inputrec*           ir,
73                  DisResRunMode         disResRunMode,
74                  DDRole                ddRole,
75                  NumRanks              numRanks,
76                  MPI_Comm              communicator,
77                  const gmx_multisim_t* ms,
78                  t_disresdata*         dd,
79                  t_state*              state,
80                  gmx_bool              bIsREMD)
81 {
82     int                  fa, nmol, npair, np;
83     history_t*           hist;
84     gmx_mtop_ilistloop_t iloop;
85     char*                ptr;
86     int                  type_min, type_max;
87
88     if (gmx_mtop_ftype_count(mtop, F_DISRES) == 0)
89     {
90         dd->nres = 0;
91
92         return;
93     }
94
95     if (fplog)
96     {
97         fprintf(fplog, "Initializing the distance restraints\n");
98     }
99
100     dd->dr_weighting = ir->eDisreWeighting;
101     dd->dr_fc        = ir->dr_fc;
102     if (EI_DYNAMICS(ir->eI))
103     {
104         dd->dr_tau = ir->dr_tau;
105     }
106     else
107     {
108         dd->dr_tau = 0.0;
109     }
110     if (dd->dr_tau == 0.0)
111     {
112         dd->dr_bMixed = FALSE;
113         dd->ETerm     = 0.0;
114     }
115     else
116     {
117         /* We store the r^-6 time averages in an array that is indexed
118          * with the local disres iatom index, so this doesn't work with DD.
119          * Note that DD is not initialized yet here, so we check that we are on multiple ranks,
120          * but there are probably also issues with e.g. NM MPI parallelization.
121          */
122         if ((disResRunMode == DisResRunMode::MDRun) && (numRanks == NumRanks::Multiple))
123         {
124             gmx_fatal(FARGS,
125                       "Time-averaged distance restraints are not supported with MPI "
126                       "parallelization. You can use OpenMP parallelization on a single node.");
127         }
128
129         dd->dr_bMixed = ir->bDisreMixed;
130         dd->ETerm     = std::exp(-(ir->delta_t / ir->dr_tau));
131     }
132     dd->ETerm1 = 1.0 - dd->ETerm;
133
134     dd->nres  = 0;
135     dd->npair = 0;
136     type_min  = INT_MAX;
137     type_max  = 0;
138     iloop     = gmx_mtop_ilistloop_init(mtop);
139     while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
140     {
141         if (nmol > 1 && !(*il)[F_DISRES].empty() && ir->eDisre != DistanceRestraintRefinement::Ensemble)
142         {
143             gmx_fatal(FARGS,
144                       "NMR distance restraints with multiple copies of the same molecule are "
145                       "currently only supported with ensemble averaging. If you just want to "
146                       "restrain distances between atom pairs using a flat-bottomed potential, use "
147                       "a restraint potential (bonds type 10) instead.");
148         }
149
150         np = 0;
151         for (fa = 0; fa < (*il)[F_DISRES].size(); fa += 3)
152         {
153             int type;
154
155             type = (*il)[F_DISRES].iatoms[fa];
156
157             np++;
158             npair = mtop->ffparams.iparams[type].disres.npair;
159             if (np == npair)
160             {
161                 dd->nres += (ir->eDisre == DistanceRestraintRefinement::Ensemble ? 1 : nmol);
162                 dd->npair += nmol * npair;
163                 np = 0;
164
165                 type_min = std::min(type_min, type);
166                 type_max = std::max(type_max, type);
167             }
168         }
169     }
170
171     if ((disResRunMode == DisResRunMode::MDRun) && (numRanks == NumRanks::Multiple) && ir->nstdisreout > 0)
172     {
173         /* With DD we currently only have local pair information available */
174         gmx_fatal(FARGS,
175                   "With MPI parallelization distance-restraint pair output is not supported. Use "
176                   "nstdisreout=0 or use OpenMP parallelization on a single node.");
177     }
178
179     /* For communicating and/or reducing (sums of) r^-6 for pairs over threads
180      * we use multiple arrays in t_disresdata. We need to have unique indices
181      * for each restraint that work over threads and MPI ranks. To this end
182      * we use the type index. These should all be in one block and there should
183      * be dd->nres types, but we check for this here.
184      * This setup currently does not allow for multiple copies of the same
185      * molecule without ensemble averaging, this is check for above.
186      */
187     GMX_RELEASE_ASSERT(
188             type_max - type_min + 1 == dd->nres,
189             "All distance restraint parameter entries in the topology should be consecutive");
190
191     dd->type_min = type_min;
192
193     snew(dd->rt, dd->npair);
194
195     if (dd->dr_tau != 0.0)
196     {
197         GMX_RELEASE_ASSERT(state != nullptr,
198                            "We need a valid state when using time-averaged distance restraints");
199
200         hist = &state->hist;
201         /* Set the "history lack" factor to 1 */
202         state->flags |= (1 << estDISRE_INITF);
203         hist->disre_initf = 1.0;
204         /* Allocate space for the r^-3 time averages */
205         state->flags |= (1 << estDISRE_RM3TAV);
206         hist->ndisrepairs = dd->npair;
207         snew(hist->disre_rm3tav, hist->ndisrepairs);
208     }
209     /* Allocate space for a copy of rm3tav,
210      * so we can call do_force without modifying the state.
211      */
212     snew(dd->rm3tav, dd->npair);
213
214     /* Allocate Rt_6 and Rtav_6 consecutively in memory so they can be
215      * averaged over the processors in one call (in calc_disre_R_6)
216      */
217     snew(dd->Rt_6, 2 * dd->nres);
218     dd->Rtav_6 = &(dd->Rt_6[dd->nres]);
219
220     ptr = getenv("GMX_DISRE_ENSEMBLE_SIZE");
221     if ((disResRunMode == DisResRunMode::MDRun) && ms != nullptr && ptr != nullptr && !bIsREMD)
222     {
223 #if GMX_MPI
224         dd->nsystems = 0;
225         sscanf(ptr, "%d", &dd->nsystems);
226         if (fplog)
227         {
228             fprintf(fplog, "Found GMX_DISRE_ENSEMBLE_SIZE set to %d systems per ensemble\n", dd->nsystems);
229         }
230         /* This check is only valid on MASTER(cr), so probably
231          * ensemble-averaged distance restraints are broken on more
232          * than one processor per simulation system. */
233         if (ddRole == DDRole::Master)
234         {
235             check_multi_int(fplog, ms, dd->nsystems, "the number of systems per ensemble", FALSE);
236         }
237         gmx_bcast(sizeof(int), &dd->nsystems, communicator);
238
239         /* We use to allow any value of nsystems which was a divisor
240          * of ms->numSimulations_. But this required an extra communicator which
241          * pulled in mpi.h in nearly all C files.
242          */
243         if (!(ms->numSimulations_ == 1 || ms->numSimulations_ == dd->nsystems))
244         {
245             gmx_fatal(FARGS,
246                       "GMX_DISRE_ENSEMBLE_SIZE (%d) is not equal to 1 or the number of systems "
247                       "(option -multidir) %d",
248                       dd->nsystems,
249                       ms->numSimulations_);
250         }
251         if (fplog)
252         {
253             fprintf(fplog, "Our ensemble consists of systems:");
254             for (int i = 0; i < dd->nsystems; i++)
255             {
256                 fprintf(fplog, " %d", (ms->simulationIndex_ / dd->nsystems) * dd->nsystems + i);
257             }
258             fprintf(fplog, "\n");
259         }
260 #else
261         GMX_UNUSED_VALUE(communicator);
262 #endif
263     }
264     else
265     {
266         dd->nsystems = 1;
267     }
268
269     if (dd->nsystems == 1)
270     {
271         dd->Rtl_6 = dd->Rt_6;
272     }
273     else
274     {
275         snew(dd->Rtl_6, dd->nres);
276     }
277
278     if (dd->npair > 0)
279     {
280         if (fplog)
281         {
282             fprintf(fplog, "There are %d distance restraints involving %d atom pairs\n", dd->nres, dd->npair);
283         }
284         /* Have to avoid g_disre de-referencing cr blindly, mdrun not
285          * doing consistency checks for ensemble-averaged distance
286          * restraints when that's not happening, and only doing those
287          * checks from appropriate processes (since check_multi_int is
288          * too broken to check whether the communication will
289          * succeed...) */
290         if ((disResRunMode == DisResRunMode::MDRun) && ms && dd->nsystems > 1 && (ddRole == DDRole::Master))
291         {
292             check_multi_int(fplog, ms, dd->nres, "the number of distance restraints", FALSE);
293         }
294         please_cite(fplog, "Tropp80a");
295         please_cite(fplog, "Torda89a");
296     }
297 }
298
299 void calc_disres_R_6(const t_commrec*      cr,
300                      const gmx_multisim_t* ms,
301                      int                   nfa,
302                      const t_iatom         forceatoms[],
303                      const rvec            x[],
304                      const t_pbc*          pbc,
305                      t_disresdata*         dd,
306                      const history_t*      hist)
307 {
308     rvec     dx;
309     real *   rt, *rm3tav, *Rtl_6, *Rt_6, *Rtav_6;
310     real     ETerm, ETerm1, cf1 = 0, cf2 = 0;
311     gmx_bool bTav;
312
313     bTav   = (dd->dr_tau != 0);
314     ETerm  = dd->ETerm;
315     ETerm1 = dd->ETerm1;
316     rt     = dd->rt;
317     rm3tav = dd->rm3tav;
318     Rtl_6  = dd->Rtl_6;
319     Rt_6   = dd->Rt_6;
320     Rtav_6 = dd->Rtav_6;
321
322     if (bTav)
323     {
324         /* scaling factor to smoothly turn on the restraint forces *
325          * when using time averaging                               */
326         dd->exp_min_t_tau = hist->disre_initf * ETerm;
327
328         cf1 = dd->exp_min_t_tau;
329         cf2 = 1.0 / (1.0 - dd->exp_min_t_tau);
330     }
331
332     for (int res = 0; res < dd->nres; res++)
333     {
334         Rtav_6[res] = 0.0;
335         Rt_6[res]   = 0.0;
336     }
337
338     /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
339      * the total number of atoms pairs is nfa/3                          */
340     for (int fa = 0; fa < nfa; fa += 3)
341     {
342         int type = forceatoms[fa];
343         int res  = type - dd->type_min;
344         int pair = fa / 3;
345         int ai   = forceatoms[fa + 1];
346         int aj   = forceatoms[fa + 2];
347
348         if (pbc)
349         {
350             pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
351         }
352         else
353         {
354             rvec_sub(x[ai], x[aj], dx);
355         }
356         real rt2  = iprod(dx, dx);
357         real rt_1 = gmx::invsqrt(rt2);
358         real rt_3 = rt_1 * rt_1 * rt_1;
359
360         rt[pair] = rt2 * rt_1;
361         if (bTav)
362         {
363             /* Here we update rm3tav in t_disresdata using the data
364              * in history_t.
365              * Thus the results stay correct when this routine
366              * is called multiple times.
367              */
368             rm3tav[pair] = cf2 * ((ETerm - cf1) * hist->disre_rm3tav[pair] + ETerm1 * rt_3);
369         }
370         else
371         {
372             rm3tav[pair] = rt_3;
373         }
374
375         /* Currently divide_bondeds_over_threads() ensures that pairs within
376          * the same restraint get assigned to the same thread, so we could
377          * run this loop thread-parallel.
378          */
379         Rt_6[res] += rt_3 * rt_3;
380         Rtav_6[res] += rm3tav[pair] * rm3tav[pair];
381     }
382
383     /* NOTE: Rt_6 and Rtav_6 are stored consecutively in memory */
384     if (cr && DOMAINDECOMP(cr))
385     {
386         gmx_sum(2 * dd->nres, dd->Rt_6, cr);
387     }
388
389     if (dd->nsystems > 1)
390     {
391         real invn = 1.0 / dd->nsystems;
392
393         for (int res = 0; res < dd->nres; res++)
394         {
395             Rtl_6[res] = Rt_6[res];
396             Rt_6[res] *= invn;
397             Rtav_6[res] *= invn;
398         }
399
400         GMX_ASSERT(cr != nullptr && ms != nullptr, "We need multisim with nsystems>1");
401         gmx_sum_sim(2 * dd->nres, dd->Rt_6, ms);
402
403         if (DOMAINDECOMP(cr))
404         {
405             gmx_bcast(2 * dd->nres, dd->Rt_6, cr->mpi_comm_mygroup);
406         }
407     }
408
409     /* Store the base forceatoms pointer, so we can re-calculate the pair
410      * index in ta_disres() for indexing pair data in t_disresdata when
411      * using thread parallelization.
412      */
413     dd->forceatomsStart = forceatoms;
414
415     dd->sumviol = 0;
416 }
417
418 real ta_disres(int              nfa,
419                const t_iatom*   forceatoms,
420                const t_iparams* ip,
421                const rvec*      x,
422                rvec4*           f,
423                rvec*            fshift,
424                const t_pbc*     pbc,
425                real gmx_unused lambda,
426                real gmx_unused* dvdlambda,
427                const t_mdatoms gmx_unused* md,
428                t_fcdata gmx_unused* fcd,
429                t_disresdata*        disresdata,
430                t_oriresdata gmx_unused* oriresdata,
431                int gmx_unused* global_atom_index)
432 {
433     const real seven_three = 7.0 / 3.0;
434
435     rvec     dx;
436     real     weight_rt_1;
437     real     smooth_fc, Rt, Rtav, rt2, *Rtl_6, *Rt_6, *Rtav_6;
438     real     k0, f_scal = 0, fmax_scal, fk_scal, fij;
439     real     tav_viol, instant_viol, mixed_viol, violtot, vtot;
440     real     tav_viol_Rtav7, instant_viol_Rtav7;
441     real     up1, up2, low;
442     gmx_bool bConservative, bMixed, bViolation;
443     gmx_bool dr_bMixed;
444
445     DistanceRestraintWeighting dr_weighting = disresdata->dr_weighting;
446     dr_bMixed                               = disresdata->dr_bMixed;
447     Rtl_6                                   = disresdata->Rtl_6;
448     Rt_6                                    = disresdata->Rt_6;
449     Rtav_6                                  = disresdata->Rtav_6;
450
451     tav_viol = instant_viol = mixed_viol = tav_viol_Rtav7 = instant_viol_Rtav7 = 0;
452
453     smooth_fc = disresdata->dr_fc;
454     if (disresdata->dr_tau != 0)
455     {
456         /* scaling factor to smoothly turn on the restraint forces *
457          * when using time averaging                               */
458         smooth_fc *= (1.0 - disresdata->exp_min_t_tau);
459     }
460
461     violtot = 0;
462     vtot    = 0;
463
464     /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
465      * the total number of atoms pairs is nfa/3                          */
466     int faOffset = static_cast<int>(forceatoms - disresdata->forceatomsStart);
467     for (int fa = 0; fa < nfa; fa += 3)
468     {
469         int type  = forceatoms[fa];
470         int npair = ip[type].disres.npair;
471         up1       = ip[type].disres.up1;
472         up2       = ip[type].disres.up2;
473         low       = ip[type].disres.low;
474         k0        = smooth_fc * ip[type].disres.kfac;
475
476         int res = type - disresdata->type_min;
477
478         /* save some flops when there is only one pair */
479         if (ip[type].disres.type != 2)
480         {
481             bConservative = (dr_weighting == DistanceRestraintWeighting::Conservative) && (npair > 1);
482             bMixed        = dr_bMixed;
483             Rt            = gmx::invsixthroot(Rt_6[res]);
484             Rtav          = gmx::invsixthroot(Rtav_6[res]);
485         }
486         else
487         {
488             /* When rtype=2 use instantaneous not ensemble averaged distance */
489             bConservative = (npair > 1);
490             bMixed        = FALSE;
491             Rt            = gmx::invsixthroot(Rtl_6[res]);
492             Rtav          = Rt;
493         }
494
495         if (Rtav > up1)
496         {
497             bViolation = TRUE;
498             tav_viol   = Rtav - up1;
499         }
500         else if (Rtav < low)
501         {
502             bViolation = TRUE;
503             tav_viol   = Rtav - low;
504         }
505         else
506         {
507             bViolation = FALSE;
508         }
509
510         if (bViolation)
511         {
512             /* Add 1/npair energy and violation for each of the npair pairs */
513             real pairFac = 1 / static_cast<real>(npair);
514
515             /* NOTE:
516              * there is no real potential when time averaging is applied
517              */
518             vtot += 0.5 * k0 * gmx::square(tav_viol) * pairFac;
519             if (!bMixed)
520             {
521                 f_scal = -k0 * tav_viol;
522                 violtot += fabs(tav_viol) * pairFac;
523             }
524             else
525             {
526                 if (Rt > up1)
527                 {
528                     if (tav_viol > 0)
529                     {
530                         instant_viol = Rt - up1;
531                     }
532                     else
533                     {
534                         bViolation = FALSE;
535                     }
536                 }
537                 else if (Rt < low)
538                 {
539                     if (tav_viol < 0)
540                     {
541                         instant_viol = Rt - low;
542                     }
543                     else
544                     {
545                         bViolation = FALSE;
546                     }
547                 }
548                 else
549                 {
550                     bViolation = FALSE;
551                 }
552                 if (bViolation)
553                 {
554                     mixed_viol = std::sqrt(tav_viol * instant_viol);
555                     f_scal     = -k0 * mixed_viol;
556                     violtot += mixed_viol * pairFac;
557                 }
558             }
559         }
560
561         if (bViolation)
562         {
563             fmax_scal = -k0 * (up2 - up1);
564             /* Correct the force for the number of restraints */
565             if (bConservative)
566             {
567                 f_scal = std::max(f_scal, fmax_scal);
568                 if (!bMixed)
569                 {
570                     f_scal *= Rtav / Rtav_6[res];
571                 }
572                 else
573                 {
574                     f_scal /= 2 * mixed_viol;
575                     tav_viol_Rtav7     = tav_viol * Rtav / Rtav_6[res];
576                     instant_viol_Rtav7 = instant_viol * Rt / Rt_6[res];
577                 }
578             }
579             else
580             {
581                 f_scal /= npair;
582                 f_scal = std::max(f_scal, fmax_scal);
583             }
584
585             /* Exert the force ... */
586
587             int pair = (faOffset + fa) / 3;
588             int ai   = forceatoms[fa + 1];
589             int aj   = forceatoms[fa + 2];
590             int ki   = CENTRAL;
591             if (pbc)
592             {
593                 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
594             }
595             else
596             {
597                 rvec_sub(x[ai], x[aj], dx);
598             }
599             rt2 = iprod(dx, dx);
600
601             weight_rt_1 = gmx::invsqrt(rt2);
602
603             if (bConservative)
604             {
605                 if (!dr_bMixed)
606                 {
607                     weight_rt_1 *= std::pow(disresdata->rm3tav[pair], seven_three);
608                 }
609                 else
610                 {
611                     weight_rt_1 *=
612                             tav_viol_Rtav7 * std::pow(disresdata->rm3tav[pair], seven_three)
613                             + instant_viol_Rtav7
614                                       / (disresdata->rt[pair] * gmx::power6(disresdata->rt[pair]));
615                 }
616             }
617
618             fk_scal = f_scal * weight_rt_1;
619
620             for (int m = 0; m < DIM; m++)
621             {
622                 fij = fk_scal * dx[m];
623
624                 f[ai][m] += fij;
625                 f[aj][m] -= fij;
626                 if (fshift)
627                 {
628                     fshift[ki][m] += fij;
629                     fshift[CENTRAL][m] -= fij;
630                 }
631             }
632         }
633     }
634
635 #pragma omp atomic
636     disresdata->sumviol += violtot;
637
638     /* Return energy */
639     return vtot;
640 }
641
642 void update_disres_history(const t_disresdata& dd, history_t* hist)
643 {
644     if (dd.dr_tau != 0)
645     {
646         /* Copy the new time averages that have been calculated
647          * in calc_disres_R_6.
648          */
649         hist->disre_initf = dd.exp_min_t_tau;
650         for (int pair = 0; pair < dd.npair; pair++)
651         {
652             hist->disre_rm3tav[pair] = dd.rm3tav[pair];
653         }
654     }
655 }