| Title: | Thorn Guide for the AEILocalInterp Thorn |
| Author(s): | Jonathan Thornburg |
| incorporating code from Thomas Radke |
| Date: | Date |
This thorn provides the CCTK_InterpLocalUniform() API for processor-local interpolation of N-dimensional data arrays. The data arrays (in general there may be many of them) must be defined on a uniformly spaced grid.
This thorn provides several variants of Lagrange and Hermite interpolation. At present the Lagrange interpolation operators support orders 1 through 6 for 1-D interpolation, and 1 through 4 for 2-D and 3-D interpolation. The Hermite interpolation operators support orders 2, 3, and 4 for 1-D interpolation, and 2 and 3 for 2-D and 3-D interpolation.
This thorn supports a number of interpolation options, including parameter-table entries to control the handling of grid boundaries and out-of-range points, non-contiguous input arrays, taking derivatives as part of the interpolation, and querying the Jacobian of the interpolation operation.
This interpolator provides:
Section A1.3.5 of this thorn guide explains the two variants in detail, but for most purposes the tensor-product variant is preferable. For convenience this is also registered under the operator name
and this is probably the best interpolation operator for general-purpose use.1
The code allows arbitrarily-shaped interpolation molecules, but at the moment only hypercube-shaped molecules are implemented. It would be easy to add additional orders and/or dimensions if desired.
This interpolator supports a number of options specified by a parameter table (a Cactus key-value table, manipulated by the Util_Table*() APIs). Note that all the interpolator options described here apply only to the current interpolator call: there is no visible state kept inside the interpolator from one call to another.
This interpolator was written by Jonathan Thornburg in winter 2001–2002. Between then and July 2003 it lived in the LocalInterp thorn in the CactusBase arrangement, but in July 2003 it was moved to the (new) AEILocalInterp thorn in the AEIThorns arrangement so it could stay GPL (Cactus policies forbid GPL code in the CactusBase arrangement).
Within Cactus, each interpolator has a character-string name; this is mapped to a Cactus interpolator handle by CCTK_InterpHandle(). For any given interpolator handle, there may be a separate interpolator defined for each of the interpolation APIs (both the processor-local ones provided by this thorn, and the global ones provided by driver-specific interpolator thorns such as PUGHInterp).
Terminology for interpolation seems to differ a bit from one author to another. Here we refer to the point-centered interpolation of a set of input arrays (defining data on a uniformly or nonuniformly spaced grid of data points) to a set of interpolation points (specified by a corresponding set of interpolation coordinates), with the results being stored in a corresponding set of output arrays. Alternatively, we may refer to cell-centered interpolation, using a grid of data cells and a set of interpolation cells. (This latter terminology is common in hydrodynamics interpolators.)
At present all the interpolators do polynomial interpolation, and allow the interpolation of multiple input arrays (to the same set of interpolation points) in a single interpolator call, using the basic algorithm:
In the future, we may add other interpolators where the choice of molecule is data-dependent (and thus may vary from one input array to the next), and/or where the entire input grid is used in each interpolation.
We define the order of the interpolation to be the order of the fitted polynomial. That is, in our terminology
linear interpolation is order 1, quadratic is order 2, cubic is order 3, etc. An order n interpolator thus has
O
(Δx)n+1)
interpolation errors for generic smooth input data. Section A1.3.5 explains how the interpolating
polynomial is defined for multi-dimensional interpolation.
Because the interpolating polynomial generally changes if the interpolation point moves from one grid cell to another, unless something special is done the interpolating function isn’t smooth, i.e. its 1st derivative is generically discontinuous, with O(Δxn) jump discontinuities each time the interpolating polynomial changes. Correspondingly, the interpolation error is generically a “bump function” which is zero at each grid point and rises to a local maximum in each grid cell.2 This is the case, for example, for tensor-product Lagrange polynomial interpolation.
[For maximum-degree Lagrange polynomial interpolation, the interpolation error may not be zero at the grid points, and the interpolant itself may not be continuous there. Section A1.3.5 explains this in detail.]
This thorn also provides Hermite polynomial interpolation, which guarantees a smooth (C1) interpolant and (for smooth input data) interpolation error. Unfortunately, this comes at the cost of a larger molecule size than the same-order Lagrange interplator, and a much larger interpolation error if the interpolation molecules are off-centered. (By default, the Hermite interpolator doesn’t allow off-centered molecules, though this can be changed via the boundary_off_centering_tolerance[] and excision_off_centering_tolerance[] parameter-table entries.)
For most purposes Lagrange polynomial interpolation is better; only use Hermite polynomial interpolation if you need a smooth interpolant.
As described in the Function Reference section of the Cactus User’s Guide, interp_coords, input_arrays, and output_arrays are actually all pointers to arrays of void * pointers, since we support a number of different Cactus data types. Internally, the interpolator casts these void * pointers to CCTK_REAL * or whatever the correct Cactus data types are. But when describing how the interpolator accesses the various arrays, for simplicity we usually gloss over this casting, i.e. we pretend that interp_coords, input_arrays, and output_arrays are pointers to arrays of CCTK_REAL * pointers. (This may become clearer once you read the next example.)
We use pt, in, and out as generic 0-origin integer subscripts into the arrays of interpolation points, input arrays, and output arrays respectively. We use (i,j,k) as a generic N_dims-vector of integer subscripts into the input array input_arrays[in]. (Thus (i,j,k) space refers to the grid of data points.) We usually only write array subscripting expressions for the 3-D case; the restrictions/generalizations to other dimensions should be obvious.
For example, for 3-D interpolation, the (x,y,z) coordinates of a typical interpolation point are given by
(Notice here that as described above, we’ve implicitly taken interp_coords to have the C type
const CCTK_REAL* interp_coords[], glossing over the casting from its actual C type of
const void* interp_coords[].)
We take axis to be a 0-origin integer specifying a coordinate axis (dimension), i.e. 0 for x, 1 for y, 2 for z, …. We take ibndry to be a 0-origin integer specifying the combination of a coordinate axis (dimension) and a minimum/maximum “end” of the grid; these are enumerated in the order xmin, xmax, ymin, ymax, zmin, zmax, ….
When describing Jacobians and domains of dependence, it’s useful to introduce the notion of molecule position, a nominal reference point for the interpolation molecule in (i,j,k) coordinate space. (For example, the molecule position might just be the (i,j,k) coordinates of the molecule’s center.) We also introduce molecule coordinates (mi,mj,jk), which are just (i,j,k) coordinates relative to the molecule position. We use m as a generic molecule coordinate. Thus (in notation which should be obvious) a generic molecule operation can be written
![]() | (A1.1) |
Note that there is no requirement that the output be semantically located at the position posn! (This may become clearer once you look at the example in section A1.3.9.) However, by convention we do assume that m = 0 is always a valid m coordinate; this simplifies pointer arithmetic.
When describing various entries in the parameter table in, we use const qualifiers on table entries to indicate that the interpolator treats them as const variables/arrays, i.e. it promises not to change them. In contrast, table entries which are not shown as const are considered read-write by the interpolator; it typically uses them to return outputs back to the caller.
The CCTK_InterpLocalUniform() interpolation API (described in detail in the “Function Reference” appendix of the Cactus Users’ Guide) includes a parameter table (a Cactus key-value table, manipulated by the Util_Table*() APIs) whose table handle is passed to the interpolator. This can be used to specify a number of options to the interpolator; it can also be used by the interpolator to return additional information back to the caller.
The only mandatory parameter for this interpolator is the interpolation order:
As noted in section A1.2.2, in our terminology linear interpolation is order 1, quadratic is order 2, cubic is order 3, etc.
If no grid boundaries or excised points are nearby, the interpolator centers the molecules around the interpolation point as much as possible. Table A1.1 gives the molecule size and details of the centering policy for each interpolation operator and order.
| Molecule | |||
| Interpolation Operator | Order | Size | MoleculeCentering Policy |
| "Lagrange polynomial interpolation" | 1 | 2-point |
|
| "Lagrange polynomial interpolation" | 2 | 3-point |
|
| "Lagrange polynomial interpolation" | 3 | 4-point |
|
| "Lagrange polynomial interpolation" | 4 | 5-point |
|
| "Lagrange polynomial interpolation" | 5 | 6-point |
|
| "Lagrange polynomial interpolation" | 6 | 7-point |
|
| "Hermite polynomial interpolation" | 2 | 4-point |
|
| "Hermite polynomial interpolation" | 3 | 6-point |
|
| "Hermite polynomial interpolation" | 4 | 6-point |
|
If an interpolation point is too near a grid boundary and/or an excised region the interpolator can either
off-center the interpolation molecule, or refuse to interpolate that point (returning a
CCTK_INTERP_ERROR_POINT_OUTSIDE3
or CCTK_INTERP_ERROR_POINT_EXCISED error code, and possibly also printing a warning message). This behavior
is controlled by the (optional) parameter-table entries
The elements of these arrays are matched up with the grid boundaries in the order
[xmin,xmax,ymin,ymax,zmin,zmax,…]. (For excision_*_tolerance[] the minimum/maximum are
interpreted as the corresponding coordinate decreasing/increasing when moving out of the active
region of the grid into the excised region.) We use ibndry as a generic 0-origin index into these
arrays.
For any given grid boundary ibndry, the interpolator doesn’t use input data from the
N_boundary_points_to_omit[ibndry] grid points closest to the boundary. In other words, these points are
effectively omitted from the input grid. For example:
More generally, this option may be useful for controlling the extent to which the interpolator uses input data from ghost and/or symmetry zones.
The (optional) parameter-table entries
control the local interpolator’s willingness to off-center (change away from the default centering) its interpolation molecules near boundaries and excised points.
To describe how these parameters work, it’s useful to define two regions in the data coordinate space:
Figure A1.1 shows an example of these regions.
| xmin (ibndry = 0) |
| xmax (ibndry = 1) |
| ymin (ibndry = 2) |
| ymax (ibndry = 3) |
Near a boundary or excised regions, the interpolator will off-center its molecules (if and) only if the interpolation point is both
where we use “*” to denote boundary or excision as appropriate.
There are four cases for these parameters:
If any of these table entries aren’t specified, the defaults are
In other words, the interpolation points may be anywhere within the valid-data region or up to 10-10 grid spacing outside it. (The 10-10 “fudge factor” helps to avoid spurious errors in case floating-point roundoff moves an interpolation point which was supposed to be just on the boundary of the valid-data region, slightly outside it.)
[In the near future these defaults may be changed: For Lagrange polynomial interpolation the new defaults would be
while for Hermite polynomial interpolation they would be
This would leave Lagrange interpolation unchanged, while for Hermite interpolation the defaults would forbid any significant off-centering of the interpolation molecules.]
By default, AEILocalInterp prints a Cactus level 1 warning message for each interpolation point which causes it to return a CCTK_INTERP_ERROR_POINT_OUTSIDE or CCTK_INTERP_ERROR_POINT_EXCISED error code. If the key suppress_warnings is present in the parameter table (it may have any datatype and value), then these messages are not printed.
By default, an interpolator just returns a single result, either 0 for success or some (negative) error code. If there are multiple interpolation points causing errors, an interpolator should return 0 for success, or the error code for the first point in error (the one with the smallest pt).
However, sometimes you would like to know the status of the interpolation for each point. The following (optional) parameter-table entries may be used to request this, and to query if the interpolator supports this:
If the interpolator supports returning per-point status, and if the key per_point_status is present in the parameter table, then the interpolator will set
in the parameter table. If per_point_status = NULL, the interpolator will set error_point_status to 0. Otherwise, the interpolator will set error_point_status to the negative of the number of points for which the per-point status is non-zero.
In multiple dimensions, there are two plausible definitions for the generic interpolating polynomial of a given order (degree) n. For convenience I’ll describe these for the 2-D case, but the generalization to any number of dimensions should be obvious:
![]() | (A1.2) |
![]() | (A1.3) |
Because it has (n + 1)2 coefficients, the tensor-product polynomial f can (and does) pass through all the (n + 1)2 input data points in a square molecule of size n + 1. This implies that the interpolation error vanishes at the input grid points, and that the overall interpolating function is continuous (up to floating-point roundoff errors). However, f does have the slightly peculiar property of having terms up to xnyn despite being of “degree” n. (For example, the “linear” interpolant for n = 1 would have xy terms, even though those are formally quadratic in the independent variables x and y.)
In contrast, the maximum-degree polynomial g has only
(n + 1)(n + 2) coefficients, so for generic input data
(i.e. input data which isn’t actually sampled from a polynomial of the form (A1.3)) g can’t pass through all the
(n + 1)2 input data points in a square molecule of size n + 1. The interpolator actually does a least-squares fit of
the polymomial g to the input data in the molecule, so in general the molecule won’t pass through any of the
data points! Moreover, each time the interpolation point crosses a grid square (for odd n) or the center lines of a
grid square (for even n), the set of points used in the molecule changes (cf. section A1.3.2), so the interpolant
generally has a jump discontinuity! For these reasons, the tensor-product choice (A1.2) is generally
preferable.
An interpolator may support different families/shapes of interpolation molecules. Hypercube-shaped molecules are the simplest and most common case, but one could also imagine (say) octagon-shaped molecules in 2-D, or some generalization of this in higher numbers of dimensions.
The (optional) parameter
may be used to set or query the molecule family.
If this key is present in the parameter table, the interpolator sets the molecule family/shape based on the value specified. If this key isn’t present in the parameter table, then the interpolator sets it to the molecule family being used.
At present only hypercubed-shaped molecules are implemented; these are specified by molecule_family set to "cube".
Sometimes the input “arrays” used by the interpolator might not be contiguous in memory. For example, we might want to do 2-D interpolation within a plane of a 3-D grid array, but the plane might or might not be contiguous in memory. (Another example would be that the input arrays might be members of a compact group.)
The following (optional) parameter-table entries specify non-contiguous input arrays:
In general, the interpolator accesses the input using the generic subscripting expression
where
and where (i,j,k) run from input_array_min_subscripts[] to input_array_max_subscripts[] inclusive (n.b. this is an inclusive range, i.e. min ≤(i,j,k) ≤max).
The defaults are that each input array is contiguous in memory, i.e. input_array_offsets[] = 0, stride determined from input_array_dims[] in the usual Fortran manner, input_array_min_subscripts[] = all 0s, and input_array_max_subscripts[] = input_array_dims[]-1. If the stride and max subscript are both specified explicitly, then the explicit input_array_dims[] argument to CCTK_InterpLocalUniform() is ignored.
If we view the input data as being samples of a smooth function, then instead of estimating values of that function at the interpolation points, the interpolator can instead or additionally estimate values of various derivatives of that function at the interpolation points. (We don’t currently implement it, but one could also imagine interpolating more general molecule operations such as Laplacians.)
The following (optional) parameter-table entries are used to specify this:
The semantics here are that
where op is an operator specified by the operation_codes[out] value as described below.
Note that operand_indices[] doesn’t directly name the inputs, but rather gives their (0-origin) subscripts in the list of inputs. This allows for a more efficient implementation in the (common) case where some of the input arrays have many different operations applied to them. (It’s most efficient to group all operations on a given input array together in the operand_indices and operation_codes arrays, as in the example in section A1.6.2.)
Negative operation_codes[out] values are reserved for future use. An operation_codes[out] value which is ≥ 0 is taken as a decimal integer encoding a coordinate partial derivative: each decimal digit means to take the coordinate partial derivative along that (1-origin) axis; the order of the digits in a number is ignored. Table A1.2 summarizes these resulting derivative codes.
| operation_codes[out] | What it means |
| 0 | interpolate the input array itself (no derivative) |
| 1 | interpolate ∂ ∂x1 of the input array |
| 2 | interpolate ∂ ∂x2 of the input array |
| 3 | interpolate ∂ ∂x3 of the input array |
| 12 or 21 | interpolate ∂2 ∂x1∂x2 of the input array |
| 13 or 31 | interpolate ∂2 ∂x1∂x3 of the input array |
| 23 or 32 | interpolate ∂2 ∂x2∂x3 of the input array |
| 11 | interpolate ∂2 ∂(x1)2 of the input array |
| 22 | interpolate ∂2 ∂(x2)2 of the input array |
| 33 | interpolate ∂2 ∂(x3)2 of the input array |
At present we do not have any #defines for the operation codes in the Cactus header files. However, you can avoid most of the software-engineering problems of having “magic numbers” for the operation codes, by using the macro
to mark all such operation_codes[] values in your code. There’s an example of this in section A1.6.2.
The Jacobian of the interpolator is defined as
![]() | (A1.4) |
We may want to know the Jacobian itself, and/or “just” the set of (i,j,k) for which this is nonzero (i.e. the Jacobian’s sparsity structure, or equivalently the domain of dependence of the output result, or equivalently the interpolation molecule size and shape) for a given out, in, and pt.5
The complexity of doing this depends (strongly!) on the structure of the Jacobian, and in particular on the answers to the following questions:
Because the different cases differ so much in the amount of information required to describe the Jacobian, it’s hard to define a single API which covers all cases without burdening the simpler cases with excessive complexity. Instead, the interpolator defines a Jacobian-structure–query API to determine which case applies, together with (conceptually) several different APIs for the different cases. (At the moment only a single case is implemented.)
The following parameter-table entries may be used to query which of the different Jacobian-structure cases applies:
The parameter molecule_family may may be used to query what molecule family is being used. This is described in detail in section A1.3.6.
If the interpolation molecule size and/or shape vary with the interpolation coordinates, the interpolator sets the parameter
to 1. Otherwise (i.e. if the interpolation molecule size and shape are independent of the interpolation coordinates) it should set this parameter to 0.
If the interpolator supports computing derivatives as described in section A1.3.8, and if the interpolation molecule’s size and/or shape varies with operation_codes[], the interpolator sets the parameter
to 1. Otherwise (i.e. if the interpolator doesn’t support computing derivatives, or if the interpolator does support computing derivatives but the interpolation molecule size and shape are independent of the operation_code[] values), the interpolator sets this parameter to 0. Note that this query tests whether the molecule size and/or shape depend on operation_codes[] in general, independent of whether there are in fact any distinct values (or even any values at all) passed in operation_codes[] in this particular interpolator call. In other words, this is a query about the basic design of the interpolator, not about this particular call.
If the interpolation molecule’s size and/or shape varies with the actual floating-point values of the input arrays, the interpolator sets the parameter
to 1. Otherwise (i.e. if the interpolation molecule size and shape are independent of the input array values; this is a necessary, but not sufficient, condition for the interpolation to be linear), the interpolator sets this parameter to 0.
If the actual floating-point values of the Jacobian (A1.4) (for a given out, in, and pt) depend on the actual floating-point values of the input arrays (i.e. if the interpolation is nonlinear), the interpolator sets the parameter
to 1. Otherwise (i.e. if the interpolation is linear) the interpolator sets this parameter to 0.
The current implementation always sets
The simplest case (and the only one for which we have defined an API at present) is when the molecules are hypercube-shaped and of (typically small) fixed size, independent of the interpolation coordinates and the actual floating-point values in the input arrays (though presumably depending on the interpolation order and on operation_code). In other words, this case applies if (and only if) the Jacobian structure information described in section A1.3.9 returns
(These are precisely the values set by the current implementation.) In the rest of this section we describe the query API for this case.
The following parameters may be used to query the molecule size:
The semantics of these are that if both of these keys are present (the values and datatypes don’t matter), then the interpolator will (re)set the values to give the (inclusive) minimum and maximum m molecule coordinates. (Note that either both of these keys should be present, or neither should be present. This simplifies the logic in the interpolator slightly.)
The following parameter may be used to query the molecule positions:
The semantics of this is that the caller should set molecule_positions[] to an array of N_dims pointers to (caller-supplied) arrays of N_interp_points CCTK_INTs each. If this key exists, then the interpolator will store the molecule positions in the pointed-to arrays.
The following parameters may be used to query the Jacobian (A1.4) itself:
If Jacobian_pointer is present in the table, then Jacobian_interp_point_stride and Jacobian_m_strides must also be present. For each out where Jacobian_pointer[out] != NULL,7 the interpolator would then store the Jacobian (A1.4) in
where
and where part is 0 for real values and the real parts of complex values, and 1 for the imaginary parts of complex values.
By appropriately setting the various stride parameters, this allows a fairly wide variety of possible storage layouts for the Jacobian.
An example may help to clarify this: Suppose we have a 1-D grid with 11 grid points, with integer subscripts 0 through 10 inclusive, and interpolation coordinates given by coord_origin = 0.0 and coord_delta = 0.1. Suppose further that we’re doing Lagrange polynomial interpolation, with order = 2 and hence (by table A1.1) using 3-point molecules. Finally, suppose that we query the Jacobian molecule positions for the N_interp_points=14 interpolation points 0.0, 0.04, 0.06, 0.10, 0.14, 0.16, 0.20, 0.80, 0.84, 0.86, 0.90, 0.94, 0.96, 1.00. Then the queries might return
The way the generalized polynomial interpolator is implemented it’s easy to also do Savitzky-Golay smoothing.8 This is best described by way of an example: Suppose we’re doing 1-D cubic interpolation, using (by table A1.1) 4-point molecules. In other words, at each interpolation point we use a cubic interpolation polynomial fitted to 4 surrounding data points. For Savitzky-Golay smoothing, we would instead least-squares fit a cubic polynomial to some larger number of surrounding data points. This combines interpolation with smoothing, so there’s less amplification of noise in the input data in the interpolation outputs.
The optional input
specifies how much (how many points) to enlarge the interpolation molecule for this. The default is 0 (no smoothing). 1 would mean to enlarge the molecule by 1 point (e.g. to use a 5-point molecule instead of the usual 4-point one for cubic interpolation). 2 would mean to enlarge by 2 points, (e.g. to use a 6-point molecule for cubic interpolation). Etc etc.
Note that in > 1 dimension, the maximum-degree Lagrange interpolation already uses more data points than the number of free parameters in the interpolating polynomials, i.e. it already does some Savitzky-Golay smoothing. For example, in 2-D a generic cubic polynomial f(x,y) = ∑ i+j≤3cijxiyj has 10 free parameters cij, which we least-squares fit to the 16 data points in the 4 × 4 molecule.
Savitzky-Golay smoothing is basically free apart from the increase in the molecule size, e.g. a smoothing=2 cubic interpolation has exactly the same cost as any other 6-point–molecule interpolation.
The current implementation has all the framework for smoothing, but only the smoothing=0 case is implemented at the moment.
Setting the optional parameter
to a positive value turns on various debug printing inside the interpolator. This is intended for debugging the interpolator itself; you need to look at the interpolator source code to interpret the output.
This thorn has a Boolean parameter log_interp_coords. If this is set to true, this thorn will write a detailed log file giving the grid information and interpolation coordinates for each interpolator call. Each processor will write a separate log file, with the file name determined via
Each file begins with a header comment describing the file format in detail. After the header comment, there is one line in the file for each call on CCTK_InterpLocalUniform() which is handled by this thorn (with log_interp_coords set to true).
Note that these log files are typically quite large, and writing them may significantly slow the simulation – you should probably only use this option for debugging purposes.
Sometimes you may want to call the interpolator, but not actually do any interpolation for some or all of the arrays:
In such cases there are several possible ways you can supress some or all of the interpolations:
Note that the combination input_arrays[in] = NULL and output_arrays[out]≠NULL is an error: conceptually you’re asking for an interpolation (the output array pointer is non-NULL) but not supplying any input data (the input array pointer is NULL).
This interpolator’s basic design is to use separate specialized code for each combination of
i.e. in practice for each distinct choice of interpolation molecule. Maple is used to generate all the interpolation coefficients. The C preprocessor is then used to generate all the specialized code from a single master “template” file. The template code uses #ifdefs to handle lower dimensions with no extra overhead, e.g. 1-D/2-D interpolation is basically just as efficient as in a purpose-built 1-D/2-D interpolator, with no extra overhead imposed by the interpolator also supporting higher-dimensional interpolation.
The Maple code which generates the interpolation coefficients is quite general-purpose, and can handle an arbitrary dimensionality and molecule size/shape. Generating new coefficients can be rather time-consuming, though, e.g. the current coefficients for 3-D for orders 1-4 take about 8 cpu minutes to generate using Maple 7 on a 1.7 GHz P4.
Note that when compiling the code in the directory src/GeneralizedPolynomial-Uniform of this thorn, you may get compiler warnings about casts discarding const qualifiers from pointers in (the #include-ed file) template.c. Don’t worry – the code is actually ok.9 You may also get compiler warnings about unused variables in this same file; again don’t worry, the code is ok.
See the README file in the source code directory AEILocalInterp/src/ and in its subdirectories for further details on the implementation.
Here’s a simple example in C, of interpolating a CCTK_REAL and a CCTK_COMPLEX 10 × 20 2-D array, at 5 interpolation points, using cubic interpolation.
Note: Since C allows arrays to be initialized only if the initializer values are compile-time constants, we have to declare the interp_coords[], input_arrays[], and output_arrays[] arrays as non-const, and set their values with ordinary (run-time) assignment statements. In C++, there’s no restriction on initializer values, so we could declare the arrays const and initialize them as part of their declarations.
In this example we interpolate the 3-metric and its 1st derivatives at a set of interpolation points on a 2-sphere.
Note: Since we’re not using C++, we again declare the pointer arrays non-const, and “initialize” them with ordinary (run-time) assignment statements. However, in this example, for greater clarity we place these assignment statements right after the declarations. Since C allows declarations only at the start of a { } block, not in the middle of a block, we nest the rest of the program in extra blocks (with the { } indented 2 spaces to distinguish them from “normal” { } pairs) to allow for further declarations. The reader will have to decide for herself whether this style is more or less ugly than separating the declarations and initializations of the pointer arrays.
Thanks to Thomas Radke for our ongoing collaboration on Cactus interpolation, to Ian Hawke and Erik Schnetter for helpful comments on the documentation, to Tom Goodale and Thomas Radke for many useful design discussions, to Erik Schnetter for bug reports, and to all the Cactus crew for a great infrastructure!
log_interp_coords | Scope: private | BOOLEAN |
Description: should we log the grid min(delta)max and the interpolation coordinates for each call
on the interpolator?
| ||
| Default: false | ||
Implements:
aeilocalinterp
This section lists all the variables which are assigned storage by thorn AEIThorns/AEILocalInterp. Storage can either last for the duration of the run (Always means that if this thorn is activated storage will be assigned, Conditional means that if this thorn is activated storage will be assigned for the duration of the run if some condition is met), or can be turned on for the duration of a schedule function.
NONE
CCTK_STARTUP
aeilocalinterp_u_startup
register cctk_interplocaluniform() interpolation operators
| After: | driver_startup | |
| Language: | c | |
| Type: | function | |
1For backwards compatability this interpolator is also registered under still another operator name, "generalized polynomial interpolation", but this is deprecated and will probably be removed at some point.
2For a further discussion of the non-smoothness of interpolation errors, see appendix F of J. Thornburg, Physical Review D 59(10), 104007.
3For backwards compatability the error code CCTK_INTERP_ERROR_POINT_X_RANGE is also defined; it’s a synonym for CCTK_INTERP_ERROR_POINT_OUTSIDE.
4In practice 999.0 is a conventional value here. The actual value doesn’t matter so long as it’s larger than any possible molecule radius.
5For something like a spline interpolator the Jacobian is generally nonzero for all (i,j,k), but exponentially small for most of them; in this case the Jacobian-query API would probably specify a cutoff to return an approximate Jacobian with reasonable sparsity.
6We can always make the interpolation molecules be the “same size and shape” by padding them with zero entries to bring them all up to the size of the largest molecule used anywhere in the grid, but this would be very inefficient if many molecules were padded in this way.
7That is, setting Jacobian_pointer[out] = NULL supresses the Jacobian query. If (as is often the case) the Jacobian is independent of out, then it’s a good idea to avoid unnecessary computations by supressing the queries in this manner for any remaining output arrays.
8See section 14.8 of Numerical Recipes (2nd edition) for a general discussion of Savitzky-Golay smoothing.
9The warnings are really a design bug in C’s const-qualifier rules when pointers point to other pointers. See question 11.10 in the (excellent!) comp.lang.c online FAQ (http://www.eskimo.com/ scs/C-faq/top.html) for further discussion. gcc 2.95.3 gives the spurious warnings, but gcc 3.2.1 doesn’t, nor does the Intel C compiler in any version I’ve tried. Also, C++ has fixed the problems in C’s const-qualifier rules, so the same code compiled as C++ shouldn’t give any warnings (and doesn’t on the compilers I’ve tried).