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 Change "cuda" suffix to "gpu"
41 * \author Artem Zhmurov <zhmurov@gmail.com>
43 * \ingroup module_mdlib
46 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_H
47 #define GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_H
49 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
50 #include "gromacs/mdtypes/group.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/classhelpers.h"
54 class GpuEventSynchronizer;
65 class UpdateConstrainCuda
69 /*! \brief Create Update-Constrain object.
71 * The constructor is given a non-nullptr \p commandStream, in which all the update and constrain
72 * routines are executed. \p xUpdatedOnDevice should mark the completion of all kernels that modify
73 * coordinates. The event is maintained outside this class and also passed to all (if any) consumers
74 * of the updated coordinates. The \p xUpdatedOnDevice also can not be a nullptr because the
75 * markEvent(...) method is called unconditionally.
77 * \param[in] ir Input record data: LINCS takes number of iterations and order of
79 * \param[in] mtop Topology of the system: SETTLE gets the masses for O and H atoms
80 * and target O-H and H-H distances from this object.
81 * \param[in] commandStream GPU stream to use. Can be nullptr.
82 * \param[in] xUpdatedOnDevice The event synchronizer to use to mark that update is done on the GPU.
84 UpdateConstrainCuda(const t_inputrec &ir,
85 const gmx_mtop_t &mtop,
86 const void *commandStream,
87 GpuEventSynchronizer *xUpdatedOnDevice);
89 ~UpdateConstrainCuda();
93 * This will extract temperature scaling factors from tcstat, transform them into the plain
94 * array and call the normal integrate method.
96 * \param[in] fReadyOnDevice Event synchronizer indicating that the forces are ready in the device memory.
97 * \param[in] dt Timestep.
98 * \param[in] updateVelocities If the velocities should be constrained.
99 * \param[in] computeVirial If virial should be updated.
100 * \param[out] virial Place to save virial tensor.
101 * \param[in] doTempCouple If the temperature coupling should be performed.
102 * \param[in] tcstat Temperature coupling data.
103 * \param[in] doPressureCouple If the temperature coupling should be applied.
104 * \param[in] dtPressureCouple Period between pressure coupling steps
105 * \param[in] velocityScalingMatrix Parrinello-Rahman velocity scaling matrix
107 void integrate(GpuEventSynchronizer *fReadyOnDevice,
109 bool updateVelocities,
113 gmx::ArrayRef<const t_grp_tcstat> tcstat,
114 bool doPressureCouple,
115 float dtPressureCouple,
116 const matrix velocityScalingMatrix);
118 /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
120 * \param[in,out] d_x Device buffer with coordinates.
121 * \param[in,out] d_v Device buffer with velocities.
122 * \param[in] d_f Device buffer with forces.
123 * \param[in] idef System topology
124 * \param[in] md Atoms data.
125 * \param[in] numTempScaleValues Number of temperature scaling groups. Zero for no temperature scaling.
127 void set(DeviceBuffer<float> d_x,
128 DeviceBuffer<float> d_v,
129 DeviceBuffer<float> d_f,
132 int numTempScaleValues);
137 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
139 * \param[in] pbc The PBC data in t_pbc format.
141 void setPbc(const t_pbc *pbc);
143 /*! \brief Return the synchronizer associated with the event indicated that the coordinates are ready on the device.
145 GpuEventSynchronizer* getCoordinatesReadySync();
149 gmx::PrivateImplPointer<Impl> impl_;