f2f768dc7b3747f38c744a25f100df6d363095b3
[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,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \libinternal \file
38  * \brief Declares interface to constraint code.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \author Mark Abraham <mark.j.abraham@gmail.com>
42  * \ingroup module_mdlib
43  * \inlibraryapi
44  */
45
46 #ifndef GMX_MDLIB_CONSTR_H
47 #define GMX_MDLIB_CONSTR_H
48
49 #include <cstdio>
50
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"
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 t_blocka;
66 struct t_commrec;
67 struct t_ilist;
68 struct t_inputrec;
69 struct t_mdatoms;
70 struct t_nrnb;
71 struct t_pbc;
72 class t_state;
73
74 namespace gmx
75 {
76
77 //! Describes supported flavours of constrained updates.
78 enum class ConstraintVariable : int
79 {
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                   */
90 };
91
92 /*! \libinternal
93  * \brief Handles constraints */
94 class Constraints
95 {
96     private:
97         /*! \brief Constructor
98          *
99          * Private to enforce use of makeConstraints() factory
100          * function. */
101         Constraints(const gmx_mtop_t     &mtop,
102                     const t_inputrec     &ir,
103                     FILE                 *log,
104                     const t_mdatoms      &md,
105                     const t_commrec      *cr,
106                     const gmx_multisim_t &ms,
107                     t_nrnb               *nrnb,
108                     gmx_wallcycle        *wcycle,
109                     bool                  pbcHandlingRequired,
110                     int                   numConstraints,
111                     int                   numSettles);
112     public:
113         /*! \brief This member type helps implement a factory
114          * function, because its objects can access the private
115          * constructor. */
116         struct CreationHelper;
117
118         ~Constraints();
119
120         /*! \brief Returns the total number of flexible constraints in the system. */
121         int numFlexibleConstraints() const;
122
123         /*! \brief Set up all the local constraints for the domain.
124          *
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);
129
130         /*! \brief Applies constraints to coordinates.
131          *
132          * When econq=ConstraintVariable::Positions constrains
133          * coordinates xprime using th directions in x, min_proj is
134          * not used.
135          *
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
140          * out.
141          *
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.
145          *
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.
150          *
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.
155          *
156          * If v!=NULL also constrain v by adding the constraint corrections / dt.
157          *
158          * If vir!=NULL calculate the constraint virial.
159          *
160          * Return whether the application of constraints succeeded without error.
161          */
162         bool apply(bool                  bLog,
163                    bool                  bEner,
164                    int64_t               step,
165                    int                   delta_step,
166                    real                  step_scaling,
167                    rvec                 *x,
168                    rvec                 *xprime,
169                    rvec                 *min_proj,
170                    matrix                box,
171                    real                  lambda,
172                    real                 *dvdlambda,
173                    rvec                 *v,
174                    tensor               *vir,
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;
182
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. */
188         real rmsd() const;
189
190     private:
191         //! Implementation type.
192         class Impl;
193         //! Implementation object.
194         PrivateImplPointer<Impl> impl_;
195 };
196
197 /*! \brief Generate a fatal error because of too many LINCS/SETTLE warnings. */
198 [[ noreturn ]] void too_many_constraint_warnings(int eConstrAlg, int warncount);
199
200 /*! \brief Returns whether constraint with parameter \p iparamsIndex is a flexible constraint */
201 static inline bool isConstraintFlexible(const t_iparams *iparams,
202                                         int              iparamsIndex)
203 {
204     GMX_ASSERT(iparams != nullptr, "Need a valid iparams array");
205
206     return (iparams[iparamsIndex].constr.dA == 0 &&
207             iparams[iparamsIndex].constr.dB == 0);
208 };
209
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.
214  */
215
216 /*! \brief Tells make_at2con how to treat flexible constraints */
217 enum class FlexibleConstraintTreatment
218 {
219     Include, //!< Include all flexible constraints
220     Exclude  //!< Exclude all flexible constraints
221 };
222
223 /*! \brief Returns the flexible constraint treatment depending on whether the integrator is dynamic */
224 FlexibleConstraintTreatment
225 flexibleConstraintTreatment(bool haveDynamicsIntegrator);
226
227 /*! \brief Returns a block struct to go from atoms to constraints
228  *
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.
232  *
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
237  */
238 t_blocka make_at2con(const gmx_moltype_t         &moltype,
239                      const t_iparams             *iparams,
240                      FlexibleConstraintTreatment  flexibleConstraintTreatment);
241
242 /*! \brief Returns a block struct to go from atoms to constraints
243  *
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.
247  *
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
253  */
254 t_blocka make_at2con(int                          numAtoms,
255                      const t_ilist               *ilist,
256                      const t_iparams             *iparams,
257                      FlexibleConstraintTreatment  flexibleConstraintTreatment);
258
259 /*! \brief Returns an array of atom to constraints lists for the moltypes */
260 const t_blocka *atom2constraints_moltype(const Constraints *constr);
261
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);
265
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. */
269 inline const int *
270 constr_iatomptr(gmx::ArrayRef<const int> iatom_constr,
271                 gmx::ArrayRef<const int> iatom_constrnc,
272                 int                      con)
273 {
274     if (con*3 < iatom_constr.size())
275     {
276         return iatom_constr.data() + con*3;
277     }
278     else
279     {
280         return iatom_constrnc.data() + con*3 - iatom_constr.size();
281     }
282 };
283
284 /*! \brief Returns whether there are inter charge group constraints */
285 bool inter_charge_group_constraints(const gmx_mtop_t &mtop);
286
287 /*! \brief Returns whether there are inter charge group settles */
288 bool inter_charge_group_settles(const gmx_mtop_t &mtop);
289
290 }  // namespace gmx
291
292 #endif