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/arrayrefwithpadding.h"
63 #include "gromacs/math/utilities.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdlib/gmx_omp_nthreads.h"
66 #include "gromacs/mdlib/lincs.h"
67 #include "gromacs/mdlib/settle.h"
68 #include "gromacs/mdlib/shake.h"
69 #include "gromacs/mdtypes/commrec.h"
70 #include "gromacs/mdtypes/inputrec.h"
71 #include "gromacs/mdtypes/md_enums.h"
72 #include "gromacs/mdtypes/state.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/pulling/pull.h"
75 #include "gromacs/timing/wallcycle.h"
76 #include "gromacs/topology/ifunc.h"
77 #include "gromacs/topology/mtop_lookup.h"
78 #include "gromacs/topology/mtop_util.h"
79 #include "gromacs/utility/arrayref.h"
80 #include "gromacs/utility/exceptions.h"
81 #include "gromacs/utility/fatalerror.h"
82 #include "gromacs/utility/gmxassert.h"
83 #include "gromacs/utility/listoflists.h"
84 #include "gromacs/utility/pleasecite.h"
85 #include "gromacs/utility/txtdump.h"
90 /* \brief Impl class for Constraints
92 * \todo Members like md, idef are valid only for the lifetime of a
93 * domain, which would be good to make clearer in the structure of the
94 * code. It should not be possible to call apply() if setConstraints()
95 * has not been called. For example, this could be achieved if
96 * setConstraints returned a valid object with such a method. That
97 * still requires that we manage the lifetime of that object
98 * correctly, however. */
99 class Constraints::Impl
102 Impl(const gmx_mtop_t& mtop_p,
103 const t_inputrec& ir_p,
106 const t_commrec* cr_p,
107 bool useUpdateGroups,
108 const gmx_multisim_t* ms,
110 gmx_wallcycle* wcycle_p,
111 bool pbcHandlingRequired,
115 void setConstraints(gmx_localtop_t* top,
118 gmx::ArrayRef<const real> masses,
119 gmx::ArrayRef<const real> inverseMasses,
120 bool hasMassPerturbedAtoms,
122 gmx::ArrayRef<const unsigned short> cFREEZE);
123 bool apply(bool bLog,
128 ArrayRefWithPadding<RVec> x,
129 ArrayRefWithPadding<RVec> xprime,
130 ArrayRef<RVec> min_proj,
134 ArrayRefWithPadding<RVec> v,
136 tensor constraintsVirial,
137 ConstraintVariable econq);
138 //! The total number of constraints.
140 //! The number of flexible constraints.
142 //! A list of atoms to constraints for each moleculetype.
143 std::vector<ListOfLists<int>> at2con_mt;
144 //! A list of atoms to settles for each moleculetype
145 std::vector<std::vector<int>> at2settle_mt;
147 Lincs* lincsd = nullptr; // TODO this should become a unique_ptr
149 std::unique_ptr<shakedata> shaked;
151 std::unique_ptr<SettleData> settled;
152 //! The maximum number of warnings.
154 //! The number of warnings for LINCS.
155 int warncount_lincs = 0;
156 //! The number of warnings for SETTLE.
157 int warncount_settle = 0;
158 //! The essential dynamics data.
159 gmx_edsam* ed = nullptr;
161 //! Thread-local virial contribution.
162 tensor* threadConstraintsVirial = { nullptr };
163 //! Did a settle error occur?
164 bool* bSettleErrorHasOccurred = nullptr;
166 //! Pointer to the global topology - only used for printing warnings.
167 const gmx_mtop_t& mtop;
168 //! Parameters for the interactions in this domain.
169 const InteractionDefinitions* idef = nullptr;
170 //! Total number of atoms.
172 //! Number of local atoms.
173 int numHomeAtoms_ = 0;
175 gmx::ArrayRef<const real> masses_;
177 gmx::ArrayRef<const real> inverseMasses_;
178 //! If there are atoms with perturbed mass (for FEP).
179 bool hasMassPerturbedAtoms_;
180 //! FEP lambda value.
182 //! Freeze groups data
183 gmx::ArrayRef<const unsigned short> cFREEZE_;
184 //! Whether we need to do pbc for handling bonds.
185 bool pbcHandlingRequired_ = false;
189 //! Communication support.
190 const t_commrec* cr = nullptr;
191 //! Multi-sim support.
192 const gmx_multisim_t* ms = nullptr;
193 //! Pulling code object, if any.
194 pull_t* pull_work = nullptr;
195 /*!\brief Input options.
197 * \todo Replace with IMdpOptions */
198 const t_inputrec& ir;
199 //! Flop counting support.
200 t_nrnb* nrnb = nullptr;
201 //! Tracks wallcycle usage.
202 gmx_wallcycle* wcycle;
205 Constraints::~Constraints() = default;
207 int Constraints::numFlexibleConstraints() const
209 return impl_->nflexcon;
212 bool Constraints::havePerturbedConstraints() const
214 const gmx_ffparams_t& ffparams = impl_->mtop.ffparams;
216 for (size_t i = 0; i < ffparams.functype.size(); i++)
218 if ((ffparams.functype[i] == F_CONSTR || ffparams.functype[i] == F_CONSTRNC)
219 && ffparams.iparams[i].constr.dA != ffparams.iparams[i].constr.dB)
228 //! Clears constraint quantities for atoms in nonlocal region.
229 static void clear_constraint_quantity_nonlocal(const gmx_domdec_t& dd, ArrayRef<RVec> q)
231 int nonlocal_at_start, nonlocal_at_end;
233 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
235 for (int at = nonlocal_at_start; at < nonlocal_at_end; at++)
241 void too_many_constraint_warnings(ConstraintAlgorithm eConstrAlg, int warncount)
244 "Too many %s warnings (%d)\n"
245 "If you know what you are doing you can %s"
246 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
247 "but normally it is better to fix the problem",
248 (eConstrAlg == ConstraintAlgorithm::Lincs) ? "LINCS" : "SETTLE",
250 (eConstrAlg == ConstraintAlgorithm::Lincs)
251 ? "adjust the lincs warning threshold in your mdp file\nor "
255 //! Writes out coordinates.
256 static void write_constr_pdb(const char* fn,
258 const gmx_mtop_t& mtop,
262 ArrayRef<const RVec> x,
267 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
269 const char * anm, *resnm;
272 if (DOMAINDECOMP(cr))
275 dd_get_constraint_range(*dd, &dd_ac0, &dd_ac1);
282 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
286 sprintf(fname, "%s.pdb", fn);
289 out = gmx_fio_fopen(fname, "w");
291 fprintf(out, "TITLE %s\n", title);
292 gmx_write_pdb_box(out, PbcType::Unset, box);
294 for (i = start; i < start + homenr; i++)
298 if (i >= dd_numHomeAtoms(*dd) && i < dd_ac0)
302 ii = dd->globalAtomIndices[i];
308 mtopGetAtomAndResidueName(mtop, ii, &molb, &anm, &resnr, &resnm, nullptr);
309 gmx_fprintf_pdb_atomline(out,
325 fprintf(out, "TER\n");
330 //! Writes out domain contents to help diagnose crashes.
331 static void dump_confs(FILE* log,
333 const gmx_mtop_t& mtop,
337 ArrayRef<const RVec> x,
338 ArrayRef<const RVec> xprime,
341 char buf[STRLEN], buf2[22];
343 char* env = getenv("GMX_SUPPRESS_DUMP");
349 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
350 write_constr_pdb(buf, "initial coordinates", mtop, start, homenr, cr, x, box);
351 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
352 write_constr_pdb(buf, "coordinates after constraining", mtop, start, homenr, cr, xprime, box);
355 fprintf(log, "Wrote pdb files with previous and current coordinates\n");
357 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
360 bool Constraints::apply(bool bLog,
365 ArrayRefWithPadding<RVec> x,
366 ArrayRefWithPadding<RVec> xprime,
367 ArrayRef<RVec> min_proj,
371 ArrayRefWithPadding<RVec> v,
373 tensor constraintsVirial,
374 ConstraintVariable econq)
376 return impl_->apply(bLog,
393 bool Constraints::Impl::apply(bool bLog,
398 ArrayRefWithPadding<RVec> x,
399 ArrayRefWithPadding<RVec> xprime,
400 ArrayRef<RVec> min_proj,
404 ArrayRefWithPadding<RVec> v,
406 tensor constraintsVirial,
407 ConstraintVariable econq)
412 real invdt, vir_fac = 0, t;
414 t_pbc pbc, *pbc_null;
418 wallcycle_start(wcycle, WallCycleCounter::Constr);
420 if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
423 "constrain called for forces displacements while not doing energy minimization, "
424 "can not do this while the LINCS and SETTLE constraint connection matrices are "
433 scaled_delta_t = step_scaling * ir.delta_t;
435 /* Prepare time step for use in constraint implementations, and
436 avoid generating inf when ir.delta_t = 0. */
443 invdt = 1.0 / scaled_delta_t;
446 if (ir.efep != FreeEnergyPerturbationType::No && EI_DYNAMICS(ir.eI))
448 /* Set the constraint lengths for the step at which this configuration
449 * is meant to be. The invmasses should not be changed.
451 lambda += delta_step * ir.fepvals->delta_lambda;
456 clear_mat(constraintsVirial);
458 const InteractionList& settle = idef->il[F_SETTLE];
459 nsettle = settle.size() / (1 + NRAL(F_SETTLE));
463 nth = gmx_omp_nthreads_get(ModuleMultiThread::Settle);
470 /* We do not need full pbc when constraints do not cross update groups
471 * i.e. when dd->constraint_comm==NULL.
472 * Note that PBC for constraints is different from PBC for bondeds.
473 * For constraints there is both forward and backward communication.
475 if (ir.pbcType != PbcType::No && (cr->dd || pbcHandlingRequired_)
476 && !(cr->dd && cr->dd->constraint_comm == nullptr))
478 /* With pbc=screw the screw has been changed to a shift
479 * by the constraint coordinate communication routine,
480 * so that here we can use normal pbc.
482 pbc_null = set_pbc_dd(&pbc, ir.pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, FALSE, box);
489 /* Communicate the coordinates required for the non-local constraints
490 * for LINCS and/or SETTLE.
494 dd_move_x_constraints(cr->dd,
496 x.unpaddedArrayRef(),
497 xprime.unpaddedArrayRef(),
498 econq == ConstraintVariable::Positions);
502 /* We need to initialize the non-local components of v.
503 * We never actually use these values, but we do increment them,
504 * so we should avoid uninitialized variables and overflows.
506 clear_constraint_quantity_nonlocal(*cr->dd, v.unpaddedArrayRef());
510 if (lincsd != nullptr)
512 bOK = constrain_lincs(bLog || bEner,
524 hasMassPerturbedAtoms_,
528 v.unpaddedArrayRef(),
535 if (!bOK && maxwarn < INT_MAX)
540 "Constraint error in algorithm %s at step %s\n",
541 enumValueToString(ConstraintAlgorithm::Lincs),
542 gmx_step_str(step, buf));
548 if (shaked != nullptr)
550 bOK = constrain_shake(log,
555 x.unpaddedArrayRef(),
556 xprime.unpaddedArrayRef(),
563 v.unpaddedArrayRef(),
569 if (!bOK && maxwarn < INT_MAX)
574 "Constraint error in algorithm %s at step %s\n",
575 enumValueToString(ConstraintAlgorithm::Shake),
576 gmx_step_str(step, buf));
584 bool bSettleErrorHasOccurred0 = false;
588 case ConstraintVariable::Positions:
589 #pragma omp parallel for num_threads(nth) schedule(static)
590 for (int th = 0; th < nth; th++)
596 clear_mat(threadConstraintsVirial[th]);
608 th == 0 ? constraintsVirial : threadConstraintsVirial[th],
609 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
611 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
613 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
616 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle * 3);
620 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle * 3);
623 case ConstraintVariable::Velocities:
624 case ConstraintVariable::Derivative:
625 case ConstraintVariable::Force:
626 case ConstraintVariable::ForceDispl:
627 #pragma omp parallel for num_threads(nth) schedule(static)
628 for (int th = 0; th < nth; th++)
632 int calcvir_atom_end;
636 calcvir_atom_end = 0;
640 calcvir_atom_end = numHomeAtoms_;
645 clear_mat(threadConstraintsVirial[th]);
648 int start_th = (nsettle * th) / nth;
649 int end_th = (nsettle * (th + 1)) / nth;
651 if (start_th >= 0 && end_th - start_th > 0)
653 settle_proj(*settled,
656 settle.iatoms.data() + start_th * (1 + NRAL(F_SETTLE)),
658 x.unpaddedArrayRef(),
659 xprime.unpaddedArrayRef(),
662 th == 0 ? constraintsVirial : threadConstraintsVirial[th]);
665 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
667 /* This is an overestimate */
668 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
670 case ConstraintVariable::Deriv_FlexCon:
671 /* Nothing to do, since the are no flexible constraints in settles */
673 default: gmx_incons("Unknown constraint quantity for settle");
678 /* Reduce the virial contributions over the threads */
679 for (int th = 1; th < nth; th++)
681 m_add(constraintsVirial, threadConstraintsVirial[th], constraintsVirial);
685 if (econq == ConstraintVariable::Positions)
687 for (int th = 1; th < nth; th++)
689 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
692 if (bSettleErrorHasOccurred0)
698 ": One or more water molecules can not be settled.\n"
699 "Check for bad contacts and/or reduce the timestep if appropriate.\n",
703 fprintf(log, "%s", buf);
705 fprintf(stderr, "%s", buf);
707 if (warncount_settle > maxwarn)
709 too_many_constraint_warnings(ConstraintAlgorithm::Count, warncount_settle);
720 /* The normal uses of constrain() pass step_scaling = 1.0.
721 * The call to constrain() for SD1 that passes step_scaling =
722 * 0.5 also passes vir = NULL, so cannot reach this
723 * assertion. This assertion should remain until someone knows
724 * that this path works for their intended purpose, and then
725 * they can use scaled_delta_t instead of ir.delta_t
727 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
730 case ConstraintVariable::Positions: vir_fac = 0.5 / (ir.delta_t * ir.delta_t); break;
731 case ConstraintVariable::Velocities: vir_fac = 0.5 / ir.delta_t; break;
732 case ConstraintVariable::Force:
733 case ConstraintVariable::ForceDispl: vir_fac = 0.5; break;
734 default: gmx_incons("Unsupported constraint quantity for virial");
739 vir_fac *= 2; /* only constraining over half the distance here */
741 for (int i = 0; i < DIM; i++)
743 for (int j = 0; j < DIM; j++)
745 constraintsVirial[i][j] *= vir_fac;
752 dump_confs(log, step, mtop, start, numHomeAtoms_, cr, x.unpaddedArrayRef(), xprime.unpaddedArrayRef(), box);
755 if (econq == ConstraintVariable::Positions)
757 if (ir.bPull && pull_have_constraint(*pull_work))
759 if (EI_DYNAMICS(ir.eI))
761 t = ir.init_t + (step + delta_step) * ir.delta_t;
767 set_pbc(&pbc, ir.pbcType, box);
768 pull_constraint(pull_work,
774 x.unpaddedArrayRef(),
775 xprime.unpaddedArrayRef(),
776 v.unpaddedArrayRef(),
779 if (ed && delta_step > 0)
781 /* apply the essential dynamics constraints here */
782 do_edsam(&ir, step, cr, xprime.unpaddedArrayRef(), v.unpaddedArrayRef(), box, ed);
785 wallcycle_stop(wcycle, WallCycleCounter::Constr);
787 const bool haveVelocities = (!v.empty() || econq == ConstraintVariable::Velocities);
788 if (haveVelocities && !cFREEZE_.empty())
790 /* Set the velocities of frozen dimensions to zero */
792 if (econq == ConstraintVariable::Velocities)
794 vRef = xprime.unpaddedArrayRef();
798 vRef = v.unpaddedArrayRef();
801 int gmx_unused numThreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
803 #pragma omp parallel for num_threads(numThreads) schedule(static)
804 for (int i = 0; i < numHomeAtoms_; i++)
806 int freezeGroup = cFREEZE_[i];
808 for (int d = 0; d < DIM; d++)
810 if (ir.opts.nFreeze[freezeGroup][d])
821 ArrayRef<real> Constraints::rmsdData() const
825 return lincs_rmsdData(impl_->lincsd);
833 real Constraints::rmsd() const
837 return lincs_rmsd(impl_->lincsd);
845 int Constraints::numConstraintsTotal()
847 return impl_->ncon_tot;
850 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
852 if (haveDynamicsIntegrator)
854 return FlexibleConstraintTreatment::Include;
858 return FlexibleConstraintTreatment::Exclude;
862 /*! \brief Returns a block struct to go from atoms to constraints
864 * The block struct will contain constraint indices with lower indices
865 * directly matching the order in F_CONSTR and higher indices matching
866 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
868 * \param[in] numAtoms The number of atoms to construct the list for
869 * \param[in] ilists The interaction lists, size F_NRE
870 * \param[in] iparams Interaction parameters, can be null when
871 * \p flexibleConstraintTreatment==Include
872 * \param[in] flexibleConstraintTreatment The flexible constraint treatment,
875 * \returns a block struct with all constraints for each atom
877 static ListOfLists<int> makeAtomsToConstraintsList(int numAtoms,
878 ArrayRef<const InteractionList> ilists,
879 ArrayRef<const t_iparams> iparams,
880 FlexibleConstraintTreatment flexibleConstraintTreatment)
882 GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || !iparams.empty(),
883 "With flexible constraint detection we need valid iparams");
885 std::vector<int> count(numAtoms);
887 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
889 const InteractionList& ilist = ilists[ftype];
890 const int stride = 1 + NRAL(ftype);
891 for (int i = 0; i < ilist.size(); i += stride)
893 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
894 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
896 for (int j = 1; j < 3; j++)
898 int a = ilist.iatoms[i + j];
905 std::vector<int> listRanges(numAtoms + 1);
906 for (int a = 0; a < numAtoms; a++)
908 listRanges[a + 1] = listRanges[a] + count[a];
911 std::vector<int> elements(listRanges[numAtoms]);
913 /* The F_CONSTRNC constraints have constraint numbers
914 * that continue after the last F_CONSTR constraint.
916 int numConstraints = 0;
917 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
919 const InteractionList& ilist = ilists[ftype];
920 const int stride = 1 + NRAL(ftype);
921 for (int i = 0; i < ilist.size(); i += stride)
923 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
924 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
926 for (int j = 1; j < 3; j++)
928 const int a = ilist.iatoms[i + j];
929 elements[listRanges[a] + count[a]++] = numConstraints;
936 return ListOfLists<int>(std::move(listRanges), std::move(elements));
939 ListOfLists<int> make_at2con(int numAtoms,
940 ArrayRef<const InteractionList> ilist,
941 ArrayRef<const t_iparams> iparams,
942 FlexibleConstraintTreatment flexibleConstraintTreatment)
944 return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
947 ListOfLists<int> make_at2con(const gmx_moltype_t& moltype,
948 gmx::ArrayRef<const t_iparams> iparams,
949 FlexibleConstraintTreatment flexibleConstraintTreatment)
951 return makeAtomsToConstraintsList(
952 moltype.atoms.nr, makeConstArrayRef(moltype.ilist), iparams, flexibleConstraintTreatment);
955 //! Return the number of flexible constraints in the \c ilist and \c iparams.
956 int countFlexibleConstraints(ArrayRef<const InteractionList> ilist, ArrayRef<const t_iparams> iparams)
959 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
961 const int numIatomsPerConstraint = 3;
962 for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
964 const int type = ilist[ftype].iatoms[i];
965 if (iparams[type].constr.dA == 0 && iparams[type].constr.dB == 0)
975 //! Returns the index of the settle to which each atom belongs.
976 static std::vector<int> make_at2settle(int natoms, const InteractionList& ilist)
978 /* Set all to no settle */
979 std::vector<int> at2s(natoms, -1);
981 const int stride = 1 + NRAL(F_SETTLE);
983 for (int s = 0; s < ilist.size(); s += stride)
985 at2s[ilist.iatoms[s + 1]] = s / stride;
986 at2s[ilist.iatoms[s + 2]] = s / stride;
987 at2s[ilist.iatoms[s + 3]] = s / stride;
993 void Constraints::Impl::setConstraints(gmx_localtop_t* top,
996 gmx::ArrayRef<const real> masses,
997 gmx::ArrayRef<const real> inverseMasses,
998 const bool hasMassPerturbedAtoms,
1000 gmx::ArrayRef<const unsigned short> cFREEZE)
1002 numAtoms_ = numAtoms;
1003 numHomeAtoms_ = numHomeAtoms;
1005 inverseMasses_ = inverseMasses;
1006 hasMassPerturbedAtoms_ = hasMassPerturbedAtoms;
1014 /* With DD we might also need to call LINCS on a domain no constraints for
1015 * communicating coordinates to other nodes that do have constraints.
1017 if (ir.eConstrAlg == ConstraintAlgorithm::Lincs)
1019 set_lincs(*idef, numAtoms_, inverseMasses_, lambda_, EI_DYNAMICS(ir.eI), cr, lincsd);
1021 if (ir.eConstrAlg == ConstraintAlgorithm::Shake)
1025 // We are using the local topology, so there are only
1026 // F_CONSTR constraints.
1027 GMX_RELEASE_ASSERT(idef->il[F_CONSTRNC].empty(),
1028 "Here we should not have no-connect constraints");
1029 make_shake_sblock_dd(shaked.get(), idef->il[F_CONSTR]);
1033 make_shake_sblock_serial(shaked.get(), &top->idef, numAtoms_);
1040 settled->setConstraints(idef->il[F_SETTLE], numHomeAtoms_, masses_, inverseMasses_);
1043 /* Make a selection of the local atoms for essential dynamics */
1046 dd_make_local_ed_indices(cr->dd, ed);
1050 void Constraints::setConstraints(gmx_localtop_t* top,
1052 const int numHomeAtoms,
1053 gmx::ArrayRef<const real> masses,
1054 gmx::ArrayRef<const real> inverseMasses,
1055 const bool hasMassPerturbedAtoms,
1057 gmx::ArrayRef<const unsigned short> cFREEZE)
1059 impl_->setConstraints(
1060 top, numAtoms, numHomeAtoms, masses, inverseMasses, hasMassPerturbedAtoms, lambda, cFREEZE);
1063 /*! \brief Makes a per-moleculetype container of mappings from atom
1064 * indices to constraint indices.
1066 * Note that flexible constraints are only enabled with a dynamical integrator. */
1067 static std::vector<ListOfLists<int>> makeAtomToConstraintMappings(const gmx_mtop_t& mtop,
1068 FlexibleConstraintTreatment flexibleConstraintTreatment)
1070 std::vector<ListOfLists<int>> mapping;
1071 mapping.reserve(mtop.moltype.size());
1072 for (const gmx_moltype_t& moltype : mtop.moltype)
1074 mapping.push_back(make_at2con(moltype, mtop.ffparams.iparams, flexibleConstraintTreatment));
1079 Constraints::Constraints(const gmx_mtop_t& mtop,
1080 const t_inputrec& ir,
1083 const t_commrec* cr,
1084 const bool useUpdateGroups,
1085 const gmx_multisim_t* ms,
1087 gmx_wallcycle* wcycle,
1088 bool pbcHandlingRequired,
1091 impl_(new Impl(mtop, ir, pull_work, log, cr, useUpdateGroups, ms, nrnb, wcycle, pbcHandlingRequired, numConstraints, numSettles))
1095 Constraints::Impl::Impl(const gmx_mtop_t& mtop_p,
1096 const t_inputrec& ir_p,
1099 const t_commrec* cr_p,
1100 const bool useUpdateGroups,
1101 const gmx_multisim_t* ms_p,
1103 gmx_wallcycle* wcycle_p,
1104 bool pbcHandlingRequired,
1107 ncon_tot(numConstraints),
1109 pbcHandlingRequired_(pbcHandlingRequired),
1113 pull_work(pull_work),
1118 if (numConstraints + numSettles > 0 && ir.epc == PressureCoupling::Mttk)
1120 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1124 if (numConstraints > 0)
1126 at2con_mt = makeAtomToConstraintMappings(mtop, flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
1128 for (const gmx_molblock_t& molblock : mtop.molblock)
1130 int count = countFlexibleConstraints(mtop.moltype[molblock.type].ilist, mtop.ffparams.iparams);
1131 nflexcon += molblock.nmol * count;
1138 fprintf(log, "There are %d flexible constraints\n", nflexcon);
1139 if (ir.fc_stepsize == 0)
1143 "WARNING: step size for flexible constraining = 0\n"
1144 " All flexible constraints will be rigid.\n"
1145 " Will try to keep all flexible constraints at their original "
1147 " but the lengths may exhibit some drift.\n\n");
1153 please_cite(log, "Hess2002");
1157 // When there are multiple PP domains and update groups are
1158 // not in use, the constraints might be split across the
1159 // domains, needing particular handling.
1160 const bool mayHaveSplitConstraints = DOMAINDECOMP(cr) && !useUpdateGroups;
1162 if (ir.eConstrAlg == ConstraintAlgorithm::Lincs)
1164 lincsd = init_lincs(
1165 log, mtop, nflexcon, at2con_mt, mayHaveSplitConstraints, ir.nLincsIter, ir.nProjOrder);
1168 if (ir.eConstrAlg == ConstraintAlgorithm::Shake)
1170 if (mayHaveSplitConstraints)
1173 "SHAKE is not supported with domain decomposition and constraints that "
1174 "cross domain boundaries, use LINCS");
1179 "For this system also velocities and/or forces need to be constrained, "
1180 "this can not be done with SHAKE, you should select LINCS");
1182 please_cite(log, "Ryckaert77a");
1185 please_cite(log, "Barth95a");
1188 shaked = std::make_unique<shakedata>();
1194 please_cite(log, "Miyamoto92a");
1196 settled = std::make_unique<SettleData>(mtop);
1198 // SETTLE with perturbed masses is not implemented. grompp now checks
1199 // for this, but old .tpr files that did this might still exist.
1200 if (haveFepPerturbedMassesInSettles(mtop))
1203 "SETTLE is not implemented for atoms whose mass is perturbed. "
1204 "You might\ninstead use normal constraints.");
1207 /* Make an atom to settle index for use in domain decomposition */
1208 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1210 at2settle_mt.emplace_back(
1211 make_at2settle(mtop.moltype[mt].atoms.nr, mtop.moltype[mt].ilist[F_SETTLE]));
1214 /* Allocate thread-local work arrays */
1215 int nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Settle);
1216 if (nthreads > 1 && threadConstraintsVirial == nullptr)
1218 snew(threadConstraintsVirial, nthreads);
1219 snew(bSettleErrorHasOccurred, nthreads);
1224 char* env = getenv("GMX_MAXCONSTRWARN");
1228 sscanf(env, "%8d", &maxwarn);
1235 fprintf(log, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1239 fprintf(stderr, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1242 warncount_lincs = 0;
1243 warncount_settle = 0;
1246 Constraints::Impl::~Impl()
1248 if (bSettleErrorHasOccurred != nullptr)
1250 sfree(bSettleErrorHasOccurred);
1252 if (threadConstraintsVirial != nullptr)
1254 sfree(threadConstraintsVirial);
1259 void Constraints::saveEdsamPointer(gmx_edsam* ed)
1264 ArrayRef<const ListOfLists<int>> Constraints::atom2constraints_moltype() const
1266 return impl_->at2con_mt;
1269 ArrayRef<const std::vector<int>> Constraints::atom2settle_moltype() const
1271 return impl_->at2settle_mt;
1274 void do_constrain_first(FILE* fplog,
1275 gmx::Constraints* constr,
1276 const t_inputrec* ir,
1279 ArrayRefWithPadding<RVec> x,
1280 ArrayRefWithPadding<RVec> v,
1284 int i, m, start, end;
1286 real dt = ir->delta_t;
1289 PaddedVector<RVec> savex(numAtoms);
1296 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n", start, numHomeAtoms, end);
1298 /* Do a first constrain to reset particles... */
1299 step = ir->init_step;
1302 char buf[STEPSTRSIZE];
1303 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n", gmx_step_str(step, buf));
1307 bool needsLogging = true;
1308 bool computeEnergy = false;
1309 bool computeVirial = false;
1310 /* constrain the current position */
1311 constr->apply(needsLogging,
1325 gmx::ConstraintVariable::Positions);
1328 /* constrain the inital velocity, and save it */
1329 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1330 constr->apply(needsLogging,
1337 v.unpaddedArrayRef(),
1344 gmx::ConstraintVariable::Velocities);
1346 /* constrain the inital velocities at t-dt/2 */
1347 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != IntegrationAlgorithm::VV)
1349 auto subX = x.paddedArrayRef().subArray(start, end);
1350 auto subV = v.paddedArrayRef().subArray(start, end);
1351 for (i = start; (i < end); i++)
1353 for (m = 0; (m < DIM); m++)
1355 /* Reverse the velocity */
1356 subV[i][m] = -subV[i][m];
1357 /* Store the position at t-dt in buf */
1358 savex[i][m] = subX[i][m] + dt * subV[i][m];
1361 /* Shake the positions at t=-dt with the positions at t=0
1362 * as reference coordinates.
1366 char buf[STEPSTRSIZE];
1367 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n", gmx_step_str(step, buf));
1370 constr->apply(needsLogging,
1376 savex.arrayRefWithPadding(),
1384 gmx::ConstraintVariable::Positions);
1386 for (i = start; i < end; i++)
1388 for (m = 0; m < DIM; m++)
1390 /* Re-reverse the velocities */
1391 subV[i][m] = -subV[i][m];
1397 void constrain_velocities(gmx::Constraints* constr,
1403 gmx_bool computeVirial,
1404 tensor constraintsVirial)
1406 if (constr != nullptr)
1408 constr->apply(do_log,
1413 state->x.arrayRefWithPadding(),
1414 state->v.arrayRefWithPadding(),
1415 state->v.arrayRefWithPadding().unpaddedArrayRef(),
1417 state->lambda[FreeEnergyPerturbationCouplingType::Bonded],
1419 ArrayRefWithPadding<RVec>(),
1422 ConstraintVariable::Velocities);
1426 void constrain_coordinates(gmx::Constraints* constr,
1431 ArrayRefWithPadding<RVec> xp,
1433 gmx_bool computeVirial,
1434 tensor constraintsVirial)
1436 if (constr != nullptr)
1438 constr->apply(do_log,
1443 state->x.arrayRefWithPadding(),
1447 state->lambda[FreeEnergyPerturbationCouplingType::Bonded],
1449 state->v.arrayRefWithPadding(),
1452 ConstraintVariable::Positions);