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