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