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