Add second LINCS atom update task
[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/real.h"
58
59 struct gmx_edsam;
60 struct gmx_localtop_t;
61 struct gmx_moltype_t;
62 struct gmx_mtop_t;
63 struct gmx_multisim_t;
64 struct gmx_wallcycle;
65 struct pull_t;
66 struct t_commrec;
67 struct t_ilist;
68 struct t_inputrec;
69 struct t_nrnb;
70 struct t_pbc;
71 class t_state;
72 enum class ConstraintAlgorithm : int;
73 namespace gmx
74 {
75 template<typename T>
76 class ArrayRefWithPadding;
77 template<typename>
78 class ListOfLists;
79 class ObservablesReducerBuilder;
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                 bool                       useUpdateGroups,
111                 const gmx_multisim_t*      ms,
112                 t_nrnb*                    nrnb,
113                 gmx_wallcycle*             wcycle,
114                 bool                       pbcHandlingRequired,
115                 ObservablesReducerBuilder* observablesReducerBuilder,
116                 int                        numConstraints,
117                 int                        numSettles);
118
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;
124
125     ~Constraints();
126
127     /*! \brief Returns the total number of flexible constraints in the system. */
128     int numFlexibleConstraints() const;
129
130     /*! \brief Returns whether the system contains perturbed constraints */
131     bool havePerturbedConstraints() const;
132
133     /*! \brief Set up all the local constraints for the domain.
134      *
135      * \todo Make this a callback that is called automatically
136      * once a new domain has been made. */
137     void setConstraints(gmx_localtop_t*                     top,
138                         int                                 numAtoms,
139                         int                                 numHomeAtoms,
140                         gmx::ArrayRef<const real>           masses,
141                         gmx::ArrayRef<const real>           inverseMasses,
142                         bool                                hasMassPerturbedAtoms,
143                         real                                lambda,
144                         gmx::ArrayRef<const unsigned short> cFREEZE);
145
146     /*! \brief Applies constraints to coordinates.
147      *
148      * When econq=ConstraintVariable::Positions constrains
149      * coordinates xprime using the directions in x, min_proj is
150      * not used.
151      *
152      * When econq=ConstraintVariable::Derivative, calculates the
153      * components xprime in the constraint directions and
154      * subtracts these components from min_proj.  So when
155      * min_proj=xprime, the constraint components are projected
156      * out.
157      *
158      * When econq=ConstraintVariable::Deriv_FlexCon, the same is
159      * done as with ConstraintVariable::Derivative, but only the
160      * components of the flexible constraints are stored.
161      *
162      * delta_step is used for determining the constraint reference lengths
163      * when lenA != lenB or will the pull code with a pulling rate.
164      * step + delta_step is the step at which the final configuration
165      * is meant to be; for update delta_step = 1.
166      *
167      * step_scaling can be used to update coordinates based on the time
168      * step multiplied by this factor. Thus, normally 1.0 is passed. The
169      * SD1 integrator uses 0.5 in one of its calls, to correct positions
170      * for half a step of changed velocities.
171      *
172      * If v!=NULL also constrain v by adding the constraint corrections / dt.
173      *
174      * If computeVirial is true, calculate the constraint virial.
175      *
176      * Return whether the application of constraints succeeded without error.
177      *
178      * /note x is non-const, because non-local atoms need to be communicated.
179      */
180     bool apply(bool                      bLog,
181                bool                      bEner,
182                int64_t                   step,
183                int                       delta_step,
184                real                      step_scaling,
185                ArrayRefWithPadding<RVec> x,
186                ArrayRefWithPadding<RVec> xprime,
187                ArrayRef<RVec>            min_proj,
188                const matrix              box,
189                real                      lambda,
190                real*                     dvdlambda,
191                ArrayRefWithPadding<RVec> v,
192                bool                      computeVirial,
193                tensor                    constraintsVirial,
194                ConstraintVariable        econq);
195     //! Links the essentialdynamics and constraint code.
196     void saveEdsamPointer(gmx_edsam* ed);
197     //! Getter for use by domain decomposition.
198     ArrayRef<const ListOfLists<int>> atom2constraints_moltype() const;
199     //! Getter for use by domain decomposition.
200     ArrayRef<const std::vector<int>> atom2settle_moltype() const;
201
202     /*! \brief Return the RMSD of the constraints when available. */
203     real rmsd() const;
204
205     /*! \brief Get the total number of constraints.
206      *
207      * \returns Total number of constraints, including SETTLE and LINCS/SHAKE constraints.
208      */
209     int numConstraintsTotal();
210
211 private:
212     //! Implementation type.
213     class Impl;
214     //! Implementation object.
215     std::unique_ptr<Impl> impl_;
216 };
217
218 /*! \brief Generate a fatal error because of too many LINCS/SETTLE warnings. */
219 [[noreturn]] void too_many_constraint_warnings(ConstraintAlgorithm eConstrAlg, int warncount);
220
221 /*! \brief Returns whether constraint with parameter \p iparamsIndex is a flexible constraint */
222 static inline bool isConstraintFlexible(ArrayRef<const t_iparams> iparams, int iparamsIndex)
223 {
224     return (iparams[iparamsIndex].constr.dA == 0 && iparams[iparamsIndex].constr.dB == 0);
225 };
226
227 /* The at2con t_blocka struct returned by the routines below
228  * contains a list of constraints per atom.
229  * The F_CONSTRNC constraints in this structure number consecutively
230  * after the F_CONSTR constraints.
231  */
232
233 /*! \brief Tells make_at2con how to treat flexible constraints */
234 enum class FlexibleConstraintTreatment
235 {
236     Include, //!< Include all flexible constraints
237     Exclude  //!< Exclude all flexible constraints
238 };
239
240 /*! \brief Returns the flexible constraint treatment depending on whether the integrator is dynamic */
241 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator);
242
243 /*! \brief Returns a ListOfLists object to go from atoms to constraints
244  *
245  * The object will contain constraint indices with lower indices
246  * directly matching the order in F_CONSTR and higher indices matching
247  * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
248  *
249  * \param[in]  moltype   The molecule data
250  * \param[in]  iparams   Interaction parameters, can be null when
251  *                       \p flexibleConstraintTreatment==Include
252  * \param[in]  flexibleConstraintTreatment  The flexible constraint treatment,
253  *                                          see enum above
254  *
255  * \returns a ListOfLists object with all constraints for each atom
256  */
257 ListOfLists<int> make_at2con(const gmx_moltype_t&           moltype,
258                              gmx::ArrayRef<const t_iparams> iparams,
259                              FlexibleConstraintTreatment    flexibleConstraintTreatment);
260
261 /*! \brief Returns a ListOfLists object to go from atoms to constraints
262  *
263  * The object will contain constraint indices with lower indices
264  * directly matching the order in F_CONSTR and higher indices matching
265  * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
266  *
267  * \param[in]  numAtoms  The number of atoms to construct the list for
268  * \param[in]  ilist     Interaction list, size F_NRE
269  * \param[in]  iparams   Interaction parameters, can be null when
270  *                       \p flexibleConstraintTreatment==Include
271  * \param[in]  flexibleConstraintTreatment  The flexible constraint treatment,
272  *                                          see enum above
273  *
274  * \returns a ListOfLists object with all constraints for each atom
275  */
276 ListOfLists<int> make_at2con(int                             numAtoms,
277                              ArrayRef<const InteractionList> ilist,
278                              ArrayRef<const t_iparams>       iparams,
279                              FlexibleConstraintTreatment     flexibleConstraintTreatment);
280
281 //! Return the number of flexible constraints in the \c ilist and \c iparams.
282 int countFlexibleConstraints(ArrayRef<const InteractionList> ilist, ArrayRef<const t_iparams> iparams);
283
284 /*! \brief Returns the constraint iatoms for a constraint number con
285  * which comes from a list where F_CONSTR and F_CONSTRNC constraints
286  * are concatenated. */
287 inline const int* constr_iatomptr(gmx::ArrayRef<const int> iatom_constr,
288                                   gmx::ArrayRef<const int> iatom_constrnc,
289                                   int                      con)
290 {
291     if (con * 3 < iatom_constr.ssize())
292     {
293         return iatom_constr.data() + con * 3;
294     }
295     else
296     {
297         return iatom_constrnc.data() + con * 3 - iatom_constr.size();
298     }
299 };
300
301 /*! \brief Returns whether there are inter charge group constraints */
302 bool inter_charge_group_constraints(const gmx_mtop_t& mtop);
303
304 /*! \brief Returns whether there are inter charge group settles */
305 bool inter_charge_group_settles(const gmx_mtop_t& mtop);
306
307 /*! \brief Constrain the initial coordinates and velocities */
308 void do_constrain_first(FILE*                     log,
309                         gmx::Constraints*         constr,
310                         const t_inputrec*         inputrec,
311                         int                       numAtoms,
312                         int                       numHomeAtoms,
313                         ArrayRefWithPadding<RVec> x,
314                         ArrayRefWithPadding<RVec> v,
315                         const matrix              box,
316                         real                      lambda);
317
318 /*! \brief Constrain velocities only.
319  *
320  * The dvdlambda contribution has to be added to the bonded interactions
321  */
322 void constrain_velocities(gmx::Constraints* constr,
323                           bool              do_log,
324                           bool              do_ene,
325                           int64_t           step,
326                           t_state*          state,
327                           real*             dvdlambda,
328                           bool              computeVirial,
329                           tensor            constraintsVirial);
330
331 /*! \brief Constraint coordinates.
332  *
333  * The dvdlambda contribution has to be added to the bonded interactions
334  */
335 void constrain_coordinates(gmx::Constraints*         constr,
336                            bool                      do_log,
337                            bool                      do_ene,
338                            int64_t                   step,
339                            t_state*                  state,
340                            ArrayRefWithPadding<RVec> xp,
341                            real*                     dvdlambda,
342                            bool                      computeVirial,
343                            tensor                    constraintsVirial);
344
345 } // namespace gmx
346
347 #endif