2052bccb7c801b96b4b4916349cbcb0338cb76bd
[alexxy/gromacs.git] / src / gromacs / mdlib / update_constrain_cuda_impl.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Implements update and constraints class using CUDA.
38  *
39  * The class combines Leap-Frog integrator with LINCS and SETTLE constraints.
40  *
41  * \todo The computational procedures in members should be integrated to improve
42  *       computational performance.
43  *
44  * \author Artem Zhmurov <zhmurov@gmail.com>
45  *
46  * \ingroup module_mdlib
47  */
48 #include "gmxpre.h"
49
50 #include "update_constrain_cuda_impl.h"
51
52 #include <assert.h>
53 #include <stdio.h>
54
55 #include <cmath>
56
57 #include <algorithm>
58
59 #include "gromacs/gpu_utils/cudautils.cuh"
60 #include "gromacs/gpu_utils/devicebuffer.h"
61 #include "gromacs/gpu_utils/gputraits.cuh"
62 #include "gromacs/gpu_utils/vectype_ops.cuh"
63 #include "gromacs/mdlib/leapfrog_cuda.cuh"
64 #include "gromacs/mdlib/lincs_cuda.cuh"
65 #include "gromacs/mdlib/settle_cuda.cuh"
66 #include "gromacs/mdlib/update_constrain_cuda.h"
67
68 namespace gmx
69 {
70
71 void UpdateConstrainCuda::Impl::integrate(GpuEventSynchronizer             *fReadyOnDevice,
72                                           const real                        dt,
73                                           const bool                        updateVelocities,
74                                           const bool                        computeVirial,
75                                           tensor                            virial,
76                                           const bool                        doTempCouple,
77                                           gmx::ArrayRef<const t_grp_tcstat> tcstat,
78                                           const bool                        doPressureCouple,
79                                           const float                       dtPressureCouple,
80                                           const matrix                      velocityScalingMatrix)
81 {
82     // Clearing virial matrix
83     // TODO There is no point in having separate virial matrix for constraints
84     clear_mat(virial);
85
86     // Make sure that the forces are ready on device before proceeding with the update.
87     fReadyOnDevice->enqueueWaitEvent(commandStream_);
88
89     // The integrate should save a copy of the current coordinates in d_xp_ and write updated once into d_x_.
90     // The d_xp_ is only needed by constraints.
91     integrator_->integrate(d_x_, d_xp_, d_v_, d_f_, dt,
92                            doTempCouple, tcstat,
93                            doPressureCouple, dtPressureCouple, velocityScalingMatrix);
94     // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
95     // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
96     // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
97     lincsCuda_->apply(d_xp_, d_x_,
98                       updateVelocities, d_v_, 1.0/dt,
99                       computeVirial, virial);
100     settleCuda_->apply(d_xp_, d_x_,
101                        updateVelocities, d_v_, 1.0/dt,
102                        computeVirial, virial);
103
104     // scaledVirial -> virial (methods above returns scaled values)
105     float scaleFactor = 0.5f/(dt*dt);
106     for (int i = 0; i < DIM; i++)
107     {
108         for (int j = 0; j < DIM; j++)
109         {
110             virial[i][j] = scaleFactor*virial[i][j];
111         }
112     }
113
114     coordinatesReady_->markEvent(commandStream_);
115
116     return;
117 }
118
119 UpdateConstrainCuda::Impl::Impl(const t_inputrec     &ir,
120                                 const gmx_mtop_t     &mtop,
121                                 const void           *commandStream,
122                                 GpuEventSynchronizer *xUpdatedOnDevice) :
123     coordinatesReady_(xUpdatedOnDevice)
124 {
125     GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
126     commandStream != nullptr ? commandStream_ = *static_cast<const CommandStream*>(commandStream) : commandStream_ = nullptr;
127
128
129     integrator_ = std::make_unique<LeapFrogCuda>(commandStream_);
130     lincsCuda_  = std::make_unique<LincsCuda>(ir.nLincsIter, ir.nProjOrder, commandStream_);
131     settleCuda_ = std::make_unique<SettleCuda>(mtop, commandStream_);
132
133 }
134
135 UpdateConstrainCuda::Impl::~Impl()
136 {
137 }
138
139 void UpdateConstrainCuda::Impl::set(DeviceBuffer<float>        d_x,
140                                     DeviceBuffer<float>        d_v,
141                                     const DeviceBuffer<float>  d_f,
142                                     const t_idef              &idef,
143                                     const t_mdatoms           &md,
144                                     const int                  numTempScaleValues)
145 {
146     GMX_ASSERT(d_x != nullptr, "Coordinates device buffer should not be null.");
147     GMX_ASSERT(d_v != nullptr, "Velocities device buffer should not be null.");
148     GMX_ASSERT(d_f != nullptr, "Forces device buffer should not be null.");
149
150     d_x_ = reinterpret_cast<float3*>(d_x);
151     d_v_ = reinterpret_cast<float3*>(d_v);
152     d_f_ = reinterpret_cast<float3*>(d_f);
153
154     numAtoms_ = md.nr;
155
156     reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, nullptr);
157
158     reallocateDeviceBuffer(&d_inverseMasses_, numAtoms_,
159                            &numInverseMasses_, &numInverseMassesAlloc_, nullptr);
160
161     // Integrator should also update something, but it does not even have a method yet
162     integrator_->set(md, numTempScaleValues, md.cTC);
163     lincsCuda_->set(idef, md);
164     settleCuda_->set(idef, md);
165 }
166
167 void UpdateConstrainCuda::Impl::setPbc(const t_pbc *pbc)
168 {
169     setPbcAiuc(pbc->ndim_ePBC, pbc->box, &pbcAiuc_);
170     integrator_->setPbc(pbc);
171     lincsCuda_->setPbc(pbc);
172     settleCuda_->setPbc(pbc);
173 }
174
175 void UpdateConstrainCuda::Impl::waitCoordinatesReadyOnDevice()
176 {
177     coordinatesReady_->waitForEvent();
178 }
179
180 GpuEventSynchronizer* UpdateConstrainCuda::Impl::getCoordinatesReadySync()
181 {
182     return coordinatesReady_;
183 }
184
185 UpdateConstrainCuda::UpdateConstrainCuda(const t_inputrec     &ir,
186                                          const gmx_mtop_t     &mtop,
187                                          const void           *commandStream,
188                                          GpuEventSynchronizer *xUpdatedOnDevice)
189     : impl_(new Impl(ir, mtop, commandStream, xUpdatedOnDevice))
190 {
191 }
192
193 UpdateConstrainCuda::~UpdateConstrainCuda() = default;
194
195 void UpdateConstrainCuda::integrate(GpuEventSynchronizer             *fReadyOnDevice,
196                                     const real                        dt,
197                                     const bool                        updateVelocities,
198                                     const bool                        computeVirial,
199                                     tensor                            virialScaled,
200                                     const bool                        doTempCouple,
201                                     gmx::ArrayRef<const t_grp_tcstat> tcstat,
202                                     const bool                        doPressureCouple,
203                                     const float                       dtPressureCouple,
204                                     const matrix                      velocityScalingMatrix)
205 {
206     impl_->integrate(fReadyOnDevice,
207                      dt, updateVelocities, computeVirial, virialScaled,
208                      doTempCouple, tcstat,
209                      doPressureCouple, dtPressureCouple, velocityScalingMatrix);
210 }
211
212 void UpdateConstrainCuda::set(DeviceBuffer<float>        d_x,
213                               DeviceBuffer<float>        d_v,
214                               const DeviceBuffer<float>  d_f,
215                               const t_idef              &idef,
216                               const t_mdatoms           &md,
217                               const int                  numTempScaleValues)
218 {
219     impl_->set(d_x, d_v, d_f, idef, md, numTempScaleValues);
220 }
221
222 void UpdateConstrainCuda::setPbc(const t_pbc *pbc)
223 {
224     impl_->setPbc(pbc);
225 }
226
227 void UpdateConstrainCuda::waitCoordinatesReadyOnDevice()
228 {
229     impl_->waitCoordinatesReadyOnDevice();
230 }
231
232 GpuEventSynchronizer* UpdateConstrainCuda::getCoordinatesReadySync()
233 {
234     return impl_->getCoordinatesReadySync();
235 }
236
237 } //namespace gmx