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/mdrun.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/mdatom.h"
72 #include "gromacs/mdtypes/state.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/pulling/pull.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 //! Clears constraint quantities for atoms in nonlocal region.
193 static void clear_constraint_quantity_nonlocal(gmx_domdec_t *dd, rvec *q)
195 int nonlocal_at_start, nonlocal_at_end, at;
197 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
199 for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
205 void too_many_constraint_warnings(int eConstrAlg, int warncount)
208 "Too many %s warnings (%d)\n"
209 "If you know what you are doing you can %s"
210 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
211 "but normally it is better to fix the problem",
212 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
213 (eConstrAlg == econtLINCS) ?
214 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
217 //! Writes out coordinates.
218 static void write_constr_pdb(const char *fn, const char *title,
219 const gmx_mtop_t &mtop,
220 int start, int homenr, const t_commrec *cr,
221 const rvec x[], matrix box)
225 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
227 const char *anm, *resnm;
230 if (DOMAINDECOMP(cr))
233 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
240 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
244 sprintf(fname, "%s.pdb", fn);
247 out = gmx_fio_fopen(fname, "w");
249 fprintf(out, "TITLE %s\n", title);
250 gmx_write_pdb_box(out, -1, box);
252 for (i = start; i < start+homenr; i++)
256 if (i >= dd_numHomeAtoms(*dd) && i < dd_ac0)
260 ii = dd->globalAtomIndices[i];
266 mtopGetAtomAndResidueName(mtop, ii, &molb, &anm, &resnr, &resnm, nullptr);
267 gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
268 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
270 fprintf(out, "TER\n");
275 //! Writes out domain contents to help diagnose crashes.
276 static void dump_confs(FILE *log, int64_t step, const gmx_mtop_t &mtop,
277 int start, int homenr, const t_commrec *cr,
278 const rvec x[], rvec xprime[], matrix box)
280 char buf[STRLEN], buf2[22];
282 char *env = getenv("GMX_SUPPRESS_DUMP");
288 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
289 write_constr_pdb(buf, "initial coordinates",
290 mtop, start, homenr, cr, x, box);
291 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
292 write_constr_pdb(buf, "coordinates after constraining",
293 mtop, start, homenr, cr, xprime, box);
296 fprintf(log, "Wrote pdb files with previous and current coordinates\n");
298 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
302 Constraints::apply(bool bLog,
315 ConstraintVariable econq)
317 return impl_->apply(bLog,
334 Constraints::Impl::apply(bool bLog,
347 ConstraintVariable econq)
353 real invdt, vir_fac = 0, t;
355 t_pbc pbc, *pbc_null;
359 wallcycle_start(wcycle, ewcCONSTR);
361 if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
363 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");
372 scaled_delta_t = step_scaling * ir.delta_t;
374 /* Prepare time step for use in constraint implementations, and
375 avoid generating inf when ir.delta_t = 0. */
382 invdt = 1.0/scaled_delta_t;
385 if (ir.efep != efepNO && EI_DYNAMICS(ir.eI))
387 /* Set the constraint lengths for the step at which this configuration
388 * is meant to be. The invmasses should not be changed.
390 lambda += delta_step*ir.fepvals->delta_lambda;
395 clear_mat(vir_r_m_dr);
397 const t_ilist *settle = &idef->il[F_SETTLE];
398 nsettle = settle->nr/(1+NRAL(F_SETTLE));
402 nth = gmx_omp_nthreads_get(emntSETTLE);
409 /* We do not need full pbc when constraints do not cross charge groups,
410 * i.e. when dd->constraint_comm==NULL.
411 * Note that PBC for constraints is different from PBC for bondeds.
412 * For constraints there is both forward and backward communication.
414 if (ir.ePBC != epbcNONE &&
415 (cr->dd || pbcHandlingRequired_) && !(cr->dd && cr->dd->constraint_comm == nullptr))
417 /* With pbc=screw the screw has been changed to a shift
418 * by the constraint coordinate communication routine,
419 * so that here we can use normal pbc.
421 pbc_null = set_pbc_dd(&pbc, ir.ePBC,
422 DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
430 /* Communicate the coordinates required for the non-local constraints
431 * for LINCS and/or SETTLE.
435 dd_move_x_constraints(cr->dd, box, x, xprime, econq == ConstraintVariable::Positions);
439 /* We need to initialize the non-local components of v.
440 * We never actually use these values, but we do increment them,
441 * so we should avoid uninitialized variables and overflows.
443 clear_constraint_quantity_nonlocal(cr->dd, v);
447 if (lincsd != nullptr)
449 bOK = constrain_lincs(bLog || bEner, ir, step, lincsd, md, cr, ms,
451 box, pbc_null, lambda, dvdlambda,
452 invdt, v, vir != nullptr, vir_r_m_dr,
454 maxwarn, &warncount_lincs);
455 if (!bOK && maxwarn < INT_MAX)
459 fprintf(log, "Constraint error in algorithm %s at step %s\n",
460 econstr_names[econtLINCS], gmx_step_str(step, buf));
466 if (shaked != nullptr)
468 bOK = constrain_shake(log, shaked,
470 *idef, ir, x, xprime, min_proj, nrnb,
472 invdt, v, vir != nullptr, vir_r_m_dr,
473 maxwarn < INT_MAX, econq);
475 if (!bOK && maxwarn < INT_MAX)
479 fprintf(log, "Constraint error in algorithm %s at step %s\n",
480 econstr_names[econtSHAKE], gmx_step_str(step, buf));
488 bool bSettleErrorHasOccurred0 = false;
492 case ConstraintVariable::Positions:
493 #pragma omp parallel for num_threads(nth) schedule(static)
494 for (th = 0; th < nth; th++)
500 clear_mat(vir_r_m_dr_th[th]);
507 invdt, v ? v[0] : nullptr,
509 th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th],
510 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
512 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
514 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
517 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
521 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
524 case ConstraintVariable::Velocities:
525 case ConstraintVariable::Derivative:
526 case ConstraintVariable::Force:
527 case ConstraintVariable::ForceDispl:
528 #pragma omp parallel for num_threads(nth) schedule(static)
529 for (th = 0; th < nth; th++)
533 int calcvir_atom_end;
537 calcvir_atom_end = 0;
541 calcvir_atom_end = md.homenr;
546 clear_mat(vir_r_m_dr_th[th]);
549 int start_th = (nsettle* th )/nth;
550 int end_th = (nsettle*(th+1))/nth;
552 if (start_th >= 0 && end_th - start_th > 0)
554 settle_proj(settled, econq,
556 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
559 xprime, min_proj, calcvir_atom_end,
560 th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th]);
563 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
565 /* This is an overestimate */
566 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
568 case ConstraintVariable::Deriv_FlexCon:
569 /* Nothing to do, since the are no flexible constraints in settles */
572 gmx_incons("Unknown constraint quantity for settle");
577 /* Reduce the virial contributions over the threads */
578 for (int th = 1; th < nth; th++)
580 m_add(vir_r_m_dr, vir_r_m_dr_th[th], vir_r_m_dr);
584 if (econq == ConstraintVariable::Positions)
586 for (int th = 1; th < nth; th++)
588 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
591 if (bSettleErrorHasOccurred0)
595 "\nstep " "%" PRId64 ": One or more water molecules can not be settled.\n"
596 "Check for bad contacts and/or reduce the timestep if appropriate.\n",
600 fprintf(log, "%s", buf);
602 fprintf(stderr, "%s", buf);
604 if (warncount_settle > maxwarn)
606 too_many_constraint_warnings(-1, warncount_settle);
617 /* The normal uses of constrain() pass step_scaling = 1.0.
618 * The call to constrain() for SD1 that passes step_scaling =
619 * 0.5 also passes vir = NULL, so cannot reach this
620 * assertion. This assertion should remain until someone knows
621 * that this path works for their intended purpose, and then
622 * they can use scaled_delta_t instead of ir.delta_t
624 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
627 case ConstraintVariable::Positions:
628 vir_fac = 0.5/(ir.delta_t*ir.delta_t);
630 case ConstraintVariable::Velocities:
631 vir_fac = 0.5/ir.delta_t;
633 case ConstraintVariable::Force:
634 case ConstraintVariable::ForceDispl:
638 gmx_incons("Unsupported constraint quantity for virial");
643 vir_fac *= 2; /* only constraining over half the distance here */
645 for (int i = 0; i < DIM; i++)
647 for (int j = 0; j < DIM; j++)
649 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
656 dump_confs(log, step, mtop, start, homenr, cr, x, xprime, box);
659 if (econq == ConstraintVariable::Positions)
661 if (ir.bPull && pull_have_constraint(ir.pull_work))
663 if (EI_DYNAMICS(ir.eI))
665 t = ir.init_t + (step + delta_step)*ir.delta_t;
671 set_pbc(&pbc, ir.ePBC, box);
672 pull_constraint(ir.pull_work, &md, &pbc, cr, ir.delta_t, t, x, xprime, v, *vir);
674 if (ed && delta_step > 0)
676 /* apply the essential dynamics constraints here */
677 do_edsam(&ir, step, cr, xprime, v, box, ed);
680 wallcycle_stop(wcycle, ewcCONSTR);
682 if (v != nullptr && md.cFREEZE)
684 /* Set the velocities of frozen dimensions to zero */
686 int gmx_unused numThreads = gmx_omp_nthreads_get(emntUpdate);
688 #pragma omp parallel for num_threads(numThreads) schedule(static)
689 for (int i = 0; i < md.homenr; i++)
691 int freezeGroup = md.cFREEZE[i];
693 for (int d = 0; d < DIM; d++)
695 if (ir.opts.nFreeze[freezeGroup][d])
706 ArrayRef<real> Constraints::rmsdData() const
710 return lincs_rmsdData(impl_->lincsd);
714 return EmptyArrayRef();
718 real Constraints::rmsd() const
722 return lincs_rmsd(impl_->lincsd);
730 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
732 if (haveDynamicsIntegrator)
734 return FlexibleConstraintTreatment::Include;
738 return FlexibleConstraintTreatment::Exclude;
742 /*! \brief Returns a block struct to go from atoms to constraints
744 * The block struct will contain constraint indices with lower indices
745 * directly matching the order in F_CONSTR and higher indices matching
746 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
748 * \param[in] numAtoms The number of atoms to construct the list for
749 * \param[in] ilists The interaction lists, size F_NRE
750 * \param[in] iparams Interaction parameters, can be null when flexibleConstraintTreatment=Include
751 * \param[in] flexibleConstraintTreatment The flexible constraint treatment, see enum above
752 * \returns a block struct with all constraints for each atom
754 template <typename T>
756 makeAtomsToConstraintsList(int numAtoms,
758 const t_iparams *iparams,
759 FlexibleConstraintTreatment flexibleConstraintTreatment)
761 GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || iparams != nullptr, "With flexible constraint detection we need valid iparams");
763 std::vector<int> count(numAtoms);
765 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
767 const T &ilist = ilists[ftype];
768 const int stride = 1 + NRAL(ftype);
769 for (int i = 0; i < ilist.size(); i += stride)
771 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include ||
772 !isConstraintFlexible(iparams, ilist.iatoms[i]))
774 for (int j = 1; j < 3; j++)
776 int a = ilist.iatoms[i + j];
784 at2con.nr = numAtoms;
785 at2con.nalloc_index = at2con.nr + 1;
786 snew(at2con.index, at2con.nalloc_index);
788 for (int a = 0; a < numAtoms; a++)
790 at2con.index[a + 1] = at2con.index[a] + count[a];
793 at2con.nra = at2con.index[at2con.nr];
794 at2con.nalloc_a = at2con.nra;
795 snew(at2con.a, at2con.nalloc_a);
797 /* The F_CONSTRNC constraints have constraint numbers
798 * that continue after the last F_CONSTR constraint.
800 int numConstraints = 0;
801 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
803 const T &ilist = ilists[ftype];
804 const int stride = 1 + NRAL(ftype);
805 for (int i = 0; i < ilist.size(); i += stride)
807 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include ||
808 !isConstraintFlexible(iparams, ilist.iatoms[i]))
810 for (int j = 1; j < 3; j++)
812 int a = ilist.iatoms[i + j];
813 at2con.a[at2con.index[a] + count[a]++] = numConstraints;
823 t_blocka make_at2con(int numAtoms,
824 const t_ilist *ilist,
825 const t_iparams *iparams,
826 FlexibleConstraintTreatment flexibleConstraintTreatment)
828 return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
831 t_blocka make_at2con(const gmx_moltype_t &moltype,
832 gmx::ArrayRef<const t_iparams> iparams,
833 FlexibleConstraintTreatment flexibleConstraintTreatment)
835 return makeAtomsToConstraintsList(moltype.atoms.nr, moltype.ilist.data(), iparams.data(), flexibleConstraintTreatment);
838 //! Return the number of flexible constraints in the \c ilist and \c iparams.
839 template <typename T>
841 countFlexibleConstraintsTemplate(const T *ilist,
842 const t_iparams *iparams)
845 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
847 const int numIatomsPerConstraint = 3;
848 for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
850 const int type = ilist[ftype].iatoms[i];
851 if (iparams[type].constr.dA == 0 &&
852 iparams[type].constr.dB == 0)
862 int countFlexibleConstraints(const t_ilist *ilist,
863 const t_iparams *iparams)
865 return countFlexibleConstraintsTemplate(ilist, iparams);
868 //! Returns the index of the settle to which each atom belongs.
869 static std::vector<int> make_at2settle(int natoms,
870 const InteractionList &ilist)
872 /* Set all to no settle */
873 std::vector<int> at2s(natoms, -1);
875 const int stride = 1 + NRAL(F_SETTLE);
877 for (int s = 0; s < ilist.size(); s += stride)
879 at2s[ilist.iatoms[s + 1]] = s/stride;
880 at2s[ilist.iatoms[s + 2]] = s/stride;
881 at2s[ilist.iatoms[s + 3]] = s/stride;
888 Constraints::Impl::setConstraints(const gmx_localtop_t &top,
895 /* With DD we might also need to call LINCS on a domain no constraints for
896 * communicating coordinates to other nodes that do have constraints.
898 if (ir.eConstrAlg == econtLINCS)
900 set_lincs(top.idef, md, EI_DYNAMICS(ir.eI), cr, lincsd);
902 if (ir.eConstrAlg == econtSHAKE)
906 // We are using the local topology, so there are only
907 // F_CONSTR constraints.
908 make_shake_sblock_dd(shaked, &idef->il[F_CONSTR], &top.cgs, cr->dd);
912 make_shake_sblock_serial(shaked, idef, md);
919 settle_set_constraints(settled,
920 &idef->il[F_SETTLE], md);
923 /* Make a selection of the local atoms for essential dynamics */
926 dd_make_local_ed_indices(cr->dd, ed);
931 Constraints::setConstraints(const gmx_localtop_t &top,
934 impl_->setConstraints(top, md);
937 /*! \brief Makes a per-moleculetype container of mappings from atom
938 * indices to constraint indices.
940 * Note that flexible constraints are only enabled with a dynamical integrator. */
941 static std::vector<t_blocka>
942 makeAtomToConstraintMappings(const gmx_mtop_t &mtop,
943 FlexibleConstraintTreatment flexibleConstraintTreatment)
945 std::vector<t_blocka> mapping;
946 mapping.reserve(mtop.moltype.size());
947 for (const gmx_moltype_t &moltype : mtop.moltype)
949 mapping.push_back(make_at2con(moltype,
950 mtop.ffparams.iparams,
951 flexibleConstraintTreatment));
956 Constraints::Constraints(const gmx_mtop_t &mtop,
957 const t_inputrec &ir,
961 const gmx_multisim_t *ms,
963 gmx_wallcycle *wcycle,
964 bool pbcHandlingRequired,
967 : impl_(new Impl(mtop,
981 Constraints::Impl::Impl(const gmx_mtop_t &mtop_p,
982 const t_inputrec &ir_p,
984 const t_mdatoms &md_p,
985 const t_commrec *cr_p,
986 const gmx_multisim_t *ms_p,
988 gmx_wallcycle *wcycle_p,
989 bool pbcHandlingRequired,
992 : ncon_tot(numConstraints),
995 pbcHandlingRequired_(pbcHandlingRequired),
1003 if (numConstraints + numSettles > 0 && ir.epc == epcMTTK)
1005 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1009 if (numConstraints > 0)
1011 at2con_mt = makeAtomToConstraintMappings(mtop,
1012 flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
1014 for (const gmx_molblock_t &molblock : mtop.molblock)
1017 countFlexibleConstraintsTemplate(mtop.moltype[molblock.type].ilist.data(),
1018 mtop.ffparams.iparams.data());
1019 nflexcon += molblock.nmol*count;
1026 fprintf(log, "There are %d flexible constraints\n",
1028 if (ir.fc_stepsize == 0)
1031 "WARNING: step size for flexible constraining = 0\n"
1032 " All flexible constraints will be rigid.\n"
1033 " Will try to keep all flexible constraints at their original length,\n"
1034 " but the lengths may exhibit some drift.\n\n");
1040 please_cite(log, "Hess2002");
1044 if (ir.eConstrAlg == econtLINCS)
1046 lincsd = init_lincs(log, mtop,
1047 nflexcon, at2con_mt,
1048 DOMAINDECOMP(cr) && cr->dd->splitConstraints,
1049 ir.nLincsIter, ir.nProjOrder);
1052 if (ir.eConstrAlg == econtSHAKE)
1054 if (DOMAINDECOMP(cr) && cr->dd->splitConstraints)
1056 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1060 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");
1062 please_cite(log, "Ryckaert77a");
1065 please_cite(log, "Barth95a");
1068 shaked = shake_init();
1074 please_cite(log, "Miyamoto92a");
1076 bInterCGsettles = inter_charge_group_settles(mtop);
1078 settled = settle_init(mtop);
1080 /* Make an atom to settle index for use in domain decomposition */
1081 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1083 at2settle_mt.emplace_back(make_at2settle(mtop.moltype[mt].atoms.nr,
1084 mtop.moltype[mt].ilist[F_SETTLE]));
1087 /* Allocate thread-local work arrays */
1088 int nthreads = gmx_omp_nthreads_get(emntSETTLE);
1089 if (nthreads > 1 && vir_r_m_dr_th == nullptr)
1091 snew(vir_r_m_dr_th, nthreads);
1092 snew(bSettleErrorHasOccurred, nthreads);
1097 char *env = getenv("GMX_MAXCONSTRWARN");
1101 sscanf(env, "%8d", &maxwarn);
1109 "Setting the maximum number of constraint warnings to %d\n",
1115 "Setting the maximum number of constraint warnings to %d\n",
1119 warncount_lincs = 0;
1120 warncount_settle = 0;
1123 Constraints::Impl::~Impl()
1125 for (auto blocka : at2con_mt)
1127 done_blocka(&blocka);
1129 if (bSettleErrorHasOccurred != nullptr)
1131 sfree(bSettleErrorHasOccurred);
1133 if (vir_r_m_dr_th != nullptr)
1135 sfree(vir_r_m_dr_th);
1137 if (settled != nullptr)
1139 settle_free(settled);
1144 void Constraints::saveEdsamPointer(gmx_edsam * ed)
1149 const ArrayRef<const t_blocka>
1150 Constraints::atom2constraints_moltype() const
1152 return impl_->at2con_mt;
1155 ArrayRef < const std::vector < int>> Constraints::atom2settle_moltype() const
1157 return impl_->at2settle_mt;
1161 bool inter_charge_group_constraints(const gmx_mtop_t &mtop)
1163 const gmx_moltype_t *molt;
1165 int *at2cg, cg, a, ftype, i;
1169 for (size_t mb = 0; mb < mtop.molblock.size() && !bInterCG; mb++)
1171 molt = &mtop.moltype[mtop.molblock[mb].type];
1173 if (molt->ilist[F_CONSTR].size() > 0 ||
1174 molt->ilist[F_CONSTRNC].size() > 0 ||
1175 molt->ilist[F_SETTLE].size() > 0)
1178 snew(at2cg, molt->atoms.nr);
1179 for (cg = 0; cg < cgs->nr; cg++)
1181 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1187 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1189 const InteractionList &il = molt->ilist[ftype];
1190 for (i = 0; i < il.size() && !bInterCG; i += 1+NRAL(ftype))
1192 if (at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 2]])
1206 bool inter_charge_group_settles(const gmx_mtop_t &mtop)
1208 const gmx_moltype_t *molt;
1210 int *at2cg, cg, a, ftype, i;
1214 for (size_t mb = 0; mb < mtop.molblock.size() && !bInterCG; mb++)
1216 molt = &mtop.moltype[mtop.molblock[mb].type];
1218 if (molt->ilist[F_SETTLE].size() > 0)
1221 snew(at2cg, molt->atoms.nr);
1222 for (cg = 0; cg < cgs->nr; cg++)
1224 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1230 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1232 const InteractionList &il = molt->ilist[ftype];
1233 for (i = 0; i < il.size() && !bInterCG; i += 1+NRAL(F_SETTLE))
1235 if (at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 2]] ||
1236 at2cg[il.iatoms[i + 1]] != at2cg[il.iatoms[i + 3]])