Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / math / densityfittingforce.cpp
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  * \brief
37  * Implements density fitting forces.
38  *
39  * \author Christian Blau <blau@kth.se>
40  * \ingroup module_math
41  */
42 #include "gmxpre.h"
43
44 #include "densityfittingforce.h"
45
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/multidimarray.h"
48
49 namespace gmx
50 {
51
52 /********************************************************************
53  * DensityFittingForce::Impl
54  */
55
56 /*! \internal \brief
57  * Private implementation class for DensityFittingForce.
58  */
59 class DensityFittingForce::Impl
60 {
61 public:
62     /*! \brief Construct densityfitting force implementation from
63      * spread of function used to generate the density and spread range.
64      * \param[in] kernelShapeParameters determine the shape of the spreading kernel
65      */
66     explicit Impl(const GaussianSpreadKernelParameters::Shape& kernelShapeParameters);
67     //! \copydoc DensityFittingForce::evaluateForce(const DensitySpreadKernelParameters::PositionAndAmplitude & localParameters, basic_mdspan<const float, dynamicExtents3D> densityDerivative)
68     RVec evaluateForce(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters,
69                        basic_mdspan<const float, dynamicExtents3D> densityDerivative);
70     //! The width of the Gaussian in lattice spacing units
71     DVec sigma_;
72     //! The spread range in lattice points
73     IVec latticeSpreadRange_;
74     //! The three one-dimensional Gaussians that are used in the force calculation
75     std::array<GaussianOn1DLattice, DIM> gauss1d_;
76     //! The outer product of a Gaussian along the z and y dimension
77     OuterProductEvaluator outerProductZY_;
78 };
79
80 DensityFittingForce::Impl::Impl(const GaussianSpreadKernelParameters::Shape& kernelShapeParameters) :
81     sigma_{ kernelShapeParameters.sigma_ },
82     latticeSpreadRange_{ kernelShapeParameters.latticeSpreadRange() },
83     gauss1d_({ GaussianOn1DLattice(latticeSpreadRange_[XX], sigma_[XX]),
84                GaussianOn1DLattice(latticeSpreadRange_[YY], sigma_[YY]),
85                GaussianOn1DLattice(latticeSpreadRange_[ZZ], sigma_[ZZ]) })
86 {
87 }
88
89 RVec DensityFittingForce::Impl::evaluateForce(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters,
90                                               basic_mdspan<const float, dynamicExtents3D> densityDerivative)
91 {
92     const IVec closestLatticePoint(roundToInt(localParameters.coordinate_[XX]),
93                                    roundToInt(localParameters.coordinate_[YY]),
94                                    roundToInt(localParameters.coordinate_[ZZ]));
95     const auto spreadRange = spreadRangeWithinLattice(
96             closestLatticePoint, densityDerivative.extents(), latticeSpreadRange_);
97
98     // do nothing if the added Gaussian will never reach the lattice
99     if (spreadRange.empty())
100     {
101         return {};
102     }
103
104     for (int dimension = XX; dimension <= ZZ; ++dimension)
105     {
106         // multiply with amplitude so that Gauss3D = (amplitude * Gauss_x) * Gauss_y * Gauss_z
107         const float gauss1DAmplitude = dimension > XX ? 1.0 : localParameters.amplitude_;
108         gauss1d_[dimension].spread(gauss1DAmplitude, localParameters.coordinate_[dimension]
109                                                              - closestLatticePoint[dimension]);
110     }
111
112     const auto spreadZY = outerProductZY_(gauss1d_[ZZ].view(), gauss1d_[YY].view());
113     const auto spreadX  = gauss1d_[XX].view();
114     const IVec spreadGridOffset(latticeSpreadRange_[XX] - closestLatticePoint[XX],
115                                 latticeSpreadRange_[YY] - closestLatticePoint[YY],
116                                 latticeSpreadRange_[ZZ] - closestLatticePoint[ZZ]);
117
118     const DVec differenceVectorScale  = { 1. / (square(sigma_[XX])), 1. / (square(sigma_[YY])),
119                                          1. / (square(sigma_[ZZ])) };
120     const DVec differenceVectorOffset = scaleByVector(
121             spreadRange.begin().toDVec() - localParameters.coordinate_.toDVec(), differenceVectorScale);
122
123     DVec differenceVector = differenceVectorOffset;
124
125     DVec force = { 0., 0., 0. };
126
127     for (int zLatticeIndex = spreadRange.begin()[ZZ]; zLatticeIndex < spreadRange.end()[ZZ];
128          ++zLatticeIndex, differenceVector[ZZ] += differenceVectorScale[ZZ])
129     {
130         auto zSliceOfDerivative = densityDerivative[zLatticeIndex];
131
132         differenceVector[YY] = differenceVectorOffset[YY];
133         for (int yLatticeIndex = spreadRange.begin()[YY]; yLatticeIndex < spreadRange.end()[YY];
134              ++yLatticeIndex, differenceVector[YY] += differenceVectorScale[YY])
135         {
136             auto       ySliceOfDerivative = zSliceOfDerivative[yLatticeIndex];
137             const auto zyPrefactor        = spreadZY(zLatticeIndex + spreadGridOffset[ZZ],
138                                               yLatticeIndex + spreadGridOffset[YY]);
139
140             differenceVector[XX] = differenceVectorOffset[XX];
141             for (int xLatticeIndex = spreadRange.begin()[XX]; xLatticeIndex < spreadRange.end()[XX];
142                  ++xLatticeIndex, differenceVector[XX] += differenceVectorScale[XX])
143             {
144                 const double preFactor = zyPrefactor * spreadX[xLatticeIndex + spreadGridOffset[XX]]
145                                          * ySliceOfDerivative[xLatticeIndex];
146                 force += preFactor * differenceVector;
147             }
148         }
149     }
150     return localParameters.amplitude_ * force.toRVec();
151 }
152
153 /********************************************************************
154  * DensityFittingForce
155  */
156
157 DensityFittingForce::DensityFittingForce(const GaussianSpreadKernelParameters::Shape& kernelShapeParameters) :
158     impl_(new Impl(kernelShapeParameters))
159 {
160 }
161
162 RVec DensityFittingForce::evaluateForce(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters,
163                                         basic_mdspan<const float, dynamicExtents3D> densityDerivative)
164 {
165     return impl_->evaluateForce(localParameters, densityDerivative);
166 }
167
168 DensityFittingForce::~DensityFittingForce() {}
169
170 DensityFittingForce::DensityFittingForce(const DensityFittingForce& other) :
171     impl_(new Impl(*other.impl_))
172 {
173 }
174
175 DensityFittingForce& DensityFittingForce::operator=(const DensityFittingForce& other)
176 {
177     *impl_ = *other.impl_;
178     return *this;
179 }
180
181 DensityFittingForce::DensityFittingForce(DensityFittingForce&&) noexcept = default;
182
183 DensityFittingForce& DensityFittingForce::operator=(DensityFittingForce&&) noexcept = default;
184
185 } // namespace gmx