2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \libinternal \file
37 * \brief Declaration of high-level functions of CUDA implementation of update and constrain class.
39 * \todo This should only list interfaces needed for libgromacs clients (e.g.
40 * management of coordinates, velocities and forces should not be here).
41 * \todo Change "cuda" suffix to "gpu"
43 * \author Artem Zhmurov <zhmurov@gmail.com>
45 * \ingroup module_mdlib
48 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_H
49 #define GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_H
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/topology/idef.h"
55 #include "gromacs/utility/classhelpers.h"
60 class UpdateConstrainCuda
64 /*! \brief Create Update-Constrain object.
66 * \param[in] numAtoms Number of atoms.
68 * \param[in] ir Input record data: LINCS takes number of iterations and order of
70 * \param[in] mtop Topology of the system: SETTLE gets the masses for O and H atoms
71 * and target O-H and H-H distances from this object.
73 UpdateConstrainCuda(int numAtoms
,
75 const gmx_mtop_t
&mtop
);
77 ~UpdateConstrainCuda();
81 * Integrates the equation of motion using Leap-Frog algorithm and applies
82 * LINCS and SETTLE constraints.
83 * Updates d_xp_ and d_v_ fields of this object.
85 * \param[in] dt Timestep
86 * \param[in] updateVelocities If the velocities should be constrained.
87 * \param[in] computeVirial If virial should be updated.
88 * \param[out] virial Place to save virial tensor.
90 void integrate(real dt
,
91 bool updateVelocities
,
96 * Update data-structures (e.g. after NB search step).
98 * \param[in] idef System topology
99 * \param[in] md Atoms data. Can be used to update masses if needed (not used now).
101 void set(const t_idef
&idef
,
102 const t_mdatoms
&md
);
107 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
109 * \param[in] pbc The PBC data in t_pbc format.
111 void setPbc(const t_pbc
*pbc
);
114 * Copy coordinates from CPU to GPU.
116 * The data are assumed to be in float3/fvec format (single precision).
118 * \param[in] h_x CPU pointer where coordinates should be copied from.
120 void copyCoordinatesToGpu(const rvec
*h_x
);
123 * Copy velocities from CPU to GPU.
125 * The data are assumed to be in float3/fvec format (single precision).
127 * \param[in] h_v CPU pointer where velocities should be copied from.
129 void copyVelocitiesToGpu(const rvec
*h_v
);
132 * Copy forces from CPU to GPU.
134 * The data are assumed to be in float3/fvec format (single precision).
136 * \param[in] h_f CPU pointer where forces should be copied from.
138 void copyForcesToGpu(const rvec
*h_f
);
141 * Copy coordinates from GPU to CPU.
143 * The data are assumed to be in float3/fvec format (single precision).
145 * \param[out] h_xp CPU pointer where coordinates should be copied to.
147 void copyCoordinatesFromGpu(rvec
*h_xp
);
150 * Copy velocities from GPU to CPU.
152 * The velocities are assumed to be in float3/fvec format (single precision).
154 * \param[in] h_v Pointer to velocities data.
156 void copyVelocitiesFromGpu(rvec
*h_v
);
159 * Copy forces from GPU to CPU.
161 * The forces are assumed to be in float3/fvec format (single precision).
163 * \param[in] h_f Pointer to forces data.
165 void copyForcesFromGpu(rvec
*h_f
);
168 * Set the internal GPU-memory x, xprime and v pointers.
170 * Data is not copied. The data are assumed to be in float3/fvec format
171 * (float3 is used internally, but the data layout should be identical).
173 * \param[in] d_x Pointer to the coordinates for the input (on GPU)
174 * \param[in] d_xp Pointer to the coordinates for the output (on GPU)
175 * \param[in] d_v Pointer to the velocities (on GPU)
176 * \param[in] d_f Pointer to the forces (on GPU)
178 void setXVFPointers(rvec
*d_x
, rvec
*d_xp
, rvec
*d_v
, rvec
*d_f
);
182 gmx::PrivateImplPointer
<Impl
> impl_
;