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