Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / disre.c
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, 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 "config.h"
41
42 #include <math.h>
43 #include <stdlib.h>
44
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/legacyheaders/types/commrec.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/legacyheaders/bondf.h"
50 #include "gromacs/legacyheaders/copyrite.h"
51 #include "gromacs/legacyheaders/disre.h"
52 #include "gromacs/legacyheaders/main.h"
53 #include "gromacs/topology/mtop_util.h"
54
55 #include "gromacs/math/vec.h"
56 #include "gromacs/pbcutil/ishift.h"
57 #include "gromacs/pbcutil/mshift.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
61
62 void init_disres(FILE *fplog, const gmx_mtop_t *mtop,
63                  t_inputrec *ir, const t_commrec *cr,
64                  t_fcdata *fcd, t_state *state, gmx_bool bIsREMD)
65 {
66     int                  fa, nmol, i, npair, np;
67     t_iparams           *ip;
68     t_disresdata        *dd;
69     history_t           *hist;
70     gmx_mtop_ilistloop_t iloop;
71     t_ilist             *il;
72     char                *ptr;
73
74     dd = &(fcd->disres);
75
76     if (gmx_mtop_ftype_count(mtop, F_DISRES) == 0)
77     {
78         dd->nres = 0;
79
80         return;
81     }
82
83     if (fplog)
84     {
85         fprintf(fplog, "Initializing the distance restraints\n");
86     }
87
88
89     if (ir->eDisre == edrEnsemble)
90     {
91         gmx_fatal(FARGS, "Sorry, distance restraints with ensemble averaging over multiple molecules in one system are not functional in this version of GROMACS");
92     }
93
94     dd->dr_weighting = ir->eDisreWeighting;
95     dd->dr_fc        = ir->dr_fc;
96     if (EI_DYNAMICS(ir->eI))
97     {
98         dd->dr_tau   = ir->dr_tau;
99     }
100     else
101     {
102         dd->dr_tau   = 0.0;
103     }
104     if (dd->dr_tau == 0.0)
105     {
106         dd->dr_bMixed = FALSE;
107         dd->ETerm     = 0.0;
108     }
109     else
110     {
111         dd->dr_bMixed = ir->bDisreMixed;
112         dd->ETerm     = exp(-(ir->delta_t/ir->dr_tau));
113     }
114     dd->ETerm1        = 1.0 - dd->ETerm;
115
116     ip = mtop->ffparams.iparams;
117
118     dd->nres  = 0;
119     dd->npair = 0;
120     iloop     = gmx_mtop_ilistloop_init(mtop);
121     while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
122     {
123         np = 0;
124         for (fa = 0; fa < il[F_DISRES].nr; fa += 3)
125         {
126             np++;
127             npair = mtop->ffparams.iparams[il[F_DISRES].iatoms[fa]].disres.npair;
128             if (np == npair)
129             {
130                 dd->nres  += (ir->eDisre == edrEnsemble ? 1 : nmol)*npair;
131                 dd->npair += nmol*npair;
132                 np         = 0;
133             }
134         }
135     }
136
137     if (cr && PAR(cr))
138     {
139         /* Temporary check, will be removed when disre is implemented with DD */
140         const char *notestr = "NOTE: atoms involved in distance restraints should be within the same domain. If this is not the case mdrun generates a fatal error. If you encounter this, use a single MPI rank (Verlet+OpenMP+GPUs work fine).";
141
142         if (MASTER(cr))
143         {
144             fprintf(stderr, "\n%s\n\n", notestr);
145         }
146         if (fplog)
147         {
148             fprintf(fplog, "%s\n", notestr);
149         }
150
151         if (dd->dr_tau != 0 || ir->eDisre == edrEnsemble || cr->ms != NULL ||
152             dd->nres != dd->npair)
153         {
154             gmx_fatal(FARGS, "Time or ensemble averaged or multiple pair distance restraints do not work (yet) with domain decomposition, use a single MPI rank%s", cr->ms ? " per simulation" : "");
155         }
156         if (ir->nstdisreout != 0)
157         {
158             if (fplog)
159             {
160                 fprintf(fplog, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n\n");
161             }
162             if (MASTER(cr))
163             {
164                 fprintf(stderr, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n");
165             }
166             ir->nstdisreout = 0;
167         }
168     }
169
170     snew(dd->rt, dd->npair);
171
172     if (dd->dr_tau != 0.0)
173     {
174         hist = &state->hist;
175         /* Set the "history lack" factor to 1 */
176         state->flags     |= (1<<estDISRE_INITF);
177         hist->disre_initf = 1.0;
178         /* Allocate space for the r^-3 time averages */
179         state->flags     |= (1<<estDISRE_RM3TAV);
180         hist->ndisrepairs = dd->npair;
181         snew(hist->disre_rm3tav, hist->ndisrepairs);
182     }
183     /* Allocate space for a copy of rm3tav,
184      * so we can call do_force without modifying the state.
185      */
186     snew(dd->rm3tav, dd->npair);
187
188     /* Allocate Rt_6 and Rtav_6 consecutively in memory so they can be
189      * averaged over the processors in one call (in calc_disre_R_6)
190      */
191     snew(dd->Rt_6, 2*dd->nres);
192     dd->Rtav_6 = &(dd->Rt_6[dd->nres]);
193
194     ptr = getenv("GMX_DISRE_ENSEMBLE_SIZE");
195     if (cr && cr->ms != NULL && ptr != NULL && !bIsREMD)
196     {
197 #ifdef GMX_MPI
198         dd->nsystems = 0;
199         sscanf(ptr, "%d", &dd->nsystems);
200         if (fplog)
201         {
202             fprintf(fplog, "Found GMX_DISRE_ENSEMBLE_SIZE set to %d systems per ensemble\n", dd->nsystems);
203         }
204         /* This check is only valid on MASTER(cr), so probably
205          * ensemble-averaged distance restraints are broken on more
206          * than one processor per simulation system. */
207         if (MASTER(cr))
208         {
209             check_multi_int(fplog, cr->ms, dd->nsystems,
210                             "the number of systems per ensemble",
211                             FALSE);
212         }
213         gmx_bcast_sim(sizeof(int), &dd->nsystems, cr);
214
215         /* We use to allow any value of nsystems which was a divisor
216          * of ms->nsim. But this required an extra communicator which
217          * was stored in t_fcdata. This pulled in mpi.h in nearly all C files.
218          */
219         if (!(cr->ms->nsim == 1 || cr->ms->nsim == dd->nsystems))
220         {
221             gmx_fatal(FARGS, "GMX_DISRE_ENSEMBLE_SIZE (%d) is not equal to 1 or the number of systems (option -multi) %d", dd->nsystems, cr->ms->nsim);
222         }
223         if (fplog)
224         {
225             fprintf(fplog, "Our ensemble consists of systems:");
226             for (i = 0; i < dd->nsystems; i++)
227             {
228                 fprintf(fplog, " %d",
229                         (cr->ms->sim/dd->nsystems)*dd->nsystems+i);
230             }
231             fprintf(fplog, "\n");
232         }
233         snew(dd->Rtl_6, dd->nres);
234 #endif
235     }
236     else
237     {
238         dd->nsystems = 1;
239         dd->Rtl_6    = dd->Rt_6;
240     }
241
242     if (dd->npair > 0)
243     {
244         if (fplog)
245         {
246             fprintf(fplog, "There are %d distance restraints involving %d atom pairs\n", dd->nres, dd->npair);
247         }
248         /* Have to avoid g_disre de-referencing cr blindly, mdrun not
249          * doing consistency checks for ensemble-averaged distance
250          * restraints when that's not happening, and only doing those
251          * checks from appropriate processes (since check_multi_int is
252          * too broken to check whether the communication will
253          * succeed...) */
254         if (cr && cr->ms && dd->nsystems > 1 && MASTER(cr))
255         {
256             check_multi_int(fplog, cr->ms, fcd->disres.nres,
257                             "the number of distance restraints",
258                             FALSE);
259         }
260         please_cite(fplog, "Tropp80a");
261         please_cite(fplog, "Torda89a");
262     }
263 }
264
265 void calc_disres_R_6(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
266                      const rvec x[], const t_pbc *pbc,
267                      t_fcdata *fcd, history_t *hist)
268 {
269     atom_id         ai, aj;
270     int             fa, res, i, pair, ki, kj, m;
271     int             type, npair, np;
272     rvec            dx;
273     real           *rt, *rm3tav, *Rtl_6, *Rt_6, *Rtav_6;
274     real            rt_1, rt_3, rt2;
275     ivec            it, jt, dt;
276     t_disresdata   *dd;
277     real            ETerm, ETerm1, cf1 = 0, cf2 = 0, invn = 0;
278     gmx_bool        bTav;
279
280     dd           = &(fcd->disres);
281     bTav         = (dd->dr_tau != 0);
282     ETerm        = dd->ETerm;
283     ETerm1       = dd->ETerm1;
284     rt           = dd->rt;
285     rm3tav       = dd->rm3tav;
286     Rtl_6        = dd->Rtl_6;
287     Rt_6         = dd->Rt_6;
288     Rtav_6       = dd->Rtav_6;
289
290     if (bTav)
291     {
292         /* scaling factor to smoothly turn on the restraint forces *
293          * when using time averaging                               */
294         dd->exp_min_t_tau = hist->disre_initf*ETerm;
295
296         cf1 = dd->exp_min_t_tau;
297         cf2 = 1.0/(1.0 - dd->exp_min_t_tau);
298     }
299
300     if (dd->nsystems > 1)
301     {
302         invn = 1.0/dd->nsystems;
303     }
304
305     /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
306      * the total number of atoms pairs is nfa/3                          */
307     res = 0;
308     fa  = 0;
309     while (fa < nfa)
310     {
311         type  = forceatoms[fa];
312         npair = ip[type].disres.npair;
313
314         Rtav_6[res] = 0.0;
315         Rt_6[res]   = 0.0;
316
317         /* Loop over the atom pairs of 'this' restraint */
318         np = 0;
319         while (fa < nfa && np < npair)
320         {
321             pair = fa/3;
322             ai   = forceatoms[fa+1];
323             aj   = forceatoms[fa+2];
324
325             if (pbc)
326             {
327                 pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
328             }
329             else
330             {
331                 rvec_sub(x[ai], x[aj], dx);
332             }
333             rt2  = iprod(dx, dx);
334             rt_1 = gmx_invsqrt(rt2);
335             rt_3 = rt_1*rt_1*rt_1;
336
337             rt[pair]         = sqrt(rt2);
338             if (bTav)
339             {
340                 /* Here we update rm3tav in t_fcdata using the data
341                  * in history_t.
342                  * Thus the results stay correct when this routine
343                  * is called multiple times.
344                  */
345                 rm3tav[pair] = cf2*((ETerm - cf1)*hist->disre_rm3tav[pair] +
346                                     ETerm1*rt_3);
347             }
348             else
349             {
350                 rm3tav[pair] = rt_3;
351             }
352
353             Rt_6[res]       += rt_3*rt_3;
354             Rtav_6[res]     += rm3tav[pair]*rm3tav[pair];
355
356             fa += 3;
357             np++;
358         }
359         if (dd->nsystems > 1)
360         {
361             Rtl_6[res]   = Rt_6[res];
362             Rt_6[res]   *= invn;
363             Rtav_6[res] *= invn;
364         }
365
366         res++;
367     }
368 }
369
370 real ta_disres(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
371                const rvec x[], rvec f[], rvec fshift[],
372                const t_pbc *pbc, const t_graph *g,
373                real gmx_unused lambda, real gmx_unused *dvdlambda,
374                const t_mdatoms gmx_unused *md, t_fcdata *fcd,
375                int gmx_unused *global_atom_index)
376 {
377     const real      sixth       = 1.0/6.0;
378     const real      seven_three = 7.0/3.0;
379
380     atom_id         ai, aj;
381     int             fa, res, npair, p, pair, ki = CENTRAL, m;
382     int             type;
383     rvec            dx;
384     real            weight_rt_1;
385     real            smooth_fc, Rt, Rtav, rt2, *Rtl_6, *Rt_6, *Rtav_6;
386     real            k0, f_scal = 0, fmax_scal, fk_scal, fij;
387     real            tav_viol, instant_viol, mixed_viol, violtot, vtot;
388     real            tav_viol_Rtav7, instant_viol_Rtav7;
389     real            up1, up2, low;
390     gmx_bool        bConservative, bMixed, bViolation;
391     ivec            it, jt, dt;
392     t_disresdata   *dd;
393     int             dr_weighting;
394     gmx_bool        dr_bMixed;
395
396     dd           = &(fcd->disres);
397     dr_weighting = dd->dr_weighting;
398     dr_bMixed    = dd->dr_bMixed;
399     Rtl_6        = dd->Rtl_6;
400     Rt_6         = dd->Rt_6;
401     Rtav_6       = dd->Rtav_6;
402
403     tav_viol = instant_viol = mixed_viol = tav_viol_Rtav7 = instant_viol_Rtav7 = 0;
404
405     smooth_fc = dd->dr_fc;
406     if (dd->dr_tau != 0)
407     {
408         /* scaling factor to smoothly turn on the restraint forces *
409          * when using time averaging                               */
410         smooth_fc *= (1.0 - dd->exp_min_t_tau);
411     }
412
413     violtot = 0;
414     vtot    = 0;
415
416     /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
417      * the total number of atoms pairs is nfa/3                          */
418     res  = 0;
419     fa   = 0;
420     while (fa < nfa)
421     {
422         type  = forceatoms[fa];
423         /* Take action depending on restraint, calculate scalar force */
424         npair = ip[type].disres.npair;
425         up1   = ip[type].disres.up1;
426         up2   = ip[type].disres.up2;
427         low   = ip[type].disres.low;
428         k0    = smooth_fc*ip[type].disres.kfac;
429
430         /* save some flops when there is only one pair */
431         if (ip[type].disres.type != 2)
432         {
433             bConservative = (dr_weighting == edrwConservative) && (npair > 1);
434             bMixed        = dr_bMixed;
435             Rt            = pow(Rt_6[res], -sixth);
436             Rtav          = pow(Rtav_6[res], -sixth);
437         }
438         else
439         {
440             /* When rtype=2 use instantaneous not ensemble avereged distance */
441             bConservative = (npair > 1);
442             bMixed        = FALSE;
443             Rt            = pow(Rtl_6[res], -sixth);
444             Rtav          = Rt;
445         }
446
447         if (Rtav > up1)
448         {
449             bViolation = TRUE;
450             tav_viol   = Rtav - up1;
451         }
452         else if (Rtav < low)
453         {
454             bViolation = TRUE;
455             tav_viol   = Rtav - low;
456         }
457         else
458         {
459             bViolation = FALSE;
460         }
461
462         if (bViolation)
463         {
464             /* NOTE:
465              * there is no real potential when time averaging is applied
466              */
467             vtot += 0.5*k0*sqr(tav_viol);
468             if (1/vtot == 0)
469             {
470                 printf("vtot is inf: %f\n", vtot);
471             }
472             if (!bMixed)
473             {
474                 f_scal   = -k0*tav_viol;
475                 violtot += fabs(tav_viol);
476             }
477             else
478             {
479                 if (Rt > up1)
480                 {
481                     if (tav_viol > 0)
482                     {
483                         instant_viol = Rt - up1;
484                     }
485                     else
486                     {
487                         bViolation = FALSE;
488                     }
489                 }
490                 else if (Rt < low)
491                 {
492                     if (tav_viol < 0)
493                     {
494                         instant_viol = Rt - low;
495                     }
496                     else
497                     {
498                         bViolation = FALSE;
499                     }
500                 }
501                 else
502                 {
503                     bViolation = FALSE;
504                 }
505                 if (bViolation)
506                 {
507                     mixed_viol = sqrt(tav_viol*instant_viol);
508                     f_scal     = -k0*mixed_viol;
509                     violtot   += mixed_viol;
510                 }
511             }
512         }
513
514         if (bViolation)
515         {
516             fmax_scal = -k0*(up2-up1);
517             /* Correct the force for the number of restraints */
518             if (bConservative)
519             {
520                 f_scal  = max(f_scal, fmax_scal);
521                 if (!bMixed)
522                 {
523                     f_scal *= Rtav/Rtav_6[res];
524                 }
525                 else
526                 {
527                     f_scal            /= 2*mixed_viol;
528                     tav_viol_Rtav7     = tav_viol*Rtav/Rtav_6[res];
529                     instant_viol_Rtav7 = instant_viol*Rt/Rt_6[res];
530                 }
531             }
532             else
533             {
534                 f_scal /= (real)npair;
535                 f_scal  = max(f_scal, fmax_scal);
536             }
537
538             /* Exert the force ... */
539
540             /* Loop over the atom pairs of 'this' restraint */
541             for (p = 0; p < npair; p++)
542             {
543                 pair = fa/3;
544                 ai   = forceatoms[fa+1];
545                 aj   = forceatoms[fa+2];
546
547                 if (pbc)
548                 {
549                     ki = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
550                 }
551                 else
552                 {
553                     rvec_sub(x[ai], x[aj], dx);
554                 }
555                 rt2 = iprod(dx, dx);
556
557                 weight_rt_1 = gmx_invsqrt(rt2);
558
559                 if (bConservative)
560                 {
561                     if (!dr_bMixed)
562                     {
563                         weight_rt_1 *= pow(dd->rm3tav[pair], seven_three);
564                     }
565                     else
566                     {
567                         weight_rt_1 *= tav_viol_Rtav7*pow(dd->rm3tav[pair], seven_three)+
568                             instant_viol_Rtav7*pow(dd->rt[pair], -7);
569                     }
570                 }
571
572                 fk_scal  = f_scal*weight_rt_1;
573
574                 if (g)
575                 {
576                     ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
577                     ki = IVEC2IS(dt);
578                 }
579
580                 for (m = 0; m < DIM; m++)
581                 {
582                     fij            = fk_scal*dx[m];
583
584                     f[ai][m]           += fij;
585                     f[aj][m]           -= fij;
586                     fshift[ki][m]      += fij;
587                     fshift[CENTRAL][m] -= fij;
588                 }
589                 fa += 3;
590             }
591         }
592         else
593         {
594             /* No violation so force and potential contributions */
595             fa += 3*npair;
596         }
597         res++;
598     }
599
600     dd->sumviol = violtot;
601
602     /* Return energy */
603     return vtot;
604 }
605
606 void update_disres_history(t_fcdata *fcd, history_t *hist)
607 {
608     t_disresdata *dd;
609     int           pair;
610
611     dd = &(fcd->disres);
612     if (dd->dr_tau != 0)
613     {
614         /* Copy the new time averages that have been calculated
615          * in calc_disres_R_6.
616          */
617         hist->disre_initf = dd->exp_min_t_tau;
618         for (pair = 0; pair < dd->npair; pair++)
619         {
620             hist->disre_rm3tav[pair] = dd->rm3tav[pair];
621         }
622     }
623 }