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, 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.
37 /*! \libinternal \file
38 * \brief Declares interface to constraint code.
40 * \author Berk Hess <hess@kth.se>
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_mdlib
46 #ifndef GMX_MDLIB_CONSTR_H
47 #define GMX_MDLIB_CONSTR_H
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/topology/idef.h"
53 #include "gromacs/utility/arrayref.h"
54 #include "gromacs/utility/basedefinitions.h"
55 #include "gromacs/utility/classhelpers.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/real.h"
60 struct gmx_localtop_t;
63 struct gmx_multisim_t;
77 //! Describes supported flavours of constrained updates.
78 enum class ConstraintVariable : int
80 Positions, /* Constrain positions (mass weighted) */
81 Velocities, /* Constrain velocities (mass weighted) */
82 Derivative, /* Constrain a derivative (mass weighted), *
83 * for instance velocity or acceleration, *
84 * constraint virial can not be calculated. */
85 Deriv_FlexCon, /* As Derivative, but only output flex. con. */
86 Force, /* Constrain forces (non mass-weighted) */
87 // TODO What does this do? Improve the comment.
88 ForceDispl /* Like Force, but free particles will have mass
89 * 1 and frozen particles mass 0 */
93 * \brief Handles constraints */
97 /*! \brief Constructor
99 * Private to enforce use of makeConstraints() factory
101 Constraints(const gmx_mtop_t &mtop,
102 const t_inputrec &ir,
106 const gmx_multisim_t *ms,
108 gmx_wallcycle *wcycle,
109 bool pbcHandlingRequired,
113 /*! \brief This member type helps implement a factory
114 * function, because its objects can access the private
116 struct CreationHelper;
120 /*! \brief Returns the total number of flexible constraints in the system. */
121 int numFlexibleConstraints() const;
123 /*! \brief Set up all the local constraints for the domain.
125 * \todo Make this a callback that is called automatically
126 * once a new domain has been made. */
127 void setConstraints(const gmx_localtop_t &top,
128 const t_mdatoms &md);
130 /*! \brief Applies constraints to coordinates.
132 * When econq=ConstraintVariable::Positions constrains
133 * coordinates xprime using th directions in x, min_proj is
136 * When econq=ConstraintVariable::Derivative, calculates the
137 * components xprime in the constraint directions and
138 * subtracts these components from min_proj. So when
139 * min_proj=xprime, the constraint components are projected
142 * When econq=ConstraintVariable::Deriv_FlexCon, the same is
143 * done as with ConstraintVariable::Derivative, but only the
144 * components of the flexible constraints are stored.
146 * delta_step is used for determining the constraint reference lengths
147 * when lenA != lenB or will the pull code with a pulling rate.
148 * step + delta_step is the step at which the final configuration
149 * is meant to be; for update delta_step = 1.
151 * step_scaling can be used to update coordinates based on the time
152 * step multiplied by this factor. Thus, normally 1.0 is passed. The
153 * SD1 integrator uses 0.5 in one of its calls, to correct positions
154 * for half a step of changed velocities.
156 * If v!=NULL also constrain v by adding the constraint corrections / dt.
158 * If vir!=NULL calculate the constraint virial.
160 * Return whether the application of constraints succeeded without error.
162 bool apply(bool bLog,
175 ConstraintVariable econq);
176 //! Links the essentialdynamics and constraint code.
177 void saveEdsamPointer(gmx_edsam *ed);
178 //! Getter for use by domain decomposition.
179 const ArrayRef<const t_blocka> atom2constraints_moltype() const;
180 //! Getter for use by domain decomposition.
181 ArrayRef < const std::vector < int>> atom2settle_moltype() const;
183 /*! \brief Return the data for reduction for determining
184 * constraint RMS relative deviations, or an empty ArrayRef
185 * when not supported for any active constraints. */
186 ArrayRef<real> rmsdData() const;
187 /*! \brief Return the RMSD of the constraints when available. */
191 //! Implementation type.
193 //! Implementation object.
194 PrivateImplPointer<Impl> impl_;
197 /*! \brief Generate a fatal error because of too many LINCS/SETTLE warnings. */
198 [[ noreturn ]] void too_many_constraint_warnings(int eConstrAlg, int warncount);
200 /*! \brief Returns whether constraint with parameter \p iparamsIndex is a flexible constraint */
201 static inline bool isConstraintFlexible(const t_iparams *iparams,
204 GMX_ASSERT(iparams != nullptr, "Need a valid iparams array");
206 return (iparams[iparamsIndex].constr.dA == 0 &&
207 iparams[iparamsIndex].constr.dB == 0);
210 /* The at2con t_blocka struct returned by the routines below
211 * contains a list of constraints per atom.
212 * The F_CONSTRNC constraints in this structure number consecutively
213 * after the F_CONSTR constraints.
216 /*! \brief Tells make_at2con how to treat flexible constraints */
217 enum class FlexibleConstraintTreatment
219 Include, //!< Include all flexible constraints
220 Exclude //!< Exclude all flexible constraints
223 /*! \brief Returns the flexible constraint treatment depending on whether the integrator is dynamic */
224 FlexibleConstraintTreatment
225 flexibleConstraintTreatment(bool haveDynamicsIntegrator);
227 /*! \brief Returns a block struct to go from atoms to constraints
229 * The block struct will contain constraint indices with lower indices
230 * directly matching the order in F_CONSTR and higher indices matching
231 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
233 * \param[in] moltype The molecule data
234 * \param[in] iparams Interaction parameters, can be null when flexibleConstraintTreatment=Include
235 * \param[in] flexibleConstraintTreatment The flexible constraint treatment, see enum above
236 * \returns a block struct with all constraints for each atom
238 t_blocka make_at2con(const gmx_moltype_t &moltype,
239 gmx::ArrayRef<const t_iparams> iparams,
240 FlexibleConstraintTreatment flexibleConstraintTreatment);
242 /*! \brief Returns a block struct to go from atoms to constraints
244 * The block struct will contain constraint indices with lower indices
245 * directly matching the order in F_CONSTR and higher indices matching
246 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
248 * \param[in] numAtoms The number of atoms to construct the list for
249 * \param[in] ilist Interaction list, size F_NRE
250 * \param[in] iparams Interaction parameters, can be null when flexibleConstraintTreatment=Include
251 * \param[in] flexibleConstraintTreatment The flexible constraint treatment, see enum above
252 * \returns a block struct with all constraints for each atom
254 t_blocka make_at2con(int numAtoms,
255 const t_ilist *ilist,
256 const t_iparams *iparams,
257 FlexibleConstraintTreatment flexibleConstraintTreatment);
259 /*! \brief Returns an array of atom to constraints lists for the moltypes */
260 const t_blocka *atom2constraints_moltype(const Constraints *constr);
262 //! Return the number of flexible constraints in the \c ilist and \c iparams.
263 int countFlexibleConstraints(const t_ilist *ilist,
264 const t_iparams *iparams);
266 /*! \brief Returns the constraint iatoms for a constraint number con
267 * which comes from a list where F_CONSTR and F_CONSTRNC constraints
268 * are concatenated. */
270 constr_iatomptr(gmx::ArrayRef<const int> iatom_constr,
271 gmx::ArrayRef<const int> iatom_constrnc,
274 if (con*3 < iatom_constr.size())
276 return iatom_constr.data() + con*3;
280 return iatom_constrnc.data() + con*3 - iatom_constr.size();
284 /*! \brief Returns whether there are inter charge group constraints */
285 bool inter_charge_group_constraints(const gmx_mtop_t &mtop);
287 /*! \brief Returns whether there are inter charge group settles */
288 bool inter_charge_group_settles(const gmx_mtop_t &mtop);