5c742161902d028b5d13166001cd29d249971517
[alexxy/gromacs.git] / src / gromacs / tables / quadraticsplinetable.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,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.
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
36 /*! \libinternal
37  * \defgroup module_tables  Classes for table interpolation
38  * \ingroup group_utilitymodules
39  *
40  * \brief Table interpolation from analytical or numerical input
41  *
42  * This module provides quadratic spline interpolation tables used
43  * both for the nonbonded kernels and listed interactions.
44  *
45  * \author Erik Lindahl <erik.lindahl@scilifelab.se>
46  */
47
48
49 /*! \libinternal \file
50  * \brief
51  * Declares classes for quadratic spline table
52  *
53  * \inlibraryapi
54  * \author Erik Lindahl <erik.lindahl@gmail.com>
55  * \ingroup module_tables
56  */
57 #ifndef GMX_TABLES_QUADRATICSPLINETABLE_H
58 #define GMX_TABLES_QUADRATICSPLINETABLE_H
59
60 #include <functional>
61 #include <initializer_list>
62 #include <vector>
63
64 #include "gromacs/simd/simd.h"
65 #include "gromacs/tables/tableinput.h"
66 #include "gromacs/utility/alignedallocator.h"
67 #include "gromacs/utility/classhelpers.h"
68 #include "gromacs/utility/exceptions.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/real.h"
71
72 namespace gmx
73 {
74
75 /*! \libinternal \brief Quadratic spline interpolation table.
76  *
77  * This class interpolates a function specified either as an analytical
78  * expression or from user-provided table data.
79  *
80  * At initialization, you provide the reference function of vectors
81  * as a list of tuples that contain a brief name, the function, and
82  * derivative for each function to tabulate. To create a table with
83  * two functions this initializer list can for instance look like
84  *
85  *     { {"LJ6", lj6Func, lj6Der}, {"LJ12", lj12Func, lj12Der} }
86  *
87  * The names are only used so exceptions during initialization can
88  * be traced to a specific table.
89  *
90  * When interpolating, there are methods to interpolate either 1, 2, or 3
91  * functions in one go. By default these interpolation routines will
92  * operate on tables with the same number of functions as specified in
93  * the interpolation method (debug builds check that this is consistent with
94  * the table). However, it is also possible to use optional template
95  * parameters that specify the total number of functions in a table, and
96  * what function index to interpolate. For instance, to interpolate the
97  * derivative of the second function (i.e., index 1) in a
98  * multi-function-table with three functions in total, you can write
99  *
100  *     table.evaluateDerivative<3,1>(x,&der);
101  *
102  * Here too, debug builds will check that the template parameters are
103  * consistent with the table.
104  *
105  * The table data is internally adjusted to guarantee that the interpolated
106  * derivative is the true derivative of the interpolated potential, which is
107  * important to avoid systematic errors for the common case when the derivative
108  * is concave/convex in the entire interval.
109  * We do this by expressing the difference in the function value
110  * at a small offset h relative to a reference value in position 0 with a forward
111  * Taylor series expanded around 0, and then doing the opposite of expressing
112  * difference in the function at position 0 relative to a reference value in
113  * position h when using a backward Taylor expansion:
114  *
115  * \f{eqnarray*}{
116  *  \Delta V & = & hV'(0) + \frac{1}{2} h^2 V''(0) + \frac{1}{6} h^3 V'''(0) + O(h^4) \\
117  *  \Delta V & = & hV'(h) - \frac{1}{2} h^2 V''(h) + \frac{1}{6} h^3 V'''(h) + O(h^4)
118  * \f}
119  *
120  * Summing the equations leads to
121  *
122  * \f[
123  *  2 \Delta V = h(V'(0) + V'(h)) + \frac{1}{2} h^2 (V''(0)-V''(h)) +
124  * \frac{1}{6}h^3(V'''(0)+V'''(h)) + O(h^4) \f]
125  *
126  * To make the second term symmetric too, we can replace it with the average of
127  * the Taylor expansion at 0 and h (i.e., using the third derivative). This gives
128  *
129  * \f[
130  *  2 \Delta V = h(V'(0) + V'(h)) - \frac{1}{12} h^3 (V'''(0)+V'''(h)) + O(h^4)
131  * \f]
132  *
133  * Thus, if we replace the derivative in the internal quadratic table data with
134  *
135  * \f[
136  *  V' - \frac{1}{12}h^2 V'''
137  * \f]
138  *
139  * we will cancel the h^3 term in the error. This will make the integral of the
140  * forces match the potential much better (The h^4 term actually disappears, so
141  * when summing over 1/h points the remaining error will be O(h^4).
142  *
143  * While it is possible to create tables only from function values
144  * (i.e., no derivatives), it is recommended to provide derivatives for higher
145  * accuracy and to avoid issues with numerical differentiation. Note that the
146  * table input should be smooth, i.e. it should not contain noise e.g. from an
147  * (iterative) Boltzmann inversion procedure - you have been warned.
148  *
149  * \note This class is responsible for fundamental interpolation of any function,
150  *       which might or might not correspond to a potential. For this reason
151  *       both input and output derivatives are proper function derivatives, and
152  *       we do not swap any signs to get forces directly from the table.
153  *
154  * \note There will be a small additional accuracy loss from the internal
155  *       operation where we calculate the epsilon offset from the nearest table
156  *       point, since the integer part we subtract can get large in those cases.
157  *       The absolute such error both in the function and derivative value will
158  *       be roughly f''*x*GMX_REAL_EPS, where x is the argument and f'' the
159  *       second derivative.
160  *       While this is technically possible to solve with extended precision
161  *       arithmetics, that would introduce extra instructions in some highly
162  *       performance-sensitive code parts. For typical GROMACS interaction
163  *       functions the derivatives will decay faster than the potential, which
164  *       means it will never play any role. For other functions it will only
165  *       cause a small increase in the relative error for arguments where the
166  *       magnitude of the function or derivative is very small.
167  *       Since we typically sum several results in GROMACS, this should never
168  *       show up in any real cases, and for this reason we choose not to do
169  *       the extended precision arithmetics.
170  *
171  * \note These routines are not suitable for table ranges starting far away
172  *       from zero, since we allocate memory and calculate indices starting from
173  *       range zero for efficiency reasons.
174  */
175 class QuadraticSplineTable
176 {
177 private:
178     /*! \brief Change that function value falls inside range when debugging
179      *
180      *  \tparam T   Lookup argument floating-point type, typically SimdReal or real.
181      *  \param  r   Lookup argument to test
182      *
183      *  \throws Debug builds will throw gmx::RangeError for values that are
184      *          larger than the upper limit of the range, or smaller than 0.
185      *          We allow the table to be called with arguments between 0 and
186      *          the lower limit of the range, since this might in theory occur
187      *          once-in-a-blue-moon with some algorithms.
188      */
189     template<typename T>
190     void rangeCheck(T gmx_unused r) const
191     {
192 #ifndef NDEBUG
193         // Check that all values fall in range when debugging
194         if (anyTrue(r < T(0.0) || T(range_.second) <= r))
195         {
196             GMX_THROW(RangeError("Interpolation input value falls outside table definition range"));
197         }
198 #endif
199     }
200
201 public:
202     /*! \brief Default tolerance for tables is 10*GMX_FLOAT_EPS
203      *
204      *  \note Even for double precision builds we set the tolerance to
205      *        one order of magnitude above the single precision epsilon.
206      */
207     static const real defaultTolerance;
208
209     /*! \brief Initialize table data from function
210      *
211      * \param analyticalInputList Initializer list with one or more functions to tabulate,
212      *                            specified as pairs containing analytical
213      *                            functions and their derivatives. The function will also
214      *                            be called for values smaller than the lower limit of the
215      *                            range, but we avoid calling it for 0.0 if that value
216      *                            is not included in the range.
217      * \param range               Range over which the function will be tabulated.
218      *                            Constructor will throw gmx::APIError for negative values.
219      *                            Due to the way the numerical derivative evaluation depends
220      *                            on machine precision internally, this range must be
221      *                            at least 0.001, or the constructor throws gmx::APIError.
222      * \param tolerance           Requested accuracy of the table. This will be used to
223      *                            calculate the required internal spacing. If this cannot
224      *                            be achieved (for instance because the table would require
225      *                            too much memory) the constructor will throw gmx::ToleranceError.
226      *
227      * \note The functions are always defined in double precision to avoid
228      *       losing accuracy when constructing tables.
229      *
230      * \note Since we fill the table for values below range.first, you can achieve
231      *       a smaller table by using a smaller range where the tolerance has to be
232      *       met, and accept that a few function calls below range.first do not
233      *       quite reach the tolerance.
234      *
235      * \warning For efficiency reasons (since this code is used in some inner
236      *       (kernels), we always allocate memory and calculate table indices
237      *       for the complete interval [0,range.second], although the data will
238      *       not be valid outside the definition range to avoid calling the
239      *       function there. This means you should \a not use this class
240      *       to tabulate functions for small ranges very far away from zero,
241      *       since you would both waste a huge amount of memory and incur
242      *       truncation errors when calculating the index.
243      *
244      * \throws gmx::ToleranceError if the requested tolerance cannot be achieved,
245      *         and gmx::APIError for other incorrect input.
246      */
247     QuadraticSplineTable(std::initializer_list<AnalyticalSplineTableInput> analyticalInputList,
248                          const std::pair<real, real>&                      range,
249                          real tolerance = defaultTolerance);
250
251     /*! \brief Initialize table data from tabulated values and derivatives
252      *
253      * \param numericalInputList  Initializer list with one or more functions to tabulate,
254      *                            specified as pairs containing containing vectors for the
255      *                            function values and their derivatives. Data points are
256      *                            separated by the spacing parameter, starting from 0.
257      *                            Values below the lower limit of the range will be used to
258      *                            attempt defining the table, but we avoid using index 0
259      *                            unless 0.0 is included in the range. Some extra points beyond
260      *                            range.second are required to re-interpolate values, so add
261      *                            some margin. The constructor will throw gmx::APIError if the
262      *                            input vectors are too short to cover the requested range
263      *                            (and they must always be at least five points).
264      * \param range               Range over which the function will be tabulated.
265      *                            Constructor will throw gmx::APIError for negative values,
266      *                            or if the value/derivative vector does not cover the
267      *                            range.
268      * \param tolerance           Requested accuracy of the table in the range. This will be
269      *                            used to calculate the required internal spacing and possibly
270      *                            re-interpolate. The constructor will throw
271      *                            gmx::ToleranceError if the input spacing is too coarse
272      *                            to achieve this accuracy.
273      *
274      * \note The input data vectors are always double precision to avoid
275      *       losing accuracy when constructing tables.
276      *
277      * \note Since we fill the table for values below range.first, you can achieve
278      *       a smaller table by using a smaller range where the tolerance has to be
279      *       met, and accept that a few function calls below range.first do not
280      *       quite reach the tolerance.
281      *
282      * \warning For efficiency reasons (since this code is used in some inner
283      *       (kernels), we always allocate memory and calculate table indices
284      *       for the complete interval [0,range.second], although the data will
285      *       not be valid outside the definition range to avoid calling the
286      *       function there. This means you should \a not use this class
287      *       to tabulate functions for small ranges very far away from zero,
288      *       since you would both waste a huge amount of memory and incur
289      *       truncation errors when calculating the index.
290      */
291     QuadraticSplineTable(std::initializer_list<NumericalSplineTableInput> numericalInputList,
292                          const std::pair<real, real>&                     range,
293                          real tolerance = defaultTolerance);
294
295
296     /************************************************************
297      *           Evaluation methods for single functions        *
298      ************************************************************/
299
300     /*! \brief Evaluate both function and derivative, single table function
301      *
302      *  This is a templated method where the template can be either real or SimdReal.
303      *
304      *  \tparam     numFuncInTable  Number of separate functions in table, default is 1
305      *  \tparam     funcIndex       Index of function to evaluate in table, default is 0
306      *  \tparam     T               Type (SimdReal or real) of lookup and result
307      *  \param      r               Points for which to evaluate function and derivative
308      *  \param[out] functionValue   Function value
309      *  \param[out] derivativeValue Function derivative
310      *
311      *  For debug builds we assert that the input values fall in the range
312      *  specified when constructing the table.
313      */
314     template<int numFuncInTable = 1, int funcIndex = 0, typename T>
315     void evaluateFunctionAndDerivative(T r, T* functionValue, T* derivativeValue) const
316     {
317         rangeCheck(r);
318         GMX_ASSERT(numFuncInTable == numFuncInTable_,
319                    "Evaluation method not matching number of functions in table");
320         GMX_ASSERT(funcIndex < numFuncInTable,
321                    "Function index not in range of the number of tables");
322
323         T    rTable   = r * T(tableScale_);
324         auto tabIndex = cvttR2I(rTable); // type is either std::int32_t or SimdInt32
325         T    eps      = rTable - trunc(rTable);
326         T    t0;
327         T    t1;
328         T    t2;
329         T t3 gmx_unused;
330
331         // Load Derivative, Delta, Function, and Zero values for each table point.
332         // The 4 refers to these four values - not any SIMD width.
333         gatherLoadBySimdIntTranspose<4 * numFuncInTable>(ddfzMultiTableData_.data() + 4 * funcIndex,
334                                                          tabIndex, &t0, &t1, &t2, &t3);
335
336         t1               = t0 + eps * t1;
337         *functionValue   = fma(eps * T(halfSpacing_), t0 + t1, t2);
338         *derivativeValue = t1;
339     }
340
341     /*! \brief Evaluate function value only, single table function
342      *
343      *  This is a templated method where the template can be either real or SimdReal.
344      *
345      *  \tparam     numFuncInTable  Number of separate functions in table, default is 1
346      *  \tparam     funcIndex       Index of function to evaluate in table, default is 0
347      *  \tparam     T               Type (SimdReal or real) of lookup and result
348      *  \param      r               Points for which to evaluate function value
349      *  \param[out] functionValue   Function value
350      *
351      *  For debug builds we assert that the input values fall in the range
352      *  specified when constructing the table.
353      */
354     template<int numFuncInTable = 1, int funcIndex = 0, typename T>
355     void evaluateFunction(T r, T* functionValue) const
356     {
357         T der gmx_unused;
358
359         evaluateFunctionAndDerivative<numFuncInTable, funcIndex>(r, functionValue, &der);
360     }
361
362     /*! \brief Evaluate function derivative only, single table function
363      *
364      *  This is a templated method where the template can be either real or SimdReal.
365      *
366      *  \tparam     numFuncInTable  Number of separate functions in table, default is 1
367      *  \tparam     funcIndex       Index of function to evaluate in table, default is 0
368      *  \tparam     T               Type (SimdReal or real) of lookup and result
369      *  \param      r               Points for which to evaluate function derivative
370      *  \param[out] derivativeValue Function derivative
371      *
372      *  For debug builds we assert that the input values fall in the range
373      *  specified when constructing the table.
374      */
375     template<int numFuncInTable = 1, int funcIndex = 0, typename T>
376     void evaluateDerivative(T r, T* derivativeValue) const
377     {
378         rangeCheck(r);
379         GMX_ASSERT(numFuncInTable == numFuncInTable_,
380                    "Evaluation method not matching number of functions in table");
381         GMX_ASSERT(funcIndex < numFuncInTable,
382                    "Function index not in range of the number of tables");
383
384         T    rTable   = r * T(tableScale_);
385         auto tabIndex = cvttR2I(rTable); // type is either std::int32_t or SimdInt32
386         T    eps      = rTable - trunc(rTable);
387         T    t0;
388         T    t1;
389         T t2 gmx_unused;
390
391         if (numFuncInTable == 1)
392         {
393             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex,
394                                                           tabIndex, &t0, &t1); // works for scalar T too
395         }
396         else
397         {
398             // This is not ideal, but we need a version of gatherLoadUBySimdIntTranspose that
399             // only loads a single value from memory to implement it better (will be written)
400             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex,
401                                                           tabIndex, &t0, &t2); // works for scalar T too
402             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex,
403                                                           tabIndex + T(1), &t1,
404                                                           &t2); // works for scalar T too
405         }
406
407         // (1-eps)*t0 + eps*t1
408         *derivativeValue = fma(t1 - t0, eps, t0);
409     }
410
411     /************************************************************
412      *             Evaluation methods for two functions         *
413      ************************************************************/
414
415     /*! \brief Evaluate both function and derivative, two table functions
416      *
417      *  This is a templated method where the template can be either real or SimdReal.
418      *
419      *  \tparam     numFuncInTable  Number of separate functions in table, default is 2
420      *  \tparam     funcIndex0       Index of 1st function to evaluate in table, default is 0
421      *  \tparam     funcIndex1       Index of 2nd function to evaluate in table, default is 1
422      *  \tparam     T                Type (SimdReal or real) of lookup and result
423      *  \param      r                Points for which to evaluate function and derivative
424      *  \param[out] functionValue1   Interpolated value for first function
425      *  \param[out] derivativeValue1 Interpolated derivative for first function
426      *  \param[out] functionValue2   Interpolated value for second function
427      *  \param[out] derivativeValue2 Interpolated derivative for second function
428      *
429      *  For debug builds we assert that the input values fall in the range
430      *  specified when constructing the table.
431      */
432     template<int numFuncInTable = 2, int funcIndex0 = 0, int funcIndex1 = 1, typename T>
433     void evaluateFunctionAndDerivative(T r, T* functionValue1, T* derivativeValue1, T* functionValue2, T* derivativeValue2) const
434     {
435         rangeCheck(r);
436         GMX_ASSERT(numFuncInTable == numFuncInTable_,
437                    "Evaluation method not matching number of functions in table");
438         GMX_ASSERT(funcIndex0 < numFuncInTable && funcIndex1 < numFuncInTable,
439                    "Function index not in range of the number of tables");
440
441         T    rTable   = r * T(tableScale_);
442         auto tabIndex = cvttR2I(rTable); // type is either std::int32_t or SimdInt32
443         T    eps      = rTable - trunc(rTable);
444         T    t0;
445         T    t1;
446         T    t2;
447         T t3 gmx_unused;
448
449         // Load Derivative, Delta, Function, and Zero values for each table point.
450         // The 4 refers to these four values - not any SIMD width.
451         gatherLoadBySimdIntTranspose<4 * numFuncInTable>(
452                 ddfzMultiTableData_.data() + 4 * funcIndex0, tabIndex, &t0, &t1, &t2, &t3);
453         t1                = t0 + eps * t1;
454         *functionValue1   = fma(eps * T(halfSpacing_), t0 + t1, t2);
455         *derivativeValue1 = t1;
456
457         gatherLoadBySimdIntTranspose<4 * numFuncInTable>(
458                 ddfzMultiTableData_.data() + 4 * funcIndex1, tabIndex, &t0, &t1, &t2, &t3);
459         t1                = t0 + eps * t1;
460         *functionValue2   = fma(eps * T(halfSpacing_), t0 + t1, t2);
461         *derivativeValue2 = t1;
462     }
463
464     /*! \brief Evaluate function value only, two table functions
465      *
466      *  This is a templated method where the template can be either real or SimdReal.
467      *
468      *  \tparam     numFuncInTable   Number of separate functions in table, default is 2
469      *  \tparam     funcIndex0       Index of 1st function to evaluate in table, default is 0
470      *  \tparam     funcIndex1       Index of 2nd function to evaluate in table, default is 1
471      *  \tparam     T                Type (SimdReal or real) of lookup and result
472      *  \param      r                Points for which to evaluate function value
473      *  \param[out] functionValue1   Interpolated value for first function
474      *  \param[out] functionValue2   Interpolated value for second function
475      *
476      *  For debug builds we assert that the input values fall in the range
477      *  specified when constructing the table.
478      */
479     template<int numFuncInTable = 2, int funcIndex0 = 0, int funcIndex1 = 1, typename T>
480     void evaluateFunction(T r, T* functionValue1, T* functionValue2) const
481     {
482         T der1 gmx_unused;
483         T der2 gmx_unused;
484
485         evaluateFunctionAndDerivative<numFuncInTable, funcIndex0, funcIndex1>(
486                 r, functionValue1, &der1, functionValue2, &der2);
487     }
488
489     /*! \brief Evaluate function derivative only, two table functions
490      *
491      *  This is a templated method where the template can be either real or SimdReal.
492      *
493      *  \tparam     numFuncInTable   Number of separate functions in table, default is 2
494      *  \tparam     funcIndex0       Index of 1st function to evaluate in table, default is 0
495      *  \tparam     funcIndex1       Index of 2nd function to evaluate in table, default is 1
496      *  \tparam     T                Type (SimdReal or real) of lookup and result
497      *  \param      r                Points for which to evaluate function derivative
498      *  \param[out] derivativeValue1 Interpolated derivative for first function
499      *  \param[out] derivativeValue2 Interpolated derivative for second function
500      *
501      *  For debug builds we assert that the input values fall in the range
502      *  specified when constructing the table.
503      */
504     template<int numFuncInTable = 2, int funcIndex0 = 0, int funcIndex1 = 1, typename T>
505     void evaluateDerivative(T r, T* derivativeValue1, T* derivativeValue2) const
506     {
507         rangeCheck(r);
508         GMX_ASSERT(numFuncInTable == numFuncInTable_,
509                    "Evaluation method not matching number of functions in table");
510         GMX_ASSERT(funcIndex0 < numFuncInTable && funcIndex1 < numFuncInTable,
511                    "Function index not in range of the number of tables");
512
513         T    rTable   = r * T(tableScale_);
514         auto tabIndex = cvttR2I(rTable); // type is either std::int32_t or SimdInt32
515         T    eps      = rTable - trunc(rTable);
516
517         if (numFuncInTable == 2 && funcIndex0 == 0 && funcIndex1 == 1)
518         {
519             T t0A, t0B, t1A, t1B;
520             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data(), tabIndex,
521                                                           &t0A, &t0B); // works for scalar T too
522             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + 2, tabIndex,
523                                                           &t1A, &t1B); // works for scalar T too
524             *derivativeValue1 = fma(t1A - t0A, eps, t0A);
525             *derivativeValue2 = fma(t1B - t0B, eps, t0B);
526         }
527         else
528         {
529             T t0, t1, t2;
530             // This is not ideal, but we need a version of gatherLoadUBySimdIntTranspose that
531             // only loads a single value from memory to implement it better (will be written)
532             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex0,
533                                                           tabIndex, &t0, &t2); // works for scalar T too
534             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex0,
535                                                           tabIndex + T(1), &t1,
536                                                           &t2); // works for scalar T too
537             *derivativeValue1 = fma(t1 - t0, eps, t0);
538
539             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex1,
540                                                           tabIndex, &t0, &t2); // works for scalar T too
541             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex1,
542                                                           tabIndex + T(1), &t1,
543                                                           &t2); // works for scalar T too
544             *derivativeValue2 = fma(t1 - t0, eps, t0);
545         }
546     }
547
548     /************************************************************
549      *            Evaluation methods for three functions        *
550      ************************************************************/
551
552
553     /*! \brief Evaluate both function and derivative, three table functions
554      *
555      *  This is a templated method where the template can be either real or SimdReal.
556      *
557      *  \tparam     numFuncInTable   Number of separate functions in table, default is 3
558      *  \tparam     funcIndex0       Index of 1st function to evaluate in table, default is 0
559      *  \tparam     funcIndex1       Index of 2nd function to evaluate in table, default is 1
560      *  \tparam     funcIndex2       Index of 3rd function to evaluate in table, default is 2
561      *  \tparam     T                Type (SimdReal or real) of lookup and result
562      *  \param      r                Points for which to evaluate function and derivative
563      *  \param[out] functionValue1   Interpolated value for first function
564      *  \param[out] derivativeValue1 Interpolated derivative for first function
565      *  \param[out] functionValue2   Interpolated value for second function
566      *  \param[out] derivativeValue2 Interpolated derivative for second function
567      *  \param[out] functionValue3   Interpolated value for third function
568      *  \param[out] derivativeValue3 Interpolated derivative for third function
569      *
570      *  For debug builds we assert that the input values fall in the range
571      *  specified when constructing the table.
572      */
573     template<int numFuncInTable = 3, int funcIndex0 = 0, int funcIndex1 = 1, int funcIndex2 = 2, typename T>
574     void evaluateFunctionAndDerivative(T  r,
575                                        T* functionValue1,
576                                        T* derivativeValue1,
577                                        T* functionValue2,
578                                        T* derivativeValue2,
579                                        T* functionValue3,
580                                        T* derivativeValue3) const
581     {
582         rangeCheck(r);
583         GMX_ASSERT(numFuncInTable == numFuncInTable,
584                    "Evaluation method not matching number of functions in table");
585         GMX_ASSERT(funcIndex0 < numFuncInTable && funcIndex1 < numFuncInTable && funcIndex2 < numFuncInTable,
586                    "Function index not in range of the number of tables");
587
588         T    rTable   = r * T(tableScale_);
589         auto tabIndex = cvttR2I(rTable); // type is either std::int32_t or SimdInt32
590         T    eps      = rTable - trunc(rTable);
591         T    t0;
592         T    t1;
593         T    t2;
594         T t3 gmx_unused;
595
596         gatherLoadBySimdIntTranspose<4 * numFuncInTable>(
597                 ddfzMultiTableData_.data() + 4 * funcIndex0, tabIndex, &t0, &t1, &t2, &t3);
598         t1                = t0 + eps * t1;
599         *functionValue1   = fma(eps * T(halfSpacing_), t0 + t1, t2);
600         *derivativeValue1 = t1;
601
602         gatherLoadBySimdIntTranspose<4 * numFuncInTable>(
603                 ddfzMultiTableData_.data() + 4 * funcIndex1, tabIndex, &t0, &t1, &t2, &t3);
604         t1                = t0 + eps * t1;
605         *functionValue2   = fma(eps * T(halfSpacing_), t0 + t1, t2);
606         *derivativeValue2 = t1;
607
608         gatherLoadBySimdIntTranspose<4 * numFuncInTable>(
609                 ddfzMultiTableData_.data() + 4 * funcIndex2, tabIndex, &t0, &t1, &t2, &t3);
610         t1                = t0 + eps * t1;
611         *functionValue3   = fma(eps * T(halfSpacing_), t0 + t1, t2);
612         *derivativeValue3 = t1;
613     }
614
615     /*! \brief Evaluate function value only, three table functions
616      *
617      *  This is a templated method where the template can be either real or SimdReal.
618      *
619      *  \tparam     numFuncInTable   Number of separate functions in table, default is 3
620      *  \tparam     funcIndex0       Index of 1st function to evaluate in table, default is 0
621      *  \tparam     funcIndex1       Index of 2nd function to evaluate in table, default is 1
622      *  \tparam     funcIndex2       Index of 3rd function to evaluate in table, default is 2
623      *  \tparam     T                Type (SimdReal or real) of lookup and result
624      *  \param      r                Points for which to evaluate function value
625      *  \param[out] functionValue1   Interpolated value for first function
626      *  \param[out] functionValue2   Interpolated value for second function
627      *  \param[out] functionValue3   Interpolated value for third function
628      *
629      *  For debug builds we assert that the input values fall in the range
630      *  specified when constructing the table.
631      */
632     template<int numFuncInTable = 3, int funcIndex0 = 0, int funcIndex1 = 1, int funcIndex2 = 2, typename T>
633     void evaluateFunction(T r, T* functionValue1, T* functionValue2, T* functionValue3) const
634     {
635         T der1 gmx_unused;
636         T der2 gmx_unused;
637         T der3 gmx_unused;
638
639         evaluateFunctionAndDerivative<numFuncInTable, funcIndex0, funcIndex1, funcIndex2>(
640                 r, functionValue1, &der1, functionValue2, &der2, functionValue3, &der3);
641     }
642
643     /*! \brief Evaluate function derivative only, three table functions
644      *
645      *  This is a templated method where the template can be either real or SimdReal.
646      *
647      *  \tparam     numFuncInTable   Number of separate functions in table, default is 3
648      *  \tparam     funcIndex0       Index of 1st function to evaluate in table, default is 0
649      *  \tparam     funcIndex1       Index of 2nd function to evaluate in table, default is 1
650      *  \tparam     funcIndex2       Index of 3rd function to evaluate in table, default is 2
651      *  \tparam     T                Type (SimdReal or real) of lookup and result
652      *  \param      r                Points for which to evaluate function derivative
653      *  \param[out] derivativeValue1 Interpolated derivative for first function
654      *  \param[out] derivativeValue2 Interpolated derivative for second function
655      *  \param[out] derivativeValue3 Interpolated derivative for third function
656      *
657      *  For debug builds we assert that the input values fall in the range
658      *  specified when constructing the table.
659      */
660     template<int numFuncInTable = 3, int funcIndex0 = 0, int funcIndex1 = 1, int funcIndex2 = 2, typename T>
661     void evaluateDerivative(T r, T* derivativeValue1, T* derivativeValue2, T* derivativeValue3) const
662     {
663         rangeCheck(r);
664         GMX_ASSERT(numFuncInTable == numFuncInTable_,
665                    "Evaluation method not matching number of functions in table");
666         GMX_ASSERT(funcIndex0 < numFuncInTable && funcIndex1 < numFuncInTable && funcIndex2 < numFuncInTable,
667                    "Function index not in range of the number of tables");
668
669         T    rTable   = r * T(tableScale_);
670         auto tabIndex = cvttR2I(rTable); // type is either std::int32_t or SimdInt32
671         T    eps      = rTable - trunc(rTable);
672
673         if (numFuncInTable == 3 && funcIndex0 == 0 && funcIndex1 == 1 && funcIndex2 == 2)
674         {
675             T t0A, t0B, t0C, t1A, t1B, t1C;
676             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data(),
677                                                           tabIndex, &t0A, &t0B);
678             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + 2,
679                                                           tabIndex, &t0C, &t1A);
680             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + 4,
681                                                           tabIndex, &t1B, &t1C);
682             *derivativeValue1 = fma(t1A - t0A, eps, t0A);
683             *derivativeValue2 = fma(t1B - t0B, eps, t0B);
684             *derivativeValue3 = fma(t1C - t0C, eps, t0C);
685         }
686         else
687         {
688             T t0, t1, t2;
689             // This is not ideal, but we need a version of gatherLoadUBySimdIntTranspose that
690             // only loads a single value from memory to implement it better (will be written)
691             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex0,
692                                                           tabIndex, &t0, &t2); // works for scalar T too
693             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex0,
694                                                           tabIndex + T(1), &t1,
695                                                           &t2); // works for scalar T too
696             *derivativeValue1 = fma(t1 - t0, eps, t0);
697
698             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex1,
699                                                           tabIndex, &t0, &t2); // works for scalar T too
700             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex1,
701                                                           tabIndex + T(1), &t1,
702                                                           &t2); // works for scalar T too
703             *derivativeValue2 = fma(t1 - t0, eps, t0);
704
705             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex2,
706                                                           tabIndex, &t0, &t2); // works for scalar T too
707             gatherLoadUBySimdIntTranspose<numFuncInTable>(derivativeMultiTableData_.data() + funcIndex2,
708                                                           tabIndex + T(1), &t1,
709                                                           &t2); // works for scalar T too
710             *derivativeValue3 = fma(t1 - t0, eps, t0);
711         }
712     }
713
714     /*! \brief Return the table spacing (distance between points)
715      *
716      *  You should never have to use this for normal code, but due to the
717      *  way tables are constructed internally we need this in the unit tests
718      *  to check relative tolerances over each interval.
719      *
720      *  \return table spacing.
721      */
722     real tableSpacing() const { return 1.0 / tableScale_; }
723
724 private:
725     std::size_t           numFuncInTable_; //!< Number of separate tabluated functions
726     std::pair<real, real> range_;          //!< Range for which table evaluation is allowed
727     real                  tableScale_;     //!< Table scale (inverse of spacing between points)
728     real                  halfSpacing_;    //!< 0.5*spacing (used for DDFZ table data)
729
730     //!< Derivative values only, with the third-derivative subtraction described in the class documentation.
731     std::vector<real> derivativeMultiTableData_;
732
733     /*! \brief Combined derivative, difference to next derivative, value, and zero.
734      *
735      *  For table point i, this vector contains the four values:
736      *  - derivative[i]
737      *  - (derivative[i+1]-derivative[i])
738      *  - value[i]
739      *  - 0.0
740      *
741      *  For the derivative terms we have subtracted the third-derivative term described
742      *  in the main class documentation.
743      *
744      *  This is typically more efficient than the individual tables, in particular
745      *  when using SIMD. The result should be identical outside the class, so this
746      *  is merely an internal implementation optimization. However, to allow
747      *  aligned SIMD loads we need to use an aligned allocator for this container.
748      *  We occasionally abbreviate this data as DDFZ.
749      */
750     std::vector<real, AlignedAllocator<real>> ddfzMultiTableData_;
751
752     // There should never be any reason to copy the table since it is read-only
753     GMX_DISALLOW_COPY_AND_ASSIGN(QuadraticSplineTable);
754 };
755
756
757 } // namespace gmx
758
759 #endif // GMX_TABLES_QUADRATICSPLINETABLE_H