Eliminate mdlib/mdrun.h
[alexxy/gromacs.git] / src / gromacs / pulling / pull.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,2019, 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
38 /*! \libinternal \file
39  *
40  *
41  * \brief
42  * This file contains datatypes and function declarations necessary
43    for mdrun to interface with the pull code.
44  *
45  * \author Berk Hess
46  *
47  * \inlibraryapi
48  */
49
50 #ifndef GMX_PULLING_PULL_H
51 #define GMX_PULLING_PULL_H
52
53 #include <cstdio>
54
55 #include "gromacs/math/vectypes.h"
56 #include "gromacs/mdtypes/pull_params.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/basedefinitions.h"
59 #include "gromacs/utility/real.h"
60
61 struct ContinuationOptions;
62 struct gmx_mtop_t;
63 struct gmx_output_env_t;
64 struct pull_coord_work_t;
65 struct pull_params_t;
66 struct t_commrec;
67 struct t_filenm;
68 struct t_inputrec;
69 struct t_mdatoms;
70 struct t_pbc;
71 class t_state;
72
73 namespace gmx
74 {
75 class ForceWithVirial;
76 class LocalAtomSetManager;
77 }
78
79 /*! \brief Returns if the pull coordinate is an angle
80  *
81  * \param[in] pcrd The pull coordinate to query the type for.
82  * \returns a boolean telling if the coordinate is of angle type.
83  */
84 bool pull_coordinate_is_angletype(const t_pull_coord *pcrd);
85
86 /*! \brief Returns the units of the pull coordinate.
87  *
88  * \param[in] pcrd The pull coordinate to query the units for.
89  * \returns a string with the units of the coordinate.
90  */
91 const char *pull_coordinate_units(const t_pull_coord *pcrd);
92
93 /*! \brief Returns the conversion factor from the pull coord init/rate unit to internal value unit.
94  *
95  * \param[in] pcrd The pull coordinate to get the conversion factor for.
96  * \returns the conversion factor.
97  */
98 double pull_conversion_factor_userinput2internal(const t_pull_coord *pcrd);
99
100 /*! \brief Returns the conversion factor from the pull coord internal value unit to the init/rate unit.
101  *
102  * \param[in] pcrd The pull coordinate to get the conversion factor for.
103  * \returns the conversion factor.
104  */
105 double pull_conversion_factor_internal2userinput(const t_pull_coord *pcrd);
106
107 /*! \brief Get the value for pull coord coord_ind.
108  *
109  * \param[in,out] pull      The pull struct.
110  * \param[in]     coord_ind Number of the pull coordinate.
111  * \param[in]     pbc       Information structure about periodicity.
112  * \returns the value of the pull coordinate.
113  */
114 double get_pull_coord_value(struct pull_t      *pull,
115                             int                 coord_ind,
116                             const struct t_pbc *pbc);
117
118 /*! \brief Registers the provider of an external potential for a coordinate.
119  *
120  * This function is only used for checking the consistency of the pull setup.
121  * For each pull coordinate of type external-potential, selected by the user
122  * in the mdp file, there has to be a module that provides this potential.
123  * The module registers itself as the provider by calling this function.
124  * The passed \p provider string has to match the string that the user
125  * passed with the potential-provider pull coordinate mdp option.
126  * This function should be called after init_pull has been called and before
127  * pull_potential is called for the first time.
128  * This function does many consistency checks and when it returns and the
129  * first call to do_potential passes, the pull setup is guaranteed to be
130  * correct (unless the module doesn't call apply_external_pull_coord_force
131  * every step or calls it with incorrect forces). This registering function
132  * will exit with a (release) assertion failure when used incorrely or
133  * with a fatal error when the user (mdp) input in inconsistent.
134  *
135  * Thread-safe for simultaneous registration from multiple threads.
136  *
137  * \param[in,out] pull         The pull struct.
138  * \param[in]     coord_index  The pull coordinate index to register the external potential for.
139  * \param[in]     provider     Provider string, should match the potential-provider pull coordinate mdp option.
140  */
141 void register_external_pull_potential(struct pull_t *pull,
142                                       int            coord_index,
143                                       const char    *provider);
144
145
146 /*! \brief Apply forces of an external potential to a pull coordinate.
147  *
148  * This function applies the external scalar force \p coord_force to
149  * the pull coordinate, distributing it over the atoms in the groups
150  * involved in the pull coordinate. The corresponding potential energy
151  * value should be added to the pull or the module's potential energy term
152  * separately by the module itself.
153  * This function should be called after pull_potential has been called and,
154  * obviously, before the coordinates are updated uses the forces.
155  *
156  * \param[in,out] pull             The pull struct.
157  * \param[in]     coord_index      The pull coordinate index to set the force for.
158  * \param[in]     coord_force      The scalar force for the pull coordinate.
159  * \param[in]     mdatoms          Atom properties, only masses are used.
160  * \param[in,out] forceWithVirial  Force and virial buffers.
161  */
162 void apply_external_pull_coord_force(struct pull_t        *pull,
163                                      int                   coord_index,
164                                      double                coord_force,
165                                      const t_mdatoms      *mdatoms,
166                                      gmx::ForceWithVirial *forceWithVirial);
167
168
169 /*! \brief Set the all the pull forces to zero.
170  *
171  * \param pull              The pull group.
172  */
173 void clear_pull_forces(struct pull_t *pull);
174
175
176 /*! \brief Determine the COM pull forces and add them to f, return the potential
177  *
178  * \param[in,out] pull   The pull struct.
179  * \param[in]     md     All atoms.
180  * \param[in]     pbc    Information struct about periodicity.
181  * \param[in]     cr     Struct for communication info.
182  * \param[in]     t      Time.
183  * \param[in]     lambda The value of lambda in FEP calculations.
184  * \param[in]     x      Positions.
185  * \param[in,out] force  Forces and virial.
186  * \param[out] dvdlambda Pull contribution to dV/d(lambda).
187  *
188  * \returns The pull potential energy.
189  */
190 real pull_potential(struct pull_t *pull, const t_mdatoms *md, struct t_pbc *pbc,
191                     const t_commrec *cr, double t, real lambda,
192                     const rvec *x, gmx::ForceWithVirial *force, real *dvdlambda);
193
194
195 /*! \brief Constrain the coordinates xp in the directions in x
196  * and also constrain v when v != NULL.
197  *
198  * \param[in,out] pull   The pull data.
199  * \param[in]     md     All atoms.
200  * \param[in]     pbc    Information struct about periodicity.
201  * \param[in]     cr     Struct for communication info.
202  * \param[in]     dt     The time step length.
203  * \param[in]     t      The time.
204  * \param[in]     x      Positions.
205  * \param[in,out] xp     Updated x, can be NULL.
206  * \param[in,out] v      Velocities, which may get a pull correction.
207  * \param[in,out] vir    The virial, which, if != NULL, gets a pull correction.
208  */
209 void pull_constraint(struct pull_t *pull, const t_mdatoms *md, struct t_pbc *pbc,
210                      const t_commrec *cr, double dt, double t,
211                      rvec *x, rvec *xp, rvec *v, tensor vir);
212
213
214 /*! \brief Make a selection of the home atoms for all pull groups.
215  * Should be called at every domain decomposition.
216  *
217  * \param cr             Structure for communication info.
218  * \param pull           The pull group.
219  */
220 void dd_make_local_pull_groups(const t_commrec *cr, struct pull_t *pull);
221
222
223 /*! \brief Allocate, initialize and return a pull work struct.
224  *
225  * \param fplog       General output file, normally md.log.
226  * \param pull_params The pull input parameters containing all pull settings.
227  * \param ir          The inputrec.
228  * \param mtop        The topology of the whole system.
229  * \param cr          Struct for communication info.
230  * \param atomSets    The manager that handles the pull atom sets
231  * \param lambda      FEP lambda.
232  */
233 struct pull_t *init_pull(FILE                      *fplog,
234                          const pull_params_t       *pull_params,
235                          const t_inputrec          *ir,
236                          const gmx_mtop_t          *mtop,
237                          const t_commrec           *cr,
238                          gmx::LocalAtomSetManager  *atomSets,
239                          real                       lambda);
240
241
242 /*! \brief Close the pull output files and delete pull.
243  *
244  * \param pull       The pull data structure.
245  */
246 void finish_pull(struct pull_t *pull);
247
248
249 /*! \brief Calculates centers of mass all pull groups.
250  *
251  * \param[in] cr       Struct for communication info.
252  * \param[in] pull     The pull data structure.
253  * \param[in] md       All atoms.
254  * \param[in] pbc      Information struct about periodicity.
255  * \param[in] t        Time, only used for cylinder ref.
256  * \param[in] x        The local positions.
257  * \param[in,out] xp   Updated x, can be NULL.
258  *
259  */
260 void pull_calc_coms(const t_commrec *cr,
261                     pull_t          *pull,
262                     const t_mdatoms *md,
263                     t_pbc           *pbc,
264                     double           t,
265                     const rvec       x[],
266                     rvec            *xp);
267
268 /*! \brief Margin for checking pull group PBC distances compared to half the box size */
269 static constexpr real c_pullGroupPbcMargin = 0.9;
270 /*! \brief Threshold (as a factor of half the box size) for accepting pull groups without explicitly set refatom */
271 static constexpr real c_pullGroupSmallGroupThreshold = 0.5;
272
273 /*! \brief Checks whether all groups that use a reference atom are within PBC restrictions
274  *
275  * Groups that use a reference atom for determining PBC should have all their
276  * atoms within half the box size from the PBC atom. The box size is used
277  * per dimension for rectangular boxes, but can be a combination of
278  * dimensions for triclinic boxes, depending on which dimensions are
279  * involved in the pull coordinates a group is involved in. A margin is specified
280  * to ensure that atoms are not too close to the maximum distance.
281  *
282  * Should be called without MPI parallelization and after pull_calc_coms()
283  * has been called at least once.
284  *
285  * \param[in] pull       The pull data structure
286  * \param[in] x          The coordinates
287  * \param[in] pbc        Information struct about periodicity
288  * \param[in] pbcMargin  The minimum margin (as a fraction) to half the box size
289  * \returns -1 when all groups obey PBC or the first group index that fails PBC
290  */
291 int pullCheckPbcWithinGroups(const pull_t &pull,
292                              const rvec   *x,
293                              const t_pbc  &pbc,
294                              real          pbcMargin);
295
296 /*! \brief Checks whether a specific group that uses a reference atom is within PBC restrictions
297  *
298  * Groups that use a reference atom for determining PBC should have all their
299  * atoms within half the box size from the PBC atom. The box size is used
300  * per dimension for rectangular boxes, but can be a combination of
301  * dimensions for triclinic boxes, depending on which dimensions are
302  * involved in the pull coordinates a group is involved in. A margin is specified
303  * to ensure that atoms are not too close to the maximum distance. Only one group is
304  * checked.
305  *
306  * Should be called without MPI parallelization and after pull_calc_coms()
307  * has been called at least once.
308  *
309  * \param[in] pull       The pull data structure
310  * \param[in] x          The coordinates
311  * \param[in] pbc        Information struct about periodicity
312  * \param[in] groupNr    The index of the group (in pull.group[]) to check
313  * \param[in] pbcMargin  The minimum margin (as a fraction) to half the box size
314  * \returns true if the group obeys PBC otherwise false
315  */
316 bool pullCheckPbcWithinGroup(const pull_t                  &pull,
317                              gmx::ArrayRef<const gmx::RVec> x,
318                              const t_pbc                   &pbc,
319                              int                            groupNr,
320                              real                           pbcMargin);
321
322 /*! \brief Returns if we have pull coordinates with potential pulling.
323  *
324  * \param[in] pull     The pull data structure.
325  */
326 gmx_bool pull_have_potential(const struct pull_t *pull);
327
328
329 /*! \brief Returns if we have pull coordinates with constraint pulling.
330  *
331  * \param[in] pull     The pull data structure.
332  */
333 gmx_bool pull_have_constraint(const struct pull_t *pull);
334
335 /*! \brief Returns the maxing distance for pulling
336  *
337  * For distance geometries, only dimensions with pcrd->params[dim]=1
338  * are included in the distance calculation.
339  * For directional geometries, only dimensions with pcrd->vec[dim]!=0
340  * are included in the distance calculation.
341  *
342  * \param[in] pcrd Pulling data structure
343  * \param[in] pbc  Information on periodic boundary conditions
344  * \returns The maximume distance
345  */
346 real max_pull_distance2(const pull_coord_work_t *pcrd,
347                         const t_pbc             *pbc);
348
349 /*! \brief Sets the previous step COM in pull to the current COM and updates the pull_com_prev_step in the state
350  *
351  * \param[in]   pull  The COM pull force calculation data structure
352  * \param[in]   state The local (to this rank) state.
353  */
354 void updatePrevStepPullCom(struct pull_t *pull, t_state *state);
355
356 /*! \brief Allocates, initializes and communicates the previous step pull COM (if that option is set to true).
357  *
358  * If ir->pull->bSetPbcRefToPrevStepCOM is not true nothing is done.
359  *
360  * \param[in] ir                     The input options/settings of the simulation.
361  * \param[in] md                     All atoms.
362  * \param[in] state                  The local (to this rank) state.
363  * \param[in] state_global           The global state.
364  * \param[in] cr                     Struct for communication info.
365  * \param[in] startingFromCheckpoint Is the simulation starting from a checkpoint?
366  */
367 void preparePrevStepPullCom(const t_inputrec *ir, const t_mdatoms *md, t_state *state, const t_state *state_global, const t_commrec *cr, bool startingFromCheckpoint);
368
369 /*! \brief Initializes the COM of the previous step (set to initial COM)
370  *
371  * \param[in] cr       Struct for communication info.
372  * \param[in] pull     The pull data structure.
373  * \param[in] md       All atoms.
374  * \param[in] pbc      Information struct about periodicity.
375  * \param[in] x        The local positions.
376  */
377 void initPullComFromPrevStep(const t_commrec *cr,
378                              pull_t          *pull,
379                              const t_mdatoms *md,
380                              t_pbc           *pbc,
381                              const rvec       x[]);
382
383 #endif