1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 /* This file is completely threadsafe - keep it that way! */
49 #include "gmx_fatal.h"
54 #include "mtop_util.h"
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, gmx_bool bIsREMD)
60 int fa, nmol, i, npair, np;
64 gmx_mtop_ilistloop_t iloop;
70 if (gmx_mtop_ftype_count(mtop, F_DISRES) == 0)
79 fprintf(fplog, "Initializing the distance restraints\n");
83 if (ir->eDisre == edrEnsemble)
85 gmx_fatal(FARGS, "Sorry, distance restraints with ensemble averaging over multiple molecules in one system are not functional in this version of GROMACS");
88 dd->dr_weighting = ir->eDisreWeighting;
89 dd->dr_fc = ir->dr_fc;
90 if (EI_DYNAMICS(ir->eI))
92 dd->dr_tau = ir->dr_tau;
98 if (dd->dr_tau == 0.0)
100 dd->dr_bMixed = FALSE;
105 dd->dr_bMixed = ir->bDisreMixed;
106 dd->ETerm = exp(-(ir->delta_t/ir->dr_tau));
108 dd->ETerm1 = 1.0 - dd->ETerm;
110 ip = mtop->ffparams.iparams;
114 iloop = gmx_mtop_ilistloop_init(mtop);
115 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
118 for (fa = 0; fa < il[F_DISRES].nr; fa += 3)
121 npair = mtop->ffparams.iparams[il[F_DISRES].iatoms[fa]].disres.npair;
124 dd->nres += (ir->eDisre == edrEnsemble ? 1 : nmol)*npair;
125 dd->npair += nmol*npair;
131 if (cr && PAR(cr) && !bPartDecomp)
133 /* Temporary check, will be removed when disre is implemented with DD */
134 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)";
138 fprintf(stderr, "\n%s\n\n", notestr);
142 fprintf(fplog, "%s\n", notestr);
145 if (dd->dr_tau != 0 || ir->eDisre == edrEnsemble || cr->ms != NULL ||
146 dd->nres != dd->npair)
148 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)");
150 if (ir->nstdisreout != 0)
154 fprintf(fplog, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n\n");
158 fprintf(stderr, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n");
164 snew(dd->rt, dd->npair);
166 if (dd->dr_tau != 0.0)
169 /* Set the "history lack" factor to 1 */
170 state->flags |= (1<<estDISRE_INITF);
171 hist->disre_initf = 1.0;
172 /* Allocate space for the r^-3 time averages */
173 state->flags |= (1<<estDISRE_RM3TAV);
174 hist->ndisrepairs = dd->npair;
175 snew(hist->disre_rm3tav, hist->ndisrepairs);
177 /* Allocate space for a copy of rm3tav,
178 * so we can call do_force without modifying the state.
180 snew(dd->rm3tav, dd->npair);
182 /* Allocate Rt_6 and Rtav_6 consecutively in memory so they can be
183 * averaged over the processors in one call (in calc_disre_R_6)
185 snew(dd->Rt_6, 2*dd->nres);
186 dd->Rtav_6 = &(dd->Rt_6[dd->nres]);
188 ptr = getenv("GMX_DISRE_ENSEMBLE_SIZE");
189 if (cr && cr->ms != NULL && ptr != NULL && !bIsREMD)
193 sscanf(ptr, "%d", &dd->nsystems);
196 fprintf(fplog, "Found GMX_DISRE_ENSEMBLE_SIZE set to %d systems per ensemble\n", dd->nsystems);
198 /* This check is only valid on MASTER(cr), so probably
199 * ensemble-averaged distance restraints are broken on more
200 * than one processor per simulation system. */
203 check_multi_int(fplog, cr->ms, dd->nsystems,
204 "the number of systems per ensemble",
207 gmx_bcast_sim(sizeof(int), &dd->nsystems, cr);
209 /* We use to allow any value of nsystems which was a divisor
210 * of ms->nsim. But this required an extra communicator which
211 * was stored in t_fcdata. This pulled in mpi.h in nearly all C files.
213 if (!(cr->ms->nsim == 1 || cr->ms->nsim == dd->nsystems))
215 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);
219 fprintf(fplog, "Our ensemble consists of systems:");
220 for (i = 0; i < dd->nsystems; i++)
222 fprintf(fplog, " %d",
223 (cr->ms->sim/dd->nsystems)*dd->nsystems+i);
225 fprintf(fplog, "\n");
227 snew(dd->Rtl_6, dd->nres);
233 dd->Rtl_6 = dd->Rt_6;
240 fprintf(fplog, "There are %d distance restraints involving %d atom pairs\n", dd->nres, dd->npair);
242 /* Have to avoid g_disre de-referencing cr blindly, mdrun not
243 * doing consistency checks for ensemble-averaged distance
244 * restraints when that's not happening, and only doing those
245 * checks from appropriate processes (since check_multi_int is
246 * too broken to check whether the communication will
248 if (cr && cr->ms && dd->nsystems > 1 && MASTER(cr))
250 check_multi_int(fplog, cr->ms, fcd->disres.nres,
251 "the number of distance restraints",
254 please_cite(fplog, "Tropp80a");
255 please_cite(fplog, "Torda89a");
259 void calc_disres_R_6(const gmx_multisim_t *ms,
260 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 if (dd->nsystems > 1)
367 gmx_sum_sim(2*dd->nres, Rt_6, ms);
372 real ta_disres(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
373 const rvec x[], rvec f[], rvec fshift[],
374 const t_pbc *pbc, const t_graph *g,
375 real lambda, real *dvdlambda,
376 const t_mdatoms *md, t_fcdata *fcd,
377 int *global_atom_index)
379 const real sixth = 1.0/6.0;
380 const real seven_three = 7.0/3.0;
383 int fa, res, npair, p, pair, ki = CENTRAL, m;
387 real smooth_fc, Rt, Rtav, rt2, *Rtl_6, *Rt_6, *Rtav_6;
388 real k0, f_scal = 0, fmax_scal, fk_scal, fij;
389 real tav_viol, instant_viol, mixed_viol, violtot, vtot;
390 real tav_viol_Rtav7, instant_viol_Rtav7;
392 gmx_bool bConservative, bMixed, bViolation;
399 dr_weighting = dd->dr_weighting;
400 dr_bMixed = dd->dr_bMixed;
405 tav_viol = instant_viol = mixed_viol = tav_viol_Rtav7 = instant_viol_Rtav7 = 0;
407 smooth_fc = dd->dr_fc;
410 /* scaling factor to smoothly turn on the restraint forces *
411 * when using time averaging */
412 smooth_fc *= (1.0 - dd->exp_min_t_tau);
418 /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
419 * the total number of atoms pairs is nfa/3 */
424 type = forceatoms[fa];
425 /* Take action depending on restraint, calculate scalar force */
426 npair = ip[type].disres.npair;
427 up1 = ip[type].disres.up1;
428 up2 = ip[type].disres.up2;
429 low = ip[type].disres.low;
430 k0 = smooth_fc*ip[type].disres.kfac;
432 /* save some flops when there is only one pair */
433 if (ip[type].disres.type != 2)
435 bConservative = (dr_weighting == edrwConservative) && (npair > 1);
437 Rt = pow(Rt_6[res], -sixth);
438 Rtav = pow(Rtav_6[res], -sixth);
442 /* When rtype=2 use instantaneous not ensemble avereged distance */
443 bConservative = (npair > 1);
445 Rt = pow(Rtl_6[res], -sixth);
452 tav_viol = Rtav - up1;
457 tav_viol = Rtav - low;
467 * there is no real potential when time averaging is applied
469 vtot += 0.5*k0*sqr(tav_viol);
472 printf("vtot is inf: %f\n", vtot);
476 f_scal = -k0*tav_viol;
477 violtot += fabs(tav_viol);
485 instant_viol = Rt - up1;
496 instant_viol = Rt - low;
509 mixed_viol = sqrt(tav_viol*instant_viol);
510 f_scal = -k0*mixed_viol;
511 violtot += mixed_viol;
518 fmax_scal = -k0*(up2-up1);
519 /* Correct the force for the number of restraints */
522 f_scal = max(f_scal, fmax_scal);
525 f_scal *= Rtav/Rtav_6[res];
529 f_scal /= 2*mixed_viol;
530 tav_viol_Rtav7 = tav_viol*Rtav/Rtav_6[res];
531 instant_viol_Rtav7 = instant_viol*Rt/Rt_6[res];
536 f_scal /= (real)npair;
537 f_scal = max(f_scal, fmax_scal);
540 /* Exert the force ... */
542 /* Loop over the atom pairs of 'this' restraint */
543 for (p = 0; p < npair; p++)
546 ai = forceatoms[fa+1];
547 aj = forceatoms[fa+2];
551 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
555 rvec_sub(x[ai], x[aj], dx);
559 weight_rt_1 = gmx_invsqrt(rt2);
565 weight_rt_1 *= pow(dd->rm3tav[pair], seven_three);
569 weight_rt_1 *= tav_viol_Rtav7*pow(dd->rm3tav[pair], seven_three)+
570 instant_viol_Rtav7*pow(dd->rt[pair], -7);
574 fk_scal = f_scal*weight_rt_1;
578 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
582 for (m = 0; m < DIM; m++)
588 fshift[ki][m] += fij;
589 fshift[CENTRAL][m] -= fij;
596 /* No violation so force and potential contributions */
602 dd->sumviol = violtot;
608 void update_disres_history(t_fcdata *fcd, history_t *hist)
616 /* Copy the new time averages that have been calculated
617 * in calc_disres_R_6.
619 hist->disre_initf = dd->exp_min_t_tau;
620 for (pair = 0; pair < dd->npair; pair++)
622 hist->disre_rm3tav[pair] = dd->rm3tav[pair];