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,
105 const t_mdatoms &md_p,
106 const t_commrec *cr_p,
107 const gmx_multisim_t *ms,
109 gmx_wallcycle *wcycle_p,
110 bool pbcHandlingRequired,
114 void setConstraints(const gmx_localtop_t &top,
115 const t_mdatoms &md);
116 bool apply(bool bLog,
129 ConstraintVariable econq);
130 //! The total number of constraints.
132 //! The number of flexible constraints.
134 //! A list of atoms to constraints for each moleculetype.
135 std::vector<t_blocka> at2con_mt;
136 //! A list of atoms to settles for each moleculetype
137 std::vector < std::vector < int>> at2settle_mt;
138 //! Whether any SETTLES cross charge-group boundaries.
139 bool bInterCGsettles = false;
141 Lincs *lincsd = nullptr; // TODO this should become a unique_ptr
143 shakedata *shaked = nullptr;
145 settledata *settled = nullptr;
146 //! The maximum number of warnings.
148 //! The number of warnings for LINCS.
149 int warncount_lincs = 0;
150 //! The number of warnings for SETTLE.
151 int warncount_settle = 0;
152 //! The essential dynamics data.
153 gmx_edsam * ed = nullptr;
155 //! Thread-local virial contribution.
156 tensor *vir_r_m_dr_th = {nullptr};
157 //! Did a settle error occur?
158 bool *bSettleErrorHasOccurred = nullptr;
160 //! Pointer to the global topology - only used for printing warnings.
161 const gmx_mtop_t &mtop;
162 //! Parameters for the interactions in this domain.
163 const t_idef *idef = nullptr;
164 //! Data about atoms in this domain.
166 //! Whether we need to do pbc for handling bonds.
167 bool pbcHandlingRequired_ = false;
171 //! Communication support.
172 const t_commrec *cr = nullptr;
173 //! Multi-sim support.
174 const gmx_multisim_t *ms = nullptr;
175 /*!\brief Input options.
177 * \todo Replace with IMdpOptions */
178 const t_inputrec &ir;
179 //! Flop counting support.
180 t_nrnb *nrnb = nullptr;
181 //! Tracks wallcycle usage.
182 gmx_wallcycle *wcycle;
185 Constraints::~Constraints() = default;
187 int Constraints::numFlexibleConstraints() const
189 return impl_->nflexcon;
192 bool Constraints::havePerturbedConstraints() const
194 const gmx_ffparams_t &ffparams = impl_->mtop.ffparams;
196 for (size_t i = 0; i < ffparams.functype.size(); i++)
198 if ((ffparams.functype[i] == F_CONSTR ||
199 ffparams.functype[i] == F_CONSTRNC) &&
200 ffparams.iparams[i].constr.dA != ffparams.iparams[i].constr.dB)
209 //! Clears constraint quantities for atoms in nonlocal region.
210 static void clear_constraint_quantity_nonlocal(gmx_domdec_t *dd, rvec *q)
212 int nonlocal_at_start, nonlocal_at_end, at;
214 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
216 for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
222 void too_many_constraint_warnings(int eConstrAlg, int warncount)
225 "Too many %s warnings (%d)\n"
226 "If you know what you are doing you can %s"
227 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
228 "but normally it is better to fix the problem",
229 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
230 (eConstrAlg == econtLINCS) ?
231 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
234 //! Writes out coordinates.
235 static void write_constr_pdb(const char *fn, const char *title,
236 const gmx_mtop_t &mtop,
237 int start, int homenr, const t_commrec *cr,
238 const rvec x[], matrix box)
242 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
244 const char *anm, *resnm;
247 if (DOMAINDECOMP(cr))
250 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
257 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
261 sprintf(fname, "%s.pdb", fn);
264 out = gmx_fio_fopen(fname, "w");
266 fprintf(out, "TITLE %s\n", title);
267 gmx_write_pdb_box(out, -1, box);
269 for (i = start; i < start+homenr; i++)
273 if (i >= dd_numHomeAtoms(*dd) && i < dd_ac0)
277 ii = dd->globalAtomIndices[i];
283 mtopGetAtomAndResidueName(mtop, ii, &molb, &anm, &resnr, &resnm, nullptr);
284 gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
285 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
287 fprintf(out, "TER\n");
292 //! Writes out domain contents to help diagnose crashes.
293 static void dump_confs(FILE *log, int64_t step, const gmx_mtop_t &mtop,
294 int start, int homenr, const t_commrec *cr,
295 const rvec x[], rvec xprime[], matrix box)
297 char buf[STRLEN], buf2[22];
299 char *env = getenv("GMX_SUPPRESS_DUMP");
305 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
306 write_constr_pdb(buf, "initial coordinates",
307 mtop, start, homenr, cr, x, box);
308 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
309 write_constr_pdb(buf, "coordinates after constraining",
310 mtop, start, homenr, cr, xprime, box);
313 fprintf(log, "Wrote pdb files with previous and current coordinates\n");
315 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
319 Constraints::apply(bool bLog,
332 ConstraintVariable econq)
334 return impl_->apply(bLog,
351 Constraints::Impl::apply(bool bLog,
364 ConstraintVariable econq)
370 real invdt, vir_fac = 0, t;
372 t_pbc pbc, *pbc_null;
376 wallcycle_start(wcycle, ewcCONSTR);
378 if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
380 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");
389 scaled_delta_t = step_scaling * ir.delta_t;
391 /* Prepare time step for use in constraint implementations, and
392 avoid generating inf when ir.delta_t = 0. */
399 invdt = 1.0/scaled_delta_t;
402 if (ir.efep != efepNO && EI_DYNAMICS(ir.eI))
404 /* Set the constraint lengths for the step at which this configuration
405 * is meant to be. The invmasses should not be changed.
407 lambda += delta_step*ir.fepvals->delta_lambda;
412 clear_mat(vir_r_m_dr);
414 const t_ilist *settle = &idef->il[F_SETTLE];
415 nsettle = settle->nr/(1+NRAL(F_SETTLE));
419 nth = gmx_omp_nthreads_get(emntSETTLE);
426 /* We do not need full pbc when constraints do not cross charge groups,
427 * i.e. when dd->constraint_comm==NULL.
428 * Note that PBC for constraints is different from PBC for bondeds.
429 * For constraints there is both forward and backward communication.
431 if (ir.ePBC != epbcNONE &&
432 (cr->dd || pbcHandlingRequired_) && !(cr->dd && cr->dd->constraint_comm == nullptr))
434 /* With pbc=screw the screw has been changed to a shift
435 * by the constraint coordinate communication routine,
436 * so that here we can use normal pbc.
438 pbc_null = set_pbc_dd(&pbc, ir.ePBC,
439 DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
447 /* Communicate the coordinates required for the non-local constraints
448 * for LINCS and/or SETTLE.
452 dd_move_x_constraints(cr->dd, box, x, xprime, econq == ConstraintVariable::Positions);
456 /* We need to initialize the non-local components of v.
457 * We never actually use these values, but we do increment them,
458 * so we should avoid uninitialized variables and overflows.
460 clear_constraint_quantity_nonlocal(cr->dd, v);
464 if (lincsd != nullptr)
466 bOK = constrain_lincs(bLog || bEner, ir, step, lincsd, md, cr, ms,
468 box, pbc_null, lambda, dvdlambda,
469 invdt, v, vir != nullptr, vir_r_m_dr,
471 maxwarn, &warncount_lincs);
472 if (!bOK && maxwarn < INT_MAX)
476 fprintf(log, "Constraint error in algorithm %s at step %s\n",
477 econstr_names[econtLINCS], gmx_step_str(step, buf));
483 if (shaked != nullptr)
485 bOK = constrain_shake(log, shaked,
487 *idef, ir, x, xprime, min_proj, nrnb,
489 invdt, v, vir != nullptr, vir_r_m_dr,
490 maxwarn < INT_MAX, econq);
492 if (!bOK && maxwarn < INT_MAX)
496 fprintf(log, "Constraint error in algorithm %s at step %s\n",
497 econstr_names[econtSHAKE], gmx_step_str(step, buf));
505 bool bSettleErrorHasOccurred0 = false;
509 case ConstraintVariable::Positions:
510 #pragma omp parallel for num_threads(nth) schedule(static)
511 for (th = 0; th < nth; th++)
517 clear_mat(vir_r_m_dr_th[th]);
524 invdt, v ? v[0] : nullptr,
526 th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th],
527 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
529 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
531 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
534 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
538 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
541 case ConstraintVariable::Velocities:
542 case ConstraintVariable::Derivative:
543 case ConstraintVariable::Force:
544 case ConstraintVariable::ForceDispl:
545 #pragma omp parallel for num_threads(nth) schedule(static)
546 for (th = 0; th < nth; th++)
550 int calcvir_atom_end;
554 calcvir_atom_end = 0;
558 calcvir_atom_end = md.homenr;
563 clear_mat(vir_r_m_dr_th[th]);
566 int start_th = (nsettle* th )/nth;
567 int end_th = (nsettle*(th+1))/nth;
569 if (start_th >= 0 && end_th - start_th > 0)
571 settle_proj(settled, econq,
573 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
576 xprime, min_proj, calcvir_atom_end,
577 th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th]);
580 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
582 /* This is an overestimate */
583 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
585 case ConstraintVariable::Deriv_FlexCon:
586 /* Nothing to do, since the are no flexible constraints in settles */
589 gmx_incons("Unknown constraint quantity for settle");
594 /* Reduce the virial contributions over the threads */
595 for (int th = 1; th < nth; th++)
597 m_add(vir_r_m_dr, vir_r_m_dr_th[th], vir_r_m_dr);
601 if (econq == ConstraintVariable::Positions)
603 for (int th = 1; th < nth; th++)
605 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
608 if (bSettleErrorHasOccurred0)
612 "\nstep " "%" PRId64 ": One or more water molecules can not be settled.\n"
613 "Check for bad contacts and/or reduce the timestep if appropriate.\n",
617 fprintf(log, "%s", buf);
619 fprintf(stderr, "%s", buf);
621 if (warncount_settle > maxwarn)
623 too_many_constraint_warnings(-1, warncount_settle);
634 /* The normal uses of constrain() pass step_scaling = 1.0.
635 * The call to constrain() for SD1 that passes step_scaling =
636 * 0.5 also passes vir = NULL, so cannot reach this
637 * assertion. This assertion should remain until someone knows
638 * that this path works for their intended purpose, and then
639 * they can use scaled_delta_t instead of ir.delta_t
641 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
644 case ConstraintVariable::Positions:
645 vir_fac = 0.5/(ir.delta_t*ir.delta_t);
647 case ConstraintVariable::Velocities:
648 vir_fac = 0.5/ir.delta_t;
650 case ConstraintVariable::Force:
651 case ConstraintVariable::ForceDispl:
655 gmx_incons("Unsupported constraint quantity for virial");
660 vir_fac *= 2; /* only constraining over half the distance here */
662 for (int i = 0; i < DIM; i++)
664 for (int j = 0; j < DIM; j++)
666 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
673 dump_confs(log, step, mtop, start, homenr, cr, x, xprime, box);
676 if (econq == ConstraintVariable::Positions)
678 if (ir.bPull && pull_have_constraint(ir.pull_work))
680 if (EI_DYNAMICS(ir.eI))
682 t = ir.init_t + (step + delta_step)*ir.delta_t;
688 set_pbc(&pbc, ir.ePBC, box);
689 pull_constraint(ir.pull_work, &md, &pbc, cr, ir.delta_t, t, x, xprime, v, *vir);
691 if (ed && delta_step > 0)
693 /* apply the essential dynamics constraints here */
694 do_edsam(&ir, step, cr, xprime, v, box, ed);
697 wallcycle_stop(wcycle, ewcCONSTR);
699 if (v != nullptr && md.cFREEZE)
701 /* Set the velocities of frozen dimensions to zero */
703 int gmx_unused numThreads = gmx_omp_nthreads_get(emntUpdate);
705 #pragma omp parallel for num_threads(numThreads) schedule(static)
706 for (int i = 0; i < md.homenr; i++)
708 int freezeGroup = md.cFREEZE[i];
710 for (int d = 0; d < DIM; d++)
712 if (ir.opts.nFreeze[freezeGroup][d])
723 ArrayRef<real> Constraints::rmsdData() const
727 return lincs_rmsdData(impl_->lincsd);
735 real Constraints::rmsd() const
739 return lincs_rmsd(impl_->lincsd);
747 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
749 if (haveDynamicsIntegrator)
751 return FlexibleConstraintTreatment::Include;
755 return FlexibleConstraintTreatment::Exclude;
759 /*! \brief Returns a block struct to go from atoms to constraints
761 * The block struct will contain constraint indices with lower indices
762 * directly matching the order in F_CONSTR and higher indices matching
763 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
765 * \param[in] numAtoms The number of atoms to construct the list for
766 * \param[in] ilists The interaction lists, size F_NRE
767 * \param[in] iparams Interaction parameters, can be null when flexibleConstraintTreatment=Include
768 * \param[in] flexibleConstraintTreatment The flexible constraint treatment, see enum above
769 * \returns a block struct with all constraints for each atom
771 template <typename T>
773 makeAtomsToConstraintsList(int numAtoms,
775 const t_iparams *iparams,
776 FlexibleConstraintTreatment flexibleConstraintTreatment)
778 GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || iparams != nullptr, "With flexible constraint detection we need valid iparams");
780 std::vector<int> count(numAtoms);
782 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
784 const T &ilist = ilists[ftype];
785 const int stride = 1 + NRAL(ftype);
786 for (int i = 0; i < ilist.size(); i += stride)
788 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include ||
789 !isConstraintFlexible(iparams, ilist.iatoms[i]))
791 for (int j = 1; j < 3; j++)
793 int a = ilist.iatoms[i + j];
801 at2con.nr = numAtoms;
802 at2con.nalloc_index = at2con.nr + 1;
803 snew(at2con.index, at2con.nalloc_index);
805 for (int a = 0; a < numAtoms; a++)
807 at2con.index[a + 1] = at2con.index[a] + count[a];
810 at2con.nra = at2con.index[at2con.nr];
811 at2con.nalloc_a = at2con.nra;
812 snew(at2con.a, at2con.nalloc_a);
814 /* The F_CONSTRNC constraints have constraint numbers
815 * that continue after the last F_CONSTR constraint.
817 int numConstraints = 0;
818 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
820 const T &ilist = ilists[ftype];
821 const int stride = 1 + NRAL(ftype);
822 for (int i = 0; i < ilist.size(); i += stride)
824 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include ||
825 !isConstraintFlexible(iparams, ilist.iatoms[i]))
827 for (int j = 1; j < 3; j++)
829 int a = ilist.iatoms[i + j];
830 at2con.a[at2con.index[a] + count[a]++] = numConstraints;
840 t_blocka make_at2con(int numAtoms,
841 const t_ilist *ilist,
842 const t_iparams *iparams,
843 FlexibleConstraintTreatment flexibleConstraintTreatment)
845 return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
848 t_blocka make_at2con(const gmx_moltype_t &moltype,
849 gmx::ArrayRef<const t_iparams> iparams,
850 FlexibleConstraintTreatment flexibleConstraintTreatment)
852 return makeAtomsToConstraintsList(moltype.atoms.nr, moltype.ilist.data(), iparams.data(), flexibleConstraintTreatment);
855 //! Return the number of flexible constraints in the \c ilist and \c iparams.
856 template <typename T>
858 countFlexibleConstraintsTemplate(const T *ilist,
859 const t_iparams *iparams)
862 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
864 const int numIatomsPerConstraint = 3;
865 for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
867 const int type = ilist[ftype].iatoms[i];
868 if (iparams[type].constr.dA == 0 &&
869 iparams[type].constr.dB == 0)
879 int countFlexibleConstraints(const t_ilist *ilist,
880 const t_iparams *iparams)
882 return countFlexibleConstraintsTemplate(ilist, iparams);
885 //! Returns the index of the settle to which each atom belongs.
886 static std::vector<int> make_at2settle(int natoms,
887 const InteractionList &ilist)
889 /* Set all to no settle */
890 std::vector<int> at2s(natoms, -1);
892 const int stride = 1 + NRAL(F_SETTLE);
894 for (int s = 0; s < ilist.size(); s += stride)
896 at2s[ilist.iatoms[s + 1]] = s/stride;
897 at2s[ilist.iatoms[s + 2]] = s/stride;
898 at2s[ilist.iatoms[s + 3]] = s/stride;
905 Constraints::Impl::setConstraints(const gmx_localtop_t &top,
912 /* With DD we might also need to call LINCS on a domain no constraints for
913 * communicating coordinates to other nodes that do have constraints.
915 if (ir.eConstrAlg == econtLINCS)
917 set_lincs(top.idef, md, EI_DYNAMICS(ir.eI), cr, lincsd);
919 if (ir.eConstrAlg == econtSHAKE)
923 // We are using the local topology, so there are only
924 // F_CONSTR constraints.
925 make_shake_sblock_dd(shaked, &idef->il[F_CONSTR], &top.cgs, cr->dd);
929 make_shake_sblock_serial(shaked, idef, md);
936 settle_set_constraints(settled,
937 &idef->il[F_SETTLE], md);
940 /* Make a selection of the local atoms for essential dynamics */
943 dd_make_local_ed_indices(cr->dd, ed);
948 Constraints::setConstraints(const gmx_localtop_t &top,
951 impl_->setConstraints(top, md);
954 /*! \brief Makes a per-moleculetype container of mappings from atom
955 * indices to constraint indices.
957 * Note that flexible constraints are only enabled with a dynamical integrator. */
958 static std::vector<t_blocka>
959 makeAtomToConstraintMappings(const gmx_mtop_t &mtop,
960 FlexibleConstraintTreatment flexibleConstraintTreatment)
962 std::vector<t_blocka> mapping;
963 mapping.reserve(mtop.moltype.size());
964 for (const gmx_moltype_t &moltype : mtop.moltype)
966 mapping.push_back(make_at2con(moltype,
967 mtop.ffparams.iparams,
968 flexibleConstraintTreatment));
973 Constraints::Constraints(const gmx_mtop_t &mtop,
974 const t_inputrec &ir,
978 const gmx_multisim_t *ms,
980 gmx_wallcycle *wcycle,
981 bool pbcHandlingRequired,
984 : impl_(new Impl(mtop,
998 Constraints::Impl::Impl(const gmx_mtop_t &mtop_p,
999 const t_inputrec &ir_p,
1001 const t_mdatoms &md_p,
1002 const t_commrec *cr_p,
1003 const gmx_multisim_t *ms_p,
1005 gmx_wallcycle *wcycle_p,
1006 bool pbcHandlingRequired,
1009 : ncon_tot(numConstraints),
1012 pbcHandlingRequired_(pbcHandlingRequired),
1020 if (numConstraints + numSettles > 0 && ir.epc == epcMTTK)
1022 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1026 if (numConstraints > 0)
1028 at2con_mt = makeAtomToConstraintMappings(mtop,
1029 flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
1031 for (const gmx_molblock_t &molblock : mtop.molblock)
1034 countFlexibleConstraintsTemplate(mtop.moltype[molblock.type].ilist.data(),
1035 mtop.ffparams.iparams.data());
1036 nflexcon += molblock.nmol*count;
1043 fprintf(log, "There are %d flexible constraints\n",
1045 if (ir.fc_stepsize == 0)
1048 "WARNING: step size for flexible constraining = 0\n"
1049 " All flexible constraints will be rigid.\n"
1050 " Will try to keep all flexible constraints at their original length,\n"
1051 " but the lengths may exhibit some drift.\n\n");
1057 please_cite(log, "Hess2002");
1061 if (ir.eConstrAlg == econtLINCS)
1063 lincsd = init_lincs(log, mtop,
1064 nflexcon, at2con_mt,
1065 DOMAINDECOMP(cr) && cr->dd->splitConstraints,
1066 ir.nLincsIter, ir.nProjOrder);
1069 if (ir.eConstrAlg == econtSHAKE)
1071 if (DOMAINDECOMP(cr) && cr->dd->splitConstraints)
1073 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1077 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");
1079 please_cite(log, "Ryckaert77a");
1082 please_cite(log, "Barth95a");
1085 shaked = shake_init();
1091 please_cite(log, "Miyamoto92a");
1093 bInterCGsettles = inter_charge_group_settles(mtop);
1095 settled = settle_init(mtop);
1097 /* Make an atom to settle index for use in domain decomposition */
1098 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1100 at2settle_mt.emplace_back(make_at2settle(mtop.moltype[mt].atoms.nr,
1101 mtop.moltype[mt].ilist[F_SETTLE]));
1104 /* Allocate thread-local work arrays */
1105 int nthreads = gmx_omp_nthreads_get(emntSETTLE);
1106 if (nthreads > 1 && vir_r_m_dr_th == nullptr)
1108 snew(vir_r_m_dr_th, nthreads);
1109 snew(bSettleErrorHasOccurred, nthreads);
1114 char *env = getenv("GMX_MAXCONSTRWARN");
1118 sscanf(env, "%8d", &maxwarn);
1126 "Setting the maximum number of constraint warnings to %d\n",
1132 "Setting the maximum number of constraint warnings to %d\n",
1136 warncount_lincs = 0;
1137 warncount_settle = 0;
1140 Constraints::Impl::~Impl()
1142 for (auto blocka : at2con_mt)
1144 done_blocka(&blocka);
1146 if (bSettleErrorHasOccurred != nullptr)
1148 sfree(bSettleErrorHasOccurred);
1150 if (vir_r_m_dr_th != nullptr)
1152 sfree(vir_r_m_dr_th);
1154 if (settled != nullptr)
1156 settle_free(settled);
1161 void Constraints::saveEdsamPointer(gmx_edsam * ed)
1166 const ArrayRef<const t_blocka>
1167 Constraints::atom2constraints_moltype() const
1169 return impl_->at2con_mt;
1172 ArrayRef < const std::vector < int>> Constraints::atom2settle_moltype() const
1174 return impl_->at2settle_mt;
1178 bool inter_charge_group_constraints(const gmx_mtop_t &mtop)
1180 const gmx_moltype_t *molt;
1182 int *at2cg, cg, a, ftype, i;
1186 for (size_t mb = 0; mb < mtop.molblock.size() && !bInterCG; mb++)
1188 molt = &mtop.moltype[mtop.molblock[mb].type];
1190 if (molt->ilist[F_CONSTR].size() > 0 ||
1191 molt->ilist[F_CONSTRNC].size() > 0 ||
1192 molt->ilist[F_SETTLE].size() > 0)
1195 snew(at2cg, molt->atoms.nr);
1196 for (cg = 0; cg < cgs->nr; cg++)
1198 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1204 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1206 const InteractionList &il = molt->ilist[ftype];
1207 for (i = 0; i < il.size() && !bInterCG; i += 1+NRAL(ftype))
1209 if (at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 2]])
1223 bool inter_charge_group_settles(const gmx_mtop_t &mtop)
1225 const gmx_moltype_t *molt;
1227 int *at2cg, cg, a, ftype, i;
1231 for (size_t mb = 0; mb < mtop.molblock.size() && !bInterCG; mb++)
1233 molt = &mtop.moltype[mtop.molblock[mb].type];
1235 if (molt->ilist[F_SETTLE].size() > 0)
1238 snew(at2cg, molt->atoms.nr);
1239 for (cg = 0; cg < cgs->nr; cg++)
1241 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1247 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1249 const InteractionList &il = molt->ilist[ftype];
1250 for (i = 0; i < il.size() && !bInterCG; i += 1+NRAL(F_SETTLE))
1252 if (at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 2]] ||
1253 at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 3]])