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, 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;
139 Lincs* lincsd = nullptr; // TODO this should become a unique_ptr
141 shakedata* shaked = nullptr;
143 settledata* settled = nullptr;
144 //! The maximum number of warnings.
146 //! The number of warnings for LINCS.
147 int warncount_lincs = 0;
148 //! The number of warnings for SETTLE.
149 int warncount_settle = 0;
150 //! The essential dynamics data.
151 gmx_edsam* ed = nullptr;
153 //! Thread-local virial contribution.
154 tensor* vir_r_m_dr_th = { nullptr };
155 //! Did a settle error occur?
156 bool* bSettleErrorHasOccurred = nullptr;
158 //! Pointer to the global topology - only used for printing warnings.
159 const gmx_mtop_t& mtop;
160 //! Parameters for the interactions in this domain.
161 const t_idef* idef = nullptr;
162 //! Data about atoms in this domain.
164 //! Whether we need to do pbc for handling bonds.
165 bool pbcHandlingRequired_ = false;
169 //! Communication support.
170 const t_commrec* cr = nullptr;
171 //! Multi-sim support.
172 const gmx_multisim_t* ms = nullptr;
173 //! Pulling code object, if any.
174 pull_t* pull_work = 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 || ffparams.functype[i] == F_CONSTRNC)
199 && ffparams.iparams[i].constr.dA != ffparams.iparams[i].constr.dB)
208 //! Clears constraint quantities for atoms in nonlocal region.
209 static void clear_constraint_quantity_nonlocal(gmx_domdec_t* dd, rvec* q)
211 int nonlocal_at_start, nonlocal_at_end, at;
213 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
215 for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
221 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) ? "adjust the lincs warning threshold in your mdp file\nor " : "\n");
233 //! Writes out coordinates.
234 static void write_constr_pdb(const char* fn,
236 const gmx_mtop_t& mtop,
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,
298 const gmx_mtop_t& mtop,
306 char buf[STRLEN], buf2[22];
308 char* env = getenv("GMX_SUPPRESS_DUMP");
314 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
315 write_constr_pdb(buf, "initial coordinates", mtop, start, homenr, cr, x, box);
316 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
317 write_constr_pdb(buf, "coordinates after constraining", mtop, start, homenr, cr, xprime, box);
320 fprintf(log, "Wrote pdb files with previous and current coordinates\n");
322 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
325 bool Constraints::apply(bool bLog,
338 ConstraintVariable econq)
340 return impl_->apply(bLog, bEner, step, delta_step, step_scaling, x, xprime, min_proj, box,
341 lambda, dvdlambda, v, vir, econq);
344 bool Constraints::Impl::apply(bool bLog,
357 ConstraintVariable econq)
363 real invdt, vir_fac = 0, t;
365 t_pbc pbc, *pbc_null;
369 wallcycle_start(wcycle, ewcCONSTR);
371 if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
374 "constrain called for forces displacements while not doing energy minimization, "
375 "can not do this while the LINCS and SETTLE constraint connection matrices are "
385 scaled_delta_t = step_scaling * ir.delta_t;
387 /* Prepare time step for use in constraint implementations, and
388 avoid generating inf when ir.delta_t = 0. */
395 invdt = 1.0 / scaled_delta_t;
398 if (ir.efep != efepNO && EI_DYNAMICS(ir.eI))
400 /* Set the constraint lengths for the step at which this configuration
401 * is meant to be. The invmasses should not be changed.
403 lambda += delta_step * ir.fepvals->delta_lambda;
408 clear_mat(vir_r_m_dr);
410 const t_ilist* settle = &idef->il[F_SETTLE];
411 nsettle = settle->nr / (1 + NRAL(F_SETTLE));
415 nth = gmx_omp_nthreads_get(emntSETTLE);
422 /* We do not need full pbc when constraints do not cross update groups
423 * i.e. when dd->constraint_comm==NULL.
424 * Note that PBC for constraints is different from PBC for bondeds.
425 * For constraints there is both forward and backward communication.
427 if (ir.ePBC != epbcNONE && (cr->dd || pbcHandlingRequired_)
428 && !(cr->dd && cr->dd->constraint_comm == nullptr))
430 /* With pbc=screw the screw has been changed to a shift
431 * by the constraint coordinate communication routine,
432 * so that here we can use normal pbc.
434 pbc_null = set_pbc_dd(&pbc, ir.ePBC, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, FALSE, box);
441 /* Communicate the coordinates required for the non-local constraints
442 * for LINCS and/or SETTLE.
446 dd_move_x_constraints(cr->dd, box, x, xprime, econq == ConstraintVariable::Positions);
450 /* We need to initialize the non-local components of v.
451 * We never actually use these values, but we do increment them,
452 * so we should avoid uninitialized variables and overflows.
454 clear_constraint_quantity_nonlocal(cr->dd, v);
458 if (lincsd != nullptr)
460 bOK = constrain_lincs(bLog || bEner, ir, step, lincsd, md, cr, ms, x, xprime, min_proj, box,
461 pbc_null, lambda, dvdlambda, invdt, v, vir != nullptr, vir_r_m_dr,
462 econq, nrnb, maxwarn, &warncount_lincs);
463 if (!bOK && maxwarn < INT_MAX)
467 fprintf(log, "Constraint error in algorithm %s at step %s\n",
468 econstr_names[econtLINCS], gmx_step_str(step, buf));
474 if (shaked != nullptr)
476 bOK = constrain_shake(log, shaked, md.invmass, *idef, ir, x, xprime, min_proj, nrnb, lambda,
477 dvdlambda, invdt, v, vir != nullptr, vir_r_m_dr, maxwarn < INT_MAX, econq);
479 if (!bOK && maxwarn < INT_MAX)
483 fprintf(log, "Constraint error in algorithm %s at step %s\n",
484 econstr_names[econtSHAKE], gmx_step_str(step, buf));
492 bool bSettleErrorHasOccurred0 = false;
496 case ConstraintVariable::Positions:
497 #pragma omp parallel for num_threads(nth) schedule(static)
498 for (int th = 0; th < nth; th++)
504 clear_mat(vir_r_m_dr_th[th]);
507 csettle(settled, nth, th, pbc_null, x[0], xprime[0], invdt, v ? v[0] : nullptr,
508 vir != nullptr, th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th],
509 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
511 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
513 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
516 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle * 3);
520 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle * 3);
523 case ConstraintVariable::Velocities:
524 case ConstraintVariable::Derivative:
525 case ConstraintVariable::Force:
526 case ConstraintVariable::ForceDispl:
527 #pragma omp parallel for num_threads(nth) schedule(static)
528 for (int th = 0; th < nth; th++)
532 int calcvir_atom_end;
536 calcvir_atom_end = 0;
540 calcvir_atom_end = md.homenr;
545 clear_mat(vir_r_m_dr_th[th]);
548 int start_th = (nsettle * th) / nth;
549 int end_th = (nsettle * (th + 1)) / nth;
551 if (start_th >= 0 && end_th - start_th > 0)
553 settle_proj(settled, econq, end_th - start_th,
554 settle->iatoms + start_th * (1 + NRAL(F_SETTLE)), pbc_null,
555 x, xprime, min_proj, calcvir_atom_end,
556 th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th]);
559 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
561 /* This is an overestimate */
562 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
564 case ConstraintVariable::Deriv_FlexCon:
565 /* Nothing to do, since the are no flexible constraints in settles */
567 default: gmx_incons("Unknown constraint quantity for settle");
572 /* Reduce the virial contributions over the threads */
573 for (int th = 1; th < nth; th++)
575 m_add(vir_r_m_dr, vir_r_m_dr_th[th], vir_r_m_dr);
579 if (econq == ConstraintVariable::Positions)
581 for (int th = 1; th < nth; th++)
583 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
586 if (bSettleErrorHasOccurred0)
592 ": One or more water molecules can not be settled.\n"
593 "Check for bad contacts and/or reduce the timestep if appropriate.\n",
597 fprintf(log, "%s", buf);
599 fprintf(stderr, "%s", buf);
601 if (warncount_settle > maxwarn)
603 too_many_constraint_warnings(-1, warncount_settle);
614 /* The normal uses of constrain() pass step_scaling = 1.0.
615 * The call to constrain() for SD1 that passes step_scaling =
616 * 0.5 also passes vir = NULL, so cannot reach this
617 * assertion. This assertion should remain until someone knows
618 * that this path works for their intended purpose, and then
619 * they can use scaled_delta_t instead of ir.delta_t
621 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
624 case ConstraintVariable::Positions: vir_fac = 0.5 / (ir.delta_t * ir.delta_t); break;
625 case ConstraintVariable::Velocities: vir_fac = 0.5 / ir.delta_t; break;
626 case ConstraintVariable::Force:
627 case ConstraintVariable::ForceDispl: vir_fac = 0.5; break;
628 default: gmx_incons("Unsupported constraint quantity for virial");
633 vir_fac *= 2; /* only constraining over half the distance here */
635 for (int i = 0; i < DIM; i++)
637 for (int j = 0; j < DIM; j++)
639 (*vir)[i][j] = vir_fac * vir_r_m_dr[i][j];
646 dump_confs(log, step, mtop, start, homenr, cr, x, xprime, box);
649 if (econq == ConstraintVariable::Positions)
651 if (ir.bPull && pull_have_constraint(pull_work))
653 if (EI_DYNAMICS(ir.eI))
655 t = ir.init_t + (step + delta_step) * ir.delta_t;
661 set_pbc(&pbc, ir.ePBC, box);
662 pull_constraint(pull_work, &md, &pbc, cr, ir.delta_t, t, x, xprime, v, *vir);
664 if (ed && delta_step > 0)
666 /* apply the essential dynamics constraints here */
667 do_edsam(&ir, step, cr, xprime, v, box, ed);
670 wallcycle_stop(wcycle, ewcCONSTR);
672 if (v != nullptr && md.cFREEZE)
674 /* Set the velocities of frozen dimensions to zero */
676 int gmx_unused numThreads = gmx_omp_nthreads_get(emntUpdate);
678 #pragma omp parallel for num_threads(numThreads) schedule(static)
679 for (int i = 0; i < md.homenr; i++)
681 int freezeGroup = md.cFREEZE[i];
683 for (int d = 0; d < DIM; d++)
685 if (ir.opts.nFreeze[freezeGroup][d])
696 ArrayRef<real> Constraints::rmsdData() const
700 return lincs_rmsdData(impl_->lincsd);
708 real Constraints::rmsd() const
712 return lincs_rmsd(impl_->lincsd);
720 int Constraints::numConstraintsTotal()
722 return impl_->ncon_tot;
725 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
727 if (haveDynamicsIntegrator)
729 return FlexibleConstraintTreatment::Include;
733 return FlexibleConstraintTreatment::Exclude;
737 /*! \brief Returns a block struct to go from atoms to constraints
739 * The block struct will contain constraint indices with lower indices
740 * directly matching the order in F_CONSTR and higher indices matching
741 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
743 * \param[in] numAtoms The number of atoms to construct the list for
744 * \param[in] ilists The interaction lists, size F_NRE
745 * \param[in] iparams Interaction parameters, can be null when
746 * flexibleConstraintTreatment=Include \param[in] flexibleConstraintTreatment The flexible
747 * constraint treatment, see enum above \returns a block struct with all constraints for each atom
750 static t_blocka makeAtomsToConstraintsList(int numAtoms,
752 const t_iparams* iparams,
753 FlexibleConstraintTreatment flexibleConstraintTreatment)
755 GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || iparams != nullptr,
756 "With flexible constraint detection we need valid iparams");
758 std::vector<int> count(numAtoms);
760 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
762 const T& ilist = ilists[ftype];
763 const int stride = 1 + NRAL(ftype);
764 for (int i = 0; i < ilist.size(); i += stride)
766 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
767 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
769 for (int j = 1; j < 3; j++)
771 int a = ilist.iatoms[i + j];
779 at2con.nr = numAtoms;
780 at2con.nalloc_index = at2con.nr + 1;
781 snew(at2con.index, at2con.nalloc_index);
783 for (int a = 0; a < numAtoms; a++)
785 at2con.index[a + 1] = at2con.index[a] + count[a];
788 at2con.nra = at2con.index[at2con.nr];
789 at2con.nalloc_a = at2con.nra;
790 snew(at2con.a, at2con.nalloc_a);
792 /* The F_CONSTRNC constraints have constraint numbers
793 * that continue after the last F_CONSTR constraint.
795 int numConstraints = 0;
796 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
798 const T& ilist = ilists[ftype];
799 const int stride = 1 + NRAL(ftype);
800 for (int i = 0; i < ilist.size(); i += stride)
802 if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
803 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
805 for (int j = 1; j < 3; j++)
807 int a = ilist.iatoms[i + j];
808 at2con.a[at2con.index[a] + count[a]++] = numConstraints;
818 t_blocka make_at2con(int numAtoms,
819 const t_ilist* ilist,
820 const t_iparams* iparams,
821 FlexibleConstraintTreatment flexibleConstraintTreatment)
823 return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
826 t_blocka make_at2con(const gmx_moltype_t& moltype,
827 gmx::ArrayRef<const t_iparams> iparams,
828 FlexibleConstraintTreatment flexibleConstraintTreatment)
830 return makeAtomsToConstraintsList(moltype.atoms.nr, moltype.ilist.data(), iparams.data(),
831 flexibleConstraintTreatment);
834 //! Return the number of flexible constraints in the \c ilist and \c iparams.
836 static int countFlexibleConstraintsTemplate(const T* ilist, const t_iparams* iparams)
839 for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
841 const int numIatomsPerConstraint = 3;
842 for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
844 const int type = ilist[ftype].iatoms[i];
845 if (iparams[type].constr.dA == 0 && iparams[type].constr.dB == 0)
855 int countFlexibleConstraints(const t_ilist* ilist, const t_iparams* iparams)
857 return countFlexibleConstraintsTemplate(ilist, iparams);
860 //! Returns the index of the settle to which each atom belongs.
861 static std::vector<int> make_at2settle(int natoms, const InteractionList& ilist)
863 /* Set all to no settle */
864 std::vector<int> at2s(natoms, -1);
866 const int stride = 1 + NRAL(F_SETTLE);
868 for (int s = 0; s < ilist.size(); s += stride)
870 at2s[ilist.iatoms[s + 1]] = s / stride;
871 at2s[ilist.iatoms[s + 2]] = s / stride;
872 at2s[ilist.iatoms[s + 3]] = s / stride;
878 void Constraints::Impl::setConstraints(const gmx_localtop_t& top, const t_mdatoms& md)
884 /* With DD we might also need to call LINCS on a domain no constraints for
885 * communicating coordinates to other nodes that do have constraints.
887 if (ir.eConstrAlg == econtLINCS)
889 set_lincs(top.idef, md, EI_DYNAMICS(ir.eI), cr, lincsd);
891 if (ir.eConstrAlg == econtSHAKE)
895 // We are using the local topology, so there are only
896 // F_CONSTR constraints.
897 make_shake_sblock_dd(shaked, &idef->il[F_CONSTR], cr->dd);
901 make_shake_sblock_serial(shaked, idef, md);
908 settle_set_constraints(settled, &idef->il[F_SETTLE], md);
911 /* Make a selection of the local atoms for essential dynamics */
914 dd_make_local_ed_indices(cr->dd, ed);
918 void Constraints::setConstraints(const gmx_localtop_t& top, const t_mdatoms& md)
920 impl_->setConstraints(top, md);
923 /*! \brief Makes a per-moleculetype container of mappings from atom
924 * indices to constraint indices.
926 * Note that flexible constraints are only enabled with a dynamical integrator. */
927 static std::vector<t_blocka> makeAtomToConstraintMappings(const gmx_mtop_t& mtop,
928 FlexibleConstraintTreatment flexibleConstraintTreatment)
930 std::vector<t_blocka> mapping;
931 mapping.reserve(mtop.moltype.size());
932 for (const gmx_moltype_t& moltype : mtop.moltype)
934 mapping.push_back(make_at2con(moltype, mtop.ffparams.iparams, flexibleConstraintTreatment));
939 Constraints::Constraints(const gmx_mtop_t& mtop,
940 const t_inputrec& ir,
945 const gmx_multisim_t* ms,
947 gmx_wallcycle* wcycle,
948 bool pbcHandlingRequired,
951 impl_(new Impl(mtop, ir, pull_work, log, md, cr, ms, nrnb, wcycle, pbcHandlingRequired, numConstraints, numSettles))
955 Constraints::Impl::Impl(const gmx_mtop_t& mtop_p,
956 const t_inputrec& ir_p,
959 const t_mdatoms& md_p,
960 const t_commrec* cr_p,
961 const gmx_multisim_t* ms_p,
963 gmx_wallcycle* wcycle_p,
964 bool pbcHandlingRequired,
967 ncon_tot(numConstraints),
970 pbcHandlingRequired_(pbcHandlingRequired),
974 pull_work(pull_work),
979 if (numConstraints + numSettles > 0 && ir.epc == epcMTTK)
981 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
985 if (numConstraints > 0)
987 at2con_mt = makeAtomToConstraintMappings(mtop, flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
989 for (const gmx_molblock_t& molblock : mtop.molblock)
991 int count = countFlexibleConstraintsTemplate(mtop.moltype[molblock.type].ilist.data(),
992 mtop.ffparams.iparams.data());
993 nflexcon += molblock.nmol * count;
1000 fprintf(log, "There are %d flexible constraints\n", nflexcon);
1001 if (ir.fc_stepsize == 0)
1005 "WARNING: step size for flexible constraining = 0\n"
1006 " All flexible constraints will be rigid.\n"
1007 " Will try to keep all flexible constraints at their original "
1009 " but the lengths may exhibit some drift.\n\n");
1015 please_cite(log, "Hess2002");
1019 if (ir.eConstrAlg == econtLINCS)
1021 lincsd = init_lincs(log, mtop, nflexcon, at2con_mt,
1022 DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd), ir.nLincsIter,
1026 if (ir.eConstrAlg == econtSHAKE)
1028 if (DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd))
1031 "SHAKE is not supported with domain decomposition and constraint that "
1032 "cross domain boundaries, use LINCS");
1037 "For this system also velocities and/or forces need to be constrained, "
1038 "this can not be done with SHAKE, you should select LINCS");
1040 please_cite(log, "Ryckaert77a");
1043 please_cite(log, "Barth95a");
1046 shaked = shake_init();
1052 please_cite(log, "Miyamoto92a");
1054 settled = settle_init(mtop);
1056 /* Make an atom to settle index for use in domain decomposition */
1057 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1059 at2settle_mt.emplace_back(
1060 make_at2settle(mtop.moltype[mt].atoms.nr, mtop.moltype[mt].ilist[F_SETTLE]));
1063 /* Allocate thread-local work arrays */
1064 int nthreads = gmx_omp_nthreads_get(emntSETTLE);
1065 if (nthreads > 1 && vir_r_m_dr_th == nullptr)
1067 snew(vir_r_m_dr_th, nthreads);
1068 snew(bSettleErrorHasOccurred, nthreads);
1073 char* env = getenv("GMX_MAXCONSTRWARN");
1077 sscanf(env, "%8d", &maxwarn);
1084 fprintf(log, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1088 fprintf(stderr, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1091 warncount_lincs = 0;
1092 warncount_settle = 0;
1095 Constraints::Impl::~Impl()
1097 for (auto blocka : at2con_mt)
1099 done_blocka(&blocka);
1101 if (bSettleErrorHasOccurred != nullptr)
1103 sfree(bSettleErrorHasOccurred);
1105 if (vir_r_m_dr_th != nullptr)
1107 sfree(vir_r_m_dr_th);
1109 if (settled != nullptr)
1111 settle_free(settled);
1116 void Constraints::saveEdsamPointer(gmx_edsam* ed)
1121 ArrayRef<const t_blocka> Constraints::atom2constraints_moltype() const
1123 return impl_->at2con_mt;
1126 ArrayRef<const std::vector<int>> Constraints::atom2settle_moltype() const
1128 return impl_->at2settle_mt;
1131 void do_constrain_first(FILE* fplog,
1132 gmx::Constraints* constr,
1133 const t_inputrec* ir,
1134 const t_mdatoms* md,
1136 ArrayRefWithPadding<RVec> x,
1137 ArrayRefWithPadding<RVec> v,
1141 int i, m, start, end;
1143 real dt = ir->delta_t;
1147 auto xRvec = as_rvec_array(x.paddedArrayRef().data());
1148 auto vRvec = as_rvec_array(v.paddedArrayRef().data());
1150 /* We need to allocate one element extra, since we might use
1151 * (unaligned) 4-wide SIMD loads to access rvec entries.
1153 snew(savex, natoms + 1);
1160 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n", start, md->homenr, end);
1162 /* Do a first constrain to reset particles... */
1163 step = ir->init_step;
1166 char buf[STEPSTRSIZE];
1167 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n", gmx_step_str(step, buf));
1171 /* constrain the current position */
1172 constr->apply(TRUE, FALSE, step, 0, 1.0, xRvec, xRvec, nullptr, box, lambda, &dvdl_dum, nullptr,
1173 nullptr, gmx::ConstraintVariable::Positions);
1176 /* constrain the inital velocity, and save it */
1177 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1178 constr->apply(TRUE, FALSE, step, 0, 1.0, xRvec, vRvec, vRvec, box, lambda, &dvdl_dum,
1179 nullptr, nullptr, gmx::ConstraintVariable::Velocities);
1181 /* constrain the inital velocities at t-dt/2 */
1182 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
1184 auto subX = x.paddedArrayRef().subArray(start, end);
1185 auto subV = v.paddedArrayRef().subArray(start, end);
1186 for (i = start; (i < end); i++)
1188 for (m = 0; (m < DIM); m++)
1190 /* Reverse the velocity */
1191 subV[i][m] = -subV[i][m];
1192 /* Store the position at t-dt in buf */
1193 savex[i][m] = subX[i][m] + dt * subV[i][m];
1196 /* Shake the positions at t=-dt with the positions at t=0
1197 * as reference coordinates.
1201 char buf[STEPSTRSIZE];
1202 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n", gmx_step_str(step, buf));
1205 constr->apply(TRUE, FALSE, step, -1, 1.0, xRvec, savex, nullptr, box, lambda, &dvdl_dum,
1206 vRvec, nullptr, gmx::ConstraintVariable::Positions);
1208 for (i = start; i < end; i++)
1210 for (m = 0; m < DIM; m++)
1212 /* Re-reverse the velocities */
1213 subV[i][m] = -subV[i][m];