Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / disre.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * GROningen Mixture of Alchemy and Childrens' Stories
35  */
36 /* This file is completely threadsafe - keep it that way! */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include "typedefs.h"
43 #include "sysstuff.h"
44 #include "smalloc.h"
45 #include "macros.h"
46 #include "vec.h"
47 #include "futil.h"
48 #include "xvgr.h"
49 #include "gmx_fatal.h"
50 #include "bondf.h"
51 #include "copyrite.h"
52 #include "disre.h"
53 #include "main.h"
54 #include "mtop_util.h"
55
56 void init_disres(FILE *fplog,const gmx_mtop_t *mtop,
57                  t_inputrec *ir,const t_commrec *cr,gmx_bool bPartDecomp,
58                  t_fcdata *fcd,t_state *state)
59 {
60     int          fa,nmol,i,npair,np;
61     t_iparams    *ip;
62     t_disresdata *dd;
63     history_t    *hist;
64     gmx_mtop_ilistloop_t iloop;
65     t_ilist      *il;
66     char         *ptr;
67     
68     dd = &(fcd->disres);
69
70     if (gmx_mtop_ftype_count(mtop,F_DISRES) == 0)
71     {
72         dd->nres = 0;
73
74         return;
75     }
76     
77     if (fplog)
78     {
79         fprintf(fplog,"Initializing the distance restraints\n");
80     }
81     
82
83     if (ir->eDisre == edrEnsemble)
84     {
85         gmx_fatal(FARGS,"Sorry, distance restraints with ensemble averaging over multiple molecules in one system are not functional in this version of GROMACS");
86     }
87
88     dd->dr_weighting = ir->eDisreWeighting;
89     dd->dr_fc        = ir->dr_fc;
90     if (EI_DYNAMICS(ir->eI))
91     {
92         dd->dr_tau   = ir->dr_tau;
93     }
94     else
95     {
96         dd->dr_tau   = 0.0;
97     }
98     if (dd->dr_tau == 0.0)
99     {
100         dd->dr_bMixed = FALSE;
101         dd->ETerm = 0.0;
102     }
103     else
104     {
105         dd->dr_bMixed = ir->bDisreMixed;
106         dd->ETerm = exp(-(ir->delta_t/ir->dr_tau));
107     }
108     dd->ETerm1        = 1.0 - dd->ETerm;
109     
110     ip = mtop->ffparams.iparams;
111
112     dd->nres  = 0;
113     dd->npair = 0;
114     iloop = gmx_mtop_ilistloop_init(mtop);
115     while (gmx_mtop_ilistloop_next(iloop,&il,&nmol)) {
116         np = 0;
117         for(fa=0; fa<il[F_DISRES].nr; fa+=3)
118         {
119             np++;
120             npair = mtop->ffparams.iparams[il[F_DISRES].iatoms[fa]].disres.npair;
121             if (np == npair)
122             {
123                 dd->nres  += (ir->eDisre==edrEnsemble ? 1 : nmol)*npair;
124                 dd->npair += nmol*npair;
125                 np = 0;
126             }
127         }
128     }
129
130     if (cr && PAR(cr) && !bPartDecomp)
131     {
132         /* Temporary check, will be removed when disre is implemented with DD */
133         const char *notestr="NOTE: atoms involved in distance restraints should be within the longest cut-off distance, if this is not the case mdrun generates a fatal error, in that case use particle decomposition (mdrun option -pd)";
134         
135         if (MASTER(cr))
136             fprintf(stderr,"\n%s\n\n",notestr);
137         if (fplog)
138             fprintf(fplog,"%s\n",notestr);
139
140         if (dd->dr_tau != 0 || ir->eDisre == edrEnsemble || cr->ms != NULL ||
141             dd->nres != dd->npair)
142         {
143             gmx_fatal(FARGS,"Time or ensemble averaged or multiple pair distance restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
144         }
145         if (ir->nstdisreout != 0)
146         {
147             if (fplog)
148             {
149                 fprintf(fplog,"\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n\n");
150             }
151             if (MASTER(cr))
152             {
153                 fprintf(stderr,"\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n");
154             }
155             ir->nstdisreout = 0;
156         }
157     }
158
159     snew(dd->rt,dd->npair);
160     
161     if (dd->dr_tau != 0.0)
162     {
163         hist = &state->hist;
164         /* Set the "history lack" factor to 1 */
165         state->flags |= (1<<estDISRE_INITF);
166         hist->disre_initf = 1.0;
167         /* Allocate space for the r^-3 time averages */
168         state->flags |= (1<<estDISRE_RM3TAV);
169         hist->ndisrepairs = dd->npair;
170         snew(hist->disre_rm3tav,hist->ndisrepairs);
171     }
172     /* Allocate space for a copy of rm3tav,
173      * so we can call do_force without modifying the state.
174      */
175     snew(dd->rm3tav,dd->npair);
176
177     /* Allocate Rt_6 and Rtav_6 consecutively in memory so they can be
178      * averaged over the processors in one call (in calc_disre_R_6)
179      */
180     snew(dd->Rt_6,2*dd->nres);
181     dd->Rtav_6 = &(dd->Rt_6[dd->nres]);
182
183     ptr = getenv("GMX_DISRE_ENSEMBLE_SIZE");
184     if (cr && cr->ms != NULL && ptr != NULL)
185     {
186 #ifdef GMX_MPI
187         dd->nsystems = 0;
188         sscanf(ptr,"%d",&dd->nsystems);
189         if (fplog)
190         {
191             fprintf(fplog,"Found GMX_DISRE_ENSEMBLE_SIZE set to %d systems per ensemble\n",dd->nsystems);
192         }
193         check_multi_int(fplog,cr->ms,dd->nsystems,
194                         "the number of systems per ensemble");
195         if (dd->nsystems <= 0 ||  cr->ms->nsim % dd->nsystems != 0)
196         {
197             gmx_fatal(FARGS,"The number of systems %d is not divisible by the number of systems per ensemble %d\n",cr->ms->nsim,dd->nsystems);
198         }
199         /* Split the inter-master communicator into different ensembles */
200         MPI_Comm_split(cr->ms->mpi_comm_masters,
201                        cr->ms->sim/dd->nsystems,
202                        cr->ms->sim,
203                        &dd->mpi_comm_ensemble);
204         if (fplog)
205         {
206             fprintf(fplog,"Our ensemble consists of systems:");
207             for(i=0; i<dd->nsystems; i++)
208             {
209                 fprintf(fplog," %d",
210                         (cr->ms->sim/dd->nsystems)*dd->nsystems+i);
211             }
212             fprintf(fplog,"\n");
213             }
214         snew(dd->Rtl_6,dd->nres);
215 #endif
216     }
217     else
218     {
219         dd->nsystems = 1;
220         dd->Rtl_6 = dd->Rt_6;
221     }
222     
223     if (dd->npair > 0)
224     {
225         if (fplog) {
226             fprintf(fplog,"There are %d distance restraints involving %d atom pairs\n",dd->nres,dd->npair);
227         }
228         if (cr && cr->ms)
229         {
230             check_multi_int(fplog,cr->ms,fcd->disres.nres,
231                             "the number of distance restraints");
232         }
233         please_cite(fplog,"Tropp80a");
234         please_cite(fplog,"Torda89a");
235     }
236 }
237
238 void calc_disres_R_6(const gmx_multisim_t *ms,
239                      int nfa,const t_iatom forceatoms[],const t_iparams ip[],
240                      const rvec x[],const t_pbc *pbc,
241                      t_fcdata *fcd,history_t *hist)
242 {
243     atom_id     ai,aj;
244     int         fa,res,i,pair,ki,kj,m;
245     int         type,npair,np;
246     rvec        dx;
247     real        *rt,*rm3tav,*Rtl_6,*Rt_6,*Rtav_6;
248     real        rt_1,rt_3,rt2;
249     ivec        it,jt,dt;
250     t_disresdata *dd;
251     real        ETerm,ETerm1,cf1=0,cf2=0,invn=0;
252     gmx_bool        bTav;
253
254     dd = &(fcd->disres);
255     bTav         = (dd->dr_tau != 0);
256     ETerm        = dd->ETerm;
257     ETerm1       = dd->ETerm1;
258     rt           = dd->rt;
259     rm3tav       = dd->rm3tav;
260     Rtl_6        = dd->Rtl_6;
261     Rt_6         = dd->Rt_6;
262     Rtav_6       = dd->Rtav_6;
263     
264     if (bTav)
265     {
266         /* scaling factor to smoothly turn on the restraint forces *
267          * when using time averaging                               */
268         dd->exp_min_t_tau = hist->disre_initf*ETerm;
269         
270         cf1 = dd->exp_min_t_tau;
271         cf2 = 1.0/(1.0 - dd->exp_min_t_tau);
272     }
273     
274     if (dd->nsystems > 1)
275     {
276         invn = 1.0/dd->nsystems;
277     }
278     
279     /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
280      * the total number of atoms pairs is nfa/3                          */
281     res = 0;
282     fa  = 0;
283     while (fa < nfa)
284     {
285         type  = forceatoms[fa];
286         npair = ip[type].disres.npair;
287         
288         Rtav_6[res] = 0.0;
289         Rt_6[res]   = 0.0;
290         
291         /* Loop over the atom pairs of 'this' restraint */
292         np = 0;
293         while (fa < nfa && np < npair)
294         {
295             pair = fa/3;
296             ai   = forceatoms[fa+1];
297             aj   = forceatoms[fa+2];
298             
299             if (pbc)
300             {
301                 pbc_dx_aiuc(pbc,x[ai],x[aj],dx);
302             }
303             else
304             {
305                 rvec_sub(x[ai],x[aj],dx);
306             }
307             rt2  = iprod(dx,dx);
308             rt_1 = gmx_invsqrt(rt2);
309             rt_3 = rt_1*rt_1*rt_1;
310             
311             rt[pair]         = sqrt(rt2);
312             if (bTav)
313             {
314                 /* Here we update rm3tav in t_fcdata using the data
315                  * in history_t.
316                  * Thus the results stay correct when this routine
317                  * is called multiple times.
318                  */
319                 rm3tav[pair] = cf2*((ETerm - cf1)*hist->disre_rm3tav[pair] +
320                                     ETerm1*rt_3);
321             }
322             else
323             {
324                 rm3tav[pair] = rt_3;
325             }
326
327             Rt_6[res]       += rt_3*rt_3;
328             Rtav_6[res]     += rm3tav[pair]*rm3tav[pair];
329
330             fa += 3;
331             np++;
332         }
333         if (dd->nsystems > 1)
334         {
335             Rtl_6[res]   = Rt_6[res];
336             Rt_6[res]   *= invn;
337             Rtav_6[res] *= invn;
338         }
339
340         res++;
341     }
342     
343 #ifdef GMX_MPI
344     if (dd->nsystems > 1)
345     {
346         gmx_sum_comm(2*dd->nres,Rt_6,dd->mpi_comm_ensemble);
347     }
348 #endif
349 }
350
351 real ta_disres(int nfa,const t_iatom forceatoms[],const t_iparams ip[],
352                const rvec x[],rvec f[],rvec fshift[],
353                const t_pbc *pbc,const t_graph *g,
354                real lambda,real *dvdlambda,
355                const t_mdatoms *md,t_fcdata *fcd,
356                int *global_atom_index)
357 {
358     const real sixth=1.0/6.0;
359     const real seven_three=7.0/3.0;
360     
361     atom_id     ai,aj;
362     int         fa,res,npair,p,pair,ki=CENTRAL,m;
363     int         type;
364     rvec        dx;
365     real        weight_rt_1;
366     real        smooth_fc,Rt,Rtav,rt2,*Rtl_6,*Rt_6,*Rtav_6;
367     real        k0,f_scal=0,fmax_scal,fk_scal,fij;
368     real        tav_viol,instant_viol,mixed_viol,violtot,vtot;
369     real        tav_viol_Rtav7,instant_viol_Rtav7;
370     real        up1,up2,low;
371     gmx_bool        bConservative,bMixed,bViolation;
372     ivec        it,jt,dt;
373     t_disresdata *dd;
374     int         dr_weighting;
375     gmx_bool        dr_bMixed;
376     
377     dd = &(fcd->disres);
378     dr_weighting = dd->dr_weighting;
379     dr_bMixed    = dd->dr_bMixed;
380     Rtl_6        = dd->Rtl_6;
381     Rt_6         = dd->Rt_6;
382     Rtav_6       = dd->Rtav_6;
383
384     tav_viol=instant_viol=mixed_viol=tav_viol_Rtav7=instant_viol_Rtav7=0;
385
386     smooth_fc = dd->dr_fc;
387     if (dd->dr_tau != 0)
388     {
389         /* scaling factor to smoothly turn on the restraint forces *
390          * when using time averaging                               */
391         smooth_fc *= (1.0 - dd->exp_min_t_tau); 
392     }
393     
394     violtot = 0;
395     vtot    = 0;
396     
397     /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
398      * the total number of atoms pairs is nfa/3                          */
399     res  = 0;
400     fa   = 0;
401     while (fa < nfa)
402     {
403         type  = forceatoms[fa];
404         /* Take action depending on restraint, calculate scalar force */
405         npair = ip[type].disres.npair;
406         up1   = ip[type].disres.up1;
407         up2   = ip[type].disres.up2;
408         low   = ip[type].disres.low;
409         k0    = smooth_fc*ip[type].disres.kfac;
410         
411         /* save some flops when there is only one pair */
412         if (ip[type].disres.type != 2)
413         {
414             bConservative = (dr_weighting == edrwConservative) && (npair > 1);
415             bMixed        = dr_bMixed;
416             Rt   = pow(Rt_6[res],-sixth);
417             Rtav = pow(Rtav_6[res],-sixth);
418         }
419         else
420         {
421             /* When rtype=2 use instantaneous not ensemble avereged distance */
422             bConservative = (npair > 1);
423             bMixed        = FALSE;
424             Rt   = pow(Rtl_6[res],-sixth);
425             Rtav = Rt;
426         }
427         
428         if (Rtav > up1)
429         {
430             bViolation = TRUE;
431             tav_viol = Rtav - up1;
432         }
433         else if (Rtav < low)
434         {
435             bViolation = TRUE;
436             tav_viol = Rtav - low;
437         }
438         else
439         {
440             bViolation = FALSE;
441         }
442         
443         if (bViolation)
444         {
445             /* NOTE:
446              * there is no real potential when time averaging is applied
447              */
448             vtot += 0.5*k0*sqr(tav_viol);
449             if (1/vtot == 0)
450             {
451                 printf("vtot is inf: %f\n",vtot);
452             }
453             if (!bMixed)
454             {
455                 f_scal   = -k0*tav_viol;
456                 violtot += fabs(tav_viol);
457             }
458             else
459             {
460                 if (Rt > up1)
461                 {
462                     if (tav_viol > 0)
463                     {
464                         instant_viol = Rt - up1;
465                     }
466                     else
467                     {
468                         bViolation = FALSE;
469                     }
470                 }
471                 else if (Rt < low)
472                 {
473                     if (tav_viol < 0)
474                     {
475                         instant_viol = Rt - low;
476                     }
477                     else
478                     {
479                         bViolation = FALSE;
480                     }
481                 }
482                 else
483                 {
484                     bViolation = FALSE;
485                 }
486                 if (bViolation)
487                 {
488                     mixed_viol = sqrt(tav_viol*instant_viol);
489                     f_scal     = -k0*mixed_viol;
490                     violtot   += mixed_viol;
491                 }
492             }
493         }
494
495         if (bViolation)
496         {
497             fmax_scal = -k0*(up2-up1);
498             /* Correct the force for the number of restraints */
499             if (bConservative)
500             {
501                 f_scal  = max(f_scal,fmax_scal);
502                 if (!bMixed)
503                 {
504                     f_scal *= Rtav/Rtav_6[res];
505                 }
506                 else
507                 {
508                     f_scal /= 2*mixed_viol;
509                     tav_viol_Rtav7     = tav_viol*Rtav/Rtav_6[res];
510                     instant_viol_Rtav7 = instant_viol*Rt/Rt_6[res];
511                 }
512             }
513             else
514             {
515                 f_scal /= (real)npair;
516                 f_scal  = max(f_scal,fmax_scal);
517             }    
518             
519             /* Exert the force ... */
520             
521             /* Loop over the atom pairs of 'this' restraint */
522             for(p=0; p<npair; p++)
523             {
524                 pair = fa/3;
525                 ai   = forceatoms[fa+1];
526                 aj   = forceatoms[fa+2];
527                 
528                 if (pbc)
529                 {
530                     ki = pbc_dx_aiuc(pbc,x[ai],x[aj],dx);
531                 }
532                 else
533                 {
534                     rvec_sub(x[ai],x[aj],dx);
535                 }
536                 rt2 = iprod(dx,dx);
537                 
538                 weight_rt_1 = gmx_invsqrt(rt2);
539                 
540                 if (bConservative)
541                 {
542                     if (!dr_bMixed)
543                     {
544                         weight_rt_1 *= pow(dd->rm3tav[pair],seven_three);
545                     }
546                     else
547                     {
548                         weight_rt_1 *= tav_viol_Rtav7*pow(dd->rm3tav[pair],seven_three)+
549                             instant_viol_Rtav7*pow(dd->rt[pair],-7);
550                     }
551                 }
552                 
553                 fk_scal  = f_scal*weight_rt_1;
554                 
555                 if (g)
556                 {
557                     ivec_sub(SHIFT_IVEC(g,ai),SHIFT_IVEC(g,aj),dt);
558                     ki=IVEC2IS(dt);
559                 }
560                 
561                 for(m=0; m<DIM; m++)
562                 {
563                     fij            = fk_scal*dx[m];
564                     
565                     f[ai][m]      += fij;
566                     f[aj][m]      -= fij;
567                     fshift[ki][m] += fij;
568                     fshift[CENTRAL][m] -= fij;
569                 }
570                 fa += 3;
571             }
572         }
573         else
574         {
575             /* No violation so force and potential contributions */
576             fa += 3*npair;
577         }
578         res++;
579     }
580     
581     dd->sumviol = violtot;
582     
583     /* Return energy */
584     return vtot;
585 }
586
587 void update_disres_history(t_fcdata *fcd,history_t *hist)
588 {
589     t_disresdata *dd;
590     int pair;
591     
592     dd = &(fcd->disres);
593     if (dd->dr_tau != 0)
594     {
595         /* Copy the new time averages that have been calculated
596          * in calc_disres_R_6.
597          */
598         hist->disre_initf = dd->exp_min_t_tau;
599         for(pair=0; pair<dd->npair; pair++)
600         {
601             hist->disre_rm3tav[pair] = dd->rm3tav[pair];
602         }
603     }
604 }