2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
39 /*! \libinternal \file
43 * This file contains datatypes and function declarations necessary
44 for mdrun to interface with the pull code.
51 #ifndef GMX_PULLING_PULL_H
52 #define GMX_PULLING_PULL_H
57 #include "gromacs/math/vectypes.h"
58 #include "gromacs/mdtypes/pull_params.h"
59 #include "gromacs/utility/arrayref.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/real.h"
64 struct gmx_output_env_t;
65 struct pull_coord_work_t;
79 class ForceWithVirial;
80 class LocalAtomSetManager;
83 /*! \brief Returns the units of the pull coordinate.
85 * \param[in] pcrd The pull coordinate to query the units for.
86 * \returns a string with the units of the coordinate.
88 const char* pull_coordinate_units(const t_pull_coord& pcrd);
90 /*! \brief Returns the conversion factor from the pull coord init/rate unit to internal value unit.
92 * \param[in] pcrd The pull coordinate to get the conversion factor for.
93 * \returns the conversion factor.
95 double pull_conversion_factor_userinput2internal(const t_pull_coord& pcrd);
97 /*! \brief Returns the conversion factor from the pull coord internal value unit to the init/rate unit.
99 * \param[in] pcrd The pull coordinate to get the conversion factor for.
100 * \returns the conversion factor.
102 double pull_conversion_factor_internal2userinput(const t_pull_coord& pcrd);
104 /*! \brief Get the value for pull coord coord_ind.
106 * \param[in,out] pull The pull struct.
107 * \param[in] coordIndex Index of the pull coordinate in the list of coordinates
108 * \param[in] pbc Information structure about periodicity.
109 * \returns the value of the pull coordinate.
111 double get_pull_coord_value(pull_t* pull, int coordIndex, const t_pbc& pbc);
113 /*! \brief Registers the provider of an external potential for a coordinate.
115 * This function is only used for checking the consistency of the pull setup.
116 * For each pull coordinate of type external-potential, selected by the user
117 * in the mdp file, there has to be a module that provides this potential.
118 * The module registers itself as the provider by calling this function.
119 * The passed \p provider string has to match the string that the user
120 * passed with the potential-provider pull coordinate mdp option.
121 * This function should be called after init_pull has been called and before
122 * pull_potential is called for the first time.
123 * This function does many consistency checks and when it returns and the
124 * first call to do_potential passes, the pull setup is guaranteed to be
125 * correct (unless the module doesn't call apply_external_pull_coord_force
126 * every step or calls it with incorrect forces). This registering function
127 * will exit with a (release) assertion failure when used incorrely or
128 * with a fatal error when the user (mdp) input in inconsistent.
130 * Thread-safe for simultaneous registration from multiple threads.
132 * \param[in,out] pull The pull struct.
133 * \param[in] coord_index The pull coordinate index to register the external potential for.
134 * \param[in] provider Provider string, should match the potential-provider pull coordinate mdp option.
136 void register_external_pull_potential(struct pull_t* pull, int coord_index, const char* provider);
139 /*! \brief Apply forces of an external potential to a pull coordinate.
141 * This function applies the external scalar force \p coord_force to
142 * the pull coordinate, distributing it over the atoms in the groups
143 * involved in the pull coordinate. The corresponding potential energy
144 * value should be added to the pull or the module's potential energy term
145 * separately by the module itself.
146 * This function should be called after pull_potential has been called and,
147 * obviously, before the coordinates are updated uses the forces.
149 * \param[in,out] pull The pull struct.
150 * \param[in] coord_index The pull coordinate index to set the force for.
151 * \param[in] coord_force The scalar force for the pull coordinate.
152 * \param[in] masses Atoms masses.
153 * \param[in,out] forceWithVirial Force and virial buffers.
155 void apply_external_pull_coord_force(struct pull_t* pull,
158 gmx::ArrayRef<const real> masses,
159 gmx::ForceWithVirial* forceWithVirial);
161 /*! \brief Set the all the pull forces to zero.
163 * \param pull The pull group.
165 void clear_pull_forces(pull_t* pull);
168 /*! \brief Determine the COM pull forces and add them to f, return the potential
170 * \param[in,out] pull The pull struct.
171 * \param[in] masses Atoms masses.
172 * \param[in] pbc Information struct about periodicity.
173 * \param[in] cr Struct for communication info.
175 * \param[in] lambda The value of lambda in FEP calculations.
176 * \param[in] x Positions.
177 * \param[in,out] force Forces and virial.
178 * \param[out] dvdlambda Pull contribution to dV/d(lambda).
180 * \returns The pull potential energy.
182 real pull_potential(pull_t* pull,
183 gmx::ArrayRef<const real> masses,
188 gmx::ArrayRef<const gmx::RVec> x,
189 gmx::ForceWithVirial* force,
193 /*! \brief Constrain the coordinates xp in the directions in x
194 * and also constrain v when v != NULL.
196 * \param[in,out] pull The pull data.
197 * \param[in] masses Atoms masses.
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.
207 void pull_constraint(struct pull_t* pull,
208 gmx::ArrayRef<const real> masses,
213 gmx::ArrayRef<gmx::RVec> x,
214 gmx::ArrayRef<gmx::RVec> xp,
215 gmx::ArrayRef<gmx::RVec> v,
219 /*! \brief Make a selection of the home atoms for all pull groups.
220 * Should be called at every domain decomposition.
222 * \param cr Structure for communication info.
223 * \param pull The pull group.
225 void dd_make_local_pull_groups(const t_commrec* cr, pull_t* pull);
228 /*! \brief Allocate, initialize and return a pull work struct.
230 * \param fplog General output file, normally md.log.
231 * \param pull_params The pull input parameters containing all pull settings.
232 * \param ir The inputrec.
233 * \param mtop The topology of the whole system.
234 * \param cr Struct for communication info.
235 * \param atomSets The manager that handles the pull atom sets
236 * \param lambda FEP lambda.
238 struct pull_t* init_pull(FILE* fplog,
239 const pull_params_t* pull_params,
240 const t_inputrec* ir,
241 const gmx_mtop_t& mtop,
243 gmx::LocalAtomSetManager* atomSets,
247 /*! \brief Close the pull output files and delete pull.
249 * \param pull The pull data structure.
251 void finish_pull(struct pull_t* pull);
254 /*! \brief Calculates centers of mass all pull groups.
256 * \param[in] cr Struct for communication info.
257 * \param[in] pull The pull data structure.
258 * \param[in] masses Atoms masses.
259 * \param[in] pbc Information struct about periodicity.
260 * \param[in] t Time, only used for cylinder ref.
261 * \param[in] x The local positions.
262 * \param[in,out] xp Updated x, can be NULL.
265 void pull_calc_coms(const t_commrec* cr,
267 gmx::ArrayRef<const real> masses,
270 gmx::ArrayRef<const gmx::RVec> x,
271 gmx::ArrayRef<gmx::RVec> xp);
273 /*! \brief Margin for checking pull group PBC distances compared to half the box size */
274 static constexpr real c_pullGroupPbcMargin = 0.9;
275 /*! \brief Threshold (as a factor of half the box size) for accepting pull groups without explicitly set refatom */
276 static constexpr real c_pullGroupSmallGroupThreshold = 0.5;
278 /*! \brief Checks whether all groups that use a reference atom are within PBC restrictions
280 * Groups that use a reference atom for determining PBC should have all their
281 * atoms within half the box size from the PBC atom. The box size is used
282 * per dimension for rectangular boxes, but can be a combination of
283 * dimensions for triclinic boxes, depending on which dimensions are
284 * involved in the pull coordinates a group is involved in. A margin is specified
285 * to ensure that atoms are not too close to the maximum distance.
287 * Should be called without MPI parallelization and after pull_calc_coms()
288 * has been called at least once.
290 * \param[in] pull The pull data structure
291 * \param[in] x The coordinates
292 * \param[in] pbc Information struct about periodicity
293 * \param[in] pbcMargin The minimum margin (as a fraction) to half the box size
294 * \returns -1 when all groups obey PBC or the first group index that fails PBC
296 int pullCheckPbcWithinGroups(const pull_t& pull, gmx::ArrayRef<const gmx::RVec> x, const t_pbc& pbc, real pbcMargin);
298 /*! \brief Checks whether a specific group that uses a reference atom is within PBC restrictions
300 * Groups that use a reference atom for determining PBC should have all their
301 * atoms within half the box size from the PBC atom. The box size is used
302 * per dimension for rectangular boxes, but can be a combination of
303 * dimensions for triclinic boxes, depending on which dimensions are
304 * involved in the pull coordinates a group is involved in. A margin is specified
305 * to ensure that atoms are not too close to the maximum distance. Only one group is
308 * Should be called without MPI parallelization and after pull_calc_coms()
309 * has been called at least once.
311 * \param[in] pull The pull data structure
312 * \param[in] x The coordinates
313 * \param[in] pbc Information struct about periodicity
314 * \param[in] groupNr The index of the group (in pull.group[]) to check
315 * \param[in] pbcMargin The minimum margin (as a fraction) to half the box size
316 * \returns true if the group obeys PBC otherwise false
318 bool pullCheckPbcWithinGroup(const pull_t& pull,
319 gmx::ArrayRef<const gmx::RVec> x,
324 /*! \brief Returns if we have pull coordinates with potential pulling.
326 * \param[in] pull The pull data structure.
328 bool pull_have_potential(const pull_t& pull);
331 /*! \brief Returns if we have pull coordinates with constraint pulling.
333 * \param[in] pull The pull data structure.
335 bool pull_have_constraint(const pull_t& pull);
337 /*! \brief Returns if inputrec has pull coordinates with constraint pulling.
339 * \param[in] pullParameters Pulling input parameters from input record.
341 bool pull_have_constraint(const pull_params_t& pullParameters);
343 /*! \brief Returns the maxing distance for pulling
345 * For distance geometries, only dimensions with pcrd->params[dim]=1
346 * are included in the distance calculation.
347 * For directional geometries, only dimensions with pcrd->vec[dim]!=0
348 * are included in the distance calculation.
350 * \param[in] pcrd Pulling data structure
351 * \param[in] pbc Information on periodic boundary conditions
352 * \returns The maximume distance
354 real max_pull_distance2(const pull_coord_work_t& pcrd, const t_pbc& pbc);
356 /*! \brief Sets the previous step COM in pull to the current COM, and optionally
357 * updates it in the provided ArrayRef
359 * \param[in] pull The COM pull force calculation data structure
360 * \param[in] comPreviousStep The COM of the previous step of each pull group
362 void updatePrevStepPullCom(pull_t* pull, std::optional<gmx::ArrayRef<double>> comPreviousStep);
364 /*! \brief Returns a copy of the previous step pull COM as flat vector
366 * Used for modular simulator checkpointing. Allows to keep the
367 * implementation details of pull_t hidden from its users.
369 * \param[in] pull The COM pull force calculation data structure
370 * \return A copy of the previous step COM
372 std::vector<double> prevStepPullCom(const pull_t* pull);
374 /*! \brief Set the previous step pull COM from a flat vector
376 * Used to restore modular simulator checkpoints. Allows to keep the
377 * implementation details of pull_t hidden from its users.
379 * \param[in] pull The COM pull force calculation data structure
380 * \param[in] prevStepPullCom The previous step COM to set
382 void setPrevStepPullCom(pull_t* pull, gmx::ArrayRef<const double> prevStepPullCom);
384 /*! \brief Allocates, initializes and communicates the previous step pull COM (if that option is set to true).
386 * If ir->pull->bSetPbcRefToPrevStepCOM is not true nothing is done.
388 * \param[in] ir The input options/settings of the simulation.
389 * \param[in] pull_work The COM pull force calculation data structure
390 * \param[in] masses Atoms masses.
391 * \param[in] state The local (to this rank) state.
392 * \param[in] state_global The global state.
393 * \param[in] cr Struct for communication info.
394 * \param[in] startingFromCheckpoint Is the simulation starting from a checkpoint?
396 void preparePrevStepPullCom(const t_inputrec* ir,
398 gmx::ArrayRef<const real> masses,
400 const t_state* state_global,
402 bool startingFromCheckpoint);
404 /*! \brief Initializes the COM of the previous step (set to initial COM)
406 * \param[in] cr Struct for communication info.
407 * \param[in] pull The pull data structure.
408 * \param[in] masses Atoms masses.
409 * \param[in] pbc Information struct about periodicity.
410 * \param[in] x The local positions.
412 void initPullComFromPrevStep(const t_commrec* cr,
414 gmx::ArrayRef<const real> masses,
416 gmx::ArrayRef<const gmx::RVec> x);
418 /*! \brief Initializes the previous step pull COM for new simulations (no reading from checkpoint).
420 * \param[in] cr Struct for communication info.
421 * \param[in] pull_work The COM pull force calculation data structure.
422 * \param[in] masses Atoms masses.
423 * \param[in] x The local positions.
424 * \param[in] box The current box matrix.
425 * \param[in] pbcType The type of periodic boundary conditions.
426 * \param[in] comPreviousStep The COM of the previous step of each pull group.
428 void preparePrevStepPullComNewSimulation(const t_commrec* cr,
430 gmx::ArrayRef<const real> masses,
431 gmx::ArrayRef<const gmx::RVec> x,
434 std::optional<gmx::ArrayRef<double>>&& comPreviousStep);