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