Moved copyrite.* to fileio from gmxlib and legacyheaders.
[gromacs.git] / src / gromacs / gmxlib / disre.cpp
blobfb075a7dc5ade450fc6393611f3ed69f207031ad
1 /*
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,2015, 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! */
38 #include "gmxpre.h"
40 #include "disre.h"
42 #include "config.h"
44 #include <cmath>
45 #include <cstdlib>
46 #include <cstring>
48 #include <algorithm>
50 #include "gromacs/fileio/copyrite.h"
51 #include "gromacs/gmxlib/main.h"
52 #include "gromacs/legacyheaders/network.h"
53 #include "gromacs/legacyheaders/types/commrec.h"
54 #include "gromacs/legacyheaders/types/fcdata.h"
55 #include "gromacs/legacyheaders/types/state.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/pbcutil/ishift.h"
58 #include "gromacs/pbcutil/mshift.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/topology/mtop_util.h"
61 #include "gromacs/topology/topology.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/smalloc.h"
66 void init_disres(FILE *fplog, const gmx_mtop_t *mtop,
67 t_inputrec *ir, const t_commrec *cr,
68 t_fcdata *fcd, t_state *state, gmx_bool bIsREMD)
70 int fa, nmol, npair, np;
71 t_disresdata *dd;
72 history_t *hist;
73 gmx_mtop_ilistloop_t iloop;
74 t_ilist *il;
75 char *ptr;
77 dd = &(fcd->disres);
79 if (gmx_mtop_ftype_count(mtop, F_DISRES) == 0)
81 dd->nres = 0;
83 return;
86 if (fplog)
88 fprintf(fplog, "Initializing the distance restraints\n");
92 if (ir->eDisre == edrEnsemble)
94 gmx_fatal(FARGS, "Sorry, distance restraints with ensemble averaging over multiple molecules in one system are not functional in this version of GROMACS");
97 dd->dr_weighting = ir->eDisreWeighting;
98 dd->dr_fc = ir->dr_fc;
99 if (EI_DYNAMICS(ir->eI))
101 dd->dr_tau = ir->dr_tau;
103 else
105 dd->dr_tau = 0.0;
107 if (dd->dr_tau == 0.0)
109 dd->dr_bMixed = FALSE;
110 dd->ETerm = 0.0;
112 else
114 dd->dr_bMixed = ir->bDisreMixed;
115 dd->ETerm = std::exp(-(ir->delta_t/ir->dr_tau));
117 dd->ETerm1 = 1.0 - dd->ETerm;
119 dd->nres = 0;
120 dd->npair = 0;
121 iloop = gmx_mtop_ilistloop_init(mtop);
122 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
124 np = 0;
125 for (fa = 0; fa < il[F_DISRES].nr; fa += 3)
127 np++;
128 npair = mtop->ffparams.iparams[il[F_DISRES].iatoms[fa]].disres.npair;
129 if (np == npair)
131 dd->nres += (ir->eDisre == edrEnsemble ? 1 : nmol)*npair;
132 dd->npair += nmol*npair;
133 np = 0;
138 if (cr && PAR(cr))
140 /* Temporary check, will be removed when disre is implemented with DD */
141 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).";
143 if (MASTER(cr))
145 fprintf(stderr, "\n%s\n\n", notestr);
147 if (fplog)
149 fprintf(fplog, "%s\n", notestr);
152 if (dd->dr_tau != 0 || ir->eDisre == edrEnsemble || cr->ms != NULL ||
153 dd->nres != dd->npair)
155 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" : "");
157 if (ir->nstdisreout != 0)
159 if (fplog)
161 fprintf(fplog, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n\n");
163 if (MASTER(cr))
165 fprintf(stderr, "\nWARNING: Can not write distance restraint data to energy file with domain decomposition\n");
167 ir->nstdisreout = 0;
171 snew(dd->rt, dd->npair);
173 if (dd->dr_tau != 0.0)
175 hist = &state->hist;
176 /* Set the "history lack" factor to 1 */
177 state->flags |= (1<<estDISRE_INITF);
178 hist->disre_initf = 1.0;
179 /* Allocate space for the r^-3 time averages */
180 state->flags |= (1<<estDISRE_RM3TAV);
181 hist->ndisrepairs = dd->npair;
182 snew(hist->disre_rm3tav, hist->ndisrepairs);
184 /* Allocate space for a copy of rm3tav,
185 * so we can call do_force without modifying the state.
187 snew(dd->rm3tav, dd->npair);
189 /* Allocate Rt_6 and Rtav_6 consecutively in memory so they can be
190 * averaged over the processors in one call (in calc_disre_R_6)
192 snew(dd->Rt_6, 2*dd->nres);
193 dd->Rtav_6 = &(dd->Rt_6[dd->nres]);
195 ptr = getenv("GMX_DISRE_ENSEMBLE_SIZE");
196 if (cr && cr->ms != NULL && ptr != NULL && !bIsREMD)
198 #ifdef GMX_MPI
199 dd->nsystems = 0;
200 sscanf(ptr, "%d", &dd->nsystems);
201 if (fplog)
203 fprintf(fplog, "Found GMX_DISRE_ENSEMBLE_SIZE set to %d systems per ensemble\n", dd->nsystems);
205 /* This check is only valid on MASTER(cr), so probably
206 * ensemble-averaged distance restraints are broken on more
207 * than one processor per simulation system. */
208 if (MASTER(cr))
210 check_multi_int(fplog, cr->ms, dd->nsystems,
211 "the number of systems per ensemble",
212 FALSE);
214 gmx_bcast_sim(sizeof(int), &dd->nsystems, cr);
216 /* We use to allow any value of nsystems which was a divisor
217 * of ms->nsim. But this required an extra communicator which
218 * was stored in t_fcdata. This pulled in mpi.h in nearly all C files.
220 if (!(cr->ms->nsim == 1 || cr->ms->nsim == dd->nsystems))
222 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);
224 if (fplog)
226 fprintf(fplog, "Our ensemble consists of systems:");
227 for (int i = 0; i < dd->nsystems; i++)
229 fprintf(fplog, " %d",
230 (cr->ms->sim/dd->nsystems)*dd->nsystems+i);
232 fprintf(fplog, "\n");
234 snew(dd->Rtl_6, dd->nres);
235 #endif
237 else
239 dd->nsystems = 1;
240 dd->Rtl_6 = dd->Rt_6;
243 if (dd->npair > 0)
245 if (fplog)
247 fprintf(fplog, "There are %d distance restraints involving %d atom pairs\n", dd->nres, dd->npair);
249 /* Have to avoid g_disre de-referencing cr blindly, mdrun not
250 * doing consistency checks for ensemble-averaged distance
251 * restraints when that's not happening, and only doing those
252 * checks from appropriate processes (since check_multi_int is
253 * too broken to check whether the communication will
254 * succeed...) */
255 if (cr && cr->ms && dd->nsystems > 1 && MASTER(cr))
257 check_multi_int(fplog, cr->ms, fcd->disres.nres,
258 "the number of distance restraints",
259 FALSE);
261 please_cite(fplog, "Tropp80a");
262 please_cite(fplog, "Torda89a");
266 void calc_disres_R_6(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
267 const rvec x[], const t_pbc *pbc,
268 t_fcdata *fcd, history_t *hist)
270 atom_id ai, aj;
271 int fa, res, pair;
272 int type, npair, np;
273 rvec dx;
274 real *rt, *rm3tav, *Rtl_6, *Rt_6, *Rtav_6;
275 real rt_1, rt_3, rt2;
276 t_disresdata *dd;
277 real ETerm, ETerm1, cf1 = 0, cf2 = 0, invn = 0;
278 gmx_bool bTav;
280 dd = &(fcd->disres);
281 bTav = (dd->dr_tau != 0);
282 ETerm = dd->ETerm;
283 ETerm1 = dd->ETerm1;
284 rt = dd->rt;
285 rm3tav = dd->rm3tav;
286 Rtl_6 = dd->Rtl_6;
287 Rt_6 = dd->Rt_6;
288 Rtav_6 = dd->Rtav_6;
290 if (bTav)
292 /* scaling factor to smoothly turn on the restraint forces *
293 * when using time averaging */
294 dd->exp_min_t_tau = hist->disre_initf*ETerm;
296 cf1 = dd->exp_min_t_tau;
297 cf2 = 1.0/(1.0 - dd->exp_min_t_tau);
300 if (dd->nsystems > 1)
302 invn = 1.0/dd->nsystems;
305 /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
306 * the total number of atoms pairs is nfa/3 */
307 res = 0;
308 fa = 0;
309 while (fa < nfa)
311 type = forceatoms[fa];
312 npair = ip[type].disres.npair;
314 Rtav_6[res] = 0.0;
315 Rt_6[res] = 0.0;
317 /* Loop over the atom pairs of 'this' restraint */
318 np = 0;
319 while (fa < nfa && np < npair)
321 pair = fa/3;
322 ai = forceatoms[fa+1];
323 aj = forceatoms[fa+2];
325 if (pbc)
327 pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
329 else
331 rvec_sub(x[ai], x[aj], dx);
333 rt2 = iprod(dx, dx);
334 rt_1 = gmx_invsqrt(rt2);
335 rt_3 = rt_1*rt_1*rt_1;
337 rt[pair] = std::sqrt(rt2);
338 if (bTav)
340 /* Here we update rm3tav in t_fcdata using the data
341 * in history_t.
342 * Thus the results stay correct when this routine
343 * is called multiple times.
345 rm3tav[pair] = cf2*((ETerm - cf1)*hist->disre_rm3tav[pair] +
346 ETerm1*rt_3);
348 else
350 rm3tav[pair] = rt_3;
353 Rt_6[res] += rt_3*rt_3;
354 Rtav_6[res] += rm3tav[pair]*rm3tav[pair];
356 fa += 3;
357 np++;
359 if (dd->nsystems > 1)
361 Rtl_6[res] = Rt_6[res];
362 Rt_6[res] *= invn;
363 Rtav_6[res] *= invn;
366 res++;
370 real ta_disres(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
371 const rvec x[], rvec f[], rvec fshift[],
372 const t_pbc *pbc, const t_graph *g,
373 real gmx_unused lambda, real gmx_unused *dvdlambda,
374 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
375 int gmx_unused *global_atom_index)
377 const real sixth = 1.0/6.0;
378 const real seven_three = 7.0/3.0;
380 atom_id ai, aj;
381 int fa, res, npair, p, pair, ki = CENTRAL, m;
382 int type;
383 rvec dx;
384 real weight_rt_1;
385 real smooth_fc, Rt, Rtav, rt2, *Rtl_6, *Rt_6, *Rtav_6;
386 real k0, f_scal = 0, fmax_scal, fk_scal, fij;
387 real tav_viol, instant_viol, mixed_viol, violtot, vtot;
388 real tav_viol_Rtav7, instant_viol_Rtav7;
389 real up1, up2, low;
390 gmx_bool bConservative, bMixed, bViolation;
391 ivec dt;
392 t_disresdata *dd;
393 int dr_weighting;
394 gmx_bool dr_bMixed;
396 dd = &(fcd->disres);
397 dr_weighting = dd->dr_weighting;
398 dr_bMixed = dd->dr_bMixed;
399 Rtl_6 = dd->Rtl_6;
400 Rt_6 = dd->Rt_6;
401 Rtav_6 = dd->Rtav_6;
403 tav_viol = instant_viol = mixed_viol = tav_viol_Rtav7 = instant_viol_Rtav7 = 0;
405 smooth_fc = dd->dr_fc;
406 if (dd->dr_tau != 0)
408 /* scaling factor to smoothly turn on the restraint forces *
409 * when using time averaging */
410 smooth_fc *= (1.0 - dd->exp_min_t_tau);
413 violtot = 0;
414 vtot = 0;
416 /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, *
417 * the total number of atoms pairs is nfa/3 */
418 res = 0;
419 fa = 0;
420 while (fa < nfa)
422 type = forceatoms[fa];
423 /* Take action depending on restraint, calculate scalar force */
424 npair = ip[type].disres.npair;
425 up1 = ip[type].disres.up1;
426 up2 = ip[type].disres.up2;
427 low = ip[type].disres.low;
428 k0 = smooth_fc*ip[type].disres.kfac;
430 /* save some flops when there is only one pair */
431 if (ip[type].disres.type != 2)
433 bConservative = (dr_weighting == edrwConservative) && (npair > 1);
434 bMixed = dr_bMixed;
435 Rt = std::pow(Rt_6[res], -sixth);
436 Rtav = std::pow(Rtav_6[res], -sixth);
438 else
440 /* When rtype=2 use instantaneous not ensemble avereged distance */
441 bConservative = (npair > 1);
442 bMixed = FALSE;
443 Rt = std::pow(Rtl_6[res], -sixth);
444 Rtav = Rt;
447 if (Rtav > up1)
449 bViolation = TRUE;
450 tav_viol = Rtav - up1;
452 else if (Rtav < low)
454 bViolation = TRUE;
455 tav_viol = Rtav - low;
457 else
459 bViolation = FALSE;
462 if (bViolation)
464 /* NOTE:
465 * there is no real potential when time averaging is applied
467 vtot += 0.5*k0*sqr(tav_viol);
468 if (1/vtot == 0)
470 printf("vtot is inf: %f\n", vtot);
472 if (!bMixed)
474 f_scal = -k0*tav_viol;
475 violtot += fabs(tav_viol);
477 else
479 if (Rt > up1)
481 if (tav_viol > 0)
483 instant_viol = Rt - up1;
485 else
487 bViolation = FALSE;
490 else if (Rt < low)
492 if (tav_viol < 0)
494 instant_viol = Rt - low;
496 else
498 bViolation = FALSE;
501 else
503 bViolation = FALSE;
505 if (bViolation)
507 mixed_viol = std::sqrt(tav_viol*instant_viol);
508 f_scal = -k0*mixed_viol;
509 violtot += mixed_viol;
514 if (bViolation)
516 fmax_scal = -k0*(up2-up1);
517 /* Correct the force for the number of restraints */
518 if (bConservative)
520 f_scal = std::max(f_scal, fmax_scal);
521 if (!bMixed)
523 f_scal *= Rtav/Rtav_6[res];
525 else
527 f_scal /= 2*mixed_viol;
528 tav_viol_Rtav7 = tav_viol*Rtav/Rtav_6[res];
529 instant_viol_Rtav7 = instant_viol*Rt/Rt_6[res];
532 else
534 f_scal /= npair;
535 f_scal = std::max(f_scal, fmax_scal);
538 /* Exert the force ... */
540 /* Loop over the atom pairs of 'this' restraint */
541 for (p = 0; p < npair; p++)
543 pair = fa/3;
544 ai = forceatoms[fa+1];
545 aj = forceatoms[fa+2];
547 if (pbc)
549 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
551 else
553 rvec_sub(x[ai], x[aj], dx);
555 rt2 = iprod(dx, dx);
557 weight_rt_1 = gmx_invsqrt(rt2);
559 if (bConservative)
561 if (!dr_bMixed)
563 weight_rt_1 *= std::pow(dd->rm3tav[pair], seven_three);
565 else
567 weight_rt_1 *= tav_viol_Rtav7*std::pow(dd->rm3tav[pair], seven_three)+
568 instant_viol_Rtav7*std::pow(dd->rt[pair], static_cast<real>(-7));
572 fk_scal = f_scal*weight_rt_1;
574 if (g)
576 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
577 ki = IVEC2IS(dt);
580 for (m = 0; m < DIM; m++)
582 fij = fk_scal*dx[m];
584 f[ai][m] += fij;
585 f[aj][m] -= fij;
586 fshift[ki][m] += fij;
587 fshift[CENTRAL][m] -= fij;
589 fa += 3;
592 else
594 /* No violation so force and potential contributions */
595 fa += 3*npair;
597 res++;
600 dd->sumviol = violtot;
602 /* Return energy */
603 return vtot;
606 void update_disres_history(t_fcdata *fcd, history_t *hist)
608 t_disresdata *dd;
609 int pair;
611 dd = &(fcd->disres);
612 if (dd->dr_tau != 0)
614 /* Copy the new time averages that have been calculated
615 * in calc_disres_R_6.
617 hist->disre_initf = dd->exp_min_t_tau;
618 for (pair = 0; pair < dd->npair; pair++)
620 hist->disre_rm3tav[pair] = dd->rm3tav[pair];