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,2021, 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/exceptions.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 t_oriresdata::t_oriresdata(FILE* fplog,
73 const gmx_mtop_t& mtop,
76 const gmx_multisim_t* ms,
77 t_state* globalState) :
78 numRestraints(gmx_mtop_ftype_count(mtop, F_ORIRES))
80 GMX_RELEASE_ASSERT(numRestraints > 0,
81 "orires() should only be called with orientation restraints present");
83 const int numFitParams = 5;
84 if (numRestraints <= numFitParams)
86 const std::string mesg = gmx::formatString(
87 "The system has %d orientation restraints, but at least %d are required, since "
88 "there are %d fitting parameters.",
92 GMX_THROW(gmx::InvalidInputError(mesg));
97 /* Since we apply fitting, we need to make molecules whole and this
98 * can not be done when periodic molecules are present.
100 GMX_THROW(gmx::InvalidInputError(
101 "Orientation restraints can not be applied when periodic molecules are present "
107 GMX_THROW(gmx::InvalidInputError(
108 "Orientation restraints do not work with MPI parallelization. Choose 1 MPI rank, "
112 GMX_RELEASE_ASSERT(globalState != nullptr, "We need a valid global state in t_oriresdata()");
117 std::vector<int> nr_ex;
120 for (const auto il : IListRange(mtop))
122 const int numOrires = il.list()[F_ORIRES].size();
123 if (il.nmol() > 1 && numOrires > 0)
125 const std::string mesg = gmx::formatString(
126 "Found %d copies of a molecule with orientation restrains while the current "
127 "code only supports a single copy. If you want to ensemble average, run "
128 "multiple copies of the system using the multi-sim feature of mdrun.",
130 GMX_THROW(gmx::InvalidInputError(mesg));
133 for (int i = 0; i < numOrires; i += 3)
135 int type = il.list()[F_ORIRES].iatoms[i];
136 int ex = mtop.ffparams.iparams[type].orires.ex;
137 if (ex >= numExperiments)
139 nr_ex.resize(ex + 1, 0);
140 numExperiments = ex + 1;
144 typeMin = std::min(typeMin, type);
145 typeMax = std::max(typeMax, type);
148 /* With domain decomposition we use the type index for indexing in global arrays */
150 typeMax - typeMin + 1 == numRestraints,
151 "All orientation restraint parameter entries in the topology should be consecutive");
153 snew(orderTensors, numExperiments);
154 /* When not doing time averaging, the instaneous and time averaged data
155 * are indentical and the pointers can point to the same memory.
157 snew(DTensors, numRestraints);
161 snew(DTensorsEnsembleAv, numRestraints);
165 DTensorsEnsembleAv = DTensors;
168 if (ir.orires_tau == 0)
170 DTensorsTimeAndEnsembleAv = DTensorsEnsembleAv;
176 snew(DTensorsTimeAndEnsembleAv, numRestraints);
177 edt = std::exp(-ir.delta_t / ir.orires_tau);
180 /* Extend the state with the orires history */
181 globalState->flags |= enumValueToBitMask(StateEntry::OrireInitF);
182 globalState->hist.orire_initf = 1;
183 globalState->flags |= enumValueToBitMask(StateEntry::OrireDtav);
184 globalState->hist.orire_Dtav.resize(numRestraints * 5);
187 orientations.resize(numRestraints);
190 orientationsEnsembleAvBuffer.resize(numRestraints);
191 orientationsEnsembleAv = orientationsEnsembleAvBuffer;
195 orientationsEnsembleAv = orientations;
197 if (ir.orires_tau == 0)
199 orientationsTimeAndEnsembleAv = orientationsEnsembleAv;
203 orientationsTimeAndEnsembleAvBuffer.resize(numRestraints);
204 orientationsTimeAndEnsembleAv = orientationsTimeAndEnsembleAvBuffer;
206 tmpEq.resize(numExperiments);
208 numReferenceAtoms = 0;
209 for (int i = 0; i < mtop.natoms; i++)
211 if (getGroupType(mtop.groups, SimulationAtomGroupType::OrientationRestraintsFit, i) == 0)
216 mref.resize(numReferenceAtoms);
217 xref.resize(numReferenceAtoms);
218 xtmp.resize(numReferenceAtoms);
220 eigenOutput.resize(numExperiments * c_numEigenRealsPerExperiment);
222 /* Determine the reference structure on the master node.
223 * Copy it to the other nodes after checking multi compatibility,
224 * so we are sure the subsystems match before copying.
226 auto x = makeArrayRef(globalState->x);
227 rvec com = { 0, 0, 0 };
230 for (const AtomProxy atomP : AtomRange(mtop))
232 const t_atom& local = atomP.atom();
233 int i = atomP.globalAtomNumber();
234 if (mtop.groups.groupNumbers[SimulationAtomGroupType::OrientationRestraintsFit].empty()
235 || mtop.groups.groupNumbers[SimulationAtomGroupType::OrientationRestraintsFit][i] == 0)
237 /* Not correct for free-energy with changing masses */
239 // Note that only one rank per sim is supported.
242 copy_rvec(x[i], xref[j]);
243 for (int d = 0; d < DIM; d++)
245 com[d] += mref[j] * x[i][d];
252 svmul(1.0 / mtot, com, com);
255 for (int j = 0; j < numReferenceAtoms; j++)
257 rvec_dec(xref[j], com);
263 fprintf(fplog, "Found %d orientation experiments\n", numExperiments);
264 for (int i = 0; i < numExperiments; i++)
266 fprintf(fplog, " experiment %d has %d restraints\n", i + 1, nr_ex[i]);
269 fprintf(fplog, " the fit group consists of %d atoms and has total mass %g\n", numReferenceAtoms, mtot);
277 " the orientation restraints are ensemble averaged over %d systems\n",
278 ms->numSimulations_);
281 check_multi_int(fplog, ms, numRestraints, "the number of orientation restraints", FALSE);
283 fplog, ms, numReferenceAtoms, "the number of fit atoms for orientation restraining", FALSE);
284 check_multi_int(fplog, ms, ir.nsteps, "nsteps", FALSE);
285 /* Copy the reference coordinates from the master to the other nodes */
286 gmx_sum_sim(DIM * numReferenceAtoms, xref[0], ms);
289 please_cite(fplog, "Hess2003");
292 t_oriresdata::~t_oriresdata()
295 if (DTensorsTimeAndEnsembleAv != DTensorsEnsembleAv)
297 sfree(DTensorsTimeAndEnsembleAv);
299 if (DTensorsEnsembleAv != DTensors)
301 sfree(DTensorsEnsembleAv);
306 void diagonalize_orires_tensors(t_oriresdata* od)
308 for (int ex = 0; ex < od->numExperiments; ex++)
310 /* Rotate the S tensor back to the reference frame */
312 mmul(od->rotationMatrix, od->orderTensors[ex], TMP);
313 mtmul(TMP, od->rotationMatrix, S);
314 for (int i = 0; i < DIM; i++)
316 for (int j = 0; j < DIM; j++)
318 od->M[i][j] = S[i][j];
322 jacobi(od->M, od->eig_diag, od->v);
325 for (int i = 0; i < DIM; i++)
329 for (int i = 0; i < DIM; i++)
331 for (int j = i + 1; j < DIM; j++)
333 if (gmx::square(od->eig_diag[ord[j]]) > gmx::square(od->eig_diag[ord[i]]))
342 for (int i = 0; i < DIM; i++)
344 od->eigenOutput[ex * t_oriresdata::c_numEigenRealsPerExperiment + i] = od->eig_diag[ord[i]];
346 for (int i = 0; i < DIM; i++)
348 for (int j = 0; j < DIM; j++)
350 od->eigenOutput[ex * t_oriresdata::c_numEigenRealsPerExperiment + 3 + 3 * i + j] =
357 void print_orires_log(FILE* log, t_oriresdata* od)
359 diagonalize_orires_tensors(od);
361 for (int ex = 0; ex < od->numExperiments; ex++)
363 const real* eig = od->eigenOutput.data() + ex * t_oriresdata::c_numEigenRealsPerExperiment;
364 fprintf(log, " Orientation experiment %d:\n", ex + 1);
365 fprintf(log, " order parameter: %g\n", eig[0]);
366 for (int i = 0; i < DIM; i++)
369 " eig: %6.3f %6.3f %6.3f %6.3f\n",
370 (eig[0] != 0) ? eig[i] / eig[0] : eig[i],
371 eig[DIM + i * DIM + XX],
372 eig[DIM + i * DIM + YY],
373 eig[DIM + i * DIM + ZZ]);
379 real calc_orires_dev(const gmx_multisim_t* ms,
381 const t_iatom forceatoms[],
382 const t_iparams ip[],
384 ArrayRef<const RVec> xWholeMolecules,
388 const history_t* hist)
390 real invn, pfac, r2, invr, corrfac, wsv2, sw, dev;
392 rvec com, r_unrot, r;
393 const real two_thr = 2.0 / 3.0;
395 const bool bTAV = (od->edt != 0);
396 const real edt = od->edt;
397 const real edt_1 = od->edt_1;
398 gmx::ArrayRef<OriresMatEq> matEq = od->tmpEq;
399 const int nref = od->numReferenceAtoms;
400 gmx::ArrayRef<real> mref = od->mref;
401 gmx::ArrayRef<const gmx::RVec> xref = od->xref;
402 gmx::ArrayRef<gmx::RVec> xtmp = od->xtmp;
406 od->exp_min_t_tau = hist->orire_initf * edt;
408 /* Correction factor to correct for the lack of history
411 corrfac = 1.0 / (1.0 - od->exp_min_t_tau);
420 invn = 1.0 / ms->numSimulations_;
430 auto* massT = md->massT;
431 auto* cORF = md->cORF;
432 for (int i = 0; i < md->nr; i++)
436 copy_rvec(xWholeMolecules[i], xtmp[j]);
438 for (int d = 0; d < DIM; d++)
440 com[d] += mref[j] * xtmp[j][d];
446 svmul(1.0 / mtot, com, com);
447 for (int j = 0; j < nref; j++)
449 rvec_dec(xtmp[j], com);
451 /* Calculate the rotation matrix to rotate x to the reference orientation */
452 calc_fit_R(DIM, nref, mref.data(), as_rvec_array(xref.data()), as_rvec_array(xtmp.data()), od->rotationMatrix);
454 for (int fa = 0; fa < nfa; fa += 3)
456 const int type = forceatoms[fa];
457 const int restraintIndex = type - od->typeMin;
460 pbc_dx_aiuc(pbc, x[forceatoms[fa + 1]], x[forceatoms[fa + 2]], r_unrot);
464 rvec_sub(x[forceatoms[fa + 1]], x[forceatoms[fa + 2]], r_unrot);
466 mvmul(od->rotationMatrix, r_unrot, r);
468 invr = gmx::invsqrt(r2);
469 /* Calculate the prefactor for the D tensor, this includes the factor 3! */
470 pfac = ip[type].orires.c * invr * invr * 3;
471 for (int i = 0; i < ip[type].orires.power; i++)
475 rvec5& Dinsl = od->DTensors[restraintIndex];
476 Dinsl[0] = pfac * (2 * r[0] * r[0] + r[1] * r[1] - r2);
477 Dinsl[1] = pfac * (2 * r[0] * r[1]);
478 Dinsl[2] = pfac * (2 * r[0] * r[2]);
479 Dinsl[3] = pfac * (2 * r[1] * r[1] + r[0] * r[0] - r2);
480 Dinsl[4] = pfac * (2 * r[1] * r[2]);
484 for (int i = 0; i < 5; i++)
486 od->DTensorsEnsembleAv[restraintIndex][i] = Dinsl[i] * invn;
493 gmx_sum_sim(5 * od->numRestraints, od->DTensorsEnsembleAv[0], ms);
496 /* Calculate the order tensor S for each experiment via optimization */
497 for (int ex = 0; ex < od->numExperiments; ex++)
499 for (int i = 0; i < 5; i++)
501 matEq[ex].rhs[i] = 0;
502 for (int j = 0; j <= i; j++)
504 matEq[ex].mat[i][j] = 0;
509 for (int fa = 0; fa < nfa; fa += 3)
511 const int type = forceatoms[fa];
512 const int restraintIndex = type - od->typeMin;
513 rvec5& Dtav = od->DTensorsTimeAndEnsembleAv[restraintIndex];
516 /* Here we update DTensorsTimeAndEnsembleAv in t_fcdata using the data in history_t.
517 * Thus the results stay correct when this routine
518 * is called multiple times.
520 for (int i = 0; i < 5; i++)
522 Dtav[i] = edt * hist->orire_Dtav[restraintIndex * 5 + i]
523 + edt_1 * od->DTensorsEnsembleAv[restraintIndex][i];
527 int ex = ip[type].orires.ex;
528 real weight = ip[type].orires.kfac;
529 /* Calculate the vector rhs and half the matrix T for the 5 equations */
530 for (int i = 0; i < 5; i++)
532 matEq[ex].rhs[i] += Dtav[i] * ip[type].orires.obs * weight;
533 for (int j = 0; j <= i; j++)
535 matEq[ex].mat[i][j] += Dtav[i] * Dtav[j] * weight;
539 /* Now we have all the data we can calculate S */
540 for (int ex = 0; ex < od->numExperiments; ex++)
542 OriresMatEq& eq = matEq[ex];
543 /* Correct corrfac and copy one half of T to the other half */
544 for (int i = 0; i < 5; i++)
546 eq.rhs[i] *= corrfac;
547 eq.mat[i][i] *= gmx::square(corrfac);
548 for (int j = 0; j < i; j++)
550 eq.mat[i][j] *= gmx::square(corrfac);
551 eq.mat[j][i] = eq.mat[i][j];
554 m_inv_gen(&eq.mat[0][0], 5, &eq.mat[0][0]);
555 /* Calculate the orientation tensor S for this experiment */
556 matrix& S = od->orderTensors[ex];
562 for (int i = 0; i < 5; i++)
564 S[0][0] += 1.5 * eq.mat[0][i] * eq.rhs[i];
565 S[0][1] += 1.5 * eq.mat[1][i] * eq.rhs[i];
566 S[0][2] += 1.5 * eq.mat[2][i] * eq.rhs[i];
567 S[1][1] += 1.5 * eq.mat[3][i] * eq.rhs[i];
568 S[1][2] += 1.5 * eq.mat[4][i] * eq.rhs[i];
573 S[2][2] = -S[0][0] - S[1][1];
576 const matrix* S = od->orderTensors;
581 for (int fa = 0; fa < nfa; fa += 3)
583 const int type = forceatoms[fa];
584 const int restraintIndex = type - od->typeMin;
585 const int ex = ip[type].orires.ex;
587 const rvec5& Dtav = od->DTensorsTimeAndEnsembleAv[restraintIndex];
588 od->orientationsTimeAndEnsembleAv[restraintIndex] =
590 * (S[ex][0][0] * Dtav[0] + S[ex][0][1] * Dtav[1] + S[ex][0][2] * Dtav[2]
591 + S[ex][1][1] * Dtav[3] + S[ex][1][2] * Dtav[4]);
594 const rvec5& Dins = od->DTensorsEnsembleAv[restraintIndex];
595 od->orientationsEnsembleAv[restraintIndex] =
597 * (S[ex][0][0] * Dins[0] + S[ex][0][1] * Dins[1] + S[ex][0][2] * Dins[2]
598 + S[ex][1][1] * Dins[3] + S[ex][1][2] * Dins[4]);
602 /* When ensemble averaging is used recalculate the local orientation
603 * for output to the energy file.
605 const rvec5& Dinsl = od->DTensors[restraintIndex];
606 od->orientations[restraintIndex] =
608 * (S[ex][0][0] * Dinsl[0] + S[ex][0][1] * Dinsl[1] + S[ex][0][2] * Dinsl[2]
609 + S[ex][1][1] * Dinsl[3] + S[ex][1][2] * Dinsl[4]);
612 dev = od->orientationsTimeAndEnsembleAv[restraintIndex] - ip[type].orires.obs;
614 wsv2 += ip[type].orires.kfac * gmx::square(dev);
615 sw += ip[type].orires.kfac;
617 od->rmsdev = std::sqrt(wsv2 / sw);
619 /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
620 for (int ex = 0; ex < od->numExperiments; ex++)
623 tmmul(od->rotationMatrix, od->orderTensors[ex], RS);
624 mmul(RS, od->rotationMatrix, od->orderTensors[ex]);
629 /* Approx. 120*nfa/3 flops */
633 const t_iatom forceatoms[],
634 const t_iparams ip[],
639 real gmx_unused lambda,
640 real gmx_unused* dvdlambda,
641 gmx::ArrayRef<const real> /*charge*/,
642 t_fcdata gmx_unused* fcd,
643 t_disresdata gmx_unused* disresdata,
644 t_oriresdata* oriresdata,
645 int gmx_unused* global_atom_index)
647 int ex, power, ki = gmx::c_centralShiftIndex;
648 real r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
655 if (oriresdata->fc != 0)
657 bTAV = (oriresdata->edt != 0);
659 smooth_fc = oriresdata->fc;
662 /* Smoothly switch on the restraining when time averaging is used */
663 smooth_fc *= (1.0 - oriresdata->exp_min_t_tau);
666 for (int fa = 0; fa < nfa; fa += 3)
668 const int type = forceatoms[fa];
669 const int ai = forceatoms[fa + 1];
670 const int aj = forceatoms[fa + 2];
671 const int restraintIndex = type - oriresdata->typeMin;
674 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
678 rvec_sub(x[ai], x[aj], r);
681 invr = gmx::invsqrt(r2);
683 ex = ip[type].orires.ex;
684 power = ip[type].orires.power;
685 fc = smooth_fc * ip[type].orires.kfac;
686 dev = oriresdata->orientationsTimeAndEnsembleAv[restraintIndex] - ip[type].orires.obs;
689 * there is no real potential when time averaging is applied
691 vtot += 0.5 * fc * gmx::square(dev);
695 /* Calculate the force as the sqrt of tav times instantaneous */
696 devins = oriresdata->orientationsEnsembleAv[restraintIndex] - ip[type].orires.obs;
697 if (dev * devins <= 0)
703 dev = std::sqrt(dev * devins);
711 pfac = fc * ip[type].orires.c * invr2;
712 for (int i = 0; i < power; i++)
716 mvmul(oriresdata->orderTensors[ex], r, Sr);
717 for (int i = 0; i < DIM; i++)
719 fij[i] = -pfac * dev * (4 * Sr[i] - 2 * (2 + power) * invr2 * iprod(Sr, r) * r[i]);
722 for (int i = 0; i < DIM; i++)
728 fshift[ki][i] += fij[i];
729 fshift[gmx::c_centralShiftIndex][i] -= fij[i];
737 /* Approx. 80*nfa/3 flops */
740 void update_orires_history(const t_oriresdata& od, history_t* hist)
744 /* Copy the new time averages that have been calculated
745 * in calc_orires_dev.
747 hist->orire_initf = od.exp_min_t_tau;
748 for (int pair = 0; pair < od.numRestraints; pair++)
750 for (int i = 0; i < 5; i++)
752 hist->orire_Dtav[pair * 5 + i] = od.DTensorsTimeAndEnsembleAv[pair][i];