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,2016,2017 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/gmxlib/network.h"
46 #include "gromacs/linearalgebra/nrjac.h"
47 #include "gromacs/math/do_fit.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdrunutility/multisim.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/mdtypes/fcdata.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/mdatom.h"
55 #include "gromacs/mdtypes/state.h"
56 #include "gromacs/pbcutil/ishift.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/topology/ifunc.h"
59 #include "gromacs/topology/mtop_util.h"
60 #include "gromacs/topology/topology.h"
61 #include "gromacs/utility/arrayref.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/pleasecite.h"
64 #include "gromacs/utility/smalloc.h"
69 // TODO This implementation of ensemble orientation restraints is nasty because
70 // a user can't just do multi-sim with single-sim orientation restraints.
72 void init_orires(FILE* fplog,
73 const gmx_mtop_t* mtop,
76 const gmx_multisim_t* ms,
80 od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
83 /* Not doing orientation restraints */
87 const int numFitParams = 5;
88 if (od->nr <= numFitParams)
91 "The system has %d orientation restraints, but at least %d are required, since "
92 "there are %d fitting parameters.",
93 od->nr, numFitParams + 1, numFitParams);
96 if (ir->bPeriodicMols)
98 /* Since we apply fitting, we need to make molecules whole and this
99 * can not be done when periodic molecules are present.
102 "Orientation restraints can not be applied when periodic molecules are present "
109 "Orientation restraints do not work with MPI parallelization. Choose 1 MPI rank, "
113 GMX_RELEASE_ASSERT(globalState != nullptr, "We need a valid global state in init_orires");
115 od->fc = ir->orires_fc;
122 int* nr_ex = nullptr;
123 int typeMin = INT_MAX;
125 gmx_mtop_ilistloop_t iloop = gmx_mtop_ilistloop_init(mtop);
127 while (const InteractionLists* il = gmx_mtop_ilistloop_next(iloop, &nmol))
129 const int numOrires = (*il)[F_ORIRES].size();
130 if (nmol > 1 && numOrires > 0)
133 "Found %d copies of a molecule with orientation restrains while the current "
134 "code only supports a single copy. If you want to ensemble average, run "
135 "multiple copies of the system using the multi-sim feature of mdrun.",
139 for (int i = 0; i < numOrires; i += 3)
141 int type = (*il)[F_ORIRES].iatoms[i];
142 int ex = mtop->ffparams.iparams[type].orires.ex;
145 srenew(nr_ex, ex + 1);
146 for (int j = od->nex; j < ex + 1; j++)
152 GMX_ASSERT(nr_ex, "Check for allocated nr_ex to keep the static analyzer happy");
155 typeMin = std::min(typeMin, type);
156 typeMax = std::max(typeMax, type);
159 /* With domain decomposition we use the type index for indexing in global arrays */
161 typeMax - typeMin + 1 == od->nr,
162 "All orientation restraint parameter entries in the topology should be consecutive");
163 /* Store typeMin so we can index array with the type offset */
164 od->typeMin = typeMin;
166 snew(od->S, od->nex);
167 /* When not doing time averaging, the instaneous and time averaged data
168 * are indentical and the pointers can point to the same memory.
170 snew(od->Dinsl, od->nr);
174 snew(od->Dins, od->nr);
178 od->Dins = od->Dinsl;
181 if (ir->orires_tau == 0)
189 snew(od->Dtav, od->nr);
190 od->edt = std::exp(-ir->delta_t / ir->orires_tau);
191 od->edt_1 = 1.0 - od->edt;
193 /* Extend the state with the orires history */
194 globalState->flags |= (1 << estORIRE_INITF);
195 globalState->hist.orire_initf = 1;
196 globalState->flags |= (1 << estORIRE_DTAV);
197 globalState->hist.norire_Dtav = od->nr * 5;
198 snew(globalState->hist.orire_Dtav, globalState->hist.norire_Dtav);
201 snew(od->oinsl, od->nr);
204 snew(od->oins, od->nr);
208 od->oins = od->oinsl;
210 if (ir->orires_tau == 0)
216 snew(od->otav, od->nr);
218 snew(od->tmpEq, od->nex);
221 for (int i = 0; i < mtop->natoms; i++)
223 if (getGroupType(mtop->groups, SimulationAtomGroupType::OrientationRestraintsFit, i) == 0)
228 snew(od->mref, od->nref);
229 snew(od->xref, od->nref);
230 snew(od->xtmp, od->nref);
232 snew(od->eig, od->nex * 12);
234 /* Determine the reference structure on the master node.
235 * Copy it to the other nodes after checking multi compatibility,
236 * so we are sure the subsystems match before copying.
238 auto x = makeArrayRef(globalState->x);
239 rvec com = { 0, 0, 0 };
242 for (const AtomProxy atomP : AtomRange(*mtop))
244 const t_atom& local = atomP.atom();
245 int i = atomP.globalAtomNumber();
246 if (mtop->groups.groupNumbers[SimulationAtomGroupType::OrientationRestraintsFit].empty()
247 || mtop->groups.groupNumbers[SimulationAtomGroupType::OrientationRestraintsFit][i] == 0)
249 /* Not correct for free-energy with changing masses */
250 od->mref[j] = local.m;
251 // Note that only one rank per sim is supported.
254 copy_rvec(x[i], od->xref[j]);
255 for (int d = 0; d < DIM; d++)
257 com[d] += od->mref[j] * x[i][d];
264 svmul(1.0 / mtot, com, com);
267 for (int j = 0; j < od->nref; j++)
269 rvec_dec(od->xref[j], com);
273 fprintf(fplog, "Found %d orientation experiments\n", od->nex);
274 for (int i = 0; i < od->nex; i++)
276 fprintf(fplog, " experiment %d has %d restraints\n", i + 1, nr_ex[i]);
281 fprintf(fplog, " the fit group consists of %d atoms and has total mass %g\n", od->nref, mtot);
285 fprintf(fplog, " the orientation restraints are ensemble averaged over %d systems\n",
286 ms->numSimulations_);
288 check_multi_int(fplog, ms, od->nr, "the number of orientation restraints", FALSE);
289 check_multi_int(fplog, ms, od->nref, "the number of fit atoms for orientation restraining", FALSE);
290 check_multi_int(fplog, ms, ir->nsteps, "nsteps", FALSE);
291 /* Copy the reference coordinates from the master to the other nodes */
292 gmx_sum_sim(DIM * od->nref, od->xref[0], ms);
295 please_cite(fplog, "Hess2003");
298 void diagonalize_orires_tensors(t_oriresdata* od)
300 if (od->M == nullptr)
303 for (int i = 0; i < DIM; i++)
307 snew(od->eig_diag, DIM);
309 for (int i = 0; i < DIM; i++)
315 for (int ex = 0; ex < od->nex; ex++)
317 /* Rotate the S tensor back to the reference frame */
319 mmul(od->R, od->S[ex], TMP);
320 mtmul(TMP, od->R, S);
321 for (int i = 0; i < DIM; i++)
323 for (int j = 0; j < DIM; j++)
325 od->M[i][j] = S[i][j];
330 jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
333 for (int i = 0; i < DIM; i++)
337 for (int i = 0; i < DIM; i++)
339 for (int j = i + 1; j < DIM; j++)
341 if (gmx::square(od->eig_diag[ord[j]]) > gmx::square(od->eig_diag[ord[i]]))
350 for (int i = 0; i < DIM; i++)
352 od->eig[ex * 12 + i] = od->eig_diag[ord[i]];
354 for (int i = 0; i < DIM; i++)
356 for (int j = 0; j < DIM; j++)
358 od->eig[ex * 12 + 3 + 3 * i + j] = od->v[j][ord[i]];
364 void print_orires_log(FILE* log, t_oriresdata* od)
368 diagonalize_orires_tensors(od);
370 for (int ex = 0; ex < od->nex; ex++)
372 eig = od->eig + ex * 12;
373 fprintf(log, " Orientation experiment %d:\n", ex + 1);
374 fprintf(log, " order parameter: %g\n", eig[0]);
375 for (int i = 0; i < DIM; i++)
377 fprintf(log, " eig: %6.3f %6.3f %6.3f %6.3f\n", (eig[0] != 0) ? eig[i] / eig[0] : eig[i],
378 eig[DIM + i * DIM + XX], eig[DIM + i * DIM + YY], eig[DIM + i * DIM + ZZ]);
384 real calc_orires_dev(const gmx_multisim_t* ms,
386 const t_iatom forceatoms[],
387 const t_iparams ip[],
389 ArrayRef<const RVec> xWholeMolecules,
396 real edt, edt_1, invn, pfac, r2, invr, corrfac, wsv2, sw, dev;
400 rvec * xref, *xtmp, com, r_unrot, r;
402 const real two_thr = 2.0 / 3.0;
406 /* This means that this is not the master node */
408 "Orientation restraints are only supported on the master rank, use fewer ranks");
411 bTAV = (od->edt != 0);
422 od->exp_min_t_tau = hist->orire_initf * edt;
424 /* Correction factor to correct for the lack of history
427 corrfac = 1.0 / (1.0 - od->exp_min_t_tau);
436 invn = 1.0 / ms->numSimulations_;
446 for (int i = 0; i < md->nr; i++)
448 if (md->cORF[i] == 0)
450 copy_rvec(xWholeMolecules[i], xtmp[j]);
451 mref[j] = md->massT[i];
452 for (int d = 0; d < DIM; d++)
454 com[d] += mref[j] * xtmp[j][d];
460 svmul(1.0 / mtot, com, com);
461 for (int j = 0; j < nref; j++)
463 rvec_dec(xtmp[j], com);
465 /* Calculate the rotation matrix to rotate x to the reference orientation */
466 calc_fit_R(DIM, nref, mref, xref, xtmp, od->R);
468 for (int fa = 0; fa < nfa; fa += 3)
470 const int type = forceatoms[fa];
471 const int restraintIndex = type - od->typeMin;
474 pbc_dx_aiuc(pbc, x[forceatoms[fa + 1]], x[forceatoms[fa + 2]], r_unrot);
478 rvec_sub(x[forceatoms[fa + 1]], x[forceatoms[fa + 2]], r_unrot);
480 mvmul(od->R, r_unrot, r);
482 invr = gmx::invsqrt(r2);
483 /* Calculate the prefactor for the D tensor, this includes the factor 3! */
484 pfac = ip[type].orires.c * invr * invr * 3;
485 for (int i = 0; i < ip[type].orires.power; i++)
489 rvec5& Dinsl = od->Dinsl[restraintIndex];
490 Dinsl[0] = pfac * (2 * r[0] * r[0] + r[1] * r[1] - r2);
491 Dinsl[1] = pfac * (2 * r[0] * r[1]);
492 Dinsl[2] = pfac * (2 * r[0] * r[2]);
493 Dinsl[3] = pfac * (2 * r[1] * r[1] + r[0] * r[0] - r2);
494 Dinsl[4] = pfac * (2 * r[1] * r[2]);
498 for (int i = 0; i < 5; i++)
500 od->Dins[restraintIndex][i] = Dinsl[i] * invn;
507 gmx_sum_sim(5 * od->nr, od->Dins[0], ms);
510 /* Calculate the order tensor S for each experiment via optimization */
511 for (int ex = 0; ex < od->nex; ex++)
513 for (int i = 0; i < 5; i++)
515 matEq[ex].rhs[i] = 0;
516 for (int j = 0; j <= i; j++)
518 matEq[ex].mat[i][j] = 0;
523 for (int fa = 0; fa < nfa; fa += 3)
525 const int type = forceatoms[fa];
526 const int restraintIndex = type - od->typeMin;
527 rvec5& Dtav = od->Dtav[restraintIndex];
530 /* Here we update Dtav in t_fcdata using the data in history_t.
531 * Thus the results stay correct when this routine
532 * is called multiple times.
534 for (int i = 0; i < 5; i++)
536 Dtav[i] = edt * hist->orire_Dtav[restraintIndex * 5 + i]
537 + edt_1 * od->Dins[restraintIndex][i];
541 int ex = ip[type].orires.ex;
542 real weight = ip[type].orires.kfac;
543 /* Calculate the vector rhs and half the matrix T for the 5 equations */
544 for (int i = 0; i < 5; i++)
546 matEq[ex].rhs[i] += Dtav[i] * ip[type].orires.obs * weight;
547 for (int j = 0; j <= i; j++)
549 matEq[ex].mat[i][j] += Dtav[i] * Dtav[j] * weight;
553 /* Now we have all the data we can calculate S */
554 for (int ex = 0; ex < od->nex; ex++)
556 OriresMatEq& eq = matEq[ex];
557 /* Correct corrfac and copy one half of T to the other half */
558 for (int i = 0; i < 5; i++)
560 eq.rhs[i] *= corrfac;
561 eq.mat[i][i] *= gmx::square(corrfac);
562 for (int j = 0; j < i; j++)
564 eq.mat[i][j] *= gmx::square(corrfac);
565 eq.mat[j][i] = eq.mat[i][j];
568 m_inv_gen(&eq.mat[0][0], 5, &eq.mat[0][0]);
569 /* Calculate the orientation tensor S for this experiment */
570 matrix& S = od->S[ex];
576 for (int i = 0; i < 5; i++)
578 S[0][0] += 1.5 * eq.mat[0][i] * eq.rhs[i];
579 S[0][1] += 1.5 * eq.mat[1][i] * eq.rhs[i];
580 S[0][2] += 1.5 * eq.mat[2][i] * eq.rhs[i];
581 S[1][1] += 1.5 * eq.mat[3][i] * eq.rhs[i];
582 S[1][2] += 1.5 * eq.mat[4][i] * eq.rhs[i];
587 S[2][2] = -S[0][0] - S[1][1];
590 const matrix* S = od->S;
595 for (int fa = 0; fa < nfa; fa += 3)
597 const int type = forceatoms[fa];
598 const int restraintIndex = type - od->typeMin;
599 const int ex = ip[type].orires.ex;
601 const rvec5& Dtav = od->Dtav[restraintIndex];
602 od->otav[restraintIndex] =
604 * (S[ex][0][0] * Dtav[0] + S[ex][0][1] * Dtav[1] + S[ex][0][2] * Dtav[2]
605 + S[ex][1][1] * Dtav[3] + S[ex][1][2] * Dtav[4]);
608 const rvec5& Dins = od->Dins[restraintIndex];
609 od->oins[restraintIndex] =
611 * (S[ex][0][0] * Dins[0] + S[ex][0][1] * Dins[1] + S[ex][0][2] * Dins[2]
612 + S[ex][1][1] * Dins[3] + S[ex][1][2] * Dins[4]);
616 /* When ensemble averaging is used recalculate the local orientation
617 * for output to the energy file.
619 const rvec5& Dinsl = od->Dinsl[restraintIndex];
620 od->oinsl[restraintIndex] =
622 * (S[ex][0][0] * Dinsl[0] + S[ex][0][1] * Dinsl[1] + S[ex][0][2] * Dinsl[2]
623 + S[ex][1][1] * Dinsl[3] + S[ex][1][2] * Dinsl[4]);
626 dev = od->otav[restraintIndex] - ip[type].orires.obs;
628 wsv2 += ip[type].orires.kfac * gmx::square(dev);
629 sw += ip[type].orires.kfac;
631 od->rmsdev = std::sqrt(wsv2 / sw);
633 /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
634 for (int ex = 0; ex < od->nex; ex++)
637 tmmul(od->R, od->S[ex], RS);
638 mmul(RS, od->R, od->S[ex]);
643 /* Approx. 120*nfa/3 flops */
647 const t_iatom forceatoms[],
648 const t_iparams ip[],
653 real gmx_unused lambda,
654 real gmx_unused* dvdlambda,
655 const t_mdatoms gmx_unused* md,
657 int gmx_unused* global_atom_index)
659 int ex, power, ki = CENTRAL;
660 real r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
663 const t_oriresdata* od;
671 bTAV = (od->edt != 0);
676 /* Smoothly switch on the restraining when time averaging is used */
677 smooth_fc *= (1.0 - od->exp_min_t_tau);
680 for (int fa = 0; fa < nfa; fa += 3)
682 const int type = forceatoms[fa];
683 const int ai = forceatoms[fa + 1];
684 const int aj = forceatoms[fa + 2];
685 const int restraintIndex = type - od->typeMin;
688 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
692 rvec_sub(x[ai], x[aj], r);
695 invr = gmx::invsqrt(r2);
697 ex = ip[type].orires.ex;
698 power = ip[type].orires.power;
699 fc = smooth_fc * ip[type].orires.kfac;
700 dev = od->otav[restraintIndex] - ip[type].orires.obs;
703 * there is no real potential when time averaging is applied
705 vtot += 0.5 * fc * gmx::square(dev);
709 /* Calculate the force as the sqrt of tav times instantaneous */
710 devins = od->oins[restraintIndex] - ip[type].orires.obs;
711 if (dev * devins <= 0)
717 dev = std::sqrt(dev * devins);
725 pfac = fc * ip[type].orires.c * invr2;
726 for (int i = 0; i < power; i++)
730 mvmul(od->S[ex], r, Sr);
731 for (int i = 0; i < DIM; i++)
733 fij[i] = -pfac * dev * (4 * Sr[i] - 2 * (2 + power) * invr2 * iprod(Sr, r) * r[i]);
736 for (int i = 0; i < DIM; i++)
742 fshift[ki][i] += fij[i];
743 fshift[CENTRAL][i] -= fij[i];
751 /* Approx. 80*nfa/3 flops */
754 void update_orires_history(const t_oriresdata& od, history_t* hist)
758 /* Copy the new time averages that have been calculated
759 * in calc_orires_dev.
761 hist->orire_initf = od.exp_min_t_tau;
762 for (int pair = 0; pair < od.nr; pair++)
764 for (int i = 0; i < 5; i++)
766 hist->orire_Dtav[pair * 5 + i] = od.Dtav[pair][i];