Add conserved quantity for Berendsen P-couple
[gromacs.git] / src / gromacs / tables / cubicsplinetable.cpp
blobf0b589ea8bb1ab7a209845aec8b41ba64bcebf9c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016, 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 cubic spline table functions
40 * \author Erik Lindahl <erik.lindahl@gmail.com>
41 * \ingroup module_tables
43 #include "gmxpre.h"
45 #include "cubicsplinetable.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 Calculate table elements from function/derivative data
71 * \param functionValue0 Function value for the present table index
72 * \param functionValue1 Function value for the next table index
73 * \param derivativeValue0 Derivative value for the present table index
74 * \param derivativeValue1 Derivative value for the next table index
75 * \param spacing Distance between table points
76 * \param Y Function value for table index
77 * \param F Component to multiply with offset eps
78 * \param G Component to multiply with eps^2
79 * \param H Component to multiply with eps^3
81 void
82 calculateCubicSplineCoefficients(double functionValue0,
83 double functionValue1,
84 double derivativeValue0,
85 double derivativeValue1,
86 double spacing,
87 double *Y,
88 double *F,
89 double *G,
90 double *H)
92 *Y = functionValue0;
93 *F = spacing * derivativeValue0;
94 *G = 3.0*( functionValue1 - functionValue0) - spacing * (derivativeValue1 + 2.0 * derivativeValue0);
95 *H = -2.0*( functionValue1 - functionValue0) + spacing * (derivativeValue1 + derivativeValue0);
98 /*! \brief Perform cubic spline interpolation in interval from function/derivative
100 * \param functionValue0 Function value for the present table index
101 * \param functionValue1 Function value for the next table index
102 * \param derivativeValue0 Derivative value for the present table index
103 * \param derivativeValue1 Derivative value for the next table index
104 * \param spacing Distance between table points
105 * \param eps Offset from lower table point for evaluation
106 * \param[out] interpolatedFunctionValue Output function value
107 * \param[out] interpolatedDerivativeValue Output derivative value
109 void
110 cubicSplineInterpolationFromFunctionAndDerivative(double functionValue0,
111 double functionValue1,
112 double derivativeValue0,
113 double derivativeValue1,
114 double spacing,
115 double eps,
116 double *interpolatedFunctionValue,
117 double *interpolatedDerivativeValue)
119 double Y, F, G, H;
121 calculateCubicSplineCoefficients(functionValue0, functionValue1,
122 derivativeValue0, derivativeValue1,
123 spacing,
124 &Y, &F, &G, &H);
126 double Fp = fma(fma(H, eps, G), eps, F);
128 *interpolatedFunctionValue = fma(Fp, eps, Y);
129 *interpolatedDerivativeValue = fma(eps, fma(2.0*eps, H, G), Fp)/spacing;
134 /*! \brief Construct the data for a single cubic table from analytical functions
136 * \param[in] function Analytical functiojn
137 * \param[in] derivative Analytical derivative
138 * \param[in] range Upper/lower limit of region to tabulate
139 * \param[in] spacing Distance between table points
140 * \param[out] yfghTableData Output cubic spline table with Y,F,G,H entries
142 void
143 fillSingleCubicSplineTableData(const std::function<double(double)> &function,
144 const std::function<double(double)> &derivative,
145 const std::pair<real, real> &range,
146 double spacing,
147 std::vector<real> *yfghTableData)
149 int endIndex = range.second / spacing + 2;
151 yfghTableData->resize(4*endIndex);
153 double maxMagnitude = 0.0001*GMX_REAL_MAX;
154 bool functionIsInRange = true;
155 std::size_t lastIndexInRange = endIndex - 1;
157 for (int i = endIndex - 1; i >= 0; i--)
159 double x = i * spacing;
160 double tmpFunctionValue;
161 double tmpDerivativeValue;
162 double nextHigherFunction;
163 double nextHigherDerivative;
164 double Y, F, G, H;
166 if (range.first > 0 && i == 0)
168 // Avoid x==0 if it is not in the range, since it can lead to
169 // singularities even if the value for i==1 was within or required magnitude
170 functionIsInRange = false;
173 if (functionIsInRange)
175 tmpFunctionValue = function(x);
176 tmpDerivativeValue = derivative(x);
177 nextHigherFunction = ((i+1) < endIndex) ? function(x+spacing) : 0.0;
178 nextHigherDerivative = ((i+1) < endIndex) ? derivative(x+spacing) : 0.0;
180 if (std::abs(tmpFunctionValue) > maxMagnitude || std::abs(tmpDerivativeValue) > maxMagnitude)
182 functionIsInRange = false; // Once this happens, it never resets to true again
186 if (functionIsInRange)
188 calculateCubicSplineCoefficients(tmpFunctionValue, nextHigherFunction,
189 tmpDerivativeValue, nextHigherDerivative,
190 spacing,
191 &Y, &F, &G, &H);
192 lastIndexInRange--;
194 else
196 double lastIndexY = (*yfghTableData)[4*lastIndexInRange];
197 double lastIndexF = (*yfghTableData)[4*lastIndexInRange + 1];
199 Y = lastIndexY + lastIndexF * (i - lastIndexInRange);
200 F = lastIndexF;
201 G = 0.0;
202 H = 0.0;
205 (*yfghTableData)[4*i ] = Y;
206 (*yfghTableData)[4*i+1] = F;
207 (*yfghTableData)[4*i+2] = G;
208 (*yfghTableData)[4*i+3] = H;
213 /*! \brief Construct the data for a single cubic table from vector data
215 * \param[in] function Input vector with function data
216 * \param[in] derivative Input vector with derivative data
217 * \param[in] inputSpacing Distance between points in input vectors
218 * \param[in] range Upper/lower limit of region to tabulate
219 * \param[in] spacing Distance between table points
220 * \param[out] yfghTableData Output cubic spline table with Y,F,G,H entries
222 void
223 fillSingleCubicSplineTableData(ConstArrayRef<double> function,
224 ConstArrayRef<double> derivative,
225 double inputSpacing,
226 const std::pair<real, real> &range,
227 double spacing,
228 std::vector<real> *yfghTableData)
230 int endIndex = range.second / spacing + 2;
232 std::vector<double> tmpFunction(endIndex);
233 std::vector<double> tmpDerivative(endIndex);
235 double maxMagnitude = 0.0001*GMX_REAL_MAX;
236 bool functionIsInRange = true;
237 std::size_t lastIndexInRange = endIndex - 1;
239 // Interpolate function and derivative values in positions needed for output
240 for (int i = endIndex - 1; i >= 0; i--)
242 double x = i * spacing;
243 double xtab = x / inputSpacing;
244 int index = xtab;
245 double eps = xtab - index;
247 if (range.first > 0 && i == 0)
249 // Avoid x==0 if it is not in the range, since it can lead to
250 // singularities even if the value for i==1 was within or required magnitude
251 functionIsInRange = false;
254 if (functionIsInRange && (std::abs(function[index]) > maxMagnitude || std::abs(derivative[index]) > maxMagnitude))
256 functionIsInRange = false; // Once this happens, it never resets to true again
259 if (functionIsInRange)
261 cubicSplineInterpolationFromFunctionAndDerivative(function[index],
262 function[index+1],
263 derivative[index],
264 derivative[index+1],
265 inputSpacing,
266 eps,
267 &(tmpFunction[i]),
268 &(tmpDerivative[i]));
269 lastIndexInRange--;
271 else
273 double lastIndexFunction = tmpFunction[lastIndexInRange];
274 double lastIndexDerivative = tmpDerivative[lastIndexInRange];
275 tmpFunction[i] = lastIndexFunction + lastIndexDerivative * (i - lastIndexInRange) * spacing;
276 tmpDerivative[i] = lastIndexDerivative;
280 yfghTableData->resize(4*endIndex);
282 for (int i = 0; i < endIndex; i++)
284 double Y, F, G, H;
286 double nextFunction = ((i+1) < endIndex) ? tmpFunction[i+1] : 0.0;
287 double nextDerivative = ((i+1) < endIndex) ? tmpDerivative[i+1] : 0.0;
289 calculateCubicSplineCoefficients(tmpFunction[i], nextFunction,
290 tmpDerivative[i], nextDerivative,
291 spacing,
292 &Y, &F, &G, &H);
293 (*yfghTableData)[4*i ] = Y;
294 (*yfghTableData)[4*i+1] = F;
295 (*yfghTableData)[4*i+2] = G;
296 (*yfghTableData)[4*i+3] = H;
301 } // namespace anonymous
305 #if GMX_DOUBLE
306 const real
307 CubicSplineTable::defaultTolerance = 1e-10;
308 #else
309 const real
310 CubicSplineTable::defaultTolerance = 10.0 * GMX_FLOAT_EPS;
311 #endif
313 CubicSplineTable::CubicSplineTable(std::initializer_list<AnalyticalSplineTableInput> analyticalInputList,
314 const std::pair<real, real> &range,
315 real tolerance)
316 : numFuncInTable_(analyticalInputList.size()), range_(range)
318 // Sanity check on input values
319 if (range_.first < 0.0 || (range_.second-range_.first) < 0.001)
321 GMX_THROW(InvalidInputError("Range to tabulate cannot include negative values and must span at least 0.001"));
324 if (tolerance < GMX_REAL_EPS)
326 GMX_THROW(ToleranceError("Table tolerance cannot be smaller than GMX_REAL_EPS"));
329 double minQuotient = GMX_REAL_MAX;
331 // loop over all functions to find smallest spacing
332 for (auto thisFuncInput : analyticalInputList)
336 internal::throwUnlessDerivativeIsConsistentWithFunction(thisFuncInput.function, thisFuncInput.derivative, range_);
338 catch (gmx::GromacsException &ex)
340 ex.prependContext("Error generating cubic spline table for function '" + thisFuncInput.desc + "'");
341 throw;
343 // Calculate the required table spacing h. The error we make with linear interpolation
344 // of the derivative will be described by the third-derivative correction term.
345 // This means we can compute the required spacing as h = sqrt(12*tolerance*min(f'/f''')),
346 // where f'/f''' is the first and third derivative of the function, respectively.
348 double thisMinQuotient = internal::findSmallestQuotientOfFunctionAndThirdDerivative(thisFuncInput.derivative, range_);
350 minQuotient = std::min(minQuotient, thisMinQuotient);
353 double spacing = 0.5 * std::cbrt(72.0 * std::sqrt(3.0) * tolerance * minQuotient);
355 tableScale_ = 1.0 / spacing;
357 if (range_.second * tableScale_ > 2e6)
359 GMX_THROW(ToleranceError("Over a million points would be required for table; decrease range or increase tolerance"));
362 // Loop over all tables again.
363 // Here we create the actual table for each function, and then
364 // combine them into a multiplexed table function.
365 std::size_t funcIndex = 0;
367 for (auto thisFuncInput : analyticalInputList)
371 std::vector<real> tmpYfghTableData;
373 fillSingleCubicSplineTableData(thisFuncInput.function,
374 thisFuncInput.derivative,
375 range_,
376 spacing,
377 &tmpYfghTableData);
379 internal::fillMultiplexedTableData(tmpYfghTableData, &yfghMultiTableData_,
380 4, numFuncInTable_, funcIndex);
382 funcIndex++;
384 catch (gmx::GromacsException &ex)
386 ex.prependContext("Error generating cubic spline table for function '" + thisFuncInput.desc + "'");
387 throw;
393 CubicSplineTable::CubicSplineTable(std::initializer_list<NumericalSplineTableInput> numericalInputList,
394 const std::pair<real, real> &range,
395 real tolerance)
396 : numFuncInTable_(numericalInputList.size()), range_(range)
398 // Sanity check on input values
399 if (range.first < 0.0 || (range.second-range.first) < 0.001)
401 GMX_THROW(InvalidInputError("Range to tabulate cannot include negative values and must span at least 0.001"));
404 if (tolerance < GMX_REAL_EPS)
406 GMX_THROW(ToleranceError("Table tolerance cannot be smaller than GMX_REAL_EPS"));
409 double minQuotient = GMX_REAL_MAX;
411 // loop over all functions to find smallest spacing
412 for (auto thisFuncInput : numericalInputList)
416 // We do not yet know what the margin is, but we need to test that we at least cover
417 // the requested range before starting to calculate derivatives
418 if (thisFuncInput.function.size() < range_.second / thisFuncInput.spacing + 1)
420 GMX_THROW(InconsistentInputError("Table input vectors must cover requested range, and a margin beyond the upper endpoint"));
423 if (thisFuncInput.function.size() != thisFuncInput.derivative.size())
425 GMX_THROW(InconsistentInputError("Function and derivative vectors have different lengths"));
428 internal::throwUnlessDerivativeIsConsistentWithFunction(thisFuncInput.function, thisFuncInput.derivative, thisFuncInput.spacing, range_);
430 catch (gmx::GromacsException &ex)
432 ex.prependContext("Error generating cubic spline table for function '" + thisFuncInput.desc + "'");
433 throw;
435 // Calculate the required table spacing h. The error we make with linear interpolation
436 // of the derivative will be described by the third-derivative correction term.
437 // This means we can compute the required spacing as h = sqrt(12*tolerance*min(f'/f''')),
438 // where f'/f''' is the first and third derivative of the function, respectively.
440 double thisMinQuotient = internal::findSmallestQuotientOfFunctionAndThirdDerivative(thisFuncInput.derivative, thisFuncInput.spacing, range_);
442 minQuotient = std::min(minQuotient, thisMinQuotient);
445 double spacing = std::cbrt(72.0 * std::sqrt(3.0) * tolerance * minQuotient);
447 tableScale_ = 1.0 / spacing;
449 if (range_.second * tableScale_ > 1e6)
451 GMX_THROW(ToleranceError("Requested tolerance would require over a million points in table"));
454 // Loop over all tables again.
455 // Here we create the actual table for each function, and then
456 // combine them into a multiplexed table function.
457 std::size_t funcIndex = 0;
459 for (auto thisFuncInput : numericalInputList)
463 if (spacing < thisFuncInput.spacing)
465 GMX_THROW(ToleranceError("Input vector spacing cannot achieve tolerance requested"));
468 std::vector<real> tmpYfghTableData;
470 fillSingleCubicSplineTableData(thisFuncInput.function,
471 thisFuncInput.derivative,
472 thisFuncInput.spacing,
473 range,
474 spacing,
475 &tmpYfghTableData);
477 internal::fillMultiplexedTableData(tmpYfghTableData, &yfghMultiTableData_,
478 4, numFuncInTable_, funcIndex);
480 funcIndex++;
482 catch (gmx::GromacsException &ex)
484 ex.prependContext("Error generating cubic spline table for function '" + thisFuncInput.desc + "'");
485 throw;
490 } // namespace gmx