0847e759a8cd2bea570a2ad25da4ca1b5dd06e37
[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,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "disre.h"
41
42 #include "config.h"
43
44 #include <cmath>
45 #include <cstdlib>
46 #include <cstring>
47
48 #include <algorithm>
49
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdrunutility/multisim.h"
54 #include "gromacs/mdtypes/commrec.h"
55 #include "gromacs/mdtypes/fcdata.h"
56 #include "gromacs/mdtypes/inputrec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/mdtypes/state.h"
59 #include "gromacs/pbcutil/ishift.h"
60 #include "gromacs/pbcutil/mshift.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                  const t_commrec*      cr,
74                  const gmx_multisim_t* ms,
75                  t_fcdata*             fcd,
76                  t_state*              state,
77                  gmx_bool              bIsREMD)
78 {
79     int                  fa, nmol, npair, np;
80     t_disresdata*        dd;
81     history_t*           hist;
82     gmx_mtop_ilistloop_t iloop;
83     char*                ptr;
84     int                  type_min, type_max;
85
86     dd = &(fcd->disres);
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 for PAR(cr),
120          * but there are probably also issues with e.g. NM MPI parallelization.
121          */
122         if (cr && PAR(cr))
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].size() > 0 && ir->eDisre != edrEnsemble)
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 == edrEnsemble ? 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 (cr && PAR(cr) && 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 (cr && 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 (MASTER(cr))
234         {
235             check_multi_int(fplog, ms, dd->nsystems, "the number of systems per ensemble", FALSE);
236         }
237         gmx_bcast_sim(sizeof(int), &dd->nsystems, cr);
238
239         /* We use to allow any value of nsystems which was a divisor
240          * of ms->nsim. But this required an extra communicator which
241          * was stored in t_fcdata. This pulled in mpi.h in nearly all C files.
242          */
243         if (!(ms->nsim == 1 || ms->nsim == 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, ms->nsim);
249         }
250         if (fplog)
251         {
252             fprintf(fplog, "Our ensemble consists of systems:");
253             for (int i = 0; i < dd->nsystems; i++)
254             {
255                 fprintf(fplog, " %d", (ms->sim / dd->nsystems) * dd->nsystems + i);
256             }
257             fprintf(fplog, "\n");
258         }
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,
280                     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 (cr && ms && dd->nsystems > 1 && MASTER(cr))
289         {
290             check_multi_int(fplog, ms, fcd->disres.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_fcdata*             fcd,
304                      history_t*            hist)
305 {
306     rvec          dx;
307     real *        rt, *rm3tav, *Rtl_6, *Rt_6, *Rtav_6;
308     t_disresdata* dd;
309     real          ETerm, ETerm1, cf1 = 0, cf2 = 0;
310     gmx_bool      bTav;
311
312     dd     = &(fcd->disres);
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_fcdata 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 (fcd->disres.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);
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                const t_graph*  g,
426                real gmx_unused lambda,
427                real gmx_unused* dvdlambda,
428                const t_mdatoms gmx_unused* md,
429                t_fcdata*                   fcd,
430                int gmx_unused* global_atom_index)
431 {
432     const real seven_three = 7.0 / 3.0;
433
434     rvec          dx;
435     real          weight_rt_1;
436     real          smooth_fc, Rt, Rtav, rt2, *Rtl_6, *Rt_6, *Rtav_6;
437     real          k0, f_scal = 0, fmax_scal, fk_scal, fij;
438     real          tav_viol, instant_viol, mixed_viol, violtot, vtot;
439     real          tav_viol_Rtav7, instant_viol_Rtav7;
440     real          up1, up2, low;
441     gmx_bool      bConservative, bMixed, bViolation;
442     ivec          dt;
443     t_disresdata* dd;
444     int           dr_weighting;
445     gmx_bool      dr_bMixed;
446
447     dd           = &(fcd->disres);
448     dr_weighting = dd->dr_weighting;
449     dr_bMixed    = dd->dr_bMixed;
450     Rtl_6        = dd->Rtl_6;
451     Rt_6         = dd->Rt_6;
452     Rtav_6       = dd->Rtav_6;
453
454     tav_viol = instant_viol = mixed_viol = tav_viol_Rtav7 = instant_viol_Rtav7 = 0;
455
456     smooth_fc = dd->dr_fc;
457     if (dd->dr_tau != 0)
458     {
459         /* scaling factor to smoothly turn on the restraint forces *
460          * when using time averaging                               */
461         smooth_fc *= (1.0 - dd->exp_min_t_tau);
462     }
463
464     violtot = 0;
465     vtot    = 0;
466
467     /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
468      * the total number of atoms pairs is nfa/3                          */
469     int faOffset = static_cast<int>(forceatoms - dd->forceatomsStart);
470     for (int fa = 0; fa < nfa; fa += 3)
471     {
472         int type  = forceatoms[fa];
473         int npair = ip[type].disres.npair;
474         up1       = ip[type].disres.up1;
475         up2       = ip[type].disres.up2;
476         low       = ip[type].disres.low;
477         k0        = smooth_fc * ip[type].disres.kfac;
478
479         int res = type - dd->type_min;
480
481         /* save some flops when there is only one pair */
482         if (ip[type].disres.type != 2)
483         {
484             bConservative = (dr_weighting == edrwConservative) && (npair > 1);
485             bMixed        = dr_bMixed;
486             Rt            = gmx::invsixthroot(Rt_6[res]);
487             Rtav          = gmx::invsixthroot(Rtav_6[res]);
488         }
489         else
490         {
491             /* When rtype=2 use instantaneous not ensemble averaged distance */
492             bConservative = (npair > 1);
493             bMixed        = FALSE;
494             Rt            = gmx::invsixthroot(Rtl_6[res]);
495             Rtav          = Rt;
496         }
497
498         if (Rtav > up1)
499         {
500             bViolation = TRUE;
501             tav_viol   = Rtav - up1;
502         }
503         else if (Rtav < low)
504         {
505             bViolation = TRUE;
506             tav_viol   = Rtav - low;
507         }
508         else
509         {
510             bViolation = FALSE;
511         }
512
513         if (bViolation)
514         {
515             /* Add 1/npair energy and violation for each of the npair pairs */
516             real pairFac = 1 / static_cast<real>(npair);
517
518             /* NOTE:
519              * there is no real potential when time averaging is applied
520              */
521             vtot += 0.5 * k0 * gmx::square(tav_viol) * pairFac;
522             if (!bMixed)
523             {
524                 f_scal = -k0 * tav_viol;
525                 violtot += fabs(tav_viol) * pairFac;
526             }
527             else
528             {
529                 if (Rt > up1)
530                 {
531                     if (tav_viol > 0)
532                     {
533                         instant_viol = Rt - up1;
534                     }
535                     else
536                     {
537                         bViolation = FALSE;
538                     }
539                 }
540                 else if (Rt < low)
541                 {
542                     if (tav_viol < 0)
543                     {
544                         instant_viol = Rt - low;
545                     }
546                     else
547                     {
548                         bViolation = FALSE;
549                     }
550                 }
551                 else
552                 {
553                     bViolation = FALSE;
554                 }
555                 if (bViolation)
556                 {
557                     mixed_viol = std::sqrt(tav_viol * instant_viol);
558                     f_scal     = -k0 * mixed_viol;
559                     violtot += mixed_viol * pairFac;
560                 }
561             }
562         }
563
564         if (bViolation)
565         {
566             fmax_scal = -k0 * (up2 - up1);
567             /* Correct the force for the number of restraints */
568             if (bConservative)
569             {
570                 f_scal = std::max(f_scal, fmax_scal);
571                 if (!bMixed)
572                 {
573                     f_scal *= Rtav / Rtav_6[res];
574                 }
575                 else
576                 {
577                     f_scal /= 2 * mixed_viol;
578                     tav_viol_Rtav7     = tav_viol * Rtav / Rtav_6[res];
579                     instant_viol_Rtav7 = instant_viol * Rt / Rt_6[res];
580                 }
581             }
582             else
583             {
584                 f_scal /= npair;
585                 f_scal = std::max(f_scal, fmax_scal);
586             }
587
588             /* Exert the force ... */
589
590             int pair = (faOffset + fa) / 3;
591             int ai   = forceatoms[fa + 1];
592             int aj   = forceatoms[fa + 2];
593             int ki   = CENTRAL;
594             if (pbc)
595             {
596                 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
597             }
598             else
599             {
600                 rvec_sub(x[ai], x[aj], dx);
601             }
602             rt2 = iprod(dx, dx);
603
604             weight_rt_1 = gmx::invsqrt(rt2);
605
606             if (bConservative)
607             {
608                 if (!dr_bMixed)
609                 {
610                     weight_rt_1 *= std::pow(dd->rm3tav[pair], seven_three);
611                 }
612                 else
613                 {
614                     weight_rt_1 *= tav_viol_Rtav7 * std::pow(dd->rm3tav[pair], seven_three)
615                                    + instant_viol_Rtav7 / (dd->rt[pair] * gmx::power6(dd->rt[pair]));
616                 }
617             }
618
619             fk_scal = f_scal * weight_rt_1;
620
621             if (g)
622             {
623                 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
624                 ki = IVEC2IS(dt);
625             }
626
627             for (int m = 0; m < DIM; m++)
628             {
629                 fij = fk_scal * dx[m];
630
631                 f[ai][m] += fij;
632                 f[aj][m] -= fij;
633                 if (fshift)
634                 {
635                     fshift[ki][m] += fij;
636                     fshift[CENTRAL][m] -= fij;
637                 }
638             }
639         }
640     }
641
642 #pragma omp atomic
643     dd->sumviol += violtot;
644
645     /* Return energy */
646     return vtot;
647 }
648
649 void update_disres_history(const t_fcdata* fcd, history_t* hist)
650 {
651     const t_disresdata* dd;
652     int                 pair;
653
654     dd = &(fcd->disres);
655     if (dd->dr_tau != 0)
656     {
657         /* Copy the new time averages that have been calculated
658          * in calc_disres_R_6.
659          */
660         hist->disre_initf = dd->exp_min_t_tau;
661         for (pair = 0; pair < dd->npair; pair++)
662         {
663             hist->disre_rm3tav[pair] = dd->rm3tav[pair];
664         }
665     }
666 }