Interpolation

TrixiParticles.interpolate_lineMethod
interpolate_line(start, end_, n_points, semi, ref_system, sol; endpoint=true,
                 smoothing_length=ref_system.smoothing_length, cut_off_bnd=true,
                 clip_negative_pressure=false)

Interpolates properties along a line in a TrixiParticles simulation. The line interpolation is accomplished by generating a series of evenly spaced points between start and end_. If endpoint is false, the line is interpolated between the start and end points, but does not include these points.

See also: interpolate_point, interpolate_plane_2d, interpolate_plane_2d_vtk, interpolate_plane_3d.

Arguments

  • start: The starting point of the line.
  • end_: The ending point of the line.
  • n_points: The number of points to interpolate along the line.
  • semi: The semidiscretization used for the simulation.
  • ref_system: The reference system for the interpolation.
  • sol: The solution state from which the properties are interpolated.

Keywords

  • endpoint=true: A boolean to include (true) or exclude (false) the end point in the interpolation.
  • smoothing_length=ref_system.smoothing_length: The smoothing length used in the interpolation.
  • cut_off_bnd=true: Boolean to indicate if quantities should be set to NaN when the point is "closer" to the boundary than to the fluid in a kernel-weighted sense. Or, in more detail, when the boundary has more influence than the fluid on the density summation in this point, i.e., when the boundary particles add more kernel-weighted mass than the fluid particles.
  • clip_negative_pressure=false: One common approach in SPH models is to clip negative pressure values, but this is unphysical. Instead we clip here during interpolation thus only impacting the local interpolated value.

Returns

  • A NamedTuple of arrays containing interpolated properties at each point along the line.
Note
  • This function is particularly useful for analyzing gradients or creating visualizations along a specified line in the SPH simulation domain.
  • The interpolation accuracy is subject to the density of particles and the chosen smoothing length.
  • With cut_off_bnd, a density-based estimation of the surface is used which is not as accurate as a real surface reconstruction.

Examples

# Interpolating along a line from [1.0, 0.0] to [1.0, 1.0] with 5 points
results = interpolate_line([1.0, 0.0], [1.0, 1.0], 5, semi, ref_system, sol)
source
TrixiParticles.interpolate_plane_2dMethod
interpolate_plane_2d(min_corner, max_corner, resolution, semi, ref_system, sol;
                     smoothing_length=ref_system.smoothing_length, cut_off_bnd=true,
                     clip_negative_pressure=false)

Interpolates properties along a plane in a TrixiParticles simulation. The region for interpolation is defined by its lower left and top right corners, with a specified resolution determining the density of the interpolation points.

The function generates a grid of points within the defined region, spaced uniformly according to the given resolution.

See also: interpolate_plane_2d_vtk, interpolate_plane_3d, interpolate_line, interpolate_point.

Arguments

  • min_corner: The lower left corner of the interpolation region.
  • max_corner: The top right corner of the interpolation region.
  • resolution: The distance between adjacent interpolation points in the grid.
  • semi: The semidiscretization used for the simulation.
  • ref_system: The reference system for the interpolation.
  • sol: The solution state from which the properties are interpolated.

Keywords

  • smoothing_length=ref_system.smoothing_length: The smoothing length used in the interpolation.
  • cut_off_bnd=true: Boolean to indicate if quantities should be set to NaN when the point is "closer" to the boundary than to the fluid in a kernel-weighted sense. Or, in more detail, when the boundary has more influence than the fluid on the density summation in this point, i.e., when the boundary particles add more kernel-weighted mass than the fluid particles.
  • clip_negative_pressure=false: One common approach in SPH models is to clip negative pressure values, but this is unphysical. Instead we clip here during interpolation thus only impacting the local interpolated value.

Returns

  • A NamedTuple of arrays containing interpolated properties at each point within the plane.
Note
  • The interpolation accuracy is subject to the density of particles and the chosen smoothing length.
  • With cut_off_bnd, a density-based estimation of the surface is used, which is not as accurate as a real surface reconstruction.

Examples

# Interpolating across a plane from [0.0, 0.0] to [1.0, 1.0] with a resolution of 0.2
results = interpolate_plane_2d([0.0, 0.0], [1.0, 1.0], 0.2, semi, ref_system, sol)
source
TrixiParticles.interpolate_plane_2d_vtkMethod
interpolate_plane_2d_vtk(min_corner, max_corner, resolution, semi, ref_system, sol;
                         smoothing_length=ref_system.smoothing_length, cut_off_bnd=true,
                         clip_negative_pressure=false, output_directory="out", filename="plane")

Interpolates properties along a plane in a TrixiParticles simulation and exports the result as a VTI file. The region for interpolation is defined by its lower left and top right corners, with a specified resolution determining the density of the interpolation points.

The function generates a grid of points within the defined region, spaced uniformly according to the given resolution.

See also: interpolate_plane_2d, interpolate_plane_3d, interpolate_line, interpolate_point.

Arguments

  • min_corner: The lower left corner of the interpolation region.
  • max_corner: The top right corner of the interpolation region.
  • resolution: The distance between adjacent interpolation points in the grid.
  • semi: The semidiscretization used for the simulation.
  • ref_system: The reference system for the interpolation.
  • sol: The solution state from which the properties are interpolated.

Keywords

  • smoothing_length=ref_system.smoothing_length: The smoothing length used in the interpolation.
  • output_directory="out": Directory to save the VTI file.
  • filename="plane": Name of the VTI file.
  • cut_off_bnd=true: Boolean to indicate if quantities should be set to NaN when the point is "closer" to the boundary than to the fluid in a kernel-weighted sense. Or, in more detail, when the boundary has more influence than the fluid on the density summation in this point, i.e., when the boundary particles add more kernel-weighted mass than the fluid particles.
  • clip_negative_pressure=false: One common approach in SPH models is to clip negative pressure values, but this is unphysical. Instead we clip here during interpolation thus only impacting the local interpolated value.
Note
  • The interpolation accuracy is subject to the density of particles and the chosen smoothing length.
  • With cut_off_bnd, a density-based estimation of the surface is used, which is not as accurate as a real surface reconstruction.

Examples

# Interpolating across a plane from [0.0, 0.0] to [1.0, 1.0] with a resolution of 0.2
results = interpolate_plane_2d([0.0, 0.0], [1.0, 1.0], 0.2, semi, ref_system, sol)
source
TrixiParticles.interpolate_plane_3dMethod
interpolate_plane_3d(point1, point2, point3, resolution, semi, ref_system, sol;
                     smoothing_length=ref_system.smoothing_length, cut_off_bnd=true,
                     clip_negative_pressure=false)

Interpolates properties along a plane in a 3D space in a TrixiParticles simulation. The plane for interpolation is defined by three points in 3D space, with a specified resolution determining the density of the interpolation points.

The function generates a grid of points on a parallelogram within the plane defined by the three points, spaced uniformly according to the given resolution.

See also: interpolate_plane_2d, interpolate_plane_2d_vtk, interpolate_line, interpolate_point.

Arguments

  • point1: The first point defining the plane.
  • point2: The second point defining the plane.
  • point3: The third point defining the plane. The points must not be collinear.
  • resolution: The distance between adjacent interpolation points in the grid.
  • semi: The semidiscretization used for the simulation.
  • ref_system: The reference system for the interpolation.
  • sol: The solution state from which the properties are interpolated.

Keywords

  • smoothing_length=ref_system.smoothing_length: The smoothing length used in the interpolation.
  • cut_off_bnd=true: Boolean to indicate if quantities should be set to NaN when the point is "closer" to the boundary than to the fluid in a kernel-weighted sense. Or, in more detail, when the boundary has more influence than the fluid on the density summation in this point, i.e., when the boundary particles add more kernel-weighted mass than the fluid particles.
  • clip_negative_pressure=false: One common approach in SPH models is to clip negative pressure values, but this is unphysical. Instead we clip here during interpolation thus only impacting the local interpolated value.

Returns

  • A NamedTuple of arrays containing interpolated properties at each point within the plane.
Note
  • The interpolation accuracy is subject to the density of particles and the chosen smoothing length.
  • With cut_off_bnd, a density-based estimation of the surface is used which is not as accurate as a real surface reconstruction.

Examples

# Interpolating across a plane defined by points [0.0, 0.0, 0.0], [1.0, 0.0, 0.0], and [0.0, 1.0, 0.0]
# with a resolution of 0.1
results = interpolate_plane_3d([0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.0, 1.0, 0.0], 0.1, semi, ref_system, sol)
source
TrixiParticles.interpolate_pointMethod
interpolate_point(points_coords::Array{Array{Float64,1},1}, semi, ref_system, sol;
                  smoothing_length=ref_system.smoothing_length, cut_off_bnd=true,
                  clip_negative_pressure=false)

interpolate_point(point_coords, semi, ref_system, sol;
                  smoothing_length=ref_system.smoothing_length, cut_off_bnd=true,
                  clip_negative_pressure=false)

Performs interpolation of properties at specified points or an array of points in a TrixiParticles simulation.

When given an array of points (points_coords), it iterates over each point and applies interpolation individually. For a single point (point_coords), it performs the interpolation at that specific location. The interpolation utilizes the same kernel function of the SPH simulation to weigh contributions from nearby particles.

See also: interpolate_line, interpolate_plane_2d, interpolate_plane_2d_vtk, interpolate_plane_3d, .

Arguments

  • points_coords: An array of point coordinates, for which to interpolate properties.
  • point_coords: The coordinates of a single point for interpolation.
  • semi: The semidiscretization used in the SPH simulation.
  • ref_system: The reference system defining the properties of the SPH particles.
  • sol: The current solution state from which properties are interpolated.

Keywords

  • smoothing_length=ref_system.smoothing_length: The smoothing length used in the interpolation.
  • cut_off_bnd=true: Boolean to indicate if quantities should be set to NaN when the point is "closer" to the boundary than to the fluid in a kernel-weighted sense. Or, in more detail, when the boundary has more influence than the fluid on the density summation in this point, i.e., when the boundary particles add more kernel-weighted mass than the fluid particles.
  • clip_negative_pressure=false: One common approach in SPH models is to clip negative pressure values, but this is unphysical. Instead we clip here during interpolation thus only impacting the local interpolated value.

Returns

  • For multiple points: A NamedTuple of arrays containing interpolated properties at each point.
  • For a single point: A NamedTuple of interpolated properties at the point.

Examples

# For a single point
result = interpolate_point([1.0, 0.5], semi, ref_system, sol)

# For multiple points
points = [[1.0, 0.5], [1.0, 0.6], [1.0, 0.7]]
results = interpolate_point(points, semi, ref_system, sol)
Note
  • This function is particularly useful for analyzing gradients or creating visualizations along a specified line in the SPH simulation domain.
  • The interpolation accuracy is subject to the density of particles and the chosen smoothing length.
  • With cut_off_bnd, a density-based estimation of the surface is used which is not as

accurate as a real surface reconstruction.

source