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