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