Interpolation¶
Description¶
This treatment interpolates a multiblock source CFD domain on another target CFD domain.
Two interpolation methods are available.
The first one is the one by default:
1st Method: Inverse Distance Weighting¶
The first method (method= idw) is the Inverse Distance Weighting one.
The complexity is N*logN thanks to the use of a fast spatial search structure (kd-tree).
The set of N known data points: \([ \ldots, (x_i, y_i), \ldots]\)
Interpolated value u at a given point x based on samples \(\displaystyle u_{i}=u(x_{i})\):
\(\displaystyle u(x) = \frac{\sum_1^N \frac{1}{d(x,x_i)^p} u_i} {\sum_1^N \frac{1}{d(x,x_i)^p}}\) if \(\displaystyle d(x,x_i) \neq 0 \forall i\)
\(\displaystyle u(x) = u_i\) if \(\displaystyle \exists i, d(x,x_i) = 0\)
x is the arbitrary interpolated point.
\(\displaystyle x_i\) is the known interpolating point.
\(\displaystyle d(x,x_i)\) is the distance from the known point
\(\displaystyle x_i\) to the unknown point x.
N is the total number of known points used in interpolation.
p is the power parameter.
2nd Method: linear¶
The second method(method= linear) is the linear one.
The complexity is N*logN too thanks to the use of a fast spatial search structure (kd-tree). But the computation is between five and a hundred times slower than idw because of more computations.
The set of N known data points: \([ \ldots, (x_i, y_i), \ldots]\)
Interpolated value u at a given point x based on samples \(\displaystyle u_{i}=u(x_{i})\):
\(\displaystyle u(x) = \sum_1^N w_i u_i\) where \(\displaystyle w_i\) is the weight of the point \(\displaystyle x_i\).
The weight \(\displaystyle w_i\) is calculated using the distance, area and volume in 1D, 2D and 3D respectively of the sub-forms created by the point in the cell in which it is located.
Parameters¶
- source:
Base
The source mesh with the variables you want to interpolate. If you want to use the linear method, the source_base have to be unstructure before. If it’s not, you can use
source_base.unstructure()
.
- source:
- target:
Base
The mesh to interpolate on. The base is modified in-place at the output.
- target:
- coordinates: list(str), default= None
The variable names that define the set of coordinates used for interpolation. The KDTree extracts closest points based on these variables. If None, use the coordinates from the source or target base.
- variables: list(str), default= None
List of variables to apply for the interpolation. If None, all variables will be interpolated, except the coordinates. Using this parameter allow to call the interpolation treatment on bases with different variable locations.
- tolerance: float, default= 1e-10
The value that tells if the distance from the closest point is smaller than this tolerance, the interpolated value is equal to the value of that closest point.
- tolerance_edge: float, default= 1e-10
It’s used to include a point in a cell if it’s on a vertex or an edge.
- duplicated_tolerance: int, default= None
Number of decimal places to round to for duplicated points detection. If decimals is negative, it specifies the number of positions to the left of the decimal point. If None, the detection will be exact.
- invdist_power: int, default= 1
The power of the inverse distance, this acts like a norm (1, 2, or infinity etc.) and can be used to smooth/refined the interpolation.
- nb_points: int, default= None
The number of closest points used for the interpolation. If not provided by the user, it will automatically be set depending on the type of cells in the mesh. For idw: Is the number of ‘close’ points used to calculate a weight. For linear: It’s the maximum number of ‘close’ points that will be used to find out in which cell the target point is located. Decrease this number to reduce the computation time but increase the risk of error in the interpolation. Increase this number to decrease the risk of missing values but increase the computation time.
- with_boundaries: bool, default= False
Whether or not to use data from the boundaries.
- method: str, default= idw
The interpolation method to use. Can be ‘idw’ or ‘linear’.
- data: tuple, default= None
These are the data, notably the ‘weights’, required to calculate the interpolation for a given mesh. This argument avoids having to recalculate the weights for a given mesh with a different solution. To retrieve this data for the first calculation, use ‘get_data’ = True. Note that this tuple is specific to the interpolation method.
- get_data: tuple, default= None
If True, treatment.execute() will return ‘target_base, data’ instead of just ‘target_base’. This data will be useful to compute another solution in the same source_base and target_base.
filter: callable, default= None
A function that returns a list of boolean values signaling where the interpolation should be performed on the target point or not. The function must accept the following variables: 1. The source base 2. The target base 3. The target zone to be interpolated 4. The target instant to be interpolated 5. A custom dictionary to be used freely by the user
The function must return a boolean array of the same shape as the target instant. If the value at the index i is True, the interpolation is performed on that point, otherwise we use the default_value
- default_value: float, default= np.nan:
The value assigned to points where the interpolation is not performed.
- create_target_instant: bool, default= False
If True, the treatment will automatically create the instants in the target base using the instants in the source base.
Returns¶
Base
The target base containing the results.
Preconditions¶
The source and target bases must have the same number of coordinates.
If the method is linear: - The source base must be unstructured; - The source base must be monozone; - The cell types that can be treat are ‘bi’, ‘tri’, ‘qua’, ‘tet’, ‘hex’, ‘pyr’.
Postconditions¶
The target base is modified in-place.
Example¶
import antares
myt = antares.Treatment('interpolation')
myt['source'] = source_base
myt['target'] = target_base
myt.execute()
Example¶
"""
This example illustrates the interpolation of a
cloud of points against another cloud of points
"""
from copy import copy
import os
if not os.path.isdir('OUTPUT'):
os.makedirs('OUTPUT')
from antares import Reader, Treatment, Writer
# ------------------
# Reading the files
# ------------------
reader = Reader('bin_tp')
reader['filename'] = os.path.join('..', 'data', 'ROTOR37', 'GENERIC', 'flow_<zone>_<instant>.dat')
source_base = reader.read()
print(source_base)
# --------------------
# Build a target base
# --------------------
# (same mesh but at cell centers)
tmp_base = copy(source_base[:, :, ('x', 'y', 'z')])
tmp_base.node_to_cell()
target_base = tmp_base.get_location('cell')
# --------------
# Interpolation
# --------------
treatment = Treatment('interpolation')
treatment['source'] = source_base
treatment['target'] = target_base
result = treatment.execute()
# Note that the result contains the same instants as the source and
# each instant contains the variable interpolated and the coordinates
# used for interpolation
print(result[0])
# >>> Zone object
# - instants : ['ite0', 'ite1', 'ite2']
print(result[0][0])
# >>> Instant object
# - variables : [('x', 'cell'), ('y', 'cell'), ('z', 'cell'), ('ro', 'cell'), ('rovx', 'cell'), ('rovy', 'cell'), ('rovz', 'cell'), ('roE', 'cell')]
# - shape : None
# -------------------
# Writing the result
# -------------------
writer = Writer('bin_tp')
writer['base'] = result
writer['filename'] = os.path.join('OUTPUT', 'ex_interpolation.plt')
writer.dump()