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