Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / math / gausstransform.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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 Gaussian function evaluations on lattices and related functionality
38  *
39  * \author Christian Blau <blau@kth.se>
40  *
41  * \ingroup module_math
42  */
43 #include "gmxpre.h"
44
45 #include "gausstransform.h"
46
47 #include <cmath>
48
49 #include <algorithm>
50 #include <array>
51
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/multidimarray.h"
54 #include "gromacs/math/utilities.h"
55
56 namespace gmx
57 {
58
59 /********************************************************************
60  * GaussianOn1DLattice::Impl
61  */
62
63 class GaussianOn1DLattice::Impl
64 {
65 public:
66     Impl(int numGridPointsForSpreadingHalfWidth, real sigma);
67     ~Impl()                 = default;
68     Impl(const Impl& other) = default;
69     Impl& operator=(const Impl& other) = default;
70
71     /*! \brief evaluate Gaussian function at all lattice points
72      * \param[in] amplitude the amplitude of the Gaussian
73      * \param[in] dx distance from the center
74      */
75     void spread(double amplitude, real dx);
76     //! Largest distance in number of gridpoints from 0
77     int numGridPointsForSpreadingHalfWidth_;
78     /*! \brief Avoid overflow for E2^offset and underflow for E3(i).
79      *
80      * Occurs when sigma is much smaller than numGridPointsForSpreadingHalfWidth_.
81      *
82      * E2^offset smaller than maximum float requires
83      * \f$exp(dx / (2*square(sigma))^numGridPointsForSpreadingHalfWidth_ \leq max_float \f$
84      * The maximum expected distance of the Gaussian center to the next lattice point is dx = 0.5,
85      * thus the maximum spread distance here is \f$4 * sigma^2 * \log(\mathrm{maxfloat})\f$ .
86      *
87      * E3(i) larger than minmium float requires
88      * exp(i^2 / 2*(sigma)^2) > min_float
89      * Thus the maximum spread distance here is \f$\sigma \sqrt(-2\log(\mathrm{minfloat}))\f$
90      */
91     int maxEvaluatedSpreadDistance_;
92     //! Width of the Gaussian function
93     double sigma_;
94     //! The result of the spreading calculation
95     std::vector<float> spreadingResult_;
96     //! Pre-calculated exp(-gridIndex^2/2 * (sigma^2)) named as in Greengard2004
97     std::vector<float> e3_;
98     /*! \brief Equal to std::floor(std::log(std::numeric_limits<float>::max())).
99      * Above expression is not constexpr and a const variable would implicitly delete default copy
100      * assignment. Therefore resorting to setting number manually.
101      */
102     static constexpr double c_logMaxFloat = 88.72284;
103     static constexpr double c_logMinFloat = -87.33654;
104 };
105
106 GaussianOn1DLattice::Impl::Impl(int numGridPointsForSpreadingHalfWidth, real sigma) :
107     numGridPointsForSpreadingHalfWidth_(numGridPointsForSpreadingHalfWidth),
108     sigma_(sigma),
109     spreadingResult_(2 * numGridPointsForSpreadingHalfWidth + 1)
110 {
111     maxEvaluatedSpreadDistance_ =
112             std::min(numGridPointsForSpreadingHalfWidth_,
113                      static_cast<int>(std::floor(4 * square(sigma) * c_logMaxFloat)) - 1);
114     maxEvaluatedSpreadDistance_ =
115             std::min(maxEvaluatedSpreadDistance_,
116                      static_cast<int>(std::floor(sigma * sqrt(-2.0 * c_logMinFloat))) - 1);
117
118     std::generate_n(
119             std::back_inserter(e3_), maxEvaluatedSpreadDistance_ + 1, [sigma, latticeIndex = 0]() mutable {
120                 return std::exp(-0.5 * square(latticeIndex++ / sigma));
121             });
122
123     std::fill(std::begin(spreadingResult_), std::end(spreadingResult_), 0.);
124 };
125
126 void GaussianOn1DLattice::Impl::spread(double amplitude, real dx)
127 {
128     /* The spreading routine implements the fast gaussian gridding as in
129      *
130      * Leslie Greengard and June-Yub Lee,
131      * "Accelerating the Nonuniform Fast Fourier Transform"
132      * SIAM REV 2004 Vol. 46, No. 3, pp. 443-454 DOI. 10.1137/S003614450343200X
133      *
134      * Following the naming conventions for e1, e2 and e3, nu = 1, m = numGridPointsForSpreadingHalfWidth_.
135      *
136      * Speed up is achieved by factorization of the exponential that is evaluted
137      * at regular lattice points i, where the distance from the
138      * Gaussian center is \f$x-i\f$:
139      *
140      * \f[
141      *      a * \exp(-(x^2-2*i*x+ i^2)/(2*\sigma^2)) =
142      *      a * \exp(-x^2/2*\sigma^2) * \exp(x/\sigma^2)^i * \exp(i/2*sigma^2) =
143      *      e_1(x) * e_2(x)^i * e_3(i)
144      * \f]
145      *
146      * Requiring only two exp evaluations per spreading operation.
147      *
148      */
149     const double e1 = amplitude * exp(-0.5 * dx * dx / square(sigma_)) / (sqrt(2 * M_PI) * sigma_);
150     spreadingResult_[numGridPointsForSpreadingHalfWidth_] = e1;
151
152     const double e2 = exp(dx / square(sigma_));
153
154     double e2pow = e2; //< powers of e2, e2^offset
155
156     // Move outwards from mid-point, using e2pow value for both points simultaneously
157     //      o    o    o<----O---->o    o    o
158     for (int offset = 1; offset < maxEvaluatedSpreadDistance_; offset++)
159     {
160         const double e1_3                                              = e1 * e3_[offset];
161         spreadingResult_[numGridPointsForSpreadingHalfWidth_ + offset] = e1_3 * e2pow;
162         spreadingResult_[numGridPointsForSpreadingHalfWidth_ - offset] = e1_3 / e2pow;
163         e2pow *= e2;
164     }
165     // separate statement for gridpoints at the end of the range avoids
166     // overflow for large sigma and saves one e2 multiplication operation
167     spreadingResult_[numGridPointsForSpreadingHalfWidth_ - maxEvaluatedSpreadDistance_] =
168             (e1 / e2pow) * e3_[maxEvaluatedSpreadDistance_];
169     spreadingResult_[numGridPointsForSpreadingHalfWidth_ + maxEvaluatedSpreadDistance_] =
170             (e1 * e2pow) * e3_[maxEvaluatedSpreadDistance_];
171 }
172
173 /********************************************************************
174  * GaussianOn1DLattice
175  */
176
177 GaussianOn1DLattice::GaussianOn1DLattice(int numGridPointsForSpreadingHalfWidth_, real sigma) :
178     impl_(new Impl(numGridPointsForSpreadingHalfWidth_, sigma))
179 {
180 }
181
182 GaussianOn1DLattice::~GaussianOn1DLattice() {}
183
184 void GaussianOn1DLattice::spread(double amplitude, real dx)
185 {
186     impl_->spread(amplitude, dx);
187 }
188
189 ArrayRef<const float> GaussianOn1DLattice::view()
190 {
191     return impl_->spreadingResult_;
192 }
193
194 GaussianOn1DLattice::GaussianOn1DLattice(const GaussianOn1DLattice& other) :
195     impl_(new Impl(*other.impl_))
196 {
197 }
198
199 GaussianOn1DLattice& GaussianOn1DLattice::operator=(const GaussianOn1DLattice& other)
200 {
201     *impl_ = *other.impl_;
202     return *this;
203 }
204
205 GaussianOn1DLattice::GaussianOn1DLattice(GaussianOn1DLattice&&) noexcept = default;
206
207 GaussianOn1DLattice& GaussianOn1DLattice::operator=(GaussianOn1DLattice&&) noexcept = default;
208
209 namespace
210 {
211
212 //! rounds real-valued coordinate to the closest integer values
213 IVec closestIntegerPoint(const RVec& coordinate)
214 {
215     return { roundToInt(coordinate[XX]), roundToInt(coordinate[YY]), roundToInt(coordinate[ZZ]) };
216 }
217
218 /*! \brief Substracts a range from a three-dimensional integer coordinate and ensures
219  * the resulting coordinate is within a lattice.
220  * \param[in] index point in lattice
221  * \param[in] range to be shifted
222  * \returns Shifted index or zero if shifted index is smaller than zero.
223  */
224 IVec rangeBeginWithinLattice(const IVec& index, const IVec& range)
225 {
226     return elementWiseMax({ 0, 0, 0 }, index - range);
227 }
228
229 /*! \brief Adds a range from a three-dimensional integer coordinate and ensures
230  * the resulting coordinate is within a lattice.
231  * \param[in] index point in lattice
232  * \param[in] extents extent of the lattice
233  * \param[in] range to be shifted
234  * \returns Shifted index or the lattice extent if shifted index is larger than the extent
235  */
236 IVec rangeEndWithinLattice(const IVec& index, const dynamicExtents3D& extents, const IVec& range)
237 {
238     IVec extentAsIvec(static_cast<int>(extents.extent(ZZ)),
239                       static_cast<int>(extents.extent(YY)),
240                       static_cast<int>(extents.extent(XX)));
241     return elementWiseMin(extentAsIvec, index + range);
242 }
243
244
245 } // namespace
246
247 /********************************************************************
248  * OuterProductEvaluator
249  */
250
251 mdspan<const float, dynamic_extent, dynamic_extent> OuterProductEvaluator::
252                                                     operator()(ArrayRef<const float> x, ArrayRef<const float> y)
253 {
254     data_.resize(ssize(x), ssize(y));
255     for (gmx::index xIndex = 0; xIndex < ssize(x); ++xIndex)
256     {
257         const auto xValue = x[xIndex];
258         std::transform(std::begin(y), std::end(y), begin(data_.asView()[xIndex]), [xValue](float yValue) {
259             return xValue * yValue;
260         });
261     }
262     return data_.asConstView();
263 }
264
265 /********************************************************************
266  * IntegerBox
267  */
268
269 IntegerBox::IntegerBox(const IVec& begin, const IVec& end) : begin_{ begin }, end_{ end } {}
270
271 const IVec& IntegerBox::begin() const
272 {
273     return begin_;
274 }
275 const IVec& IntegerBox::end() const
276 {
277     return end_;
278 }
279
280 bool IntegerBox::empty() const
281 {
282     return !((begin_[XX] < end_[XX]) && (begin_[YY] < end_[YY]) && (begin_[ZZ] < end_[ZZ]));
283 }
284
285 IntegerBox spreadRangeWithinLattice(const IVec& center, dynamicExtents3D extent, IVec range)
286 {
287     const IVec begin = rangeBeginWithinLattice(center, range);
288     const IVec end   = rangeEndWithinLattice(center, extent, range);
289     return { begin, end };
290 }
291 /********************************************************************
292  * GaussianSpreadKernel
293  */
294
295 IVec GaussianSpreadKernelParameters::Shape::latticeSpreadRange() const
296 {
297     DVec range(std::ceil(sigma_[XX] * spreadWidthMultiplesOfSigma_),
298                std::ceil(sigma_[YY] * spreadWidthMultiplesOfSigma_),
299                std::ceil(sigma_[ZZ] * spreadWidthMultiplesOfSigma_));
300     return range.toIVec();
301 }
302
303 /********************************************************************
304  * GaussTransform3D::Impl
305  */
306
307 /*! \internal \brief
308  * Private implementation class for GaussTransform3D.
309  */
310 class GaussTransform3D::Impl
311 {
312 public:
313     //! Construct from extent and spreading width and range
314     Impl(const dynamicExtents3D& extent, const GaussianSpreadKernelParameters::Shape& kernelShapeParameters);
315     ~Impl() = default;
316     //! Copy constructor
317     Impl(const Impl& other) = default;
318     //! Copy assignment
319     Impl& operator=(const Impl& other) = default;
320     //! Add another gaussian
321     void add(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParamters);
322     //! The width of the Gaussian in lattice spacing units
323     BasicVector<double> sigma_;
324     //! The spread range in lattice points
325     IVec spreadRange_;
326     //! The result of the Gauss transform
327     MultiDimArray<std::vector<float>, dynamicExtents3D> data_;
328     //! The outer product of a Gaussian along the z and y dimension
329     OuterProductEvaluator outerProductZY_;
330     //! The three one-dimensional Gaussians, whose outer product is added to the Gauss transform
331     std::array<GaussianOn1DLattice, DIM> gauss1d_;
332 };
333
334 GaussTransform3D::Impl::Impl(const dynamicExtents3D&                      extent,
335                              const GaussianSpreadKernelParameters::Shape& kernelShapeParameters) :
336     sigma_{ kernelShapeParameters.sigma_ },
337     spreadRange_{ kernelShapeParameters.latticeSpreadRange() },
338     data_{ extent },
339     gauss1d_({ GaussianOn1DLattice(spreadRange_[XX], sigma_[XX]),
340                GaussianOn1DLattice(spreadRange_[YY], sigma_[YY]),
341                GaussianOn1DLattice(spreadRange_[ZZ], sigma_[ZZ]) })
342 {
343 }
344
345 void GaussTransform3D::Impl::add(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters)
346 {
347     const IVec closestLatticePoint = closestIntegerPoint(localParameters.coordinate_);
348     const auto spreadRange =
349             spreadRangeWithinLattice(closestLatticePoint, data_.asView().extents(), spreadRange_);
350
351     // do nothing if the added Gaussian will never reach the lattice
352     if (spreadRange.empty())
353     {
354         return;
355     }
356
357     for (int dimension = XX; dimension <= ZZ; ++dimension)
358     {
359         // multiply with amplitude so that Gauss3D = (amplitude * Gauss_x) * Gauss_y * Gauss_z
360         const float gauss1DAmplitude = dimension > XX ? 1.0 : localParameters.amplitude_;
361         gauss1d_[dimension].spread(
362                 gauss1DAmplitude, localParameters.coordinate_[dimension] - closestLatticePoint[dimension]);
363     }
364
365     const auto spreadZY         = outerProductZY_(gauss1d_[ZZ].view(), gauss1d_[YY].view());
366     const auto spreadX          = gauss1d_[XX].view();
367     const IVec spreadGridOffset = spreadRange_ - closestLatticePoint;
368
369     // \todo optimize these loops if performance critical
370     // The looping strategy uses that the last, x-dimension is contiguous in the memory layout
371     for (int zLatticeIndex = spreadRange.begin()[ZZ]; zLatticeIndex < spreadRange.end()[ZZ]; ++zLatticeIndex)
372     {
373         const auto zSlice = data_.asView()[zLatticeIndex];
374
375         for (int yLatticeIndex = spreadRange.begin()[YY]; yLatticeIndex < spreadRange.end()[YY]; ++yLatticeIndex)
376         {
377             const auto  ySlice      = zSlice[yLatticeIndex];
378             const float zyPrefactor = spreadZY(zLatticeIndex + spreadGridOffset[ZZ],
379                                                yLatticeIndex + spreadGridOffset[YY]);
380
381             for (int xLatticeIndex = spreadRange.begin()[XX]; xLatticeIndex < spreadRange.end()[XX];
382                  ++xLatticeIndex)
383             {
384                 const float xPrefactor = spreadX[xLatticeIndex + spreadGridOffset[XX]];
385                 ySlice[xLatticeIndex] += zyPrefactor * xPrefactor;
386             }
387         }
388     }
389 }
390
391 /********************************************************************
392  * GaussTransform3D
393  */
394
395 GaussTransform3D::GaussTransform3D(const dynamicExtents3D&                      extent,
396                                    const GaussianSpreadKernelParameters::Shape& kernelShapeParameters) :
397     impl_(new Impl(extent, kernelShapeParameters))
398 {
399 }
400
401 void GaussTransform3D::add(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters)
402 {
403     impl_->add(localParameters);
404 }
405
406 void GaussTransform3D::setZero()
407 {
408     std::fill(begin(impl_->data_), end(impl_->data_), 0.);
409 }
410
411 basic_mdspan<float, dynamicExtents3D> GaussTransform3D::view()
412 {
413     return impl_->data_.asView();
414 }
415
416 basic_mdspan<const float, dynamicExtents3D> GaussTransform3D::constView() const
417 {
418     return impl_->data_.asConstView();
419 }
420
421 GaussTransform3D::~GaussTransform3D() {}
422
423 GaussTransform3D::GaussTransform3D(const GaussTransform3D& other) : impl_(new Impl(*other.impl_)) {}
424
425 GaussTransform3D& GaussTransform3D::operator=(const GaussTransform3D& other)
426 {
427     *impl_ = *other.impl_;
428     return *this;
429 }
430
431 GaussTransform3D::GaussTransform3D(GaussTransform3D&&) noexcept = default;
432
433 GaussTransform3D& GaussTransform3D::operator=(GaussTransform3D&&) noexcept = default;
434
435 } // namespace gmx