Use ArrayRef in signatures of constraints and some pulling
[alexxy/gromacs.git] / src / gromacs / mdlib / constr.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2021, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /*! \libinternal \file
39  * \brief Declares interface to constraint code.
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \author Mark Abraham <mark.j.abraham@gmail.com>
43  * \ingroup module_mdlib
44  * \inlibraryapi
45  */
46
47 #ifndef GMX_MDLIB_CONSTR_H
48 #define GMX_MDLIB_CONSTR_H
49
50 #include <cstdio>
51
52 #include <memory>
53
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/topology/idef.h"
56 #include "gromacs/utility/arrayref.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/real.h"
59
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_nrnb;
71 struct t_pbc;
72 class t_state;
73 enum class ConstraintAlgorithm : int;
74 namespace gmx
75 {
76 template<typename T>
77 class ArrayRefWithPadding;
78 template<typename>
79 class ListOfLists;
80
81 //! Describes supported flavours of constrained updates.
82 enum class ConstraintVariable : int
83 {
84     Positions,     /* Constrain positions (mass weighted)             */
85     Velocities,    /* Constrain velocities (mass weighted)            */
86     Derivative,    /* Constrain a derivative (mass weighted),         *
87                     * for instance velocity or acceleration,          *
88                     * constraint virial can not be calculated.        */
89     Deriv_FlexCon, /* As Derivative, but only output flex. con.       */
90     Force,         /* Constrain forces (non mass-weighted)            */
91     // TODO What does this do? Improve the comment.
92     ForceDispl /* Like Force, but free particles will have mass
93                 * 1 and frozen particles mass 0                   */
94 };
95
96 /*! \libinternal
97  * \brief Handles constraints */
98 class Constraints
99 {
100 private:
101     /*! \brief Constructor
102      *
103      * Private to enforce use of makeConstraints() factory
104      * function. */
105     Constraints(const gmx_mtop_t&     mtop,
106                 const t_inputrec&     ir,
107                 pull_t*               pull_work,
108                 FILE*                 log,
109                 const t_commrec*      cr,
110                 const gmx_multisim_t* ms,
111                 t_nrnb*               nrnb,
112                 gmx_wallcycle*        wcycle,
113                 bool                  pbcHandlingRequired,
114                 int                   numConstraints,
115                 int                   numSettles);
116
117 public:
118     /*! \brief This member type helps implement a factory
119      * function, because its objects can access the private
120      * constructor. */
121     struct CreationHelper;
122
123     ~Constraints();
124
125     /*! \brief Returns the total number of flexible constraints in the system. */
126     int numFlexibleConstraints() const;
127
128     /*! \brief Returns whether the system contains perturbed constraints */
129     bool havePerturbedConstraints() const;
130
131     /*! \brief Set up all the local constraints for the domain.
132      *
133      * \todo Make this a callback that is called automatically
134      * once a new domain has been made. */
135     void setConstraints(gmx_localtop_t*                     top,
136                         int                                 numAtoms,
137                         int                                 numHomeAtoms,
138                         gmx::ArrayRef<const real>           masses,
139                         gmx::ArrayRef<const real>           inverseMasses,
140                         bool                                hasMassPerturbedAtoms,
141                         real                                lambda,
142                         gmx::ArrayRef<const unsigned short> cFREEZE);
143
144     /*! \brief Applies constraints to coordinates.
145      *
146      * When econq=ConstraintVariable::Positions constrains
147      * coordinates xprime using the directions in x, min_proj is
148      * not used.
149      *
150      * When econq=ConstraintVariable::Derivative, calculates the
151      * components xprime in the constraint directions and
152      * subtracts these components from min_proj.  So when
153      * min_proj=xprime, the constraint components are projected
154      * out.
155      *
156      * When econq=ConstraintVariable::Deriv_FlexCon, the same is
157      * done as with ConstraintVariable::Derivative, but only the
158      * components of the flexible constraints are stored.
159      *
160      * delta_step is used for determining the constraint reference lengths
161      * when lenA != lenB or will the pull code with a pulling rate.
162      * step + delta_step is the step at which the final configuration
163      * is meant to be; for update delta_step = 1.
164      *
165      * step_scaling can be used to update coordinates based on the time
166      * step multiplied by this factor. Thus, normally 1.0 is passed. The
167      * SD1 integrator uses 0.5 in one of its calls, to correct positions
168      * for half a step of changed velocities.
169      *
170      * If v!=NULL also constrain v by adding the constraint corrections / dt.
171      *
172      * If computeVirial is true, calculate the constraint virial.
173      *
174      * Return whether the application of constraints succeeded without error.
175      *
176      * /note x is non-const, because non-local atoms need to be communicated.
177      */
178     bool apply(bool                      bLog,
179                bool                      bEner,
180                int64_t                   step,
181                int                       delta_step,
182                real                      step_scaling,
183                ArrayRefWithPadding<RVec> x,
184                ArrayRefWithPadding<RVec> xprime,
185                ArrayRef<RVec>            min_proj,
186                const matrix              box,
187                real                      lambda,
188                real*                     dvdlambda,
189                ArrayRefWithPadding<RVec> v,
190                bool                      computeVirial,
191                tensor                    constraintsVirial,
192                ConstraintVariable        econq);
193     //! Links the essentialdynamics and constraint code.
194     void saveEdsamPointer(gmx_edsam* ed);
195     //! Getter for use by domain decomposition.
196     ArrayRef<const ListOfLists<int>> atom2constraints_moltype() const;
197     //! Getter for use by domain decomposition.
198     ArrayRef<const std::vector<int>> atom2settle_moltype() const;
199
200     /*! \brief Return the data for reduction for determining
201      * constraint RMS relative deviations, or an empty ArrayRef
202      * when not supported for any active constraints. */
203     ArrayRef<real> rmsdData() const;
204     /*! \brief Return the RMSD of the constraints when available. */
205     real rmsd() const;
206
207     /*! \brief Get the total number of constraints.
208      *
209      * \returns Total number of constraints, including SETTLE and LINCS/SHAKE constraints.
210      */
211     int numConstraintsTotal();
212
213 private:
214     //! Implementation type.
215     class Impl;
216     //! Implementation object.
217     std::unique_ptr<Impl> impl_;
218 };
219
220 /*! \brief Generate a fatal error because of too many LINCS/SETTLE warnings. */
221 [[noreturn]] void too_many_constraint_warnings(ConstraintAlgorithm eConstrAlg, int warncount);
222
223 /*! \brief Returns whether constraint with parameter \p iparamsIndex is a flexible constraint */
224 static inline bool isConstraintFlexible(ArrayRef<const t_iparams> iparams, int iparamsIndex)
225 {
226     return (iparams[iparamsIndex].constr.dA == 0 && iparams[iparamsIndex].constr.dB == 0);
227 };
228
229 /* The at2con t_blocka struct returned by the routines below
230  * contains a list of constraints per atom.
231  * The F_CONSTRNC constraints in this structure number consecutively
232  * after the F_CONSTR constraints.
233  */
234
235 /*! \brief Tells make_at2con how to treat flexible constraints */
236 enum class FlexibleConstraintTreatment
237 {
238     Include, //!< Include all flexible constraints
239     Exclude  //!< Exclude all flexible constraints
240 };
241
242 /*! \brief Returns the flexible constraint treatment depending on whether the integrator is dynamic */
243 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator);
244
245 /*! \brief Returns a ListOfLists object to go from atoms to constraints
246  *
247  * The object will contain constraint indices with lower indices
248  * directly matching the order in F_CONSTR and higher indices matching
249  * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
250  *
251  * \param[in]  moltype   The molecule data
252  * \param[in]  iparams   Interaction parameters, can be null when
253  *                       \p flexibleConstraintTreatment==Include
254  * \param[in]  flexibleConstraintTreatment  The flexible constraint treatment,
255  *                                          see enum above
256  *
257  * \returns a ListOfLists object with all constraints for each atom
258  */
259 ListOfLists<int> make_at2con(const gmx_moltype_t&           moltype,
260                              gmx::ArrayRef<const t_iparams> iparams,
261                              FlexibleConstraintTreatment    flexibleConstraintTreatment);
262
263 /*! \brief Returns a ListOfLists object to go from atoms to constraints
264  *
265  * The object will contain constraint indices with lower indices
266  * directly matching the order in F_CONSTR and higher indices matching
267  * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
268  *
269  * \param[in]  numAtoms  The number of atoms to construct the list for
270  * \param[in]  ilist     Interaction list, size F_NRE
271  * \param[in]  iparams   Interaction parameters, can be null when
272  *                       \p flexibleConstraintTreatment==Include
273  * \param[in]  flexibleConstraintTreatment  The flexible constraint treatment,
274  *                                          see enum above
275  *
276  * \returns a ListOfLists object with all constraints for each atom
277  */
278 ListOfLists<int> make_at2con(int                             numAtoms,
279                              ArrayRef<const InteractionList> ilist,
280                              ArrayRef<const t_iparams>       iparams,
281                              FlexibleConstraintTreatment     flexibleConstraintTreatment);
282
283 //! Return the number of flexible constraints in the \c ilist and \c iparams.
284 int countFlexibleConstraints(ArrayRef<const InteractionList> ilist, ArrayRef<const t_iparams> iparams);
285
286 /*! \brief Returns the constraint iatoms for a constraint number con
287  * which comes from a list where F_CONSTR and F_CONSTRNC constraints
288  * are concatenated. */
289 inline const int* constr_iatomptr(gmx::ArrayRef<const int> iatom_constr,
290                                   gmx::ArrayRef<const int> iatom_constrnc,
291                                   int                      con)
292 {
293     if (con * 3 < iatom_constr.ssize())
294     {
295         return iatom_constr.data() + con * 3;
296     }
297     else
298     {
299         return iatom_constrnc.data() + con * 3 - iatom_constr.size();
300     }
301 };
302
303 /*! \brief Returns whether there are inter charge group constraints */
304 bool inter_charge_group_constraints(const gmx_mtop_t& mtop);
305
306 /*! \brief Returns whether there are inter charge group settles */
307 bool inter_charge_group_settles(const gmx_mtop_t& mtop);
308
309 /*! \brief Constrain the initial coordinates and velocities */
310 void do_constrain_first(FILE*                     log,
311                         gmx::Constraints*         constr,
312                         const t_inputrec*         inputrec,
313                         int                       numAtoms,
314                         int                       numHomeAtoms,
315                         ArrayRefWithPadding<RVec> x,
316                         ArrayRefWithPadding<RVec> v,
317                         const matrix              box,
318                         real                      lambda);
319
320 /*! \brief Constrain velocities only.
321  *
322  * The dvdlambda contribution has to be added to the bonded interactions
323  */
324 void constrain_velocities(gmx::Constraints* constr,
325                           bool              do_log,
326                           bool              do_ene,
327                           int64_t           step,
328                           t_state*          state,
329                           real*             dvdlambda,
330                           bool              computeVirial,
331                           tensor            constraintsVirial);
332
333 /*! \brief Constraint coordinates.
334  *
335  * The dvdlambda contribution has to be added to the bonded interactions
336  */
337 void constrain_coordinates(gmx::Constraints*         constr,
338                            bool                      do_log,
339                            bool                      do_ene,
340                            int64_t                   step,
341                            t_state*                  state,
342                            ArrayRefWithPadding<RVec> xp,
343                            real*                     dvdlambda,
344                            bool                      computeVirial,
345                            tensor                    constraintsVirial);
346
347 } // namespace gmx
348
349 #endif