b22f365358da3a593ae10f53b233dd89e7fb4ea1
[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     return;
111 }
112
113 UpdateConstrainCuda::Impl::Impl(const t_inputrec  &ir,
114                                 const gmx_mtop_t  &mtop,
115                                 const void        *commandStream)
116 {
117     commandStream != nullptr ? commandStream_ = *static_cast<const CommandStream*>(commandStream) : commandStream_ = nullptr;
118
119     integrator_ = std::make_unique<LeapFrogCuda>(commandStream_);
120     lincsCuda_  = std::make_unique<LincsCuda>(ir.nLincsIter, ir.nProjOrder, commandStream_);
121     settleCuda_ = std::make_unique<SettleCuda>(mtop, commandStream_);
122
123 }
124
125 UpdateConstrainCuda::Impl::~Impl()
126 {
127 }
128
129 void UpdateConstrainCuda::Impl::set(DeviceBuffer<float>        d_x,
130                                     DeviceBuffer<float>        d_v,
131                                     const DeviceBuffer<float>  d_f,
132                                     const t_idef              &idef,
133                                     const t_mdatoms           &md,
134                                     const int                  numTempScaleValues)
135 {
136     GMX_ASSERT(d_x != nullptr, "Coordinates device buffer should not be null.");
137     GMX_ASSERT(d_v != nullptr, "Velocities device buffer should not be null.");
138     GMX_ASSERT(d_f != nullptr, "Forces device buffer should not be null.");
139
140     d_x_ = reinterpret_cast<float3*>(d_x);
141     d_v_ = reinterpret_cast<float3*>(d_v);
142     d_f_ = reinterpret_cast<float3*>(d_f);
143
144     numAtoms_ = md.nr;
145
146     reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, nullptr);
147
148     reallocateDeviceBuffer(&d_inverseMasses_, numAtoms_,
149                            &numInverseMasses_, &numInverseMassesAlloc_, nullptr);
150
151     // Integrator should also update something, but it does not even have a method yet
152     integrator_->set(md, numTempScaleValues, md.cTC);
153     lincsCuda_->set(idef, md);
154     settleCuda_->set(idef, md);
155 }
156
157 void UpdateConstrainCuda::Impl::setPbc(const t_pbc *pbc)
158 {
159     setPbcAiuc(pbc->ndim_ePBC, pbc->box, &pbcAiuc_);
160     integrator_->setPbc(pbc);
161     lincsCuda_->setPbc(pbc);
162     settleCuda_->setPbc(pbc);
163 }
164
165 UpdateConstrainCuda::UpdateConstrainCuda(const t_inputrec  &ir,
166                                          const gmx_mtop_t  &mtop,
167                                          const void        *commandStream)
168     : impl_(new Impl(ir, mtop, commandStream))
169 {
170 }
171
172 UpdateConstrainCuda::~UpdateConstrainCuda() = default;
173
174 void UpdateConstrainCuda::integrate(const real                        dt,
175                                     const bool                        updateVelocities,
176                                     const bool                        computeVirial,
177                                     tensor                            virialScaled,
178                                     const bool                        doTempCouple,
179                                     gmx::ArrayRef<const t_grp_tcstat> tcstat,
180                                     const bool                        doPressureCouple,
181                                     const float                       dtPressureCouple,
182                                     const matrix                      velocityScalingMatrix)
183 {
184     impl_->integrate(dt, updateVelocities, computeVirial, virialScaled,
185                      doTempCouple, tcstat,
186                      doPressureCouple, dtPressureCouple, velocityScalingMatrix);
187 }
188
189 void UpdateConstrainCuda::set(DeviceBuffer<float>        d_x,
190                               DeviceBuffer<float>        d_v,
191                               const DeviceBuffer<float>  d_f,
192                               const t_idef              &idef,
193                               const t_mdatoms           &md,
194                               const int                  numTempScaleValues)
195 {
196     impl_->set(d_x, d_v, d_f, idef, md, numTempScaleValues);
197 }
198
199 void UpdateConstrainCuda::setPbc(const t_pbc *pbc)
200 {
201     impl_->setPbc(pbc);
202 }
203
204 } //namespace gmx