Merge branch release-2018
[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, 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/basedefinitions.h"
58 #include "gromacs/utility/real.h"
59
60 struct ContinuationOptions;
61 struct gmx_mtop_t;
62 struct gmx_output_env_t;
63 struct pull_coord_work_t;
64 struct pull_params_t;
65 struct t_commrec;
66 struct t_filenm;
67 struct t_inputrec;
68 struct t_mdatoms;
69 struct t_pbc;
70
71 namespace gmx
72 {
73 class ForceWithVirial;
74 class LocalAtomSetManager;
75 }
76
77 /*! \brief Returns if the pull coordinate is an angle
78  *
79  * \param[in] pcrd The pull coordinate to query the type for.
80  * \returns a boolean telling if the coordinate is of angle type.
81  */
82 bool pull_coordinate_is_angletype(const t_pull_coord *pcrd);
83
84 /*! \brief Returns the units of the pull coordinate.
85  *
86  * \param[in] pcrd The pull coordinate to query the units for.
87  * \returns a string with the units of the coordinate.
88  */
89 const char *pull_coordinate_units(const t_pull_coord *pcrd);
90
91 /*! \brief Returns the conversion factor from the pull coord init/rate unit to internal value unit.
92  *
93  * \param[in] pcrd The pull coordinate to get the conversion factor for.
94  * \returns the conversion factor.
95  */
96 double pull_conversion_factor_userinput2internal(const t_pull_coord *pcrd);
97
98 /*! \brief Returns the conversion factor from the pull coord internal value unit to the init/rate unit.
99  *
100  * \param[in] pcrd The pull coordinate to get the conversion factor for.
101  * \returns the conversion factor.
102  */
103 double pull_conversion_factor_internal2userinput(const t_pull_coord *pcrd);
104
105 /*! \brief Get the value for pull coord coord_ind.
106  *
107  * \param[in,out] pull      The pull struct.
108  * \param[in]     coord_ind Number of the pull coordinate.
109  * \param[in]     pbc       Information structure about periodicity.
110  * \returns the value of the pull coordinate.
111  */
112 double get_pull_coord_value(struct pull_t      *pull,
113                             int                 coord_ind,
114                             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,
140                                       int            coord_index,
141                                       const char    *provider);
142
143
144 /*! \brief Apply forces of an external potential to a pull coordinate.
145  *
146  * This function applies the external scalar force \p coord_force to
147  * the pull coordinate, distributing it over the atoms in the groups
148  * involved in the pull coordinate. The corresponding potential energy
149  * value should be added to the pull or the module's potential energy term
150  * separately by the module itself.
151  * This function should be called after pull_potential has been called and,
152  * obviously, before the coordinates are updated uses the forces.
153  *
154  * \param[in,out] pull             The pull struct.
155  * \param[in]     coord_index      The pull coordinate index to set the force for.
156  * \param[in]     coord_force      The scalar force for the pull coordinate.
157  * \param[in]     mdatoms          Atom properties, only masses are used.
158  * \param[in,out] forceWithVirial  Force and virial buffers.
159  */
160 void apply_external_pull_coord_force(struct pull_t        *pull,
161                                      int                   coord_index,
162                                      double                coord_force,
163                                      const t_mdatoms      *mdatoms,
164                                      gmx::ForceWithVirial *forceWithVirial);
165
166
167 /*! \brief Set the all the pull forces to zero.
168  *
169  * \param pull              The pull group.
170  */
171 void clear_pull_forces(struct pull_t *pull);
172
173
174 /*! \brief Determine the COM pull forces and add them to f, return the potential
175  *
176  * \param[in,out] pull   The pull struct.
177  * \param[in]     md     All atoms.
178  * \param[in]     pbc    Information struct about periodicity.
179  * \param[in]     cr     Struct for communication info.
180  * \param[in]     t      Time.
181  * \param[in]     lambda The value of lambda in FEP calculations.
182  * \param[in]     x      Positions.
183  * \param[in,out] force  Forces and virial.
184  * \param[out] dvdlambda Pull contribution to dV/d(lambda).
185  *
186  * \returns The pull potential energy.
187  */
188 real pull_potential(struct pull_t *pull, const t_mdatoms *md, struct t_pbc *pbc,
189                     const t_commrec *cr, double t, real lambda,
190                     const rvec *x, gmx::ForceWithVirial *force, real *dvdlambda);
191
192
193 /*! \brief Constrain the coordinates xp in the directions in x
194  * and also constrain v when v != NULL.
195  *
196  * \param[in,out] pull   The pull data.
197  * \param[in]     md     All atoms.
198  * \param[in]     pbc    Information struct about periodicity.
199  * \param[in]     cr     Struct for communication info.
200  * \param[in]     dt     The time step length.
201  * \param[in]     t      The time.
202  * \param[in]     x      Positions.
203  * \param[in,out] xp     Updated x, can be NULL.
204  * \param[in,out] v      Velocities, which may get a pull correction.
205  * \param[in,out] vir    The virial, which, if != NULL, gets a pull correction.
206  */
207 void pull_constraint(struct pull_t *pull, const t_mdatoms *md, struct t_pbc *pbc,
208                      const t_commrec *cr, double dt, double t,
209                      rvec *x, rvec *xp, rvec *v, tensor vir);
210
211
212 /*! \brief Make a selection of the home atoms for all pull groups.
213  * Should be called at every domain decomposition.
214  *
215  * \param cr             Structure for communication info.
216  * \param pull           The pull group.
217  */
218 void dd_make_local_pull_groups(const t_commrec *cr, struct pull_t *pull);
219
220
221 /*! \brief Allocate, initialize and return a pull work struct.
222  *
223  * \param fplog       General output file, normally md.log.
224  * \param pull_params The pull input parameters containing all pull settings.
225  * \param ir          The inputrec.
226  * \param mtop        The topology of the whole system.
227  * \param cr          Struct for communication info.
228  * \param atomSets    The manager that handles the pull atom sets
229  * \param lambda      FEP lambda.
230  */
231 struct pull_t *init_pull(FILE                      *fplog,
232                          const pull_params_t       *pull_params,
233                          const t_inputrec          *ir,
234                          const gmx_mtop_t          *mtop,
235                          const t_commrec           *cr,
236                          gmx::LocalAtomSetManager  *atomSets,
237                          real                       lambda);
238
239 /*! \brief Set up and open the pull output files, when requested.
240  *
241  * NOTE: This should only be called on the master rank and only when
242  *       doing dynamics (e.g. not with energy minimization).
243  *
244  * \param pull        The pull work data struct
245  * \param nfile       Number of files.
246  * \param fnm         Standard filename struct.
247  * \param oenv        Output options.
248  * \param continuationOptions  Options for continuing from checkpoint file
249  */
250 void init_pull_output_files(pull_t                    *pull,
251                             int                        nfile,
252                             const t_filenm             fnm[],
253                             const gmx_output_env_t    *oenv,
254                             const ContinuationOptions &continuationOptions);
255
256 /*! \brief Close the pull output files.
257  *
258  * \param pull       The pull group.
259  */
260 void finish_pull(struct pull_t *pull);
261
262
263 /*! \brief Print the pull output (x and/or f)
264  *
265  * \param pull     The pull data structure.
266  * \param step     Time step number.
267  * \param time     Time.
268  */
269 void pull_print_output(struct pull_t *pull, int64_t step, double time);
270
271
272 /*! \brief Calculates centers of mass all pull groups.
273  *
274  * \param[in] cr       Struct for communication info.
275  * \param[in] pull     The pull data structure.
276  * \param[in] md       All atoms.
277  * \param[in] pbc      Information struct about periodicity.
278  * \param[in] t        Time, only used for cylinder ref.
279  * \param[in] x        The local positions.
280  * \param[in,out] xp   Updated x, can be NULL.
281  *
282  */
283 void pull_calc_coms(const t_commrec *cr,
284                     pull_t          *pull,
285                     const t_mdatoms *md,
286                     t_pbc           *pbc,
287                     double           t,
288                     const rvec       x[],
289                     rvec            *xp);
290
291 /*! \brief Margin for checking PBC distances compared to half the box size in pullCheckPbcWithinGroups() */
292 static constexpr real c_pullGroupPbcMargin = 0.9;
293
294 /*! \brief Checks whether all groups that use a reference atom are within PBC restrictions
295  *
296  * Groups that use a reference atom for determining PBC should have all their
297  * atoms within half the box size from the PBC atom. The box size is used
298  * per dimension for rectangular boxes, but can be a combination of
299  * dimensions for triclinic boxes, depending on which dimensions are
300  * involved in the pull coordinates a group is involved in.
301  *
302  * Should be called without MPI parallelization and after pull_calc_coms()
303  * has been called at least once.
304  *
305  * \param[in] pull  The pull data structure
306  * \param[in] x     The coordinates
307  * \param[in] pbc   Information struct about periodicity
308  * \returns -1 when all groups obey PBC or the first group index that fails PBC
309  */
310 int pullCheckPbcWithinGroups(const pull_t &pull,
311                              const rvec   *x,
312                              const t_pbc  &pbc);
313
314 /*! \brief Returns if we have pull coordinates with potential pulling.
315  *
316  * \param[in] pull     The pull data structure.
317  */
318 gmx_bool pull_have_potential(const struct pull_t *pull);
319
320
321 /*! \brief Returns if we have pull coordinates with constraint pulling.
322  *
323  * \param[in] pull     The pull data structure.
324  */
325 gmx_bool pull_have_constraint(const struct pull_t *pull);
326
327 /*! \brief Returns the maxing distance for pulling
328  *
329  * For distance geometries, only dimensions with pcrd->params[dim]=1
330  * are included in the distance calculation.
331  * For directional geometries, only dimensions with pcrd->vec[dim]!=0
332  * are included in the distance calculation.
333  *
334  * \param[in] pcrd Pulling data structure
335  * \param[in] pbc  Information on periodic boundary conditions
336  * \returns The maximume distance
337  */
338 real max_pull_distance2(const pull_coord_work_t *pcrd,
339                         const t_pbc             *pbc);
340
341 #endif