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