Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / tables / quadraticsplinetable.cpp
blob04b32c171ad3318e6bf51583eb79c32522c48763
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018, 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.
36 /*! \internal \file
37 * \brief
38 * Implements classes for quadratic spline table functions
40 * \author Erik Lindahl <erik.lindahl@gmail.com>
41 * \ingroup module_tables
43 #include "gmxpre.h"
45 #include "quadraticsplinetable.h"
47 #include <cmath>
49 #include <algorithm>
50 #include <functional>
51 #include <initializer_list>
52 #include <utility>
53 #include <vector>
55 #include "gromacs/tables/tableinput.h"
56 #include "gromacs/utility/alignedallocator.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/exceptions.h"
59 #include "gromacs/utility/real.h"
61 #include "splineutil.h"
63 namespace gmx
66 namespace
69 /*! \brief Construct the data for a single quadratic table from analytical functions
71 * \param[in] function Analytical functiojn
72 * \param[in] derivative Analytical derivative
73 * \param[in] range Upper/lower limit of region to tabulate
74 * \param[in] spacing Distance between table points
75 * \param[out] functionTableData Output table with function data
76 * \param[out] derivativeTableData OUtput table with (adjusted) derivative data
78 void
79 fillSingleQuadraticSplineTableData(const std::function<double(double)> &function,
80 const std::function<double(double)> &derivative,
81 const std::pair<real, real> &range,
82 double spacing,
83 std::vector<real> *functionTableData,
84 std::vector<real> *derivativeTableData)
86 std::size_t endIndex = static_cast<std::size_t>(range.second / spacing + 2);
88 functionTableData->resize(endIndex);
89 derivativeTableData->resize(endIndex);
91 double maxMagnitude = 0.0001*GMX_REAL_MAX;
92 bool functionIsInRange = true;
93 std::size_t lastIndexInRange = endIndex - 1;
95 for (int i = endIndex - 1; i >= 0; i--)
97 double x = i * spacing;
98 double tmpFunctionValue;
99 double tmpDerivativeValue;
101 if (range.first > 0 && i == 0)
103 // Avoid x==0 if it is not in the range, since it can lead to
104 // singularities even if the value for i==1 was within or required magnitude
105 functionIsInRange = false;
108 if (functionIsInRange)
110 tmpFunctionValue = function(x);
112 // Calculate third derivative term (2nd derivative of the derivative)
113 // Make sure we stay in range. In practice this means we use one-sided
114 // interpolation at the interval endpoints (indentical to an offset for 3-point formula)
115 const double h = std::pow( GMX_DOUBLE_EPS, 0.25 );
116 double y = std::min( std::max(x, range.first + h), range.second - h);
117 double thirdDerivativeValue = ( derivative(y+h) - 2.0 * derivative(y) + derivative(y-h) ) / ( h * h );
119 tmpDerivativeValue = derivative(x) - spacing * spacing * thirdDerivativeValue / 12.0;
121 if (std::abs(tmpFunctionValue) > maxMagnitude || std::abs(tmpDerivativeValue) > maxMagnitude)
123 functionIsInRange = false; // Once this happens, it never resets to true again
127 if (functionIsInRange)
129 (*functionTableData)[i] = tmpFunctionValue;
130 (*derivativeTableData)[i] = tmpDerivativeValue;
131 lastIndexInRange--;
133 else
135 // Once the function or derivative (more likely) has reached very large values,
136 // we simply make a linear function from the last in-range value of the derivative.
137 double lastIndexFunction = (*functionTableData)[lastIndexInRange];
138 double lastIndexDerivative = (*derivativeTableData)[lastIndexInRange];
139 (*functionTableData)[i] = lastIndexFunction + lastIndexDerivative * (i - lastIndexInRange) * spacing;
140 (*derivativeTableData)[i] = lastIndexDerivative;
146 /*! \brief Construct the data for a single quadratic table from vector data
148 * \param[in] function Input vector with function data
149 * \param[in] derivative Input vector with derivative data
150 * \param[in] inputSpacing Distance between points in input vectors
151 * \param[in] range Upper/lower limit of region to tabulate
152 * \param[in] spacing Distance between table points
153 * \param[out] functionTableData Output table with function data
154 * \param[out] derivativeTableData OUtput table with (adjusted) derivative data
156 void
157 fillSingleQuadraticSplineTableData(ArrayRef<const double> function,
158 ArrayRef<const double> derivative,
159 double inputSpacing,
160 const std::pair<real, real> &range,
161 double spacing,
162 std::vector<real> *functionTableData,
163 std::vector<real> *derivativeTableData)
165 std::size_t endIndex = static_cast<std::size_t>(range.second / spacing + 2);
167 functionTableData->resize(endIndex);
168 derivativeTableData->resize(endIndex);
170 std::vector<double> thirdDerivative(internal::vectorSecondDerivative(derivative, inputSpacing));
172 double maxMagnitude = 0.0001*GMX_REAL_MAX;
173 bool functionIsInRange = true;
174 std::size_t lastIndexInRange = endIndex - 1;
176 for (int i = endIndex - 1; i >= 0; i--)
178 double x = i * spacing;
179 double tmpFunctionValue;
180 double tmpDerivativeValue;
182 if (range.first > 0 && i == 0)
184 // Avoid x==0 if it is not in the range, since it can lead to
185 // singularities even if the value for i==1 was within or required magnitude
186 functionIsInRange = false;
189 if (functionIsInRange)
191 // Step 1: Interpolate the function value at x from input table.
192 double inputXTab = x / inputSpacing;
193 int inputIndex = static_cast<std::size_t>(inputXTab);
194 double inputEps = inputXTab - inputIndex;
196 // Linear interpolation of input derivative and third derivative
197 double thirdDerivativeValue = (1.0 - inputEps) * thirdDerivative[inputIndex] + inputEps * thirdDerivative[inputIndex+1];
198 double derivativeValue = (1.0 - inputEps) * derivative[inputIndex] + inputEps * derivative[inputIndex+1];
200 // Quadratic interpolation for function value
201 tmpFunctionValue = function[inputIndex] + 0.5 * (derivative[inputIndex] + derivativeValue) * inputEps * inputSpacing;
202 tmpDerivativeValue = derivativeValue - spacing * spacing * thirdDerivativeValue / 12.0;
204 if (std::abs(tmpFunctionValue) > maxMagnitude || std::abs(tmpDerivativeValue) > maxMagnitude)
206 functionIsInRange = false; // Once this happens, it never resets to true again
210 if (functionIsInRange)
212 (*functionTableData)[i] = tmpFunctionValue;
213 (*derivativeTableData)[i] = tmpDerivativeValue;
214 lastIndexInRange--;
216 else
218 // Once the function or derivative (more likely) has reached very large values,
219 // we simply make a linear function from the last in-range value of the derivative.
220 double lastIndexFunction = (*functionTableData)[lastIndexInRange];
221 double lastIndexDerivative = (*derivativeTableData)[lastIndexInRange];
222 (*functionTableData)[i] = lastIndexFunction + lastIndexDerivative * (i - lastIndexInRange) * spacing;
223 (*derivativeTableData)[i] = lastIndexDerivative;
228 /*! \brief Create merged DDFZ vector from function & derivative data
230 * \param functionTableData Function values
231 * \param derivativeTableData Derivative values. We have already subtracted the
232 * small third derivative component when calling this
233 * function, but in practice it is just an arbitrary
234 * vector here.
235 * \param ddfzTableData Vector four times longer, filled with
236 * the derivative, the difference to the next derivative
237 * point, the function value, and zero.
239 * \throws If the vector lengths do not match.
241 void
242 fillDdfzTableData(const std::vector<real> &functionTableData,
243 const std::vector<real> &derivativeTableData,
244 std::vector<real> *ddfzTableData)
246 GMX_ASSERT(functionTableData.size() == derivativeTableData.size(), "Mismatching vector lengths");
248 std::size_t points = functionTableData.size();
250 ddfzTableData->resize(4 * points);
252 for (std::size_t i = 0; i < points; i++)
254 (*ddfzTableData)[4*i] = derivativeTableData[i];
256 double nextDerivative = ( i < functionTableData.size() - 1 ) ? derivativeTableData[i+1] : 0.0;
258 (*ddfzTableData)[4*i + 1] = nextDerivative - derivativeTableData[i];
259 (*ddfzTableData)[4*i + 2] = functionTableData[i];
260 (*ddfzTableData)[4*i + 3] = 0.0;
264 } // namespace
268 const real
269 QuadraticSplineTable::defaultTolerance = 10.0 * GMX_FLOAT_EPS;
272 QuadraticSplineTable::QuadraticSplineTable(std::initializer_list<AnalyticalSplineTableInput> analyticalInputList,
273 const std::pair<real, real> &range,
274 real tolerance)
275 : numFuncInTable_(analyticalInputList.size()), range_(range)
277 // Sanity check on input values
278 if (range_.first < 0.0 || (range_.second-range_.first) < 0.001)
280 GMX_THROW(InvalidInputError("Range to tabulate cannot include negative values and must span at least 0.001"));
283 if (tolerance < GMX_REAL_EPS)
285 GMX_THROW(ToleranceError("Table tolerance cannot be smaller than GMX_REAL_EPS"));
288 double minQuotient = GMX_REAL_MAX;
290 // loop over all functions to find smallest spacing
291 for (const auto &thisFuncInput : analyticalInputList)
295 internal::throwUnlessDerivativeIsConsistentWithFunction(thisFuncInput.function, thisFuncInput.derivative, range_);
297 catch (gmx::GromacsException &ex)
299 ex.prependContext("Error generating quadratic spline table for function '" + thisFuncInput.desc + "'");
300 throw;
302 // Calculate the required table spacing h. The error we make with linear interpolation
303 // of the derivative will be described by the third-derivative correction term.
304 // This means we can compute the required spacing as h = sqrt(12*tolerance*min(f'/f''')),
305 // where f'/f''' is the first and third derivative of the function, respectively.
307 double thisMinQuotient = internal::findSmallestQuotientOfFunctionAndSecondDerivative(thisFuncInput.derivative, range_);
309 minQuotient = std::min(minQuotient, thisMinQuotient);
312 double spacing = std::sqrt(12.0 * tolerance * minQuotient);
314 halfSpacing_ = 0.5 * spacing;
315 tableScale_ = 1.0 / spacing;
317 if (range_.second * tableScale_ > 1e6)
319 GMX_THROW(ToleranceError("Over a million points would be required for table; decrease range or increase tolerance"));
322 // Loop over all tables again.
323 // Here we create the actual table for each function, and then
324 // combine them into a multiplexed table function.
325 std::size_t funcIndex = 0;
327 for (const auto &thisFuncInput : analyticalInputList)
331 std::vector<real> tmpFuncTableData;
332 std::vector<real> tmpDerTableData;
333 std::vector<real> tmpDdfzTableData;
335 fillSingleQuadraticSplineTableData(thisFuncInput.function,
336 thisFuncInput.derivative,
337 range_,
338 spacing,
339 &tmpFuncTableData,
340 &tmpDerTableData);
342 fillDdfzTableData(tmpFuncTableData, tmpDerTableData, &tmpDdfzTableData);
344 internal::fillMultiplexedTableData(tmpDerTableData, &derivativeMultiTableData_,
345 1, numFuncInTable_, funcIndex);
347 internal::fillMultiplexedTableData(tmpDdfzTableData, &ddfzMultiTableData_,
348 4, numFuncInTable_, funcIndex);
350 funcIndex++;
352 catch (gmx::GromacsException &ex)
354 ex.prependContext("Error generating quadratic spline table for function '" + thisFuncInput.desc + "'");
355 throw;
361 QuadraticSplineTable::QuadraticSplineTable(std::initializer_list<NumericalSplineTableInput> numericalInputList,
362 const std::pair<real, real> &range,
363 real tolerance)
364 : numFuncInTable_(numericalInputList.size()), range_(range)
366 // Sanity check on input values
367 if (range.first < 0.0 || (range.second-range.first) < 0.001)
369 GMX_THROW(InvalidInputError("Range to tabulate cannot include negative values and must span at least 0.001"));
372 if (tolerance < GMX_REAL_EPS)
374 GMX_THROW(ToleranceError("Table tolerance cannot be smaller than GMX_REAL_EPS"));
377 double minQuotient = GMX_REAL_MAX;
379 // loop over all functions to find smallest spacing
380 for (auto thisFuncInput : numericalInputList)
384 // We do not yet know what the margin is, but we need to test that we at least cover
385 // the requested range before starting to calculate derivatives
386 if (thisFuncInput.function.size() < range_.second / thisFuncInput.spacing + 1)
388 GMX_THROW(InconsistentInputError("Table input vectors must cover requested range, and a margin beyond the upper endpoint"));
391 if (thisFuncInput.function.size() != thisFuncInput.derivative.size())
393 GMX_THROW(InconsistentInputError("Function and derivative vectors have different lengths"));
396 internal::throwUnlessDerivativeIsConsistentWithFunction(thisFuncInput.function, thisFuncInput.derivative, thisFuncInput.spacing, range_);
398 catch (gmx::GromacsException &ex)
400 ex.prependContext("Error generating quadratic spline table for function '" + thisFuncInput.desc + "'");
401 throw;
403 // Calculate the required table spacing h. The error we make with linear interpolation
404 // of the derivative will be described by the third-derivative correction term.
405 // This means we can compute the required spacing as h = sqrt(12*tolerance*min(f'/f''')),
406 // where f'/f''' is the first and third derivative of the function, respectively.
407 // Since we already have an analytical form of the derivative, we reduce the numerical
408 // errors by calculating the quotient of the function and second derivative of the
409 // input-derivative-analytical function instead.
411 double thisMinQuotient = internal::findSmallestQuotientOfFunctionAndSecondDerivative(thisFuncInput.derivative, thisFuncInput.spacing, range_);
413 minQuotient = std::min(minQuotient, thisMinQuotient);
416 double spacing = std::sqrt(12.0 * tolerance * minQuotient);
418 halfSpacing_ = 0.5 * spacing;
419 tableScale_ = 1.0 / spacing;
421 if (range_.second * tableScale_ > 1e6)
423 GMX_THROW(ToleranceError("Requested tolerance would require over a million points in table"));
426 // Loop over all tables again.
427 // Here we create the actual table for each function, and then
428 // combine them into a multiplexed table function.
429 std::size_t funcIndex = 0;
431 for (auto thisFuncInput : numericalInputList)
435 if (spacing < thisFuncInput.spacing)
437 GMX_THROW(ToleranceError("Input vector spacing cannot achieve tolerance requested"));
440 std::vector<real> tmpFuncTableData;
441 std::vector<real> tmpDerTableData;
442 std::vector<real> tmpDdfzTableData;
444 fillSingleQuadraticSplineTableData(thisFuncInput.function,
445 thisFuncInput.derivative,
446 thisFuncInput.spacing,
447 range,
448 spacing,
449 &tmpFuncTableData,
450 &tmpDerTableData);
452 fillDdfzTableData(tmpFuncTableData, tmpDerTableData, &tmpDdfzTableData);
454 internal::fillMultiplexedTableData(tmpDerTableData, &derivativeMultiTableData_,
455 1, numFuncInTable_, funcIndex);
457 internal::fillMultiplexedTableData(tmpDdfzTableData, &ddfzMultiTableData_,
458 4, numFuncInTable_, funcIndex);
460 funcIndex++;
462 catch (gmx::GromacsException &ex)
464 ex.prependContext("Error generating quadratic spline table for function '" + thisFuncInput.desc + "'");
465 throw;
470 } // namespace gmx