EarthBox.ModelDataContainer.InterpolationContainer.Interpolation — Type
Interpolation <: CollectionContainerData structure containing parameter and array objects for marker-grid interpolation.
Fields
parameters::Parameters: Parameter groups for interpolation optionsarrays::Arrays: Array groups for grid and marker interpolation weights
Constructor
Interpolation(ynum::Int, xnum::Int, marknum::Int)::InterpolationCreate a new Interpolation collection with the given grid and marker dimensions.
Arguments
ynum::Int: Number of grid points in y-directionxnum::Int: Number of grid points in x-directionmarknum::Int: Number of markers in the model
EarthBox.ModelDataContainer.InterpolationContainer.ParameterCollection.Parameters — Type
Parameters <: AbstractParameterCollectionParameter collection for marker-grid interpolation.
Fields
interp_options::MarkerAndGridInterpOptions: Interpolation method options
Constructor
Parameters()Create a new Parameters collection with default interpolation options.
EarthBox.ModelDataContainer.InterpolationContainer.ParameterCollection.MarkerAndGridInterpOptionsGroup.MarkerAndGridInterpOptions — Type
MarkerAndGridInterpOptions <: AbstractParameterGroupParameter group for marker-grid interpolation options.
Fields
iuse_initial_temp_interp::ParameterInt: Interpolate nodal temperatures back to markers at start: 0 off; 1 oniuse_harmonic_avg_normal_viscosity::ParameterInt: Harmonic averaging of shear viscosity for normal viscosity: 0 off; 1 onobj_list::Vector{ParameterInt}: List of parameter objects
Nested Dot Access
iuse_initial_temp_interp = model.interpolation.parameters.interp_options.iuse_initial_temp_interp.valueiuse_harmonic_avg_normal_viscosity = model.interpolation.parameters.interp_options.iuse_harmonic_avg_normal_viscosity.value
Constructor
MarkerAndGridInterpOptions()Create a new MarkerAndGridInterpOptions parameter group with default values.
Returns
MarkerAndGridInterpOptions: New parameter group with initialized values
EarthBox.ModelDataContainer.InterpolationContainer.ArrayCollection.Arrays — Type
Arrays <: AbstractArrayCollectionData structure containing array groups for marker-grid interpolation weights.
Fields
grid_weights::GridWeights: Summed weight arrays on grid nodesmarker_weights::MarkerWeights: Individual marker weight arrays
Constructor
Arrays(ynum::Int, xnum::Int, marknum::Int)::ArraysCreate a new Arrays collection with the given grid and marker dimensions.
Arguments
ynum::Int: Number of grid points in y-directionxnum::Int: Number of grid points in x-directionmarknum::Int: Number of markers in the model
EarthBox.ModelDataContainer.InterpolationContainer.ArrayCollection.GridWeightsGroup.GridWeights — Type
GridWeights <: AbstractArrayGroupArray group containing summed interpolation weight arrays on grid nodes.
Fields
wtnodes::ScalarArray2DState:(ynum, xnum): Summed marker weight factors on basic grid nodes using inclusive method with maximum search radius for markers.wtetas::ScalarArray2DState:(ynum, xnum): Summed marker weight factors on basic grid nodes using exclusive method with shorter search radius for markers.wtetan::ScalarArray2DState:(ynum - 1, xnum - 1): Summed marker weight factors for normal viscosity on pressure nodes.wtnodes_vy::ScalarArray2DState:(ynum + 1, xnum): Summed marker weight factors on Vy staggered grid nodes using inclusive method with maximum search radius.
Nested Dot Access
wtnodes = model.interpolation.arrays.grid_weights.wtnodes.arraywtetas = model.interpolation.arrays.grid_weights.wtetas.arraywtetan = model.interpolation.arrays.grid_weights.wtetan.arraywtnodes_vy = model.interpolation.arrays.grid_weights.wtnodes_vy.array
Constructor
GridWeights(ynum::Int, xnum::Int)Arguments
ynum::Int: Number of grid points in y-directionxnum::Int: Number of grid points in x-direction
Returns
GridWeights: New GridWeights array group with initialized arrays
EarthBox.ModelDataContainer.InterpolationContainer.ArrayCollection.MarkerWeightsGroup.MarkerWeights — Type
MarkerWeights <: AbstractArrayGroupArray group containing individual marker weight arrays for interpolation.
Fields
marker_wtforULnode::MarkerArrayFloat1DState:(marknum): Marker weight for upper-left basic grid node.marker_wtforLLnode::MarkerArrayFloat1DState:(marknum): Marker weight for lower-left basic grid node.marker_wtforURnode::MarkerArrayFloat1DState:(marknum): Marker weight for upper-right basic grid node.marker_wtforLRnode::MarkerArrayFloat1DState:(marknum): Marker weight for lower-right basic grid node.marker_wtforCnode::MarkerArrayFloat1DState:(marknum): Marker weight for central pressure grid node.marker_wtforULnodeVy::MarkerArrayFloat1DState:(marknum): Marker weight for upper-left Vy staggered grid node.marker_wtforLLnodeVy::MarkerArrayFloat1DState:(marknum): Marker weight for lower-left Vy staggered grid node.marker_wtforURnodeVy::MarkerArrayFloat1DState:(marknum): Marker weight for upper-right Vy staggered grid node.marker_wtforLRnodeVy::MarkerArrayFloat1DState:(marknum): Marker weight for lower-right Vy staggered grid node.
Nested Dot Access
marker_wtforULnode = model.interpolation.arrays.marker_weights.marker_wtforULnode.arraymarker_wtforLLnode = model.interpolation.arrays.marker_weights.marker_wtforLLnode.arraymarker_wtforURnode = model.interpolation.arrays.marker_weights.marker_wtforURnode.arraymarker_wtforLRnode = model.interpolation.arrays.marker_weights.marker_wtforLRnode.arraymarker_wtforCnode = model.interpolation.arrays.marker_weights.marker_wtforCnode.arraymarker_wtforULnodeVy = model.interpolation.arrays.marker_weights.marker_wtforULnodeVy.arraymarker_wtforLLnodeVy = model.interpolation.arrays.marker_weights.marker_wtforLLnodeVy.arraymarker_wtforURnodeVy = model.interpolation.arrays.marker_weights.marker_wtforURnodeVy.arraymarker_wtforLRnodeVy = model.interpolation.arrays.marker_weights.marker_wtforLRnodeVy.array
Constructor
MarkerWeights(marknum::Int)Arguments
marknum::Int: Number of markers in the model
Returns
MarkerWeights: New MarkerWeights array group with initialized arrays