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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief Defines the high-level constraint code.
40 * \author Berk Hess <hess@kth.se>
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_mdlib
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/essentialdynamics/edsam.h"
57 #include "gromacs/fileio/confio.h"
58 #include "gromacs/fileio/gmxfio.h"
59 #include "gromacs/fileio/pdbio.h"
60 #include "gromacs/gmxlib/nrnb.h"
61 #include "gromacs/math/utilities.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdlib/lincs.h"
65 #include "gromacs/mdlib/settle.h"
66 #include "gromacs/mdlib/shake.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/mdatom.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/block.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/pleasecite.h"
84 #include "gromacs/utility/smalloc.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_mdatoms &md_p,
107 const t_commrec *cr_p,
108 const gmx_multisim_t *ms,
110 gmx_wallcycle *wcycle_p,
111 bool pbcHandlingRequired,
115 void setConstraints(const gmx_localtop_t &top,
116 const t_mdatoms &md);
117 bool apply(bool bLog,
130 ConstraintVariable econq);
131 //! The total number of constraints.
133 //! The number of flexible constraints.
135 //! A list of atoms to constraints for each moleculetype.
136 std::vector<t_blocka> at2con_mt;
137 //! A list of atoms to settles for each moleculetype
138 std::vector < std::vector < int>> at2settle_mt;
139 //! Whether any SETTLES cross charge-group boundaries.
140 bool bInterCGsettles = false;
142 Lincs *lincsd = nullptr; // TODO this should become a unique_ptr
144 shakedata *shaked = nullptr;
146 settledata *settled = nullptr;
147 //! The maximum number of warnings.
149 //! The number of warnings for LINCS.
150 int warncount_lincs = 0;
151 //! The number of warnings for SETTLE.
152 int warncount_settle = 0;
153 //! The essential dynamics data.
154 gmx_edsam * ed = nullptr;
156 //! Thread-local virial contribution.
157 tensor *vir_r_m_dr_th = {nullptr};
158 //! Did a settle error occur?
159 bool *bSettleErrorHasOccurred = nullptr;
161 //! Pointer to the global topology - only used for printing warnings.
162 const gmx_mtop_t &mtop;
163 //! Parameters for the interactions in this domain.
164 const t_idef *idef = nullptr;
165 //! Data about atoms in this domain.
167 //! Whether we need to do pbc for handling bonds.
168 bool pbcHandlingRequired_ = false;
172 //! Communication support.
173 const t_commrec *cr = nullptr;
174 //! Multi-sim support.
175 const gmx_multisim_t *ms = nullptr;
176 //! Pulling code object, if any.
177 pull_t *pull_work = nullptr;
178 /*!\brief Input options.
180 * \todo Replace with IMdpOptions */
181 const t_inputrec &ir;
182 //! Flop counting support.
183 t_nrnb *nrnb = nullptr;
184 //! Tracks wallcycle usage.
185 gmx_wallcycle *wcycle;
188 Constraints::~Constraints() = default;
190 int Constraints::numFlexibleConstraints() const
192 return impl_->nflexcon;
195 bool Constraints::havePerturbedConstraints() const
197 const gmx_ffparams_t &ffparams = impl_->mtop.ffparams;
199 for (size_t i = 0; i < ffparams.functype.size(); i++)
201 if ((ffparams.functype[i] == F_CONSTR ||
202 ffparams.functype[i] == F_CONSTRNC) &&
203 ffparams.iparams[i].constr.dA != ffparams.iparams[i].constr.dB)
212 //! Clears constraint quantities for atoms in nonlocal region.
213 static void clear_constraint_quantity_nonlocal(gmx_domdec_t *dd, rvec *q)
215 int nonlocal_at_start, nonlocal_at_end, at;
217 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
219 for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
225 void too_many_constraint_warnings(int eConstrAlg, int warncount)
228 "Too many %s warnings (%d)\n"
229 "If you know what you are doing you can %s"
230 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
231 "but normally it is better to fix the problem",
232 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
233 (eConstrAlg == econtLINCS) ?
234 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
237 //! Writes out coordinates.
238 static void write_constr_pdb(const char *fn, const char *title,
239 const gmx_mtop_t &mtop,
240 int start, int homenr, const t_commrec *cr,
241 const rvec x[], matrix box)
245 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
247 const char *anm, *resnm;
250 if (DOMAINDECOMP(cr))
253 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
260 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
264 sprintf(fname, "%s.pdb", fn);
267 out = gmx_fio_fopen(fname, "w");
269 fprintf(out, "TITLE %s\n", title);
270 gmx_write_pdb_box(out, -1, box);
272 for (i = start; i < start+homenr; i++)
276 if (i >= dd_numHomeAtoms(*dd) && i < dd_ac0)
280 ii = dd->globalAtomIndices[i];
286 mtopGetAtomAndResidueName(mtop, ii, &molb, &anm, &resnr, &resnm, nullptr);
287 gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
288 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
290 fprintf(out, "TER\n");
295 //! Writes out domain contents to help diagnose crashes.
296 static void dump_confs(FILE *log, int64_t step, const gmx_mtop_t &mtop,
297 int start, int homenr, const t_commrec *cr,
298 const rvec x[], rvec xprime[], matrix box)
300 char buf[STRLEN], buf2[22];
302 char *env = getenv("GMX_SUPPRESS_DUMP");
308 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
309 write_constr_pdb(buf, "initial coordinates",
310 mtop, start, homenr, cr, x, box);
311 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
312 write_constr_pdb(buf, "coordinates after constraining",
313 mtop, start, homenr, cr, xprime, box);
316 fprintf(log, "Wrote pdb files with previous and current coordinates\n");
318 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
322 Constraints::apply(bool bLog,
335 ConstraintVariable econq)
337 return impl_->apply(bLog,
354 Constraints::Impl::apply(bool bLog,
367 ConstraintVariable econq)
373 real invdt, vir_fac = 0, t;
375 t_pbc pbc, *pbc_null;
379 wallcycle_start(wcycle, ewcCONSTR);
381 if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
383 gmx_incons("constrain called for forces displacements while not doing energy minimization, can not do this while the LINCS and SETTLE constraint connection matrices are mass weighted");
392 scaled_delta_t = step_scaling * ir.delta_t;
394 /* Prepare time step for use in constraint implementations, and
395 avoid generating inf when ir.delta_t = 0. */
402 invdt = 1.0/scaled_delta_t;
405 if (ir.efep != efepNO && EI_DYNAMICS(ir.eI))
407 /* Set the constraint lengths for the step at which this configuration
408 * is meant to be. The invmasses should not be changed.
410 lambda += delta_step*ir.fepvals->delta_lambda;
415 clear_mat(vir_r_m_dr);
417 const t_ilist *settle = &idef->il[F_SETTLE];
418 nsettle = settle->nr/(1+NRAL(F_SETTLE));
422 nth = gmx_omp_nthreads_get(emntSETTLE);
429 /* We do not need full pbc when constraints do not cross charge groups,
430 * i.e. when dd->constraint_comm==NULL.
431 * Note that PBC for constraints is different from PBC for bondeds.
432 * For constraints there is both forward and backward communication.
434 if (ir.ePBC != epbcNONE &&
435 (cr->dd || pbcHandlingRequired_) && !(cr->dd && cr->dd->constraint_comm == nullptr))
437 /* With pbc=screw the screw has been changed to a shift
438 * by the constraint coordinate communication routine,
439 * so that here we can use normal pbc.
441 pbc_null = set_pbc_dd(&pbc, ir.ePBC,
442 DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
450 /* Communicate the coordinates required for the non-local constraints
451 * for LINCS and/or SETTLE.
455 dd_move_x_constraints(cr->dd, box, x, xprime, econq == ConstraintVariable::Positions);
459 /* We need to initialize the non-local components of v.
460 * We never actually use these values, but we do increment them,
461 * so we should avoid uninitialized variables and overflows.
463 clear_constraint_quantity_nonlocal(cr->dd, v);
467 if (lincsd != nullptr)
469 bOK = constrain_lincs(bLog || bEner, ir, step, lincsd, md, cr, ms,
471 box, pbc_null, lambda, dvdlambda,
472 invdt, v, vir != nullptr, vir_r_m_dr,
474 maxwarn, &warncount_lincs);
475 if (!bOK && maxwarn < INT_MAX)
479 fprintf(log, "Constraint error in algorithm %s at step %s\n",
480 econstr_names[econtLINCS], gmx_step_str(step, buf));
486 if (shaked != nullptr)
488 bOK = constrain_shake(log, shaked,
490 *idef, ir, x, xprime, min_proj, nrnb,
492 invdt, v, vir != nullptr, vir_r_m_dr,
493 maxwarn < INT_MAX, econq);
495 if (!bOK && maxwarn < INT_MAX)
499 fprintf(log, "Constraint error in algorithm %s at step %s\n",
500 econstr_names[econtSHAKE], gmx_step_str(step, buf));
508 bool bSettleErrorHasOccurred0 = false;
512 case ConstraintVariable::Positions:
513 #pragma omp parallel for num_threads(nth) schedule(static)
514 for (int th = 0; th < nth; th++)
520 clear_mat(vir_r_m_dr_th[th]);
527 invdt, v ? v[0] : nullptr,
529 th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th],
530 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
532 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
534 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
537 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
541 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
544 case ConstraintVariable::Velocities:
545 case ConstraintVariable::Derivative:
546 case ConstraintVariable::Force:
547 case ConstraintVariable::ForceDispl:
548 #pragma omp parallel for num_threads(nth) schedule(static)
549 for (int th = 0; th < nth; th++)
553 int calcvir_atom_end;
557 calcvir_atom_end = 0;
561 calcvir_atom_end = md.homenr;
566 clear_mat(vir_r_m_dr_th[th]);
569 int start_th = (nsettle* th )/nth;
570 int end_th = (nsettle*(th+1))/nth;
572 if (start_th >= 0 && end_th - start_th > 0)
574 settle_proj(settled, econq,
576 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
579 xprime, min_proj, calcvir_atom_end,
580 th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th]);
583 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
585 /* This is an overestimate */
586 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
588 case ConstraintVariable::Deriv_FlexCon:
589 /* Nothing to do, since the are no flexible constraints in settles */
592 gmx_incons("Unknown constraint quantity for settle");
597 /* Reduce the virial contributions over the threads */
598 for (int th = 1; th < nth; th++)
600 m_add(vir_r_m_dr, vir_r_m_dr_th[th], vir_r_m_dr);
604 if (econq == ConstraintVariable::Positions)
606 for (int th = 1; th < nth; th++)
608 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
611 if (bSettleErrorHasOccurred0)
615 "\nstep " "%" PRId64 ": One or more water molecules can not be settled.\n"
616 "Check for bad contacts and/or reduce the timestep if appropriate.\n",
620 fprintf(log, "%s", buf);
622 fprintf(stderr, "%s", buf);
624 if (warncount_settle > maxwarn)
626 too_many_constraint_warnings(-1, warncount_settle);
637 /* The normal uses of constrain() pass step_scaling = 1.0.
638 * The call to constrain() for SD1 that passes step_scaling =
639 * 0.5 also passes vir = NULL, so cannot reach this
640 * assertion. This assertion should remain until someone knows
641 * that this path works for their intended purpose, and then
642 * they can use scaled_delta_t instead of ir.delta_t
644 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
647 case ConstraintVariable::Positions:
648 vir_fac = 0.5/(ir.delta_t*ir.delta_t);
650 case ConstraintVariable::Velocities:
651 vir_fac = 0.5/ir.delta_t;
653 case ConstraintVariable::Force:
654 case ConstraintVariable::ForceDispl:
658 gmx_incons("Unsupported constraint quantity for virial");
663 vir_fac *= 2; /* only constraining over half the distance here */
665 for (int i = 0; i < DIM; i++)
667 for (int j = 0; j < DIM; j++)
669 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
676 dump_confs(log, step, mtop, start, homenr, cr, x, xprime, box);
679 if (econq == ConstraintVariable::Positions)
681 if (ir.bPull && pull_have_constraint(pull_work))
683 if (EI_DYNAMICS(ir.eI))
685 t = ir.init_t + (step + delta_step)*ir.delta_t;
691 set_pbc(&pbc, ir.ePBC, box);
692 pull_constraint(pull_work, &md, &pbc, cr, ir.delta_t, t, x, xprime, v, *vir);
694 if (ed && delta_step > 0)
696 /* apply the essential dynamics constraints here */
697 do_edsam(&ir, step, cr, xprime, v, box, ed);
700 wallcycle_stop(wcycle, ewcCONSTR);
702 if (v != nullptr && md.cFREEZE)
704 /* Set the velocities of frozen dimensions to zero */
706 int gmx_unused numThreads = gmx_omp_nthreads_get(emntUpdate);
708 #pragma omp parallel for num_threads(numThreads) schedule(static)
709 for (int i = 0; i < md.homenr; i++)
711 int freezeGroup = md.cFREEZE[i];
713 for (int d = 0; d < DIM; d++)
715 if (ir.opts.nFreeze[freezeGroup][d])
726 ArrayRef<real> Constraints::rmsdData() const
730 return lincs_rmsdData(impl_->lincsd);
738 real Constraints::rmsd() const
742 return lincs_rmsd(impl_->lincsd);
750 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
752 if (haveDynamicsIntegrator)
754 return FlexibleConstraintTreatment::Include;
758 return FlexibleConstraintTreatment::Exclude;
762 /*! \brief Returns a block struct to go from atoms to constraints
764 * The block struct will contain constraint indices with lower indices
765 * directly matching the order in F_CONSTR and higher indices matching
766 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
768 * \param[in] numAtoms The number of atoms to construct the list for
769 * \param[in] ilists The interaction lists, size F_NRE
770 * \param[in] iparams Interaction parameters, can be null when flexibleConstraintTreatment=Include
771 * \param[in] flexibleConstraintTreatment The flexible constraint treatment, see enum above
772 * \returns a block struct with all constraints for each atom
774 template <typename T>
776 makeAtomsToConstraintsList(int numAtoms,
778 const t_iparams *iparams,
779 FlexibleConstraintTreatment flexibleConstraintTreatment)
781 GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || iparams != nullptr, "With flexible constraint detection we need valid iparams");
783 std::vector<int> count(numAtoms);
785 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
787 const T &ilist = ilists[ftype];
788 const int stride = 1 + NRAL(ftype);
789 for (int i = 0; i < ilist.size(); i += stride)
791 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include ||
792 !isConstraintFlexible(iparams, ilist.iatoms[i]))
794 for (int j = 1; j < 3; j++)
796 int a = ilist.iatoms[i + j];
804 at2con.nr = numAtoms;
805 at2con.nalloc_index = at2con.nr + 1;
806 snew(at2con.index, at2con.nalloc_index);
808 for (int a = 0; a < numAtoms; a++)
810 at2con.index[a + 1] = at2con.index[a] + count[a];
813 at2con.nra = at2con.index[at2con.nr];
814 at2con.nalloc_a = at2con.nra;
815 snew(at2con.a, at2con.nalloc_a);
817 /* The F_CONSTRNC constraints have constraint numbers
818 * that continue after the last F_CONSTR constraint.
820 int numConstraints = 0;
821 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
823 const T &ilist = ilists[ftype];
824 const int stride = 1 + NRAL(ftype);
825 for (int i = 0; i < ilist.size(); i += stride)
827 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include ||
828 !isConstraintFlexible(iparams, ilist.iatoms[i]))
830 for (int j = 1; j < 3; j++)
832 int a = ilist.iatoms[i + j];
833 at2con.a[at2con.index[a] + count[a]++] = numConstraints;
843 t_blocka make_at2con(int numAtoms,
844 const t_ilist *ilist,
845 const t_iparams *iparams,
846 FlexibleConstraintTreatment flexibleConstraintTreatment)
848 return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
851 t_blocka make_at2con(const gmx_moltype_t &moltype,
852 gmx::ArrayRef<const t_iparams> iparams,
853 FlexibleConstraintTreatment flexibleConstraintTreatment)
855 return makeAtomsToConstraintsList(moltype.atoms.nr, moltype.ilist.data(), iparams.data(), flexibleConstraintTreatment);
858 //! Return the number of flexible constraints in the \c ilist and \c iparams.
859 template <typename T>
861 countFlexibleConstraintsTemplate(const T *ilist,
862 const t_iparams *iparams)
865 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
867 const int numIatomsPerConstraint = 3;
868 for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
870 const int type = ilist[ftype].iatoms[i];
871 if (iparams[type].constr.dA == 0 &&
872 iparams[type].constr.dB == 0)
882 int countFlexibleConstraints(const t_ilist *ilist,
883 const t_iparams *iparams)
885 return countFlexibleConstraintsTemplate(ilist, iparams);
888 //! Returns the index of the settle to which each atom belongs.
889 static std::vector<int> make_at2settle(int natoms,
890 const InteractionList &ilist)
892 /* Set all to no settle */
893 std::vector<int> at2s(natoms, -1);
895 const int stride = 1 + NRAL(F_SETTLE);
897 for (int s = 0; s < ilist.size(); s += stride)
899 at2s[ilist.iatoms[s + 1]] = s/stride;
900 at2s[ilist.iatoms[s + 2]] = s/stride;
901 at2s[ilist.iatoms[s + 3]] = s/stride;
908 Constraints::Impl::setConstraints(const gmx_localtop_t &top,
915 /* With DD we might also need to call LINCS on a domain no constraints for
916 * communicating coordinates to other nodes that do have constraints.
918 if (ir.eConstrAlg == econtLINCS)
920 set_lincs(top.idef, md, EI_DYNAMICS(ir.eI), cr, lincsd);
922 if (ir.eConstrAlg == econtSHAKE)
926 // We are using the local topology, so there are only
927 // F_CONSTR constraints.
928 make_shake_sblock_dd(shaked, &idef->il[F_CONSTR], cr->dd);
932 make_shake_sblock_serial(shaked, idef, md);
939 settle_set_constraints(settled,
940 &idef->il[F_SETTLE], md);
943 /* Make a selection of the local atoms for essential dynamics */
946 dd_make_local_ed_indices(cr->dd, ed);
951 Constraints::setConstraints(const gmx_localtop_t &top,
954 impl_->setConstraints(top, md);
957 /*! \brief Makes a per-moleculetype container of mappings from atom
958 * indices to constraint indices.
960 * Note that flexible constraints are only enabled with a dynamical integrator. */
961 static std::vector<t_blocka>
962 makeAtomToConstraintMappings(const gmx_mtop_t &mtop,
963 FlexibleConstraintTreatment flexibleConstraintTreatment)
965 std::vector<t_blocka> mapping;
966 mapping.reserve(mtop.moltype.size());
967 for (const gmx_moltype_t &moltype : mtop.moltype)
969 mapping.push_back(make_at2con(moltype,
970 mtop.ffparams.iparams,
971 flexibleConstraintTreatment));
976 Constraints::Constraints(const gmx_mtop_t &mtop,
977 const t_inputrec &ir,
982 const gmx_multisim_t *ms,
984 gmx_wallcycle *wcycle,
985 bool pbcHandlingRequired,
988 : impl_(new Impl(mtop,
1003 Constraints::Impl::Impl(const gmx_mtop_t &mtop_p,
1004 const t_inputrec &ir_p,
1007 const t_mdatoms &md_p,
1008 const t_commrec *cr_p,
1009 const gmx_multisim_t *ms_p,
1011 gmx_wallcycle *wcycle_p,
1012 bool pbcHandlingRequired,
1015 : ncon_tot(numConstraints),
1018 pbcHandlingRequired_(pbcHandlingRequired),
1022 pull_work(pull_work),
1027 if (numConstraints + numSettles > 0 && ir.epc == epcMTTK)
1029 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1033 if (numConstraints > 0)
1035 at2con_mt = makeAtomToConstraintMappings(mtop,
1036 flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
1038 for (const gmx_molblock_t &molblock : mtop.molblock)
1041 countFlexibleConstraintsTemplate(mtop.moltype[molblock.type].ilist.data(),
1042 mtop.ffparams.iparams.data());
1043 nflexcon += molblock.nmol*count;
1050 fprintf(log, "There are %d flexible constraints\n",
1052 if (ir.fc_stepsize == 0)
1055 "WARNING: step size for flexible constraining = 0\n"
1056 " All flexible constraints will be rigid.\n"
1057 " Will try to keep all flexible constraints at their original length,\n"
1058 " but the lengths may exhibit some drift.\n\n");
1064 please_cite(log, "Hess2002");
1068 if (ir.eConstrAlg == econtLINCS)
1070 lincsd = init_lincs(log, mtop,
1071 nflexcon, at2con_mt,
1072 DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd),
1073 ir.nLincsIter, ir.nProjOrder);
1076 if (ir.eConstrAlg == econtSHAKE)
1078 if (DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd))
1080 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1084 gmx_fatal(FARGS, "For this system also velocities and/or forces need to be constrained, this can not be done with SHAKE, you should select LINCS");
1086 please_cite(log, "Ryckaert77a");
1089 please_cite(log, "Barth95a");
1092 shaked = shake_init();
1098 please_cite(log, "Miyamoto92a");
1100 bInterCGsettles = inter_charge_group_settles(mtop);
1102 settled = settle_init(mtop);
1104 /* Make an atom to settle index for use in domain decomposition */
1105 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1107 at2settle_mt.emplace_back(make_at2settle(mtop.moltype[mt].atoms.nr,
1108 mtop.moltype[mt].ilist[F_SETTLE]));
1111 /* Allocate thread-local work arrays */
1112 int nthreads = gmx_omp_nthreads_get(emntSETTLE);
1113 if (nthreads > 1 && vir_r_m_dr_th == nullptr)
1115 snew(vir_r_m_dr_th, nthreads);
1116 snew(bSettleErrorHasOccurred, nthreads);
1121 char *env = getenv("GMX_MAXCONSTRWARN");
1125 sscanf(env, "%8d", &maxwarn);
1133 "Setting the maximum number of constraint warnings to %d\n",
1139 "Setting the maximum number of constraint warnings to %d\n",
1143 warncount_lincs = 0;
1144 warncount_settle = 0;
1147 Constraints::Impl::~Impl()
1149 for (auto blocka : at2con_mt)
1151 done_blocka(&blocka);
1153 if (bSettleErrorHasOccurred != nullptr)
1155 sfree(bSettleErrorHasOccurred);
1157 if (vir_r_m_dr_th != nullptr)
1159 sfree(vir_r_m_dr_th);
1161 if (settled != nullptr)
1163 settle_free(settled);
1168 void Constraints::saveEdsamPointer(gmx_edsam * ed)
1173 ArrayRef<const t_blocka>
1174 Constraints::atom2constraints_moltype() const
1176 return impl_->at2con_mt;
1179 ArrayRef < const std::vector < int>> Constraints::atom2settle_moltype() const
1181 return impl_->at2settle_mt;
1185 bool inter_charge_group_constraints(const gmx_mtop_t &mtop)
1187 const gmx_moltype_t *molt;
1189 int *at2cg, cg, a, ftype, i;
1193 for (size_t mb = 0; mb < mtop.molblock.size() && !bInterCG; mb++)
1195 molt = &mtop.moltype[mtop.molblock[mb].type];
1197 if (molt->ilist[F_CONSTR].size() > 0 ||
1198 molt->ilist[F_CONSTRNC].size() > 0 ||
1199 molt->ilist[F_SETTLE].size() > 0)
1202 snew(at2cg, molt->atoms.nr);
1203 for (cg = 0; cg < cgs->nr; cg++)
1205 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1211 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1213 const InteractionList &il = molt->ilist[ftype];
1214 for (i = 0; i < il.size() && !bInterCG; i += 1+NRAL(ftype))
1216 if (at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 2]])
1230 bool inter_charge_group_settles(const gmx_mtop_t &mtop)
1232 const gmx_moltype_t *molt;
1234 int *at2cg, cg, a, ftype, i;
1238 for (size_t mb = 0; mb < mtop.molblock.size() && !bInterCG; mb++)
1240 molt = &mtop.moltype[mtop.molblock[mb].type];
1242 if (molt->ilist[F_SETTLE].size() > 0)
1245 snew(at2cg, molt->atoms.nr);
1246 for (cg = 0; cg < cgs->nr; cg++)
1248 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1254 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1256 const InteractionList &il = molt->ilist[ftype];
1257 for (i = 0; i < il.size() && !bInterCG; i += 1+NRAL(F_SETTLE))
1259 if (at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 2]] ||
1260 at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 3]])
1274 void do_constrain_first(FILE *fplog, gmx::Constraints *constr,
1275 const t_inputrec *ir, const t_mdatoms *md,
1278 int i, m, start, end;
1280 real dt = ir->delta_t;
1284 /* We need to allocate one element extra, since we might use
1285 * (unaligned) 4-wide SIMD loads to access rvec entries.
1287 snew(savex, state->natoms + 1);
1294 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
1295 start, md->homenr, end);
1297 /* Do a first constrain to reset particles... */
1298 step = ir->init_step;
1301 char buf[STEPSTRSIZE];
1302 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
1303 gmx_step_str(step, buf));
1307 /* constrain the current position */
1308 constr->apply(TRUE, FALSE,
1310 state->x.rvec_array(), state->x.rvec_array(), nullptr,
1312 state->lambda[efptBONDED], &dvdl_dum,
1313 nullptr, nullptr, gmx::ConstraintVariable::Positions);
1316 /* constrain the inital velocity, and save it */
1317 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1318 constr->apply(TRUE, FALSE,
1320 state->x.rvec_array(), state->v.rvec_array(), state->v.rvec_array(),
1322 state->lambda[efptBONDED], &dvdl_dum,
1323 nullptr, nullptr, gmx::ConstraintVariable::Velocities);
1325 /* constrain the inital velocities at t-dt/2 */
1326 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
1328 auto x = makeArrayRef(state->x).subArray(start, end);
1329 auto v = makeArrayRef(state->v).subArray(start, end);
1330 for (i = start; (i < end); i++)
1332 for (m = 0; (m < DIM); m++)
1334 /* Reverse the velocity */
1336 /* Store the position at t-dt in buf */
1337 savex[i][m] = x[i][m] + dt*v[i][m];
1340 /* Shake the positions at t=-dt with the positions at t=0
1341 * as reference coordinates.
1345 char buf[STEPSTRSIZE];
1346 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
1347 gmx_step_str(step, buf));
1350 constr->apply(TRUE, FALSE,
1352 state->x.rvec_array(), savex, nullptr,
1354 state->lambda[efptBONDED], &dvdl_dum,
1355 state->v.rvec_array(), nullptr, gmx::ConstraintVariable::Positions);
1357 for (i = start; i < end; i++)
1359 for (m = 0; m < DIM; m++)
1361 /* Re-reverse the velocities */