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