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 by 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.
39 * \brief Defines the high-level constraint code.
41 * \author Berk Hess <hess@kth.se>
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
43 * \ingroup module_mdlib
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/essentialdynamics/edsam.h"
58 #include "gromacs/fileio/confio.h"
59 #include "gromacs/fileio/gmxfio.h"
60 #include "gromacs/fileio/pdbio.h"
61 #include "gromacs/gmxlib/nrnb.h"
62 #include "gromacs/math/utilities.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdlib/lincs.h"
66 #include "gromacs/mdlib/settle.h"
67 #include "gromacs/mdlib/shake.h"
68 #include "gromacs/mdtypes/commrec.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/pulling/pull.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/topology/ifunc.h"
76 #include "gromacs/topology/mtop_lookup.h"
77 #include "gromacs/topology/mtop_util.h"
78 #include "gromacs/utility/exceptions.h"
79 #include "gromacs/utility/fatalerror.h"
80 #include "gromacs/utility/gmxassert.h"
81 #include "gromacs/utility/listoflists.h"
82 #include "gromacs/utility/pleasecite.h"
83 #include "gromacs/utility/txtdump.h"
88 /* \brief Impl class for Constraints
90 * \todo Members like md, idef are valid only for the lifetime of a
91 * domain, which would be good to make clearer in the structure of the
92 * code. It should not be possible to call apply() if setConstraints()
93 * has not been called. For example, this could be achieved if
94 * setConstraints returned a valid object with such a method. That
95 * still requires that we manage the lifetime of that object
96 * correctly, however. */
97 class Constraints::Impl
100 Impl(const gmx_mtop_t& mtop_p,
101 const t_inputrec& ir_p,
104 const t_commrec* cr_p,
105 const gmx_multisim_t* ms,
107 gmx_wallcycle* wcycle_p,
108 bool pbcHandlingRequired,
112 void setConstraints(gmx_localtop_t* top,
117 bool hasMassPerturbedAtoms,
119 unsigned short* cFREEZE);
120 bool apply(bool bLog,
125 ArrayRefWithPadding<RVec> x,
126 ArrayRefWithPadding<RVec> xprime,
127 ArrayRef<RVec> min_proj,
131 ArrayRefWithPadding<RVec> v,
133 tensor constraintsVirial,
134 ConstraintVariable econq);
135 //! The total number of constraints.
137 //! The number of flexible constraints.
139 //! A list of atoms to constraints for each moleculetype.
140 std::vector<ListOfLists<int>> at2con_mt;
141 //! A list of atoms to settles for each moleculetype
142 std::vector<std::vector<int>> at2settle_mt;
144 Lincs* lincsd = nullptr; // TODO this should become a unique_ptr
146 std::unique_ptr<shakedata> shaked;
148 std::unique_ptr<SettleData> settled;
149 //! The maximum number of warnings.
151 //! The number of warnings for LINCS.
152 int warncount_lincs = 0;
153 //! The number of warnings for SETTLE.
154 int warncount_settle = 0;
155 //! The essential dynamics data.
156 gmx_edsam* ed = nullptr;
158 //! Thread-local virial contribution.
159 tensor* threadConstraintsVirial = { nullptr };
160 //! Did a settle error occur?
161 bool* bSettleErrorHasOccurred = nullptr;
163 //! Pointer to the global topology - only used for printing warnings.
164 const gmx_mtop_t& mtop;
165 //! Parameters for the interactions in this domain.
166 const InteractionDefinitions* idef = nullptr;
167 //! Total number of atoms.
169 //! Number of local atoms.
170 int numHomeAtoms_ = 0;
174 real* inverseMasses_;
175 //! If there are atoms with perturbed mass (for FEP).
176 bool hasMassPerturbedAtoms_;
177 //! FEP lambda value.
179 //! Freeze groups data
180 unsigned short* cFREEZE_;
181 //! Whether we need to do pbc for handling bonds.
182 bool pbcHandlingRequired_ = false;
186 //! Communication support.
187 const t_commrec* cr = nullptr;
188 //! Multi-sim support.
189 const gmx_multisim_t* ms = nullptr;
190 //! Pulling code object, if any.
191 pull_t* pull_work = nullptr;
192 /*!\brief Input options.
194 * \todo Replace with IMdpOptions */
195 const t_inputrec& ir;
196 //! Flop counting support.
197 t_nrnb* nrnb = nullptr;
198 //! Tracks wallcycle usage.
199 gmx_wallcycle* wcycle;
202 Constraints::~Constraints() = default;
204 int Constraints::numFlexibleConstraints() const
206 return impl_->nflexcon;
209 bool Constraints::havePerturbedConstraints() const
211 const gmx_ffparams_t& ffparams = impl_->mtop.ffparams;
213 for (size_t i = 0; i < ffparams.functype.size(); i++)
215 if ((ffparams.functype[i] == F_CONSTR || ffparams.functype[i] == F_CONSTRNC)
216 && ffparams.iparams[i].constr.dA != ffparams.iparams[i].constr.dB)
225 //! Clears constraint quantities for atoms in nonlocal region.
226 static void clear_constraint_quantity_nonlocal(gmx_domdec_t* dd, ArrayRef<RVec> q)
228 int nonlocal_at_start, nonlocal_at_end, at;
230 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
232 for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
238 void too_many_constraint_warnings(ConstraintAlgorithm eConstrAlg, int warncount)
241 "Too many %s warnings (%d)\n"
242 "If you know what you are doing you can %s"
243 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
244 "but normally it is better to fix the problem",
245 (eConstrAlg == ConstraintAlgorithm::Lincs) ? "LINCS" : "SETTLE",
247 (eConstrAlg == ConstraintAlgorithm::Lincs)
248 ? "adjust the lincs warning threshold in your mdp file\nor "
252 //! Writes out coordinates.
253 static void write_constr_pdb(const char* fn,
255 const gmx_mtop_t& mtop,
259 ArrayRef<const RVec> x,
264 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
266 const char * anm, *resnm;
269 if (DOMAINDECOMP(cr))
272 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
279 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
283 sprintf(fname, "%s.pdb", fn);
286 out = gmx_fio_fopen(fname, "w");
288 fprintf(out, "TITLE %s\n", title);
289 gmx_write_pdb_box(out, PbcType::Unset, box);
291 for (i = start; i < start + homenr; i++)
295 if (i >= dd_numHomeAtoms(*dd) && i < dd_ac0)
299 ii = dd->globalAtomIndices[i];
305 mtopGetAtomAndResidueName(mtop, ii, &molb, &anm, &resnr, &resnm, nullptr);
306 gmx_fprintf_pdb_atomline(out,
322 fprintf(out, "TER\n");
327 //! Writes out domain contents to help diagnose crashes.
328 static void dump_confs(FILE* log,
330 const gmx_mtop_t& mtop,
334 ArrayRef<const RVec> x,
335 ArrayRef<const RVec> xprime,
338 char buf[STRLEN], buf2[22];
340 char* env = getenv("GMX_SUPPRESS_DUMP");
346 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
347 write_constr_pdb(buf, "initial coordinates", mtop, start, homenr, cr, x, box);
348 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
349 write_constr_pdb(buf, "coordinates after constraining", mtop, start, homenr, cr, xprime, box);
352 fprintf(log, "Wrote pdb files with previous and current coordinates\n");
354 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
357 bool Constraints::apply(bool bLog,
362 ArrayRefWithPadding<RVec> x,
363 ArrayRefWithPadding<RVec> xprime,
364 ArrayRef<RVec> min_proj,
368 ArrayRefWithPadding<RVec> v,
370 tensor constraintsVirial,
371 ConstraintVariable econq)
373 return impl_->apply(bLog,
390 bool Constraints::Impl::apply(bool bLog,
395 ArrayRefWithPadding<RVec> x,
396 ArrayRefWithPadding<RVec> xprime,
397 ArrayRef<RVec> min_proj,
401 ArrayRefWithPadding<RVec> v,
403 tensor constraintsVirial,
404 ConstraintVariable econq)
409 real invdt, vir_fac = 0, t;
411 t_pbc pbc, *pbc_null;
415 wallcycle_start(wcycle, ewcCONSTR);
417 if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
420 "constrain called for forces displacements while not doing energy minimization, "
421 "can not do this while the LINCS and SETTLE constraint connection matrices are "
430 scaled_delta_t = step_scaling * ir.delta_t;
432 /* Prepare time step for use in constraint implementations, and
433 avoid generating inf when ir.delta_t = 0. */
440 invdt = 1.0 / scaled_delta_t;
443 if (ir.efep != FreeEnergyPerturbationType::No && EI_DYNAMICS(ir.eI))
445 /* Set the constraint lengths for the step at which this configuration
446 * is meant to be. The invmasses should not be changed.
448 lambda += delta_step * ir.fepvals->delta_lambda;
453 clear_mat(constraintsVirial);
455 const InteractionList& settle = idef->il[F_SETTLE];
456 nsettle = settle.size() / (1 + NRAL(F_SETTLE));
460 nth = gmx_omp_nthreads_get(emntSETTLE);
467 /* We do not need full pbc when constraints do not cross update groups
468 * i.e. when dd->constraint_comm==NULL.
469 * Note that PBC for constraints is different from PBC for bondeds.
470 * For constraints there is both forward and backward communication.
472 if (ir.pbcType != PbcType::No && (cr->dd || pbcHandlingRequired_)
473 && !(cr->dd && cr->dd->constraint_comm == nullptr))
475 /* With pbc=screw the screw has been changed to a shift
476 * by the constraint coordinate communication routine,
477 * so that here we can use normal pbc.
479 pbc_null = set_pbc_dd(&pbc, ir.pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, FALSE, box);
486 /* Communicate the coordinates required for the non-local constraints
487 * for LINCS and/or SETTLE.
491 dd_move_x_constraints(cr->dd,
493 x.unpaddedArrayRef(),
494 xprime.unpaddedArrayRef(),
495 econq == ConstraintVariable::Positions);
499 /* We need to initialize the non-local components of v.
500 * We never actually use these values, but we do increment them,
501 * so we should avoid uninitialized variables and overflows.
503 clear_constraint_quantity_nonlocal(cr->dd, v.unpaddedArrayRef());
507 if (lincsd != nullptr)
509 bOK = constrain_lincs(bLog || bEner,
521 hasMassPerturbedAtoms_,
525 v.unpaddedArrayRef(),
532 if (!bOK && maxwarn < INT_MAX)
537 "Constraint error in algorithm %s at step %s\n",
538 enumValueToString(ConstraintAlgorithm::Lincs),
539 gmx_step_str(step, buf));
545 if (shaked != nullptr)
547 bOK = constrain_shake(log,
552 x.unpaddedArrayRef(),
553 xprime.unpaddedArrayRef(),
560 v.unpaddedArrayRef(),
566 if (!bOK && maxwarn < INT_MAX)
571 "Constraint error in algorithm %s at step %s\n",
572 enumValueToString(ConstraintAlgorithm::Shake),
573 gmx_step_str(step, buf));
581 bool bSettleErrorHasOccurred0 = false;
585 case ConstraintVariable::Positions:
586 #pragma omp parallel for num_threads(nth) schedule(static)
587 for (int th = 0; th < nth; th++)
593 clear_mat(threadConstraintsVirial[th]);
605 th == 0 ? constraintsVirial : threadConstraintsVirial[th],
606 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
608 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
610 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
613 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle * 3);
617 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle * 3);
620 case ConstraintVariable::Velocities:
621 case ConstraintVariable::Derivative:
622 case ConstraintVariable::Force:
623 case ConstraintVariable::ForceDispl:
624 #pragma omp parallel for num_threads(nth) schedule(static)
625 for (int th = 0; th < nth; th++)
629 int calcvir_atom_end;
633 calcvir_atom_end = 0;
637 calcvir_atom_end = numHomeAtoms_;
642 clear_mat(threadConstraintsVirial[th]);
645 int start_th = (nsettle * th) / nth;
646 int end_th = (nsettle * (th + 1)) / nth;
648 if (start_th >= 0 && end_th - start_th > 0)
650 settle_proj(*settled,
653 settle.iatoms.data() + start_th * (1 + NRAL(F_SETTLE)),
655 x.unpaddedArrayRef(),
656 xprime.unpaddedArrayRef(),
659 th == 0 ? constraintsVirial : threadConstraintsVirial[th]);
662 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
664 /* This is an overestimate */
665 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
667 case ConstraintVariable::Deriv_FlexCon:
668 /* Nothing to do, since the are no flexible constraints in settles */
670 default: gmx_incons("Unknown constraint quantity for settle");
675 /* Reduce the virial contributions over the threads */
676 for (int th = 1; th < nth; th++)
678 m_add(constraintsVirial, threadConstraintsVirial[th], constraintsVirial);
682 if (econq == ConstraintVariable::Positions)
684 for (int th = 1; th < nth; th++)
686 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
689 if (bSettleErrorHasOccurred0)
695 ": One or more water molecules can not be settled.\n"
696 "Check for bad contacts and/or reduce the timestep if appropriate.\n",
700 fprintf(log, "%s", buf);
702 fprintf(stderr, "%s", buf);
704 if (warncount_settle > maxwarn)
706 too_many_constraint_warnings(ConstraintAlgorithm::Count, warncount_settle);
717 /* The normal uses of constrain() pass step_scaling = 1.0.
718 * The call to constrain() for SD1 that passes step_scaling =
719 * 0.5 also passes vir = NULL, so cannot reach this
720 * assertion. This assertion should remain until someone knows
721 * that this path works for their intended purpose, and then
722 * they can use scaled_delta_t instead of ir.delta_t
724 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
727 case ConstraintVariable::Positions: vir_fac = 0.5 / (ir.delta_t * ir.delta_t); break;
728 case ConstraintVariable::Velocities: vir_fac = 0.5 / ir.delta_t; break;
729 case ConstraintVariable::Force:
730 case ConstraintVariable::ForceDispl: vir_fac = 0.5; break;
731 default: gmx_incons("Unsupported constraint quantity for virial");
736 vir_fac *= 2; /* only constraining over half the distance here */
738 for (int i = 0; i < DIM; i++)
740 for (int j = 0; j < DIM; j++)
742 constraintsVirial[i][j] *= vir_fac;
749 dump_confs(log, step, mtop, start, numHomeAtoms_, cr, x.unpaddedArrayRef(), xprime.unpaddedArrayRef(), box);
752 if (econq == ConstraintVariable::Positions)
754 if (ir.bPull && pull_have_constraint(*pull_work))
756 if (EI_DYNAMICS(ir.eI))
758 t = ir.init_t + (step + delta_step) * ir.delta_t;
764 set_pbc(&pbc, ir.pbcType, box);
765 pull_constraint(pull_work,
771 as_rvec_array(x.unpaddedArrayRef().data()),
772 as_rvec_array(xprime.unpaddedArrayRef().data()),
773 as_rvec_array(v.unpaddedArrayRef().data()),
776 if (ed && delta_step > 0)
778 /* apply the essential dynamics constraints here */
782 as_rvec_array(xprime.unpaddedArrayRef().data()),
783 as_rvec_array(v.unpaddedArrayRef().data()),
788 wallcycle_stop(wcycle, ewcCONSTR);
790 const bool haveVelocities = (!v.empty() || econq == ConstraintVariable::Velocities);
791 if (haveVelocities && cFREEZE_)
793 /* Set the velocities of frozen dimensions to zero */
795 if (econq == ConstraintVariable::Velocities)
797 vRef = xprime.unpaddedArrayRef();
801 vRef = v.unpaddedArrayRef();
804 int gmx_unused numThreads = gmx_omp_nthreads_get(emntUpdate);
806 #pragma omp parallel for num_threads(numThreads) schedule(static)
807 for (int i = 0; i < numHomeAtoms_; i++)
809 int freezeGroup = cFREEZE_[i];
811 for (int d = 0; d < DIM; d++)
813 if (ir.opts.nFreeze[freezeGroup][d])
824 ArrayRef<real> Constraints::rmsdData() const
828 return lincs_rmsdData(impl_->lincsd);
836 real Constraints::rmsd() const
840 return lincs_rmsd(impl_->lincsd);
848 int Constraints::numConstraintsTotal()
850 return impl_->ncon_tot;
853 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
855 if (haveDynamicsIntegrator)
857 return FlexibleConstraintTreatment::Include;
861 return FlexibleConstraintTreatment::Exclude;
865 /*! \brief Returns a block struct to go from atoms to constraints
867 * The block struct will contain constraint indices with lower indices
868 * directly matching the order in F_CONSTR and higher indices matching
869 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
871 * \param[in] numAtoms The number of atoms to construct the list for
872 * \param[in] ilists The interaction lists, size F_NRE
873 * \param[in] iparams Interaction parameters, can be null when
874 * \p flexibleConstraintTreatment==Include
875 * \param[in] flexibleConstraintTreatment The flexible constraint treatment,
878 * \returns a block struct with all constraints for each atom
880 static ListOfLists<int> makeAtomsToConstraintsList(int numAtoms,
881 ArrayRef<const InteractionList> ilists,
882 ArrayRef<const t_iparams> iparams,
883 FlexibleConstraintTreatment flexibleConstraintTreatment)
885 GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || !iparams.empty(),
886 "With flexible constraint detection we need valid iparams");
888 std::vector<int> count(numAtoms);
890 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
892 const InteractionList& ilist = ilists[ftype];
893 const int stride = 1 + NRAL(ftype);
894 for (int i = 0; i < ilist.size(); i += stride)
896 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
897 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
899 for (int j = 1; j < 3; j++)
901 int a = ilist.iatoms[i + j];
908 std::vector<int> listRanges(numAtoms + 1);
909 for (int a = 0; a < numAtoms; a++)
911 listRanges[a + 1] = listRanges[a] + count[a];
914 std::vector<int> elements(listRanges[numAtoms]);
916 /* The F_CONSTRNC constraints have constraint numbers
917 * that continue after the last F_CONSTR constraint.
919 int numConstraints = 0;
920 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
922 const InteractionList& ilist = ilists[ftype];
923 const int stride = 1 + NRAL(ftype);
924 for (int i = 0; i < ilist.size(); i += stride)
926 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
927 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
929 for (int j = 1; j < 3; j++)
931 const int a = ilist.iatoms[i + j];
932 elements[listRanges[a] + count[a]++] = numConstraints;
939 return ListOfLists<int>(std::move(listRanges), std::move(elements));
942 ListOfLists<int> make_at2con(int numAtoms,
943 ArrayRef<const InteractionList> ilist,
944 ArrayRef<const t_iparams> iparams,
945 FlexibleConstraintTreatment flexibleConstraintTreatment)
947 return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
950 ListOfLists<int> make_at2con(const gmx_moltype_t& moltype,
951 gmx::ArrayRef<const t_iparams> iparams,
952 FlexibleConstraintTreatment flexibleConstraintTreatment)
954 return makeAtomsToConstraintsList(
955 moltype.atoms.nr, makeConstArrayRef(moltype.ilist), iparams, flexibleConstraintTreatment);
958 //! Return the number of flexible constraints in the \c ilist and \c iparams.
959 int countFlexibleConstraints(ArrayRef<const InteractionList> ilist, ArrayRef<const t_iparams> iparams)
962 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
964 const int numIatomsPerConstraint = 3;
965 for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
967 const int type = ilist[ftype].iatoms[i];
968 if (iparams[type].constr.dA == 0 && iparams[type].constr.dB == 0)
978 //! Returns the index of the settle to which each atom belongs.
979 static std::vector<int> make_at2settle(int natoms, const InteractionList& ilist)
981 /* Set all to no settle */
982 std::vector<int> at2s(natoms, -1);
984 const int stride = 1 + NRAL(F_SETTLE);
986 for (int s = 0; s < ilist.size(); s += stride)
988 at2s[ilist.iatoms[s + 1]] = s / stride;
989 at2s[ilist.iatoms[s + 2]] = s / stride;
990 at2s[ilist.iatoms[s + 3]] = s / stride;
996 void Constraints::Impl::setConstraints(gmx_localtop_t* top,
1000 real* inverseMasses,
1001 const bool hasMassPerturbedAtoms,
1003 unsigned short* cFREEZE)
1005 numAtoms_ = numAtoms;
1006 numHomeAtoms_ = numHomeAtoms;
1008 inverseMasses_ = inverseMasses;
1009 hasMassPerturbedAtoms_ = hasMassPerturbedAtoms;
1017 /* With DD we might also need to call LINCS on a domain no constraints for
1018 * communicating coordinates to other nodes that do have constraints.
1020 if (ir.eConstrAlg == ConstraintAlgorithm::Lincs)
1022 set_lincs(*idef, numAtoms_, inverseMasses_, lambda_, EI_DYNAMICS(ir.eI), cr, lincsd);
1024 if (ir.eConstrAlg == ConstraintAlgorithm::Shake)
1028 // We are using the local topology, so there are only
1029 // F_CONSTR constraints.
1030 GMX_RELEASE_ASSERT(idef->il[F_CONSTRNC].empty(),
1031 "Here we should not have no-connect constraints");
1032 make_shake_sblock_dd(shaked.get(), idef->il[F_CONSTR]);
1036 make_shake_sblock_serial(shaked.get(), &top->idef, numAtoms_);
1043 settled->setConstraints(idef->il[F_SETTLE], numHomeAtoms_, masses_, inverseMasses_);
1046 /* Make a selection of the local atoms for essential dynamics */
1049 dd_make_local_ed_indices(cr->dd, ed);
1053 void Constraints::setConstraints(gmx_localtop_t* top,
1055 const int numHomeAtoms,
1057 real* inverseMasses,
1058 const bool hasMassPerturbedAtoms,
1060 unsigned short* cFREEZE)
1062 impl_->setConstraints(
1063 top, numAtoms, numHomeAtoms, masses, inverseMasses, hasMassPerturbedAtoms, lambda, cFREEZE);
1066 /*! \brief Makes a per-moleculetype container of mappings from atom
1067 * indices to constraint indices.
1069 * Note that flexible constraints are only enabled with a dynamical integrator. */
1070 static std::vector<ListOfLists<int>> makeAtomToConstraintMappings(const gmx_mtop_t& mtop,
1071 FlexibleConstraintTreatment flexibleConstraintTreatment)
1073 std::vector<ListOfLists<int>> mapping;
1074 mapping.reserve(mtop.moltype.size());
1075 for (const gmx_moltype_t& moltype : mtop.moltype)
1077 mapping.push_back(make_at2con(moltype, mtop.ffparams.iparams, flexibleConstraintTreatment));
1082 Constraints::Constraints(const gmx_mtop_t& mtop,
1083 const t_inputrec& ir,
1086 const t_commrec* cr,
1087 const gmx_multisim_t* ms,
1089 gmx_wallcycle* wcycle,
1090 bool pbcHandlingRequired,
1093 impl_(new Impl(mtop, ir, pull_work, log, cr, ms, nrnb, wcycle, pbcHandlingRequired, numConstraints, numSettles))
1097 Constraints::Impl::Impl(const gmx_mtop_t& mtop_p,
1098 const t_inputrec& ir_p,
1101 const t_commrec* cr_p,
1102 const gmx_multisim_t* ms_p,
1104 gmx_wallcycle* wcycle_p,
1105 bool pbcHandlingRequired,
1108 ncon_tot(numConstraints),
1110 pbcHandlingRequired_(pbcHandlingRequired),
1114 pull_work(pull_work),
1119 if (numConstraints + numSettles > 0 && ir.epc == PressureCoupling::Mttk)
1121 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1125 if (numConstraints > 0)
1127 at2con_mt = makeAtomToConstraintMappings(mtop, flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
1129 for (const gmx_molblock_t& molblock : mtop.molblock)
1131 int count = countFlexibleConstraints(mtop.moltype[molblock.type].ilist, mtop.ffparams.iparams);
1132 nflexcon += molblock.nmol * count;
1139 fprintf(log, "There are %d flexible constraints\n", nflexcon);
1140 if (ir.fc_stepsize == 0)
1144 "WARNING: step size for flexible constraining = 0\n"
1145 " All flexible constraints will be rigid.\n"
1146 " Will try to keep all flexible constraints at their original "
1148 " but the lengths may exhibit some drift.\n\n");
1154 please_cite(log, "Hess2002");
1158 if (ir.eConstrAlg == ConstraintAlgorithm::Lincs)
1160 lincsd = init_lincs(log,
1164 DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd),
1169 if (ir.eConstrAlg == ConstraintAlgorithm::Shake)
1171 if (DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd))
1174 "SHAKE is not supported with domain decomposition and constraint that "
1175 "cross domain boundaries, use LINCS");
1180 "For this system also velocities and/or forces need to be constrained, "
1181 "this can not be done with SHAKE, you should select LINCS");
1183 please_cite(log, "Ryckaert77a");
1186 please_cite(log, "Barth95a");
1189 shaked = std::make_unique<shakedata>();
1195 please_cite(log, "Miyamoto92a");
1197 settled = std::make_unique<SettleData>(mtop);
1199 /* Make an atom to settle index for use in domain decomposition */
1200 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1202 at2settle_mt.emplace_back(
1203 make_at2settle(mtop.moltype[mt].atoms.nr, mtop.moltype[mt].ilist[F_SETTLE]));
1206 /* Allocate thread-local work arrays */
1207 int nthreads = gmx_omp_nthreads_get(emntSETTLE);
1208 if (nthreads > 1 && threadConstraintsVirial == nullptr)
1210 snew(threadConstraintsVirial, nthreads);
1211 snew(bSettleErrorHasOccurred, nthreads);
1216 char* env = getenv("GMX_MAXCONSTRWARN");
1220 sscanf(env, "%8d", &maxwarn);
1227 fprintf(log, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1231 fprintf(stderr, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1234 warncount_lincs = 0;
1235 warncount_settle = 0;
1238 Constraints::Impl::~Impl()
1240 if (bSettleErrorHasOccurred != nullptr)
1242 sfree(bSettleErrorHasOccurred);
1244 if (threadConstraintsVirial != nullptr)
1246 sfree(threadConstraintsVirial);
1251 void Constraints::saveEdsamPointer(gmx_edsam* ed)
1256 ArrayRef<const ListOfLists<int>> Constraints::atom2constraints_moltype() const
1258 return impl_->at2con_mt;
1261 ArrayRef<const std::vector<int>> Constraints::atom2settle_moltype() const
1263 return impl_->at2settle_mt;
1266 void do_constrain_first(FILE* fplog,
1267 gmx::Constraints* constr,
1268 const t_inputrec* ir,
1271 ArrayRefWithPadding<RVec> x,
1272 ArrayRefWithPadding<RVec> v,
1276 int i, m, start, end;
1278 real dt = ir->delta_t;
1281 PaddedVector<RVec> savex(numAtoms);
1288 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n", start, numHomeAtoms, end);
1290 /* Do a first constrain to reset particles... */
1291 step = ir->init_step;
1294 char buf[STEPSTRSIZE];
1295 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n", gmx_step_str(step, buf));
1299 bool needsLogging = true;
1300 bool computeEnergy = false;
1301 bool computeVirial = false;
1302 /* constrain the current position */
1303 constr->apply(needsLogging,
1317 gmx::ConstraintVariable::Positions);
1320 /* constrain the inital velocity, and save it */
1321 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1322 constr->apply(needsLogging,
1329 v.unpaddedArrayRef(),
1336 gmx::ConstraintVariable::Velocities);
1338 /* constrain the inital velocities at t-dt/2 */
1339 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != IntegrationAlgorithm::VV)
1341 auto subX = x.paddedArrayRef().subArray(start, end);
1342 auto subV = v.paddedArrayRef().subArray(start, end);
1343 for (i = start; (i < end); i++)
1345 for (m = 0; (m < DIM); m++)
1347 /* Reverse the velocity */
1348 subV[i][m] = -subV[i][m];
1349 /* Store the position at t-dt in buf */
1350 savex[i][m] = subX[i][m] + dt * subV[i][m];
1353 /* Shake the positions at t=-dt with the positions at t=0
1354 * as reference coordinates.
1358 char buf[STEPSTRSIZE];
1359 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n", gmx_step_str(step, buf));
1362 constr->apply(needsLogging,
1368 savex.arrayRefWithPadding(),
1376 gmx::ConstraintVariable::Positions);
1378 for (i = start; i < end; i++)
1380 for (m = 0; m < DIM; m++)
1382 /* Re-reverse the velocities */
1383 subV[i][m] = -subV[i][m];
1389 void constrain_velocities(gmx::Constraints* constr,
1395 gmx_bool computeVirial,
1396 tensor constraintsVirial)
1398 if (constr != nullptr)
1400 constr->apply(do_log,
1405 state->x.arrayRefWithPadding(),
1406 state->v.arrayRefWithPadding(),
1407 state->v.arrayRefWithPadding().unpaddedArrayRef(),
1409 state->lambda[FreeEnergyPerturbationCouplingType::Bonded],
1411 ArrayRefWithPadding<RVec>(),
1414 ConstraintVariable::Velocities);
1418 void constrain_coordinates(gmx::Constraints* constr,
1423 ArrayRefWithPadding<RVec> xp,
1425 gmx_bool computeVirial,
1426 tensor constraintsVirial)
1428 if (constr != nullptr)
1430 constr->apply(do_log,
1435 state->x.arrayRefWithPadding(),
1439 state->lambda[FreeEnergyPerturbationCouplingType::Bonded],
1441 state->v.arrayRefWithPadding(),
1444 ConstraintVariable::Positions);