Minor cleanup in constraint and shellfc code
[gromacs.git] / src / gromacs / mdlib / constr.h
blob216ea479b64c47c867103b6ff627c9045c8a6d41
1 /*
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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /*! \libinternal \file
39 * \brief Declares interface to constraint code.
41 * \author Berk Hess <hess@kth.se>
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
43 * \ingroup module_mdlib
44 * \inlibraryapi
47 #ifndef GMX_MDLIB_CONSTR_H
48 #define GMX_MDLIB_CONSTR_H
50 #include <cstdio>
52 #include "gromacs/math/vectypes.h"
53 #include "gromacs/topology/idef.h"
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/classhelpers.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/real.h"
60 struct gmx_edsam;
61 struct gmx_localtop_t;
62 struct gmx_moltype_t;
63 struct gmx_mtop_t;
64 struct gmx_multisim_t;
65 struct gmx_wallcycle;
66 struct pull_t;
67 struct t_commrec;
68 struct t_ilist;
69 struct t_inputrec;
70 struct t_mdatoms;
71 struct t_nrnb;
72 struct t_pbc;
73 class t_state;
75 namespace gmx
77 template<typename T>
78 class ArrayRefWithPadding;
79 template<typename>
80 class ListOfLists;
82 //! Describes supported flavours of constrained updates.
83 enum class ConstraintVariable : int
85 Positions, /* Constrain positions (mass weighted) */
86 Velocities, /* Constrain velocities (mass weighted) */
87 Derivative, /* Constrain a derivative (mass weighted), *
88 * for instance velocity or acceleration, *
89 * constraint virial can not be calculated. */
90 Deriv_FlexCon, /* As Derivative, but only output flex. con. */
91 Force, /* Constrain forces (non mass-weighted) */
92 // TODO What does this do? Improve the comment.
93 ForceDispl /* Like Force, but free particles will have mass
94 * 1 and frozen particles mass 0 */
97 /*! \libinternal
98 * \brief Handles constraints */
99 class Constraints
101 private:
102 /*! \brief Constructor
104 * Private to enforce use of makeConstraints() factory
105 * function. */
106 Constraints(const gmx_mtop_t& mtop,
107 const t_inputrec& ir,
108 pull_t* pull_work,
109 FILE* log,
110 const t_mdatoms& md,
111 const t_commrec* cr,
112 const gmx_multisim_t* ms,
113 t_nrnb* nrnb,
114 gmx_wallcycle* wcycle,
115 bool pbcHandlingRequired,
116 int numConstraints,
117 int numSettles);
119 public:
120 /*! \brief This member type helps implement a factory
121 * function, because its objects can access the private
122 * constructor. */
123 struct CreationHelper;
125 ~Constraints();
127 /*! \brief Returns the total number of flexible constraints in the system. */
128 int numFlexibleConstraints() const;
130 /*! \brief Returns whether the system contains perturbed constraints */
131 bool havePerturbedConstraints() const;
133 /*! \brief Set up all the local constraints for the domain.
135 * \todo Make this a callback that is called automatically
136 * once a new domain has been made. */
137 void setConstraints(const gmx_localtop_t& top, const t_mdatoms& md);
139 /*! \brief Applies constraints to coordinates.
141 * When econq=ConstraintVariable::Positions constrains
142 * coordinates xprime using the directions in x, min_proj is
143 * not used.
145 * When econq=ConstraintVariable::Derivative, calculates the
146 * components xprime in the constraint directions and
147 * subtracts these components from min_proj. So when
148 * min_proj=xprime, the constraint components are projected
149 * out.
151 * When econq=ConstraintVariable::Deriv_FlexCon, the same is
152 * done as with ConstraintVariable::Derivative, but only the
153 * components of the flexible constraints are stored.
155 * delta_step is used for determining the constraint reference lengths
156 * when lenA != lenB or will the pull code with a pulling rate.
157 * step + delta_step is the step at which the final configuration
158 * is meant to be; for update delta_step = 1.
160 * step_scaling can be used to update coordinates based on the time
161 * step multiplied by this factor. Thus, normally 1.0 is passed. The
162 * SD1 integrator uses 0.5 in one of its calls, to correct positions
163 * for half a step of changed velocities.
165 * If v!=NULL also constrain v by adding the constraint corrections / dt.
167 * If vir!=NULL calculate the constraint virial.
169 * Return whether the application of constraints succeeded without error.
171 * /note x is non-const, because non-local atoms need to be communicated.
173 bool apply(bool bLog,
174 bool bEner,
175 int64_t step,
176 int delta_step,
177 real step_scaling,
178 ArrayRefWithPadding<RVec> x,
179 ArrayRefWithPadding<RVec> xprime,
180 ArrayRef<RVec> min_proj,
181 const matrix box,
182 real lambda,
183 real* dvdlambda,
184 ArrayRefWithPadding<RVec> v,
185 tensor* vir,
186 ConstraintVariable econq);
187 //! Links the essentialdynamics and constraint code.
188 void saveEdsamPointer(gmx_edsam* ed);
189 //! Getter for use by domain decomposition.
190 ArrayRef<const ListOfLists<int>> atom2constraints_moltype() const;
191 //! Getter for use by domain decomposition.
192 ArrayRef<const std::vector<int>> atom2settle_moltype() const;
194 /*! \brief Return the data for reduction for determining
195 * constraint RMS relative deviations, or an empty ArrayRef
196 * when not supported for any active constraints. */
197 ArrayRef<real> rmsdData() const;
198 /*! \brief Return the RMSD of the constraints when available. */
199 real rmsd() const;
201 /*! \brief Get the total number of constraints.
203 * \returns Total number of constraints, including SETTLE and LINCS/SHAKE constraints.
205 int numConstraintsTotal();
207 private:
208 //! Implementation type.
209 class Impl;
210 //! Implementation object.
211 PrivateImplPointer<Impl> impl_;
214 /*! \brief Generate a fatal error because of too many LINCS/SETTLE warnings. */
215 [[noreturn]] void too_many_constraint_warnings(int eConstrAlg, int warncount);
217 /*! \brief Returns whether constraint with parameter \p iparamsIndex is a flexible constraint */
218 static inline bool isConstraintFlexible(const t_iparams* iparams, int iparamsIndex)
220 GMX_ASSERT(iparams != nullptr, "Need a valid iparams array");
222 return (iparams[iparamsIndex].constr.dA == 0 && iparams[iparamsIndex].constr.dB == 0);
225 /* The at2con t_blocka struct returned by the routines below
226 * contains a list of constraints per atom.
227 * The F_CONSTRNC constraints in this structure number consecutively
228 * after the F_CONSTR constraints.
231 /*! \brief Tells make_at2con how to treat flexible constraints */
232 enum class FlexibleConstraintTreatment
234 Include, //!< Include all flexible constraints
235 Exclude //!< Exclude all flexible constraints
238 /*! \brief Returns the flexible constraint treatment depending on whether the integrator is dynamic */
239 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator);
241 /*! \brief Returns a ListOfLists object to go from atoms to constraints
243 * The object will contain constraint indices with lower indices
244 * directly matching the order in F_CONSTR and higher indices matching
245 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
247 * \param[in] moltype The molecule data
248 * \param[in] iparams Interaction parameters, can be null when
249 * \p flexibleConstraintTreatment==Include
250 * \param[in] flexibleConstraintTreatment The flexible constraint treatment,
251 * see enum above
253 * \returns a ListOfLists object with all constraints for each atom
255 ListOfLists<int> make_at2con(const gmx_moltype_t& moltype,
256 gmx::ArrayRef<const t_iparams> iparams,
257 FlexibleConstraintTreatment flexibleConstraintTreatment);
259 /*! \brief Returns a ListOfLists object to go from atoms to constraints
261 * The object will contain constraint indices with lower indices
262 * directly matching the order in F_CONSTR and higher indices matching
263 * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
265 * \param[in] numAtoms The number of atoms to construct the list for
266 * \param[in] ilist Interaction list, size F_NRE
267 * \param[in] iparams Interaction parameters, can be null when
268 * \p flexibleConstraintTreatment==Include
269 * \param[in] flexibleConstraintTreatment The flexible constraint treatment,
270 * see enum above
272 * \returns a ListOfLists object with all constraints for each atom
274 ListOfLists<int> make_at2con(int numAtoms,
275 const t_ilist* ilist,
276 const t_iparams* iparams,
277 FlexibleConstraintTreatment flexibleConstraintTreatment);
279 //! Return the number of flexible constraints in the \c ilist and \c iparams.
280 int countFlexibleConstraints(const t_ilist* ilist, const t_iparams* iparams);
282 /*! \brief Returns the constraint iatoms for a constraint number con
283 * which comes from a list where F_CONSTR and F_CONSTRNC constraints
284 * are concatenated. */
285 inline const int* constr_iatomptr(gmx::ArrayRef<const int> iatom_constr,
286 gmx::ArrayRef<const int> iatom_constrnc,
287 int con)
289 if (con * 3 < iatom_constr.ssize())
291 return iatom_constr.data() + con * 3;
293 else
295 return iatom_constrnc.data() + con * 3 - iatom_constr.size();
299 /*! \brief Returns whether there are inter charge group constraints */
300 bool inter_charge_group_constraints(const gmx_mtop_t& mtop);
302 /*! \brief Returns whether there are inter charge group settles */
303 bool inter_charge_group_settles(const gmx_mtop_t& mtop);
305 /*! \brief Constrain the initial coordinates and velocities */
306 void do_constrain_first(FILE* log,
307 gmx::Constraints* constr,
308 const t_inputrec* inputrec,
309 const t_mdatoms* md,
310 int natoms,
311 ArrayRefWithPadding<RVec> x,
312 ArrayRefWithPadding<RVec> v,
313 const matrix box,
314 real lambda);
316 } // namespace gmx
318 #endif