2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
48 #include "gromacs/fileio/futil.h"
50 #include "gmx_fatal.h"
55 #include "mtop_util.h"
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)
61 int fa, nmol, i, npair, np;
65 gmx_mtop_ilistloop_t iloop;
71 if (gmx_mtop_ftype_count(mtop, F_DISRES) == 0)
80 fprintf(fplog, "Initializing the distance restraints\n");
84 if (ir->eDisre == edrEnsemble)
86 gmx_fatal(FARGS, "Sorry, distance restraints with ensemble averaging over multiple molecules in one system are not functional in this version of GROMACS");
89 dd->dr_weighting = ir->eDisreWeighting;
90 dd->dr_fc = ir->dr_fc;
91 if (EI_DYNAMICS(ir->eI))
93 dd->dr_tau = ir->dr_tau;
99 if (dd->dr_tau == 0.0)
101 dd->dr_bMixed = FALSE;
106 dd->dr_bMixed = ir->bDisreMixed;
107 dd->ETerm = exp(-(ir->delta_t/ir->dr_tau));
109 dd->ETerm1 = 1.0 - dd->ETerm;
111 ip = mtop->ffparams.iparams;
115 iloop = gmx_mtop_ilistloop_init(mtop);
116 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
119 for (fa = 0; fa < il[F_DISRES].nr; fa += 3)
122 npair = mtop->ffparams.iparams[il[F_DISRES].iatoms[fa]].disres.npair;
125 dd->nres += (ir->eDisre == edrEnsemble ? 1 : nmol)*npair;
126 dd->npair += nmol*npair;
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).";
139 fprintf(stderr, "\n%s\n\n", notestr);
143 fprintf(fplog, "%s\n", notestr);
146 if (dd->dr_tau != 0 || ir->eDisre == edrEnsemble || cr->ms != NULL ||
147 dd->nres != dd->npair)
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" : "");
151 if (ir->nstdisreout != 0)
155 fprintf(fplog, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n\n");
159 fprintf(stderr, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n");
165 snew(dd->rt, dd->npair);
167 if (dd->dr_tau != 0.0)
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);
178 /* Allocate space for a copy of rm3tav,
179 * so we can call do_force without modifying the state.
181 snew(dd->rm3tav, dd->npair);
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)
186 snew(dd->Rt_6, 2*dd->nres);
187 dd->Rtav_6 = &(dd->Rt_6[dd->nres]);
189 ptr = getenv("GMX_DISRE_ENSEMBLE_SIZE");
190 if (cr && cr->ms != NULL && ptr != NULL && !bIsREMD)
194 sscanf(ptr, "%d", &dd->nsystems);
197 fprintf(fplog, "Found GMX_DISRE_ENSEMBLE_SIZE set to %d systems per ensemble\n", dd->nsystems);
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. */
204 check_multi_int(fplog, cr->ms, dd->nsystems,
205 "the number of systems per ensemble",
208 gmx_bcast_sim(sizeof(int), &dd->nsystems, cr);
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.
214 if (!(cr->ms->nsim == 1 || cr->ms->nsim == dd->nsystems))
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);
220 fprintf(fplog, "Our ensemble consists of systems:");
221 for (i = 0; i < dd->nsystems; i++)
223 fprintf(fplog, " %d",
224 (cr->ms->sim/dd->nsystems)*dd->nsystems+i);
226 fprintf(fplog, "\n");
228 snew(dd->Rtl_6, dd->nres);
234 dd->Rtl_6 = dd->Rt_6;
241 fprintf(fplog, "There are %d distance restraints involving %d atom pairs\n", dd->nres, dd->npair);
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
249 if (cr && cr->ms && dd->nsystems > 1 && MASTER(cr))
251 check_multi_int(fplog, cr->ms, fcd->disres.nres,
252 "the number of distance restraints",
255 please_cite(fplog, "Tropp80a");
256 please_cite(fplog, "Torda89a");
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)
265 int fa, res, i, pair, ki, kj, m;
268 real *rt, *rm3tav, *Rtl_6, *Rt_6, *Rtav_6;
269 real rt_1, rt_3, rt2;
272 real ETerm, ETerm1, cf1 = 0, cf2 = 0, invn = 0;
276 bTav = (dd->dr_tau != 0);
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;
291 cf1 = dd->exp_min_t_tau;
292 cf2 = 1.0/(1.0 - dd->exp_min_t_tau);
295 if (dd->nsystems > 1)
297 invn = 1.0/dd->nsystems;
300 /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
301 * the total number of atoms pairs is nfa/3 */
306 type = forceatoms[fa];
307 npair = ip[type].disres.npair;
312 /* Loop over the atom pairs of 'this' restraint */
314 while (fa < nfa && np < npair)
317 ai = forceatoms[fa+1];
318 aj = forceatoms[fa+2];
322 pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
326 rvec_sub(x[ai], x[aj], dx);
329 rt_1 = gmx_invsqrt(rt2);
330 rt_3 = rt_1*rt_1*rt_1;
332 rt[pair] = sqrt(rt2);
335 /* Here we update rm3tav in t_fcdata using the data
337 * Thus the results stay correct when this routine
338 * is called multiple times.
340 rm3tav[pair] = cf2*((ETerm - cf1)*hist->disre_rm3tav[pair] +
348 Rt_6[res] += rt_3*rt_3;
349 Rtav_6[res] += rm3tav[pair]*rm3tav[pair];
354 if (dd->nsystems > 1)
356 Rtl_6[res] = Rt_6[res];
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)
372 const real sixth = 1.0/6.0;
373 const real seven_three = 7.0/3.0;
376 int fa, res, npair, p, pair, ki = CENTRAL, m;
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;
385 gmx_bool bConservative, bMixed, bViolation;
392 dr_weighting = dd->dr_weighting;
393 dr_bMixed = dd->dr_bMixed;
398 tav_viol = instant_viol = mixed_viol = tav_viol_Rtav7 = instant_viol_Rtav7 = 0;
400 smooth_fc = dd->dr_fc;
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);
411 /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
412 * the total number of atoms pairs is nfa/3 */
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;
425 /* save some flops when there is only one pair */
426 if (ip[type].disres.type != 2)
428 bConservative = (dr_weighting == edrwConservative) && (npair > 1);
430 Rt = pow(Rt_6[res], -sixth);
431 Rtav = pow(Rtav_6[res], -sixth);
435 /* When rtype=2 use instantaneous not ensemble avereged distance */
436 bConservative = (npair > 1);
438 Rt = pow(Rtl_6[res], -sixth);
445 tav_viol = Rtav - up1;
450 tav_viol = Rtav - low;
460 * there is no real potential when time averaging is applied
462 vtot += 0.5*k0*sqr(tav_viol);
465 printf("vtot is inf: %f\n", vtot);
469 f_scal = -k0*tav_viol;
470 violtot += fabs(tav_viol);
478 instant_viol = Rt - up1;
489 instant_viol = Rt - low;
502 mixed_viol = sqrt(tav_viol*instant_viol);
503 f_scal = -k0*mixed_viol;
504 violtot += mixed_viol;
511 fmax_scal = -k0*(up2-up1);
512 /* Correct the force for the number of restraints */
515 f_scal = max(f_scal, fmax_scal);
518 f_scal *= Rtav/Rtav_6[res];
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];
529 f_scal /= (real)npair;
530 f_scal = max(f_scal, fmax_scal);
533 /* Exert the force ... */
535 /* Loop over the atom pairs of 'this' restraint */
536 for (p = 0; p < npair; p++)
539 ai = forceatoms[fa+1];
540 aj = forceatoms[fa+2];
544 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
548 rvec_sub(x[ai], x[aj], dx);
552 weight_rt_1 = gmx_invsqrt(rt2);
558 weight_rt_1 *= pow(dd->rm3tav[pair], seven_three);
562 weight_rt_1 *= tav_viol_Rtav7*pow(dd->rm3tav[pair], seven_three)+
563 instant_viol_Rtav7*pow(dd->rt[pair], -7);
567 fk_scal = f_scal*weight_rt_1;
571 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
575 for (m = 0; m < DIM; m++)
581 fshift[ki][m] += fij;
582 fshift[CENTRAL][m] -= fij;
589 /* No violation so force and potential contributions */
595 dd->sumviol = violtot;
601 void update_disres_history(t_fcdata *fcd, history_t *hist)
609 /* Copy the new time averages that have been calculated
610 * in calc_disres_R_6.
612 hist->disre_initf = dd->exp_min_t_tau;
613 for (pair = 0; pair < dd->npair; pair++)
615 hist->disre_rm3tav[pair] = dd->rm3tav[pair];