Apply clang-format to source tree
[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
90     // into d_x_. The d_xp_ is only needed by constraints.
91     integrator_->integrate(d_x_, d_xp_, d_v_, d_f_, dt, doTempCouple, tcstat, doPressureCouple,
92                            dtPressureCouple, velocityScalingMatrix);
93     // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
94     // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
95     // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
96     lincsCuda_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial);
97     settleCuda_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial);
98
99     // scaledVirial -> virial (methods above returns scaled values)
100     float scaleFactor = 0.5f / (dt * dt);
101     for (int i = 0; i < DIM; i++)
102     {
103         for (int j = 0; j < DIM; j++)
104         {
105             virial[i][j] = scaleFactor * virial[i][j];
106         }
107     }
108
109     coordinatesReady_->markEvent(commandStream_);
110
111     return;
112 }
113
114 UpdateConstrainCuda::Impl::Impl(const t_inputrec&     ir,
115                                 const gmx_mtop_t&     mtop,
116                                 const void*           commandStream,
117                                 GpuEventSynchronizer* xUpdatedOnDevice) :
118     coordinatesReady_(xUpdatedOnDevice)
119 {
120     GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
121     commandStream != nullptr ? commandStream_ = *static_cast<const CommandStream*>(commandStream)
122                              : commandStream_ = nullptr;
123
124
125     integrator_ = std::make_unique<LeapFrogCuda>(commandStream_);
126     lincsCuda_  = std::make_unique<LincsCuda>(ir.nLincsIter, ir.nProjOrder, commandStream_);
127     settleCuda_ = std::make_unique<SettleCuda>(mtop, commandStream_);
128 }
129
130 UpdateConstrainCuda::Impl::~Impl() {}
131
132 void UpdateConstrainCuda::Impl::set(DeviceBuffer<float>       d_x,
133                                     DeviceBuffer<float>       d_v,
134                                     const DeviceBuffer<float> d_f,
135                                     const t_idef&             idef,
136                                     const t_mdatoms&          md,
137                                     const int                 numTempScaleValues)
138 {
139     GMX_ASSERT(d_x != nullptr, "Coordinates device buffer should not be null.");
140     GMX_ASSERT(d_v != nullptr, "Velocities device buffer should not be null.");
141     GMX_ASSERT(d_f != nullptr, "Forces device buffer should not be null.");
142
143     d_x_ = reinterpret_cast<float3*>(d_x);
144     d_v_ = reinterpret_cast<float3*>(d_v);
145     d_f_ = reinterpret_cast<float3*>(d_f);
146
147     numAtoms_ = md.nr;
148
149     reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, nullptr);
150
151     reallocateDeviceBuffer(&d_inverseMasses_, numAtoms_, &numInverseMasses_,
152                            &numInverseMassesAlloc_, nullptr);
153
154     // Integrator should also update something, but it does not even have a method yet
155     integrator_->set(md, numTempScaleValues, md.cTC);
156     lincsCuda_->set(idef, md);
157     settleCuda_->set(idef, md);
158 }
159
160 void UpdateConstrainCuda::Impl::setPbc(const t_pbc* pbc)
161 {
162     setPbcAiuc(pbc->ndim_ePBC, pbc->box, &pbcAiuc_);
163     integrator_->setPbc(pbc);
164     lincsCuda_->setPbc(pbc);
165     settleCuda_->setPbc(pbc);
166 }
167
168 GpuEventSynchronizer* UpdateConstrainCuda::Impl::getCoordinatesReadySync()
169 {
170     return coordinatesReady_;
171 }
172
173 UpdateConstrainCuda::UpdateConstrainCuda(const t_inputrec&     ir,
174                                          const gmx_mtop_t&     mtop,
175                                          const void*           commandStream,
176                                          GpuEventSynchronizer* xUpdatedOnDevice) :
177     impl_(new Impl(ir, mtop, commandStream, xUpdatedOnDevice))
178 {
179 }
180
181 UpdateConstrainCuda::~UpdateConstrainCuda() = default;
182
183 void UpdateConstrainCuda::integrate(GpuEventSynchronizer*             fReadyOnDevice,
184                                     const real                        dt,
185                                     const bool                        updateVelocities,
186                                     const bool                        computeVirial,
187                                     tensor                            virialScaled,
188                                     const bool                        doTempCouple,
189                                     gmx::ArrayRef<const t_grp_tcstat> tcstat,
190                                     const bool                        doPressureCouple,
191                                     const float                       dtPressureCouple,
192                                     const matrix                      velocityScalingMatrix)
193 {
194     impl_->integrate(fReadyOnDevice, dt, updateVelocities, computeVirial, virialScaled, doTempCouple,
195                      tcstat, doPressureCouple, dtPressureCouple, velocityScalingMatrix);
196 }
197
198 void UpdateConstrainCuda::set(DeviceBuffer<float>       d_x,
199                               DeviceBuffer<float>       d_v,
200                               const DeviceBuffer<float> d_f,
201                               const t_idef&             idef,
202                               const t_mdatoms&          md,
203                               const int                 numTempScaleValues)
204 {
205     impl_->set(d_x, d_v, d_f, idef, md, numTempScaleValues);
206 }
207
208 void UpdateConstrainCuda::setPbc(const t_pbc* pbc)
209 {
210     impl_->setPbc(pbc);
211 }
212
213 GpuEventSynchronizer* UpdateConstrainCuda::getCoordinatesReadySync()
214 {
215     return impl_->getCoordinatesReadySync();
216 }
217
218 } // namespace gmx