GeMA
The GeMA main application
Interpolation functions

Data interpolation in GeMA can be done either locally, affecting a single cell / node, or globally, affecting the whole mesh. Available local interpolation functions are (Figure 1):

  • N2P - Nodes to Point: Calculates the value of a point inside a cell/element from nodal points;
  • G2P - Gauss to Point: Calculates the value of a point inside an element from Gauss points;
  • G2N - Gauss to Nodes: Extrapolates values from Gauss points to all element nodes;
  • GN2N - Gauss Neighbors to Nodes: Calculates the value of a node from Gauss points on neighbor elements.
  • NC2P - Node Cloud to Point: Calculates the value of a point from a cloud of mesh nodes.
  • GC2P - Gauss Cloud to Point: Calculates the value of a point from a cloud of Gauss points.

interpolationFunctionsLocal.png
Figure 1 - Local interpolation functions

Global interpolation functions are used for transfering data to every node or Gauss points in a destination mesh from data in node or Gauss points, in the same source mesh or in another one, using one of the following strategies (Figure 2):

  • EB_G2N - Element Based Gauss to Nodes: Extrapolates inside each element using the Gauss to Nodes strategy and averages each element's contribution to a node. Gauss points and nodes should belong to the same mesh;
  • NB_GN2N - Node Based Gauss Neighbors to Nodes: Calculates nodal values directly from the contribution of its neighbor Gauss nodes using the Gauss Neighbors to Nodes strategy. Gauss points and nodes should belong to the same mesh;
  • MESH_N2N - Mesh Nodes to mesh Nodes: Calculates nodal values on a destination mesh from the nodal values on a distinct source mesh using the given local interpolation method.
  • MESH_N2G - Mesh Nodes to mesh Gauss: Calculates Gauss point values on a destination mesh from the nodal values on a source mesh using the given local interpolation method. Source and destination meshes might be the same;
  • MESH_G2N - Mesh Gauss to mesh Nodes: Calculates nodal values on a destination mesh from the Gauss point values on a source mesh using the given local interpolation method. Source and destination meshes might be the same;
  • MESH_G2G - Mesh Gauss to mesh Gauss: Calculates gauss point values on a destination mesh from the Gauss point values on a source mesh using the global spline interpolation method. Source and destination meshes might be the same if the source and destination Gauss sets use differente integration rules;
  • SPLINE_N2N - Spline mesh Nodes to mesh Nodes: Calculates nodal values on a destination mesh from the nodal values on a distinct source mesh using the global spline interpolation method.
  • SPLINE_N2G - Spline mesh Nodes to mesh Gauss: Calculates Gauss point values on a destination mesh from the nodal values on a source mesh using the global spline interpolation method. Source and destination meshes might be the same;
  • SPLINE_G2N - Spline mesh Gauss to mesh Nodes: Calculates nodal values on a destination mesh from the Gauss point values on a source mesh using the global spline interpolation method. Source and destination meshes might be the same;
  • SPLINE_G2G - Spline mesh Gauss to mesh Gauss: Calculates gauss point values on a destination mesh from the Gauss point values on a source mesh using the global spline interpolation method. Source and destination meshes might be the same if the source and destination Gauss sets use differente integration rules;

The set of source nodes or Gauss points and/or the set of destination nodes or gauss points can be restricted by using node sets or cell group objects.

interpolationFunctionsGlobal.png
Figure 2 - Global interpolation functions

Several interpolation algorithms are available, and each of the above operations can be done with many of them, as sumarized at the tables below. A detailed description of the algorithms can be found on the Interpolation methods page.

Table 1: Summary of available interpolation methods

Method Description Parameter
shape Shape function interpolation. Available only for element meshes. -
lshape Shape function interpolation using a linear shape function. Available only for element meshes.
Usefull when dealing with super-parametric scenarios.
-
nn Nearest neighbor interpolation. k (1)
idw Inverse distance weighted interpolation. p (2)
mls Moving least squares interpolation. -
default or "" Default interpolator. Equal to "shape" for element meshes and to "idw" with p = 2 for cell meshes.
(1) If an integer parameter 'k' is given, the 'k' closest points will be averaged to return the interpolated value. The default value for 'k' is 1.
(2) The interpolated value is the average of the known values weighted by inverse distance from the interpolation point raised to a power 'p' provided as a real number parameter. If 'p' is absent, a power of 2 is used.

Table 2: Summary of allowed interpolation methods per type of interpolation function

Interpolation function shape lshape nn idw mls
N2P - Nodes to Point Yes (1) Yes (1) Yes Yes Yes
G2P - Gauss to Point Yes (2) Yes (2) Yes Yes Yes
G2N - Gauss to Nodes Yes Yes No No No
GN2N - Gauss Neighbors to Nodes No No Yes (3) Yes (3) Yes (3)
NC2P - Node Cloud to Point No No Yes Yes Yes
GC2P - Gauss Cloud to Point No No Yes Yes Yes
EB_G2N - Element Based Gauss to Nodes Yes Yes No No No
NB_GN2N - Node Based Gauss Neighbors to Nodes No No Yes (3) Yes (3) Yes (3)
MESH_N2N - Mesh Nodes to Mesh Nodes No No Yes Yes Yes
MESH_N2G - Mesh Nodes to Mesh Gauss No No Yes Yes Yes
MESH_G2N - Mesh Gauss to Mesh Nodes No No Yes Yes Yes
MESH_G2G - Mesh Gauss to Mesh Gauss No No Yes Yes Yes
SPLINE_N2N - Spline Mesh Nodes to Mesh Nodes (4) No No No No No
SPLINE_N2G - Spline Mesh Nodes to Mesh Gauss (4) No No No No No
SPLINE_G2N - Spline Mesh Gauss to Mesh Nodes (4) No No No No No
SPLINE_G2G - Spline Mesh Gauss to Mesh Gauss (4) No No No No No
(1) Can be used for element meshes only.
(2) Performed as a composed Gauss to Nodes + Nodes to Point operation.
(3) Requires a mesh with support for topological queries.
(4) Spline methods have its independent interpolation method.

The same interpolation algorithms are also used when interpolating nodal parameter values inside a cell user function or when transfering data from one mesh to another.

At the orchestration script, local interpolations are available by creating interpolator objects that receive as parameters the interpolation method, the required accessor objects for the known data to be interpolated and any other necessary interpolation function options. More than one data set can be interpolated together, including mixed sets with scalar and multi-dimensional accessors. The resulting object provides an interpolate() method that can then be called for interpolating cell values multiple times. Global interpolation is available through a set of global functions.

Index:

Local interpolation functions and methods

NodeToCellPointInterpolator(type, param, nodeAcList, coordAc, mesh)
NodeToElementPointInterpolator(type, param, nodeAcList, coordAc, mesh, gaussAcList)
Description: Creates a new "Nodes to Point" interpolator, capable of calculating values for a point inside a cell/element from node values. The called function should match the mesh type, so the NodeToCellPointInterpolator() version should be used together with cell meshes, while the NodeToElementPointINterpolator() with element meshes. Keep in mind that some interpolation methods are not available for cell meshes (most notably the shape function class of interpolators). For efficiency, multiple value sets can be interpolated at the same time. Also, the created object can be setup for saving interpolation results to a set of Gauss accessors instead of returning them.
Important: This call only creates and initializes the interpolator. The interpolation operation itself will take place when calling the interpolate() method over the returned object. This method can be called multiple times for interpolating values over different cells/elements.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. In particular, shape function methods are valid only for element meshes. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
nodeAcList Either a single node value accessor for the node attribute or state variable to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations).
coordAc An accessor for retrieving mesh node coordinates (see mesh:nodeCoordAccessor()).
mesh The mesh owning the given accessors. This parameter is required only if the desired interpolation algorithm requires this information. For most algorithms it can be set to nil.
gaussAcList Either a single Gauss accessor or a list with multiple accessors. When the desired interpolated points are integration points, providing a list with Gauss accessors ensures that interpolation results are automatically saved to the given accessors when interpolate() is supplied with an integration rule + integration point (the number of integration points in the rule must be compatible with the rule used by the accessors). The given list must have exactly the same number of accessors as in the nodeAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the Gauss list. Dimensions on the source and destination accessors must be equal. Alternativelly, destination accessors can be all scalars if the dimension parameter is used when calling interpolate().
Returns: Returns the new "Nodes to Point" interpolator object.

Example: See the examples below for the n2p:interpolate() method.

n2p:interpolate(cell, coord, dim)
n2p:interpolate(elem, coord, cartesian, dim)
n2p:interpolate(elem, intRule, intPoint, dim)
Description: Interpolates values on the given coordinate, belonging to the given cell, for the data sets configured when creating the interpolator object. The first signature should be used with cell meshes while the other two with element meshes. When using the third signature, the iterpolation coordinate is given implicitly by the coordinate of the given integration point. In this case, if the object was created with a set of Gauss accessors, results will be written directly to the accessors and won't be returned by the function. Otherwise, the function will return either a scalar number or a vector object, depending on the interpolated data dimension, if the interpolator was created with a single node data accessor. For objects created with multiple accessors, either a table with multiple numbers or a table with multiple vectors will be returned (the first only if all the interpolated values are scalars or if the dim parameter was given).
Parameters: cell / elem The cell / element object containing the interpolation point.
coord The interpolation coordinate. Can be either a Lua table or a vector object. When working with cell meshes, provided coordinates must be Cartesian coordinates. For element meshes, they can be either Cartesian or natural coordinates, depending on the interpolation method. If needed, natural coordinates are automatically converted to Cartesian coordinates, but the reverse is not done, so, if the interpolation method requires natural coordinates (like shape function based methods), the provided coordinates must be natural coordinates.
cartesian Required flag when using the second method signature (for element mehes). Indicates the type of values given in the coord parameter (true for Cartesian, false for natural).
intRule The integration rule object when using the third method signature. Together with intPoint defines the integration point over which the interpolation will be performed. If the interpolator object was created with a set of Gauss accessors, results will be written directly to those accessors. In that case, the integration rule must be compatible with the rule used to configure the Gauss accessors (check the ruleSet attribute used when creating the Gauss accessor).
intPoint The integration point index (a value between 1 and intRule:numPoints()) used together with intRule for providing the integration point over which the interpolation will be performed.
dim An optional parameter that can be used when dealing with multi-dimensional data to request that the interpolation be done over this dimension only. The interpolation will then return a single scalar value or a list of scalars when dealing with multiple data sets. In the later case, the given dimension must be a valid dimension for every one of the input data sets. For matrix data, the dimension is taken in relation to its linearized representation (in column major format).
Returns: If the object was created with Gauss accessors and the third signature was used, returns the number of source nodes used to interpolate. If the interpolator was created with a single node data accessor, returns either a scalar number or a vector object, depending on the interpolated data dimension, and the number of source nodes used to interpolate. For objects created with multiple accessors, returns either a table with multiple numbers (all accessors are scalars or dim was used) or a table with multiple vectors (each with its own particular dimension) and the number of source nodes used to interpolate. This method can also return only nil when mls interpolation method fails (see mls method).

Examples:

-- Example 1: Interpolates temperature data (scalar) using the shape function interpolation method
local m = modelData:mesh('myMesh')
local coordAc = m:nodeCoordAccessor()
local TAc = m:nodeValueAccessor('T')
local interpolator = NodeToElementPointInterpolator('shape', nil, TAc, coordAc)
-- Get the target elements. On this example, element 1 is a quad4 and element 2 a tri3
local quadElem = m:cell(1)
local triElem = m:cell(2)
-- Interpolate data on the cell centroid ({0, 0} for a quad and {1/3, 1/3, 1/3} for a tri)
-- Results are scalar values
local T1 = interpolator:interpolate(quadElem, {0, 0}, false)
local T2 = interpolator:interpolate(triElem, {1/3, 1/3, 1/3}, false)
-- Example 2: Interpolates temperature AND displacement data (vector(3)) using the moving least squares method
local uAc = m:nodeValueAccessor('u')
local interpolator = NodeToElementPointInterpolator('mls', nil, {TAc, uAc}, coordAc)
-- Interpolates data on the quad element Gauss points.
-- Results are a list with vectors, one for the temperature (dim = 1) and
-- another for the displacement (dim = 3)
local ir = m:elementIntegrationRule(quadElem:type(), 1)
for i = 1, ir:numPoints() do
local results = interpolator:interpolate(quadElem, ir, i)
local T = results[1](1) -- results[1] is a vector with dimension 1
local ux, uy = results[2](1), results[2](2) -- results[2] is a vector with dimension 3
end
-- Example 3: Same as example 2 but saving results directly to Gauss accessors
local gaussTAc = m:gaussAttributeAccessor('gaussT')
local gaussUAc = m:gaussAttributeAccessor('gaussU')
local interpolator = NodeToElementPointInterpolator('mls', nil, {TAc, uAc}, coordAc, m, {gaussTAc, gaussUAc})
for i = 1, ir:numPoints() do
interpolator:interpolate(quadElem, ir, i)
end


GaussToElementPointInterpolator(type, param, gaussAcList, coordAc, mesh)
Description: Creates a new "Gauss to Point" interpolator, capable of calculating values for a point inside an element from Gauss point values. For efficiency, multiple value sets can be interpolated at the same time.
Important: This call only creates and initializes the interpolator. The interpolation operation itself will take place when calling the interpolate() method over the returned object. This method can be called multiple times for interpolating values over different elements.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
gaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations). If a table is given, all the accessors must share the same integration rule (check the ruleSet attribute used when creating the Gauss accessor).
coordAc An accessor for retrieving mesh node coordinates (see mesh:nodeCoordAccessor()).
mesh The mesh owning the given accessors. Required.
Returns: Returns the new "Gauss to Point" interpolator object.

Example: See the examples below for the g2p:interpolate() method.

g2p:interpolate(elem, coord, cartesian, dim)
Description: Interpolates values on the given coordinate, belonging to the given element, for the data sets configured when creating the interpolator object. This method returns either a scalar number or a vector object, depending on the interpolated data dimension, if the interpolator was created with a single Gauss accessor. For objects created with multiple accessors, either a table with multiple numbers or a table with multiple vectors will be returned (the first only if all the interpolated values are scalars or if the dim parameter was given).
Parameters: elem The element object containing the interpolation point.
coord The interpolation coordinate. Can be either a Lua table or a vector object. Values can be either Cartesian or natural coordinates, depending on the interpolation method. If needed, natural coordinates are automatically converted to Cartesian coordinates, but the reverse is not done, so, if the interpolation method requires natural coordinates (like shape function based methods), the provided coordinates must be natural coordinates.
cartesian Indicates the type of values given in the coord parameter (true for Cartesian, false for natural).
dim An optional parameter that can be used when dealing with multi-dimensional data to request that the interpolation be done over this dimension only. The interpolation will then return a single scalar value or a list of scalars when dealing with multiple data sets. In the later case, the given dimension must be a valid dimension for every one of the input data sets. For matrix data, the dimension is taken in relation to its linearized representation (in column major format).
Returns: If the interpolator was created with a single Gauss accessor, returns either a scalar number or a vector object, depending on the interpolated data dimension. For objects created with multiple accessors, returns either a table with multiple numbers (all accessors are scalars or dim was used) or a table with multiple vectors (each with its own particular dimension). Also returns the number of source gauss points used to interpolate. This method can also return only nil when mls interpolation method fails (see mls method).

Example:

-- Interpolates Stress data (vector) using the shape function interpolation method
local m = modelData:mesh('myMesh')
local coordAc = m:nodeCoordAccessor()
local SAc = m:gaussAttributeAccessor('S')
local interpolator = GaussToElementPointInterpolator('shape', nil, SAc, coordAc, m)
-- Get the target elements. On this example, element 1 is a quad4 and element 2 a tri3
local quadElem = m:cell(1)
local triElem = m:cell(2)
-- Interpolate data on the cell centroid ({0, 0} for a quad and {1/3, 1/3, 1/3} for a tri)
-- Results are vector values
local S1 = interpolator:interpolate(quadElem, {0, 0}, false)
local S2 = interpolator:interpolate(triElem, {1/3, 1/3, 1/3}, false)


GaussToElementNodeInterpolator(type, param, gaussAcList, coordAc, mesh)
Description: Creates a new "Gauss to Nodes" interpolator, capable of extrapolating values from Gauss points to the element nodes. For efficiency, multiple value sets can be interpolated at the same time.
Important: This call only creates and initializes the interpolator. The interpolation operation itself will take place when calling the interpolate() method over the returned object. This method can be called multiple times for extrapolating values over different elements.
For transfering values from Gauss points to nodes over the entire mesh, consider using the elementBasedMeshGaussToNodes() function.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
gaussAcList Either a single Gauss accessor for the Gauss attribute to be extrapolated or a table with multiple accessors for extrapolating multiple data simultaneously (much more efficient than doing multiple single extrapolations). If a table is given, all the accessors must share the same integration rule (check the ruleSet attribute used when creating the Gauss accessor).
coordAc An accessor for retrieving mesh node coordinates (see mesh:nodeCoordAccessor()).
mesh The mesh owning the given accessors. Required.
Returns: Returns the new "Gauss to Nodes" interpolator object.

Example: See the examples below for the g2n:interpolate() method.

g2n:interpolate(elem, dim)
Description: Extrapolates values on the Gauss points of the given element to the element's nodes for the data sets configured when creating the interpolator object. This method returns either a matrix object or a list of matrix objects, depending if the interpolator was created with a single Gauss accessor or multiple ones. Either way, the returned matrices will have a number of lines equal to the extrapolated data dimension and a number of columns equal to the number of element nodes.
Parameters: elem The element object.
dim An optional parameter that can be used when dealing with multi-dimensional data to request that the extrapolation be done over this dimension only. The extrapolation will then return a matrix with just one line or a list of line matrices when dealing with multiple data sets. In the later case, the given dimension must be a valid dimension for every one of the input data sets. For matrix data, the dimension is taken in relation to its linearized representation (in column major format).
Returns: If the interpolator was created with a single Gauss accessor, returns a matrix object with a number of lines equal to the extrapolated data dimension and a number of columns equal to the number of element nodes. For objects created with multiple accessors, returns a table with one matrix per extrapolated data. Also returns the number of source gauss points used to extrapolate.

Examples:

-- Example 1: Extrapolates Stress data using the shape function interpolation method
local m = modelData:mesh('myMesh')
local coordAc = m:nodeCoordAccessor()
local SAc = m:gaussAttributeAccessor('S')
local interpolator = GaussToElementNodeInterpolator('shape', nil, SAc, coordAc, m)
-- Get the target element.
local e = m:cell(1)
-- Extrapolates data to the element nodes. Result is a matrix with number of lines equal to
-- the stress data dimension and columns to the number of element nodes
local nodeS = interpolator:interpolate(e)
-- Example 2: Extrapolates Stress and Strain data together using the shape function interpolation method
local eAc = m:gaussAttributeAccessor('e')
local interpolator = GaussToElementNodeInterpolator('shape', nil, {SAc, eAc}, coordAc, m)
-- Extrapolates data to the element nodes. Result is a list with two matrices, each one with a number
-- of lines equal to the stress and strain data dimension and columns to the number of element nodes
local results = interpolator:interpolate(e)
local nodeS = results[1]
local nodeE = results[2]


GaussNeighborsToNodeInterpolator(type, param, closest, gaussAcList, coordAc, mesh, activeOnly, nodeAcList)
Description: Creates a new "Gauss Neighbors to Nodes" interpolator, capable of calculating node values from the Gauss point data from neighboring elements. For efficiency, multiple value sets can be interpolated at the same time. Also, the created object can be setup for saving interpolation results to a set of node accessors instead of returning them.
Important: This call only creates and initializes the interpolator. The interpolation operation itself will take place when calling the interpolate() method over the returned object. This method can be called multiple times for interpolating values for different nodes.
For transfering values from Gauss points to nodes over the entire mesh, consider using the nodeBasedMeshGaussToNodes() function.
Requires a mesh with support for topological queries.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
closest A boolean value defining if only the closest Gauss point from each element that contains the destination point will be used on the interpolation (true - Figure 1, GN2N closest) or if all element Gauss points will be used (false - Figure 1, GN2N all). When closest is true and the interpolator is applied to an intermediate (quadratic) edge node, depending on the integration rule, there can be more than one Gauss point at the same distance from the node. In that case both of them will be used (side image).
elementClosestRule.png
gaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations). If a table is given, all the accessors must share the same integration rule (check the ruleSet attribute used when creating the Gauss accessor).
coordAc An accessor for retrieving mesh node coordinates (see mesh:nodeCoordAccessor()).
mesh The mesh owning the given accessors. Required.
activeOnly A boolean value defining if only Gauss points from active cells will be considered (true) or not (false) while interpolating.
nodeAcList An optional parameter. Either a single node value accessor or a list with multiple accessors. When given, ensures that interpolation results are automatically written to the node accessors when interpolate() is called. The given list must have exactly the same number of accessors as in the gaussAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the node list. Dimensions on the source and destination accessors must be equal. Alternativelly, destination accessors can be all scalars if the dimension parameter is used when calling interpolate().
Returns: Returns the new "Gauss Neighbors to Nodes" interpolator object.

Example: See the examples below for the gn2n:interpolate() method.

gn2n:interpolate(node, dim)
Description: Interpolates values on the given node, based on the values on Gauss points at neighbor elements, for the data sets configured when creating the interpolator object. If the object was created with a set of node accessors, results will be written directly to the accessors and won't be returned by the function. Otherwise, the function will return either a scalar number or a vector object, depending on the interpolated data dimension, if the interpolator was created with a single Gauss accessor. For objects created with multiple accessors, either a table with multiple numbers or a table with multiple vectors will be returned (the first only if all the interpolated values are scalars or if the dim parameter was given).
Parameters: node The mesh index for the node to be interpolated.
dim An optional parameter that can be used when dealing with multi-dimensional data to request that the interpolation be done over this dimension only. The interpolation will then return a single scalar value or a list of scalars when dealing with multiple data sets. In the later case, the given dimension must be a valid dimension for every one of the input data sets. For matrix data, the dimension is taken in relation to its linearized representation (in column major format).
Returns: If the object was created with node accessors, returns the number of source gauss points used to interpolate. If the interpolator was created with a single Gauss data accessor, returns either a scalar number or a vector object, depending on the interpolated data dimension, and the number of source gauss points used to interpolate. For objects created with multiple accessors, returns either a table with multiple numbers (all accessors are scalars or dim was used) or a table with multiple vectors (each with its own particular dimension) and the number of source gauss points used to interpolate. This method can also return only nil if there is no (active) cell around the given node index or when mls interpolation method fails (see mls method).

Examples:

-- Example 1: Interpolates Stress data using the shape function interpolation method
local m = modelData:mesh('myMesh')
local coordAc = m:nodeCoordAccessor()
local SAc = m:gaussAttributeAccessor('S')
local interpolator = GaussNeighborsToNodeInterpolator('shape', nil, false, SAc, coordAc, m, true)
-- Interpolates data for node 15. Result is a vector with the interpolated stress
local nodeS = interpolator:interpolate(15)
-- Example 2: Interpolates Stress and Strain data together using the moving least squares
-- interpolation method and only the closests Gauss points
local eAc = m:gaussAttributeAccessor('e')
local interpolator = GaussNeighborsToNodeInterpolator('mls', nil, true, {SAc, eAc}, coordAc, m, true)
-- Interpolates data for node 15. Result is a list with two vectors with the interpolated stress
-- and strain data
local results = interpolator:interpolate(15)
local nodeS = results[1]
local nodeE = results[2]
-- Example 3: Same as example 2 but saving results directly to node accessors
local nodeSAc = m:nodeValueAccessor('nodeS')
local nodeEAc = m:nodeValueAccessor('nodeE')
local interpolator = GaussNeighborsToNodeInterpolator('mls', nil, true, {SAc, eAc},
coordAc, m, true, {nodeSAc, nodeEAc})
-- Interpolates data for node 15. Results are saved directly to the node accessors
interpolator:interpolate(15)


NodeCloudToPointInterpolator(type, param, searchDomain, searchMode, bucketIndex, nodeAcList)
Description: Creates a new "Node Cloud to Point" interpolator, capable of calculating point values from a cloud of mesh nodes. The nodes belonging to the input cloud (its compact domain) can be defined by a radius around the interpolation point or by defining the number of clostest nodes to include. For efficiency, multiple value sets can be interpolated at the same time.
Important: This call only creates and initializes the interpolator. The interpolation operation itself will take place when calling the interpolate() method over the returned object. This method can be called multiple times for interpolating values for different points.
For transfering values from nodes over the entire mesh to nodes (on a different mesh) or Gauss points, consider using the meshNodesToMeshNodes() or the meshNodesToMeshGauss() functions.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
searchDomain A number providing either the radius of the point cloud (expressed in the same unit as the query coordinate) or the number of closest points to include on it.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest nodes. Allowed values are "radius" or "closest".
bucketIndex The spatial (bucket) index object containing the mesh nodes. Must implement the nearestNode search capability. Keep in mind that if the bucket index was created with a node set filter, only nodes in the node set will be used as sources for the interpolation.
nodeAcList Either a single node value accessor for the Node attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations).
Returns: Returns the new "Node Cloud to Point" interpolator object.

Example: See the examples below for the nc2p:interpolate() method.

nc2p:interpolate(coord, dim)
Description: Interpolates values on the given cartesian coordinate of the target point, based on the values on neighboring nodes, for the data sets configured when creating the interpolator object. The function returns either a scalar number or a vector object, depending on the interpolated data dimension, if the interpolator was created with a single node accessor. For objects created with multiple accessors, either a table with multiple numbers or a table with multiple vectors will be returned (the first only if all the interpolated values are scalars or if the dim parameter was given).
Parameters: coord The interpolation cartesian coordinate of the target point. Can be either a Lua table or a vector object.
dim An optional parameter that can be used when dealing with multi-dimensional data to request that the interpolation be done over this dimension only. The interpolation will then return a single scalar value or a list of scalars when dealing with multiple data sets. In the later case, the given dimension must be a valid dimension for every one of the input data sets. For matrix data, the dimension is taken in relation to its linearized representation (in column major format).
Returns: If the interpolator was created with a single node data accessor, returns either a scalar number or a vector object, depending on the interpolated data dimension. For objects created with multiple accessors, returns either a table with multiple numbers (all accessors are scalars or dim was used) or a table with multiple vectors (each with its own particular dimension). This method can also return nil if there are no nodes inside the given radius r around the target point or when mls interpolation method fails (see mls method). If the function does not return nil, it also returns the number of source nodes used to interpolate.

Examples:

-- Example 1: Interpolates node data at origin using the mls function interpolation method
-- and considering the 8 closest nodes around the origin
local m = modelData:mesh('myMesh')
local naAc = m:nodeValueAccessor('na')
local bucketOptions = {} -- default options will be used
local bucketIndex = mm.buildBucketIndex('myMesh', 'node', bucketOptions)
local interp1 = NodeCloudToPointInterpolator('mls', nil, 8, 'closest', bucketIndex, naAc)
local r, n = interp1:interpolate(Vector({0, 0}))
assert(n == 8)
-- Example 2: Interpolates node data at origin using the inverse distance function interpolation method
-- and considering the closest nodes inside a radius of 8.5 around the origin
local interp2 = NodeCloudToPointInterpolator('idw', 1, 8.5, 'radius', bucketIndex, naAc)
r, n = interp2:interpolate(Vector({0, 0}))


GaussCloudToPointInterpolator(type, param, searchDomain, searchMode, bucketIndex, gaussAcList)
Description: Creates a new "Gauss Cloud to Point" interpolator, capable of calculating point values from a cloud of Gauss points. The Gauss points belonging to the input cloud (its compact domain) can be defined by a radius around the interpolation point or by defining the number of clostest Gauss points to include. For efficiency, multiple value sets can be interpolated at the same time.
Important: This call only creates and initializes the interpolator. The interpolation operation itself will take place when calling the interpolate() method over the returned object. This method can be called multiple times for interpolating values for different points.
For transfering values from Gauss points over the entire mesh to nodes or Gauss points, consider using the meshGaussToMeshNodes() or the meshGaussToMeshGauss() functions.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
searchDomain A number providing either the radius of the point cloud (expressed in the same unit as the query coordinate) or the number of closest points to include on it.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest Gauss points. Allowed values are "radius" or "closest".
bucketIndex The spatial (bucket) index object containing the mesh nodes. Must implement the nearestGauss search capability.
gaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations).
Returns: Returns the new "Gauss Cloud to Point" interpolator object.

Example: See the examples below for the gc2p:interpolate() method.

gc2p:interpolate(coord, dim)
Description: Interpolates values on the given cartesian coordinate of the target point, based on the values on neighboring Gauss points, for the data sets configured when creating the interpolator object. The function returns either a scalar number or a vector object, depending on the interpolated data dimension, if the interpolator was created with a single Gauss accessor. For objects created with multiple accessors, either a table with multiple numbers or a table with multiple vectors will be returned (the first only if all the interpolated values are scalars or if the dim parameter was given).
Parameters: coord The interpolation cartesian coordinate of the target point. Can be either a Lua table or a vector object.
dim An optional parameter that can be used when dealing with multi-dimensional data to request that the interpolation be done over this dimension only. The interpolation will then return a single scalar value or a list of scalars when dealing with multiple data sets. In the later case, the given dimension must be a valid dimension for every one of the input data sets. For matrix data, the dimension is taken in relation to its linearized representation (in column major format).
Returns: If the interpolator was created with a single Gauss data accessor, returns either a scalar number or a vector object, depending on the interpolated data dimension. For objects created with multiple accessors, returns either a table with multiple numbers (all accessors are scalars or dim was used) or a table with multiple vectors (each with its own particular dimension). This method can also return nil if there are no Gauss points inside the given radius r around the target point or when mls interpolation method fails (see mls method). If the function does not return nil, it also returns the number of source Gauss points used to interpolate.

Examples:

-- Example 1: Interpolates Gauss data at origin using the mls function interpolation method
-- and considering the 8 closest Gauss points around the origin
local m = modelData:mesh('myMesh')
local gaAc = m:gaussAttributeAccessor('ga')
local bucketOptions = {} -- default options will be used
local bucketIndex = mm.buildBucketIndex('myMesh', 'gauss', bucketOptions)
local interp1 = GaussCloudToPointInterpolator('mls', nil, 8, 'closest', bucketIndex, gaAc)
local r, n = interp1:interpolate(Vector({0, 0}))
assert(n == 8)
-- Example 2: Interpolates Gauss data at origin using the inverse distance function interpolation method
-- and considering the closest Gauss points inside a radius of 8.5 around the origin
local interp2 = GaussCloudToPointInterpolator('idw', 1, 8.5, 'radius', bucketIndex, gaAc)
r, n = interp2:interpolate(Vector({0, 0}))


Global interpolation functions

elementBasedMeshGaussToNodes(type, param, avgMode, gaussAcList, nodeAcList, coordAc, gsMesh)
Description: A global function that can be used to transfer data on Gauss points to element nodes, for the entire mesh, by traversing mesh elements and extrapolating values in each element's Gauss points to their nodes. Since values calcuated for the same node by adjacent elements are in general different, an aggregation procedure is done fore averaging the contribution from each cell. Roughly equivalent to creating a "Gauss to Nodes" interpolator, applying it to every element on the mesh and averaging resulting node contributions.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
avgMode A string specifying the averaging method that will be applied to aggregate the interpolated values from each element sharing the same node. Available methods are:
mean - Node value is the arithmetic mean of the values calculated for that node by each element.
area - Node value is the arithmetic mean of the values calculated for that node by each element, weighted by the element characteristic dimension (length in 1D, area in 2D and volume in 3D).
gaussAcList Either a single Gauss accessor for the Gauss attribute to be extrapolated or a table with multiple accessors for extrapolating multiple data simultaneously (much more efficient than doing multiple single extrapolations). If a table is given, all the accessors must share the same integration rule (check the ruleSet attribute used when creating the Gauss accessor).
nodeAcList Either a single node value accessor or a list with multiple accessors for the node attributes / state variables where the results will be written. The given list must have exactly the same number of accessors as in the gaussAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the node list. Dimensions on the source and destination accessors must be equal.
coordAc An accessor for retrieving mesh node coordinates (see mesh:nodeCoordAccessor()).
gsMesh Either the mesh owning the given accessors or a cell group set object to limit the interpolation to a specified set of mesh elements.
Returns: Returns true on success, false if the given interpolation type is not supported by this operation.

Examples:

-- Interpolates Stress and Strain data together using the shape function interpolation method
-- to transfer values from Gauss points to nodes
local m = modelData:mesh('myMesh')
local coordAc = m:nodeCoordAccessor()
-- Source Gauss accessors
local SAc = m:gaussAttributeAccessor('S')
local eAc = m:gaussAttributeAccessor('e')
-- Destination node accessors
local nodeSAc = m:nodeValueAccessor('nodeS')
local nodeEAc = m:nodeValueAccessor('nodeE')
elementBasedMeshGaussToNodes('shape', nil, 'mean', {SAc, eAc}, {nodeSAc, nodeEAc}, coordAc, m)


nodeBasedMeshGaussToNodes(type, param, closest, gaussAcList, nodeAcList, coordAc, gsMesh)
Description: A global function that can be used to transfer data on Gauss points to element nodes, for the entire mesh, by traversing mesh nodes and interpolating values from data on surrounding Gauss points from element's that include the node on their definition. Roughly equivalent to creating a "Gauss Neighbors to Nodes" interpolator and applying it to every node on the mesh. Requires a mesh with support for topological queries.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
closest A boolean value defining if only the closest Gauss point from each element that contains the destination point will be used on the interpolation (true - Figure 1, GN2N closest) or if all element Gauss points will be used (false - Figure 1, GN2N all). When closest is true and the interpolator is applied to an intermediate (quadratic) edge node, depending on the integration rule, there can be more than one Gauss point at the same distance from the node. In that case both of them will be used (side image).
elementClosestRule.png
gaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations). If a table is given, all the accessors must share the same integration rule (check the ruleSet attribute used when creating the Gauss accessor).
nodeAcList Either a single node value accessor or a list with multiple accessors for the node attributes / state variables where the results will be written. The given list must have exactly the same number of accessors as in the gaussAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the node list. Dimensions on the source and destination accessors must be equal.
coordAc An accessor for retrieving mesh node coordinates (see mesh:nodeCoordAccessor()).
gsMesh Either the mesh owning the given accessors or a cell group set object to limit the interpolation to a specified set of mesh elements.
Returns: Returns true on success, false if the given interpolation type is not supported by this operation.

Examples:

-- Interpolates Stress and Strain data together using the moving least squares interpolation method
-- to transfer values from Gauss points to nodes
local m = modelData:mesh('myMesh')
local coordAc = m:nodeCoordAccessor()
-- Source Gauss accessors
local SAc = m:gaussAttributeAccessor('S')
local eAc = m:gaussAttributeAccessor('e')
-- Destination node accessors
local nodeSAc = m:nodeValueAccessor('nodeS')
local nodeEAc = m:nodeValueAccessor('nodeE')
nodeBasedMeshGaussToNodes('shape', nil, true, {SAc, eAc}, {nodeSAc, nodeEAc}, coordAc, m)


meshNodesToMeshNodes(type, param, searchDomain, searchMode, srcBucketIndex, srcNodeAcList, dstNsMesh, dstNodeAcList)
Description: A global function that can be used to transfer node values from a source mesh to node values of a distinct destination mesh using the given local interpolation method. Each destination node value is calculated independently by a local interpolation using a cloud based approach (NC2P) where the source data nodes are given by either a radius around the destination node or by a constant number of closest nodes.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
searchDomain A number providing either the radius or the number of closest nodes in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest nodes. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh nodes. Must implement the nearestNode search capability. Keep in mind that if the bucket index was created with a node set filter, only nodes in the node set will be used as sources for the interpolation.
srcNodeAcList Either a single node value accessor for the node attribute / state variable to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations).
dstNsMesh Either the mesh owning the given destination accessors or a node set object to limit the interpolation to a specified set of mesh nodes.
dstNodeAcList Either a single node value accessor or a list with multiple accessors for the node attributes / state variables where the results will be written. The given list must have exactly the same number of accessors as in the srcNodeAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal.
Returns: Returns true on success, false if the interpolation could not be performed.

Examples:

-- Transfers 'sv' node data from 'mesh1' to nodes of 'mesh2' using the idw interpolation method
-- The results are writen to 'svSplineInterp1' and 'svSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Bucket index with node references
local bucketIndex = mm.buildBucketIndex('mesh1', 'node', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source node accessor
local srcNodeAc = mesh1:nodeValueAccessor('sv')
-- Destination node accessors
local dstNodeAc1 = mesh2:nodeValueAccessor('svSplineInterp1')
local dstNodeAc2 = mesh2:nodeValueAccessor('svSplineInterp2')
meshNodesToMeshNodes('idw', 2, 7, 'closest', bucketIndex, srcNodeAc, mesh2, dstNodeAc1)
meshNodesToMeshNodes('idw', 2, 10.0, 'radius', bucketIndex, srcNodeAc, mesh2, dstNodeAc2)


meshNodesToMeshGauss(type, param, searchDomain, searchMode, srcBucketIndex, srcNodeAcList, dstGsMesh, activeOnly, dstGaussAcList)
Description: A global function that can be used to transfer node values from a source mesh to Gauss point values of a destination mesh using the given local interpolation method. Source and destination meshes might be the same. Each destination Gauss point value is calculated independently by a local interpolation using a cloud based approach (NC2P) where the source data nodes are given by either a radius around the destination Gauss point or by a constant number of closest nodes.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
searchDomain A number providing either the radius or the number of closest nodes in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest nodes. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh nodes. Must implement the nearestNode search capability. Keep in mind that if the bucket index was created with a node set filter, only nodes in the node set will be used as sources for the interpolation.
srcNodeAcList Either a single node value accessor for the node attribute / state variable to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations).
dstGsMesh Either the mesh owning the given destination accessors or a cell group set object to limit the interpolation to a specified set of mesh elements.
activeOnly If set to true, only active cells in the destination mesh (or in the destination cell group set) will have their Gauss point values calculated.
dstGaussAcList Either a single Gauss accessor or a list with multiple accessors for the Gauss attributes where the results will be written. The given list must have exactly the same number of accessors as in the srcNodeAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal. The rule set of every accessor in the list must be the same.
Returns: Returns true on success, false if the interpolation could not be performed.

Examples:

-- Transfers 'sv' node data from 'mesh1' to Gauss points of 'mesh2' using the idw interpolation method
-- The results are writen to 'gaSplineInterp1' and 'gaSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Bucket index with node references
local bucketIndex = mm.buildBucketIndex('mesh1', 'node', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source node accessor
local srcNodeAc = mesh1:nodeValueAccessor('sv')
-- Destination Gauss accessors
local dstGaussAc1 = mesh2:gaussAttributeAccessor('gaSplineInterp1')
local dstGaussAc2 = mesh2:gaussAttributeAccessor('gaSplineInterp2')
meshNodesToMeshGauss('idw', 2, 7, 'closest', bucketIndex, srcNodeAc, mesh2, false, dstGaussAc1)
meshNodesToMeshGauss('idw', 2, 10.0, 'radius', bucketIndex, srcNodeAc, mesh2, false, dstGaussAc2)


meshGaussToMeshNodes(type, param, searchDomain, searchMode, srcBucketIndex, srcGaussAcList, dstNsMesh, dstNodeAcList)
Description: A global function that can be used to transfer Gauss point values from a source mesh to node values of a destination mesh using the given local interpolation method. Source and destination meshes might be the same. Each destination node value is calculated independently by a local interpolation using a cloud based approach (GC2P) where the source Gauss points are given by either a radius around the destination node or by a constant number of closest Gauss points.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
searchDomain A number providing either the radius or the number of closest points in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest points. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh Gauss points. Must implement the nearestGauss search capability. Keep in mind that if the bucket index was created with a cell group set filter, only Gauss points in the cell group set will be used as sources for the interpolation. Also, the bucket index rule set must be compatible with the source Gauss accessor rule sets. The active only option for the bucket index also affects which source Gauss points will be used for the interpolation.
srcGaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations). The rule set of every accessor in the list must be the same and equal to the bucket index rule set.
dstNsMesh Either the mesh owning the given destination accessors or a node set object to limit the interpolation to a specified set of mesh nodes.
dstNodeAcList Either a single node value accessor or a list with multiple accessors for the node attributes / state variables where the results will be written. The given list must have exactly the same number of accessors as in the srcGaussAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal.
Returns: Returns true on success, false if the interpolation could not be performed.

Examples:

-- Transfers 'ga' Gauss data from 'mesh1' to nodes of 'mesh2' using the idw interpolation method
-- The results are writen to 'svSplineInterp1' and 'svSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Bucket index with Gauss references
local bucketIndex = mm.buildBucketIndex('mesh1', 'gauss', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source Gauss accessor
local srcGaussAc = mesh1:gaussAttributeAccessor('ga')
-- Destination node accessors
local dstNodeAc1 = mesh2:nodeValueAccessor('svSplineInterp1')
local dstNodeAc2 = mesh2:nodeValueAccessor('svSplineInterp2')
meshGaussToMeshNodes('idw', 2, 7, 'closest', bucketIndex, srcGaussAc, mesh2, dstNodeAc1)
meshGaussToMeshNodes('idw', 2, 10.0, 'radius', bucketIndex, srcGaussAc, mesh2, dstNodeAc2)


meshGaussToMeshGauss(type, param, searchDomain, searchMode, srcBucketIndex, srcGaussAcList, dstGsMesh, activeOnly, dstGaussAcList)
Description: A global function that can be used to transfer Gauss point values from a source mesh to Gauss point values of a destination mesh using the given local interpolation method. Source and destination meshes might be the same if the source and destination Gauss sets use differente integration rules. Each destination Gauss point value is calculated independently by a local interpolation using a cloud based approach (GC2P) where the source Gauss points are given by either a radius around the destination Gauss point or by a constant number of closest Gauss points.
Parameters: type A string with the name of the desired interpolation method. Accepts the method names summarized on Table 2. Check the method documentation for algorithm details. An empty string or the "default" string will select the standard algorithm for the mesh type.
param An interpolation method parameter whose type and valid values depends on the chosen interpolator. Can be nil but must be explicitly provided. See the summary on the table above and the method documentation for further details.
searchDomain A number providing either the radius or the number of closest points in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest points. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh Gauss points. Must implement the nearestGauss search capability. Keep in mind that if the bucket index was created with a cell group set filter, only Gauss points in the cell group set will be used as sources for the interpolation. Also, the bucket index rule set must be compatible with the source Gauss accessor rule sets. The active only option for the bucket index also affects which source Gauss points will be used for the interpolation.
srcGaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations). The rule set of every accessor in the list must be the same and equal to the bucket index rule set.
dstGsMesh Either the mesh owning the given destination accessors or a cell group set object to limit the interpolation to a specified set of mesh elements.
activeOnly If set to true, only active cells in the destination mesh (or in the destination cell group set) will have their Gauss point values calculated.
dstGaussAcList Either a single Gauss accessor or a list with multiple accessors for the Gauss attributes where the results will be written. The given list must have exactly the same number of accessors as in the srcGaussAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal. The rule set of every accessor in the list must be the same.
Returns: Returns true on success, false if the interpolation could not be performed.

Examples:

-- Transfers 'ga' Gauss data from 'mesh1' to Gauss points of 'mesh2' using the idw interpolation method
-- The results are writen to 'gaSplineInterp1' and 'gaSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Bucket index with Gauss references
local bucketIndex = mm.buildBucketIndex('mesh1', 'gauss', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source Gauss accessor
local srcGaussAc = mesh1:gaussAttributeAccessor('ga')
-- Destination node accessors
local dstGaussAc1 = mesh2:gaussAttributeAccessor('gaSplineInterp1')
local dstGaussAc2 = mesh2:gaussAttributeAccessor('gaSplineInterp2')
meshGaussToMeshGauss('idw', 2, 7, 'closest', bucketIndex, srcGaussAc, mesh2, false, dstGaussAc1)
meshGaussToMeshGauss('idw', 2, 10.0, 'radius', bucketIndex, srcGaussAc, mesh2, false, dstGaussAc2)


splineMeshNodesToMeshNodes(searchDomain, searchMode, srcBucketIndex, srcNodeAcList, dstNsMesh, dstNodeAcList, solver)
Description: A global function that can be used to transfer node values from a source mesh to node values of a distinct destination mesh using the spline interpolation method. Each destination node value is calculated via a spline interpolation. The CDR (compact domain radius) of the method can be both constant or defined using a constant number of closest nodes.
Parameters: searchDomain A number providing either the radius or the number of closest nodes in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest nodes. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh nodes. Must implement the nearestNode search capability. Keep in mind that if the bucket index was created with a node set filter, only nodes in the node set will be used as sources for the interpolation.
srcNodeAcList Either a single node value accessor for the node attribute / state variable to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations).
dstNsMesh Either the mesh owning the given destination accessors or a node set object to limit the interpolation to a specified set of mesh nodes.
dstNodeAcList Either a single node value accessor or a list with multiple accessors for the node attributes / state variables where the results will be written. The given list must have exactly the same number of accessors as in the srcNodeAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal.
solver A string with the name of the desired numerical solver used to compute the spline radial coefficients.
Returns: Returns true on success, false if the interpolation could not be performed (solver error, for example).

Examples:

-- Transfers 'sv' node data from 'mesh1' to nodes of 'mesh2' using the spline interpolation method
-- The results are writen to 'svSplineInterp1' and 'svSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Expects that a NumericalSolver named 'solver' was included in the solution definition. A minimal solver definition could be:
--
-- NumericalSolver
-- {
-- id = 'solver',
-- typeName = 'ArmadilloSolver',
-- }
-- Bucket index with node references
local bucketIndex = mm.buildBucketIndex('mesh1', 'node', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source node accessor
local srcNodeAc = mesh1:nodeValueAccessor('sv')
-- Destination node accessors
local dstNodeAc1 = mesh2:nodeValueAccessor('svSplineInterp1')
local dstNodeAc2 = mesh2:nodeValueAccessor('svSplineInterp2')
splineMeshNodesToMeshNodes( 7, 'closest', bucketIndex, srcNodeAc, mesh2, dstNodeAc1, 'solver')
splineMeshNodesToMeshNodes(10.0, 'radius', bucketIndex, srcNodeAc, mesh2, dstNodeAc2, 'solver')


splineMeshNodesToMeshGauss(searchDomain, searchMode, srcBucketIndex, srcNodeAcList, dstGsMesh, activeOnly, dstGaussAcList, solver)
Description: A global function that can be used to transfer node values from a source mesh to Gauss point values of a destination mesh using the spline interpolation method. Source and destination meshes might be the same. Each destination Gauss point value is calculated via a spline interpolation. The CDR (compact domain radius) of the method can be both constant or defined using a constant number of closest nodes.
Parameters: searchDomain A number providing either the radius or the number of closest nodes in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest nodes. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh nodes. Must implement the nearestNode search capability. Keep in mind that if the bucket index was created with a node set filter, only nodes in the node set will be used as sources for the interpolation.
srcNodeAcList Either a single node value accessor for the node attribute / state variable to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations).
dstGsMesh Either the mesh owning the given destination accessors or a cell group set object to limit the interpolation to a specified set of mesh elements.
activeOnly If set to true, only active cells in the destination mesh (or in the destination cell group set) will have their Gauss point values calculated.
dstGaussAcList Either a single Gauss accessor or a list with multiple accessors for the Gauss attributes where the results will be written. The given list must have exactly the same number of accessors as in the srcNodeAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal. The rule set of every accessor in the list must be the same.
solver A string with the name of the desired numerical solver used to compute the spline radial coefficients.
Returns: Returns true on success, false if the interpolation could not be performed (solver error, for example).

Examples:

-- Transfers 'sv' node data from 'mesh1' to Gauss points of 'mesh2' using the spline interpolation method
-- The results are writen to 'gaSplineInterp1' and 'gaSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Expects that a NumericalSolver named 'solver' was included in the solution definition. A minimal solver definition could be:
--
-- NumericalSolver
-- {
-- id = 'solver',
-- typeName = 'ArmadilloSolver',
-- }
-- Bucket index with node references
local bucketIndex = mm.buildBucketIndex('mesh1', 'node', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source node accessor
local srcNodeAc = mesh1:nodeValueAccessor('sv')
-- Destination Gauss accessors
local dstGaussAc1 = mesh2:gaussAttributeAccessor('gaSplineInterp1')
local dstGaussAc2 = mesh2:gaussAttributeAccessor('gaSplineInterp2')
splineMeshNodesToMeshGauss( 7, 'closest', bucketIndex, srcNodeAc, mesh2, false, dstGaussAc1, 'solver')
splineMeshNodesToMeshGauss(10.0, 'radius', bucketIndex, srcNodeAc, mesh2, false, dstGaussAc2, 'solver')


splineMeshGaussToMeshNodes(searchDomain, searchMode, srcBucketIndex, srcGaussAcList, dstNsMesh, dstNodeAcList, solver)
Description: A global function that can be used to transfer Gauss point values from a source mesh to node values of a destination mesh using the spline interpolation method. Source and destination meshes might be the same. Each destination node value is calculated via a spline interpolation. The CDR (compact domain radius) of the method can be both constant or defined using a constant number of closest points.
Parameters: searchDomain A number providing either the radius or the number of closest points in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest points. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh Gauss points. Must implement the nearestGauss search capability. Keep in mind that if the bucket index was created with a cell group set filter, only Gauss points in the cell group set will be used as sources for the interpolation. Also, the bucket index rule set must be compatible with the source Gauss accessor rule sets. The active only option for the bucket index also affects which source Gauss points will be used for the interpolation.
srcGaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations). The rule set of every accessor in the list must be the same and equal to the bucket index rule set.
dstNsMesh Either the mesh owning the given destination accessors or a node set object to limit the interpolation to a specified set of mesh nodes.
dstNodeAcList Either a single node value accessor or a list with multiple accessors for the node attributes / state variables where the results will be written. The given list must have exactly the same number of accessors as in the srcGaussAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal.
solver A string with the name of the desired numerical solver used to compute the spline radial coefficients.
Returns: Returns true on success, false if the interpolation could not be performed (solver error, for example).

Examples:

-- Transfers 'ga' Gauss data from 'mesh1' to nodes of 'mesh2' using the spline interpolation method
-- The results are writen to 'svSplineInterp1' and 'svSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Expects that a NumericalSolver named 'solver' was included in the solution definition. A minimal solver definition could be:
--
-- NumericalSolver
-- {
-- id = 'solver',
-- typeName = 'ArmadilloSolver',
-- }
-- Bucket index with Gauss references
local bucketIndex = mm.buildBucketIndex('mesh1', 'gauss', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source Gauss accessor
local srcGaussAc = mesh1:gaussAttributeAccessor('ga')
-- Destination node accessors
local dstNodeAc1 = mesh2:nodeValueAccessor('svSplineInterp1')
local dstNodeAc2 = mesh2:nodeValueAccessor('svSplineInterp2')
splineMeshGaussToMeshNodes( 7, 'closest', bucketIndex, srcGaussAc, mesh2, dstNodeAc1, 'solver')
splineMeshGaussToMeshNodes(10.0, 'radius', bucketIndex, srcGaussAc, mesh2, dstNodeAc2, 'solver')


splineMeshGaussToMeshGauss(searchDomain, searchMode, srcBucketIndex, srcGaussAcList, dstGsMesh, activeOnly, dstGaussAcList, solver)
Description: A global function that can be used to transfer Gauss point values from a source mesh to Gauss point values of a destination mesh using the spline interpolation method. Source and destination meshes might be the same if the source and destination Gauss sets use differente integration rules. Each destination Gauss point value is calculated via a spline interpolation. The CDR (compact domain radius) of the method can be both constant or defined using a constant number of closest points.
Parameters: searchDomain A number providing either the radius or the number of closest points in the method's CDR.
searchMode A string defining if the searchDomain parameter contains a radius or a number of closest points. Allowed values are "radius" or "closest".
srcBucketIndex The spatial (bucket) index object containing the mesh Gauss points. Must implement the nearestGauss search capability. Keep in mind that if the bucket index was created with a cell group set filter, only Gauss points in the cell group set will be used as sources for the interpolation. Also, the bucket index rule set must be compatible with the source Gauss accessor rule sets. The active only option for the bucket index also affects which source Gauss points will be used for the interpolation.
srcGaussAcList Either a single Gauss accessor for the Gauss attribute to be interpolated or a table with multiple accessors for interpolating multiple data simultaneously (much more efficient than doing multiple single interpolations). The rule set of every accessor in the list must be the same and equal to the bucket index rule set.
dstGsMesh Either the mesh owning the given destination accessors or a cell group set object to limit the interpolation to a specified set of mesh elements.
activeOnly If set to true, only active cells in the destination mesh (or in the destination cell group set) will have their Gauss point values calculated.
dstGaussAcList Either a single Gauss accessor or a list with multiple accessors for the Gauss attributes where the results will be written. The given list must have exactly the same number of accessors as in the srcGaussAcList and ordered such as the interpolation result from the i'th accessor on the source list will be saved to the i'th accessor on the destination list. Dimensions on the source and destination accessors must be equal. The rule set of every accessor in the list must be the same.
solver A string with the name of the desired numerical solver used to compute the spline radial coefficients.
Returns: Returns true on success, false if the interpolation could not be performed (solver error, for example).

Examples:

-- Transfers 'ga' Gauss data from 'mesh1' to Gauss points of 'mesh2' using the spline interpolation method
-- The results are writen to 'gaSplineInterp1' and 'gaSplineInterp2', respectively considering 'closest' and 'radius' search modes
-- Expects that a NumericalSolver named 'solver' was included in the solution definition. A minimal solver definition could be:
--
-- NumericalSolver
-- {
-- id = 'solver',
-- typeName = 'ArmadilloSolver',
-- }
-- Bucket index with Gauss references
local bucketIndex = mm.buildBucketIndex('mesh1', 'gauss', {})
local mesh1 = modelData:mesh('mesh1')
local mesh2 = modelData:mesh('mesh2')
-- Source Gauss accessor
local srcGaussAc = mesh1:gaussAttributeAccessor('ga')
-- Destination node accessors
local dstGaussAc1 = mesh2:gaussAttributeAccessor('gaSplineInterp1')
local dstGaussAc2 = mesh2:gaussAttributeAccessor('gaSplineInterp2')
splineMeshGaussToMeshGauss( 7, 'closest', bucketIndex, srcGaussAc, mesh2, false, dstGaussAc1, 'solver')
splineMeshGaussToMeshGauss(10.0, 'radius', bucketIndex, srcGaussAc, mesh2, false, dstGaussAc2, 'solver')