Update instructions in containers.rst
[gromacs.git] / src / gromacs / math / gausstransform.cpp
blob111e2c7ccd0510f14459a4eae22dc820a1d51c3b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements Gaussian function evaluations on lattices and related functionality
39 * \author Christian Blau <blau@kth.se>
41 * \ingroup module_math
43 #include "gmxpre.h"
45 #include "gausstransform.h"
47 #include <cmath>
49 #include <algorithm>
50 #include <array>
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/multidimarray.h"
54 #include "gromacs/math/utilities.h"
56 namespace gmx
59 /********************************************************************
60 * GaussianOn1DLattice::Impl
63 class GaussianOn1DLattice::Impl
65 public:
66 Impl(int numGridPointsForSpreadingHalfWidth, real sigma);
67 ~Impl() = default;
68 Impl(const Impl& other) = default;
69 Impl& operator=(const Impl& other) = default;
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
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).
80 * Occurs when sigma is much smaller than numGridPointsForSpreadingHalfWidth_.
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$ .
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$
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.
102 static constexpr double c_logMaxFloat = 88.72284;
103 static constexpr double c_logMinFloat = -87.33654;
106 GaussianOn1DLattice::Impl::Impl(int numGridPointsForSpreadingHalfWidth, real sigma) :
107 numGridPointsForSpreadingHalfWidth_(numGridPointsForSpreadingHalfWidth),
108 sigma_(sigma),
109 spreadingResult_(2 * numGridPointsForSpreadingHalfWidth + 1)
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);
118 std::generate_n(std::back_inserter(e3_), maxEvaluatedSpreadDistance_ + 1,
119 [sigma, latticeIndex = 0]() mutable {
120 return std::exp(-0.5 * square(latticeIndex++ / sigma));
123 std::fill(std::begin(spreadingResult_), std::end(spreadingResult_), 0.);
126 void GaussianOn1DLattice::Impl::spread(double amplitude, real dx)
128 /* The spreading routine implements the fast gaussian gridding as in
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
134 * Following the naming conventions for e1, e2 and e3, nu = 1, m = numGridPointsForSpreadingHalfWidth_.
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$:
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]
146 * Requiring only two exp evaluations per spreading operation.
149 const double e1 = amplitude * exp(-0.5 * dx * dx / square(sigma_)) / (sqrt(2 * M_PI) * sigma_);
150 spreadingResult_[numGridPointsForSpreadingHalfWidth_] = e1;
152 const double e2 = exp(dx / square(sigma_));
154 double e2pow = e2; //< powers of e2, e2^offset
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++)
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;
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_];
173 /********************************************************************
174 * GaussianOn1DLattice
177 GaussianOn1DLattice::GaussianOn1DLattice(int numGridPointsForSpreadingHalfWidth_, real sigma) :
178 impl_(new Impl(numGridPointsForSpreadingHalfWidth_, sigma))
182 GaussianOn1DLattice::~GaussianOn1DLattice() {}
184 void GaussianOn1DLattice::spread(double amplitude, real dx)
186 impl_->spread(amplitude, dx);
189 ArrayRef<const float> GaussianOn1DLattice::view()
191 return impl_->spreadingResult_;
194 GaussianOn1DLattice::GaussianOn1DLattice(const GaussianOn1DLattice& other) :
195 impl_(new Impl(*other.impl_))
199 GaussianOn1DLattice& GaussianOn1DLattice::operator=(const GaussianOn1DLattice& other)
201 *impl_ = *other.impl_;
202 return *this;
205 GaussianOn1DLattice::GaussianOn1DLattice(GaussianOn1DLattice&&) noexcept = default;
207 GaussianOn1DLattice& GaussianOn1DLattice::operator=(GaussianOn1DLattice&&) noexcept = default;
209 namespace
212 //! rounds real-valued coordinate to the closest integer values
213 IVec closestIntegerPoint(const RVec& coordinate)
215 return { roundToInt(coordinate[XX]), roundToInt(coordinate[YY]), roundToInt(coordinate[ZZ]) };
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.
224 IVec rangeBeginWithinLattice(const IVec& index, const IVec& range)
226 return elementWiseMax({ 0, 0, 0 }, index - range);
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
236 IVec rangeEndWithinLattice(const IVec& index, const dynamicExtents3D& extents, const IVec& range)
238 IVec extentAsIvec(static_cast<int>(extents.extent(ZZ)), static_cast<int>(extents.extent(YY)),
239 static_cast<int>(extents.extent(XX)));
240 return elementWiseMin(extentAsIvec, index + range);
244 } // namespace
246 /********************************************************************
247 * OuterProductEvaluator
250 mdspan<const float, dynamic_extent, dynamic_extent> OuterProductEvaluator::
251 operator()(ArrayRef<const float> x, ArrayRef<const float> y)
253 data_.resize(ssize(x), ssize(y));
254 for (gmx::index xIndex = 0; xIndex < ssize(x); ++xIndex)
256 const auto xValue = x[xIndex];
257 std::transform(std::begin(y), std::end(y), begin(data_.asView()[xIndex]),
258 [xValue](float yValue) { return xValue * yValue; });
260 return data_.asConstView();
263 /********************************************************************
264 * IntegerBox
267 IntegerBox::IntegerBox(const IVec& begin, const IVec& end) : begin_{ begin }, end_{ end } {}
269 const IVec& IntegerBox::begin() const
271 return begin_;
273 const IVec& IntegerBox::end() const
275 return end_;
278 bool IntegerBox::empty() const
280 return !((begin_[XX] < end_[XX]) && (begin_[YY] < end_[YY]) && (begin_[ZZ] < end_[ZZ]));
283 IntegerBox spreadRangeWithinLattice(const IVec& center, dynamicExtents3D extent, IVec range)
285 const IVec begin = rangeBeginWithinLattice(center, range);
286 const IVec end = rangeEndWithinLattice(center, extent, range);
287 return { begin, end };
289 /********************************************************************
290 * GaussianSpreadKernel
293 IVec GaussianSpreadKernelParameters::Shape::latticeSpreadRange() const
295 DVec range(std::ceil(sigma_[XX] * spreadWidthMultiplesOfSigma_),
296 std::ceil(sigma_[YY] * spreadWidthMultiplesOfSigma_),
297 std::ceil(sigma_[ZZ] * spreadWidthMultiplesOfSigma_));
298 return range.toIVec();
301 /********************************************************************
302 * GaussTransform3D::Impl
305 /*! \internal \brief
306 * Private implementation class for GaussTransform3D.
308 class GaussTransform3D::Impl
310 public:
311 //! Construct from extent and spreading width and range
312 Impl(const dynamicExtents3D& extent, const GaussianSpreadKernelParameters::Shape& kernelShapeParameters);
313 ~Impl() = default;
314 //! Copy constructor
315 Impl(const Impl& other) = default;
316 //! Copy assignment
317 Impl& operator=(const Impl& other) = default;
318 //! Add another gaussian
319 void add(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParamters);
320 //! The width of the Gaussian in lattice spacing units
321 BasicVector<double> sigma_;
322 //! The spread range in lattice points
323 IVec spreadRange_;
324 //! The result of the Gauss transform
325 MultiDimArray<std::vector<float>, dynamicExtents3D> data_;
326 //! The outer product of a Gaussian along the z and y dimension
327 OuterProductEvaluator outerProductZY_;
328 //! The three one-dimensional Gaussians, whose outer product is added to the Gauss transform
329 std::array<GaussianOn1DLattice, DIM> gauss1d_;
332 GaussTransform3D::Impl::Impl(const dynamicExtents3D& extent,
333 const GaussianSpreadKernelParameters::Shape& kernelShapeParameters) :
334 sigma_{ kernelShapeParameters.sigma_ },
335 spreadRange_{ kernelShapeParameters.latticeSpreadRange() },
336 data_{ extent },
337 gauss1d_({ GaussianOn1DLattice(spreadRange_[XX], sigma_[XX]),
338 GaussianOn1DLattice(spreadRange_[YY], sigma_[YY]),
339 GaussianOn1DLattice(spreadRange_[ZZ], sigma_[ZZ]) })
343 void GaussTransform3D::Impl::add(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters)
345 const IVec closestLatticePoint = closestIntegerPoint(localParameters.coordinate_);
346 const auto spreadRange =
347 spreadRangeWithinLattice(closestLatticePoint, data_.asView().extents(), spreadRange_);
349 // do nothing if the added Gaussian will never reach the lattice
350 if (spreadRange.empty())
352 return;
355 for (int dimension = XX; dimension <= ZZ; ++dimension)
357 // multiply with amplitude so that Gauss3D = (amplitude * Gauss_x) * Gauss_y * Gauss_z
358 const float gauss1DAmplitude = dimension > XX ? 1.0 : localParameters.amplitude_;
359 gauss1d_[dimension].spread(gauss1DAmplitude, localParameters.coordinate_[dimension]
360 - closestLatticePoint[dimension]);
363 const auto spreadZY = outerProductZY_(gauss1d_[ZZ].view(), gauss1d_[YY].view());
364 const auto spreadX = gauss1d_[XX].view();
365 const IVec spreadGridOffset = spreadRange_ - closestLatticePoint;
367 // \todo optimize these loops if performance critical
368 // The looping strategy uses that the last, x-dimension is contiguous in the memory layout
369 for (int zLatticeIndex = spreadRange.begin()[ZZ]; zLatticeIndex < spreadRange.end()[ZZ]; ++zLatticeIndex)
371 const auto zSlice = data_.asView()[zLatticeIndex];
373 for (int yLatticeIndex = spreadRange.begin()[YY]; yLatticeIndex < spreadRange.end()[YY]; ++yLatticeIndex)
375 const auto ySlice = zSlice[yLatticeIndex];
376 const float zyPrefactor = spreadZY(zLatticeIndex + spreadGridOffset[ZZ],
377 yLatticeIndex + spreadGridOffset[YY]);
379 for (int xLatticeIndex = spreadRange.begin()[XX]; xLatticeIndex < spreadRange.end()[XX];
380 ++xLatticeIndex)
382 const float xPrefactor = spreadX[xLatticeIndex + spreadGridOffset[XX]];
383 ySlice[xLatticeIndex] += zyPrefactor * xPrefactor;
389 /********************************************************************
390 * GaussTransform3D
393 GaussTransform3D::GaussTransform3D(const dynamicExtents3D& extent,
394 const GaussianSpreadKernelParameters::Shape& kernelShapeParameters) :
395 impl_(new Impl(extent, kernelShapeParameters))
399 void GaussTransform3D::add(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters)
401 impl_->add(localParameters);
404 void GaussTransform3D::setZero()
406 std::fill(begin(impl_->data_), end(impl_->data_), 0.);
409 basic_mdspan<float, dynamicExtents3D> GaussTransform3D::view()
411 return impl_->data_.asView();
414 basic_mdspan<const float, dynamicExtents3D> GaussTransform3D::constView() const
416 return impl_->data_.asConstView();
419 GaussTransform3D::~GaussTransform3D() {}
421 GaussTransform3D::GaussTransform3D(const GaussTransform3D& other) : impl_(new Impl(*other.impl_)) {}
423 GaussTransform3D& GaussTransform3D::operator=(const GaussTransform3D& other)
425 *impl_ = *other.impl_;
426 return *this;
429 GaussTransform3D::GaussTransform3D(GaussTransform3D&&) noexcept = default;
431 GaussTransform3D& GaussTransform3D::operator=(GaussTransform3D&&) noexcept = default;
433 } // namespace gmx