AXOM
Axom provides a robust, flexible software infrastructure for the development of multi-physics applications and computational tools.
axom::quest Namespace Reference

Namespaces

 detail
 
 internal
 
 shaping
 

Classes

class  Delaunay
 A class for incremental generation of a 2D or 3D Delaunay triangulation. More...
 
class  DistributedClosestPoint
 Encapsulated the Distributed closest point query for a collection of query points over an "object mesh". More...
 
class  InOutOctree
 
class  GridFunctionView
 Provides a view over an MFEM grid function. MFEM grid functions are assumed to live in host memory. This class performs data movement needed to access the grid function data within a GPU device lambda. This view is limited in scope, though could be expanded in the future. More...
 
class  IntersectionShaper
 
class  PointInCellTraits
 A traits class for the mesh associated with a PointInCell query. More...
 
class  PointInCell
 A class to accelerate Point-In-Cell queries on a computational mesh. More...
 
class  SamplingShaper
 Concrete class for sample based shaping. More...
 
class  ScatteredInterpolation
 A class to perform scattered data interpolation at arbitrary points over an input point set. More...
 
class  Shaper
 
class  SignedDistance
 

Enumerations

enum  SearchStatus { NEIGHBOR_NOT_FOUND = -1 }
 
enum class  WatertightStatus : signed char { WATERTIGHT = 0 , NOT_WATERTIGHT , CHECK_FAILED }
 
enum class  SignedDistExec { CPU = 0 , OpenMP = 1 , GPU = 2 }
 

Functions

Nearest Neighbor query
void all_nearest_neighbors (const double *x, const double *y, const double *z, const int *region, int n, double limit, int *neighbor, double *sqdistance)
 Given a list of point locations and regions, for each point, find the closest point in a different region within a given search radius. More...
 
Visualize octahedra as tet mesh
int mesh_from_discretized_polyline (axom::ArrayView< OctType > &octs, int octcount, int segcount, mint::Mesh *&mesh)
 Produces a mesh of tets from an array of Octahedra. More...
 
Mesh test and repair
template<typename ExecSpace , typename FloatType >
void findTriMeshIntersectionsBVH (mint::UnstructuredMesh< mint::SINGLE_SHAPE > *surface_mesh, std::vector< std::pair< int, int >> &intersections, std::vector< int > &degenerateIndices, double intersectionThreshold=1E-8)
 Find self-intersections and degenerate triangles in a surface mesh utilizing a Bounding Volume Hierarchy. More...
 
template<typename ExecSpace , typename FloatType >
void findTriMeshIntersectionsImplicitGrid (mint::UnstructuredMesh< mint::SINGLE_SHAPE > *surface_mesh, std::vector< std::pair< int, int >> &intersections, std::vector< int > &degenerateIndices, int spatialIndexResolution=0, double intersectionThreshold=1E-8)
 Find self-intersections and degenerate triangles in a surface mesh utilizing an implicit grid spatial index. More...
 
template<typename ExecSpace , typename FloatType >
void findTriMeshIntersectionsUniformGrid (mint::UnstructuredMesh< mint::SINGLE_SHAPE > *surface_mesh, std::vector< std::pair< int, int >> &intersections, std::vector< int > &degenerateIndices, int spatialIndexResolution=0, double intersectionThreshold=1E-8)
 Find self-intersections and degenerate triangles in a surface mesh utilizing an uniform grid spatial index. More...
 
void findTriMeshIntersections (mint::UnstructuredMesh< mint::SINGLE_SHAPE > *surface_mesh, std::vector< std::pair< int, int >> &intersections, std::vector< int > &degenerateIndices, int spatialIndexResolution=0, double intersectionThreshold=1E-8)
 Find self-intersections and degenerate triangles in a surface mesh utilizing a Uniform Grid. More...
 
WatertightStatus isSurfaceMeshWatertight (mint::UnstructuredMesh< mint::SINGLE_SHAPE > *surface_mesh)
 Check a surface mesh for holes using its face relation. More...
 
void weldTriMeshVertices (mint::UnstructuredMesh< mint::SINGLE_SHAPE > **surface_mesh, double eps)
 Mesh repair function to weld vertices that are closer than eps. More...
 
InOut query – initialization and finalization functions
int inout_init (const std::string &file, MPI_Comm comm=MPI_COMM_SELF)
 Initializes the quest inout query from a mesh file. More...
 
int inout_init (mint::Mesh *&mesh, MPI_Comm comm=MPI_COMM_SELF)
 Initialize the inout query using a pre-loaded mesh. More...
 
int inout_finalize ()
 Finalizes the inout query. More...
 
bool inout_initialized ()
 Predicate to test whether the inout query has been initialized. More...
 
InOut query – properties and querying functions
Note
These must be called after inout_init()
bool inout_evaluate (double x, double y, double z=0.)
 Tests if the point (x, y, z) is inside the contained volume. More...
 
int inout_evaluate (const double *x, const double *y, const double *z, int npoints, int *res)
 Tests an array of points for containment. More...
 
int inout_mesh_min_bounds (double *coords)
 Returns the lower coordinates of the mesh's bounding box. More...
 
int inout_mesh_max_bounds (double *coords)
 Returns the upper coordinates of the mesh's bounding box. More...
 
int inout_mesh_center_of_mass (double *coords)
 Returns the center of mass of the mesh. More...
 
int inout_get_dimension ()
 Gets the spatial dimension of the query. More...
 
InOut query – setup options and parameters
Note
These must be called before inout_init()
int inout_set_dimension (int dim)
 Sets the spatial dimension of the query. More...
 
int inout_set_verbose (bool verbosity)
 Enables/disables verbose logging output. More...
 
int inout_set_vertex_weld_threshold (double thresh)
 Sets the cutoff distance for welding vertices during initialization. More...
 
int inout_set_segments_per_knot_span (int segmentsPerKnotSpan)
 Sets the number of samples for each knot span (2D only) More...
 
Signed Distance Query Initialization Methods
int signed_distance_init (const std::string &file, MPI_Comm comm=MPI_COMM_SELF)
 Initializes the Signed Distance Query with a surface given in an STL formatted file. More...
 
int signed_distance_init (const mint::Mesh *m, MPI_Comm comm=MPI_COMM_SELF)
 Initializes the Signed Distance Query with the given surface mesh. More...
 
bool signed_distance_initialized ()
 Checks if the Signed Distance Query has been initialized. More...
 
Signed Distance Query Options
void signed_distance_set_dimension (int dim)
 Sets the dimension for the Signed Distance Query. More...
 
void signed_distance_set_closed_surface (bool status)
 Indicates whether the input to the signed distance consists of a water-tight surface mesh, or not. More...
 
void signed_distance_set_compute_signs (bool computeSign)
 Sets whether the distance query should compute or ignore the sign. More...
 
void signed_distance_set_allocator (int allocatorID)
 Sets the allocator to use for creating internal signed distance query data structures. More...
 
void signed_distance_set_verbose (bool status)
 Enables/Disables verbose output for the Signed Distance Query. More...
 
void signed_distance_use_shared_memory (bool status)
 Enable/Disable the use of MPI-3 on-node shared memory for storing the surface mesh. By default this option is disabled. More...
 
void signed_distance_set_execution_space (SignedDistExec execSpace)
 Set the execution space in which to run signed distance queries. By default this option is set to SIGNED_DIST_EVAL_CPU. More...
 
Signed Distance Query Evaluation Methods
double signed_distance_evaluate (double x, double y, double z=0.0)
 Evaluates the signed distance function at the given point. More...
 
double signed_distance_evaluate (double x, double y, double z, double &cp_x, double &cp_y, double &cp_z, double &n_x, double &n_y, double &n_z)
 Evaluates the signed distance function at the given 3D point. More...
 
void signed_distance_evaluate (const double *x, const double *y, const double *z, int npoints, double *phi)
 Evaluates the signed distance function at the given set of points. More...
 
void signed_distance_get_mesh_bounds (double *lo, double *hi)
 Computes the bounds of the specified input mesh supplied to the Signed Distance Query. More...
 
Signed Distance Query Finalization Methods
void signed_distance_finalize ()
 Finalizes the SignedDistance query. More...
 

Variables

constexpr int QUEST_INOUT_SUCCESS = 0
 
constexpr int QUEST_INOUT_FAILED = -1
 

Discretize primitive shapes to linear shapes

using SphereType = primal::Sphere< double, 3 >
 
using OctType = primal::Octahedron< double, 3 >
 
using Point2D = primal::Point< double, 2 >
 
bool discretize (const SphereType &s, int levels, axom::Array< OctType > &out, int &octcount)
 Given a primitive sphere and a refinement level, allocate and return a list of Octahedra approximating the shape. More...
 
template<typename ExecSpace >
bool discretize (axom::Array< Point2D > &polyline, int len, int levels, axom::Array< OctType > &out, int &octcount)
 Given a 2D polyline revolved around the positive X-axis, allocate and return a list of Octahedra approximating the shape. More...
 

Typedef Documentation

◆ SphereType

using axom::quest::SphereType = typedef primal::Sphere<double, 3>

◆ OctType

using axom::quest::OctType = typedef primal::Octahedron<double, 3>

◆ Point2D

using axom::quest::Point2D = typedef primal::Point<double, 2>

Enumeration Type Documentation

◆ SearchStatus

Enumerator
NEIGHBOR_NOT_FOUND 

◆ WatertightStatus

enum axom::quest::WatertightStatus : signed char
strong

Enumeration indicating mesh watertightness

Enumerator
WATERTIGHT 

Each edge in a surface mesh is incident in two cells.

NOT_WATERTIGHT 

Each edge is incident in one or two cells.

CHECK_FAILED 

Calculation failed (possibly a non-manifold mesh)

◆ SignedDistExec

Enumerator
CPU 
OpenMP 
GPU 

Function Documentation

◆ all_nearest_neighbors()

void axom::quest::all_nearest_neighbors ( const double *  x,
const double *  y,
const double *  z,
const int *  region,
int  n,
double  limit,
int *  neighbor,
double *  sqdistance 
)

Given a list of point locations and regions, for each point, find the closest point in a different region within a given search radius.

Parameters
[in]xX-coordinates of input points
[in]yY-coordinates of input points
[in]zZ-coordinates of input points
[in]regionRegion of each point
[in]nNumber of points
[in]limitMax distance for all-nearest-neighbors query
[out]neighborIndex of nearest neighbor not in the same class (or NEIGHBOR_NOT_FOUND)
[out]sqdistanceSquared distance to nearest neighbor
Precondition
x, y, z, and region have n entries
neighbor is allocated with room for n entries

This method inserts all points p at (x[i], y[i], z[i]) into a UniformGrid index. Then for each point p, it gets the UniformGrid bins that overlap the box (p - (limit, limit, limit), p + (limit, limit, limit). The method compares p to each point in this list of bins and returns the index of the closest point.

We expect the use of the UniformGrid will result in a substantial time savings over a brute-force all-to-all algorithm, but the query's run time is dependent on the point distribution.

◆ discretize() [1/2]

bool axom::quest::discretize ( const SphereType s,
int  levels,
axom::Array< OctType > &  out,
int &  octcount 
)

Given a primitive sphere and a refinement level, allocate and return a list of Octahedra approximating the shape.

Parameters
[in]sThe sphere to approximate
[in]levelsThe number of refinements to perform, in addition to a central level-zero octahedron
[out]outThe newly-initialized Array of octahedra representing s
[out]octcountThe number of elements in out
Returns
false for invalid input or error in computation; true otherwise

This routine generates O(4^level) octahedra. That's exponential growth. Use appropriate caution.

This routine initializes an Array, out.

◆ discretize() [2/2]

template<typename ExecSpace >
bool axom::quest::discretize ( axom::Array< Point2D > &  polyline,
int  len,
int  levels,
axom::Array< OctType > &  out,
int &  octcount 
)

Given a 2D polyline revolved around the positive X-axis, allocate and return a list of Octahedra approximating the shape.

Parameters
[in]polylineThe polyline to revolve around the X-axis
[in]lenThe number of points in polyline
[in]levelsThe number of refinements to perform, in addition to a central level-zero octahedron in each segment
[out]outThe newly-initialized Array of octahedra representing the revolved polyline
[out]octcountThe number of elements in out
Returns
false for invalid input or error in computation; true otherwise

This routine generates n*O(2^level) octahedra, where n is the number of segments in polyline (one less than the length). That's exponential growth. Use appropriate caution.

This routine initializes an Array, out.

◆ mesh_from_discretized_polyline()

int axom::quest::mesh_from_discretized_polyline ( axom::ArrayView< OctType > &  octs,
int  octcount,
int  segcount,
mint::Mesh *&  mesh 
)

Produces a mesh of tets from an array of Octahedra.

Parameters
[in]octsArray of Octahedron objects
[in]octcountLength of octs
[in]segcountNumber of segments in the polyline originating the octs
[out]meshPointer to a mesh object where the mesh will be generated.
Returns
status error code, zero on success

This function creates a new Mint mesh out of an array of Octahedron objects. The mesh will be valid for any oct array, but the fields it includes make sense for the output of discretizing a revolved polyline:

  • octahedron_volume: octahedron volume calculated by summing tetrahedron volumes.
  • oct_as_polyhedron_volume: octahedron volume calculated by changing the oct into a Polyhedron and reporting the Polyhedron's volume. Any discrepancy between this field and octahedron_volume is a bug.
  • segment_index: the zero-based index of the segment to which the parent octahedron belongs.
  • octahedron_index: the zero-based index of the parent octahedron within its segment.
  • level_of_refinement: the level of refinement to which parent octahedron belongs.
Note
Ownership of the mesh object is passed to the caller. Consequently, the caller is responsible for properly deallocating the mesh object that the return mesh pointer points to.

◆ findTriMeshIntersectionsBVH()

template<typename ExecSpace , typename FloatType >
void axom::quest::findTriMeshIntersectionsBVH ( mint::UnstructuredMesh< mint::SINGLE_SHAPE > *  surface_mesh,
std::vector< std::pair< int, int >> &  intersections,
std::vector< int > &  degenerateIndices,
double  intersectionThreshold = 1E-8 
)

Find self-intersections and degenerate triangles in a surface mesh utilizing a Bounding Volume Hierarchy.

Parameters
[in]surface_meshA triangle mesh in three dimensions
[out]intersectionPairs of indices of intersecting mesh triangles
[out]degenerateIndicesindices of degenerate mesh triangles
[in]intersectionThresholdTolerance threshold for triangle intersection tests (default: 1E-8) After running this function over a surface mesh, intersection will be filled with pairs of indices of intersecting triangles and degenerateIndices will be filled with the indices of the degenerate triangles in the mesh. Triangles that share vertex pairs (adjacent triangles in a watertight surface mesh) are not reported as intersecting. Degenerate triangles are not reported as intersecting other triangles.

References AXOM_ANNOTATE_SCOPE, axom::Array< T, DIM, SPACE >::begin(), axom::Array< T, DIM, SPACE >::end(), axom::Array< T, DIM, SPACE >::resize(), axom::Array< T, DIM, SPACE >::size(), and SLIC_INFO.

◆ findTriMeshIntersectionsImplicitGrid()

template<typename ExecSpace , typename FloatType >
void axom::quest::findTriMeshIntersectionsImplicitGrid ( mint::UnstructuredMesh< mint::SINGLE_SHAPE > *  surface_mesh,
std::vector< std::pair< int, int >> &  intersections,
std::vector< int > &  degenerateIndices,
int  spatialIndexResolution = 0,
double  intersectionThreshold = 1E-8 
)

Find self-intersections and degenerate triangles in a surface mesh utilizing an implicit grid spatial index.

Parameters
[in]surface_meshA triangle mesh in three dimensions
[out]intersectionPairs of indices of intersecting mesh triangles
[out]degenerateIndicesindices of degenerate mesh triangles
[in]spatialIndexResolutionThe grid resolution for the index structure (default: 0)
[in]intersectionThresholdTolerance threshold for triangle intersection tests (default: 1E-8) After running this function over a surface mesh, intersection will be filled with pairs of indices of intersecting triangles and degenerateIndices will be filled with the indices of the degenerate triangles in the mesh. Triangles that share vertex pairs (adjacent triangles in a watertight surface mesh) are not reported as intersecting. Degenerate triangles are not reported as intersecting other triangles.

This function uses a quest::ImplicitGrid spatial index. Input spatialIndexResolution specifies the bin size for the UniformGrid. The default value of 0 causes this routine to calculate a heuristic bin size based on the cube root of the number of cells in the mesh.

References AXOM_ANNOTATE_SCOPE, axom::Array< T, DIM, SPACE >::begin(), axom::Array< T, DIM, SPACE >::end(), axom::Array< T, DIM, SPACE >::resize(), axom::Array< T, DIM, SPACE >::size(), and SLIC_INFO.

◆ findTriMeshIntersectionsUniformGrid()

template<typename ExecSpace , typename FloatType >
void axom::quest::findTriMeshIntersectionsUniformGrid ( mint::UnstructuredMesh< mint::SINGLE_SHAPE > *  surface_mesh,
std::vector< std::pair< int, int >> &  intersections,
std::vector< int > &  degenerateIndices,
int  spatialIndexResolution = 0,
double  intersectionThreshold = 1E-8 
)

Find self-intersections and degenerate triangles in a surface mesh utilizing an uniform grid spatial index.

Parameters
[in]surface_meshA triangle mesh in three dimensions
[out]intersectionPairs of indices of intersecting mesh triangles
[out]degenerateIndicesindices of degenerate mesh triangles
[in]spatialIndexResolutionThe grid resolution for the index structure (default: 0)
[in]intersectionThresholdTolerance threshold for triangle intersection tests (default: 1E-8) After running this function over a surface mesh, intersection will be filled with pairs of indices of intersecting triangles and degenerateIndices will be filled with the indices of the degenerate triangles in the mesh. Triangles that share vertex pairs (adjacent triangles in a watertight surface mesh) are not reported as intersecting. Degenerate triangles are not reported as intersecting other triangles.

This function uses a quest::UniformGrid spatial index. Input spatialIndexResolution specifies the bin size for the UniformGrid. The default value of 0 causes this routine to calculate a heuristic bin size based on the cube root of the number of cells in the mesh.

References AXOM_ANNOTATE_SCOPE, axom::Array< T, DIM, SPACE >::begin(), axom::Array< T, DIM, SPACE >::end(), axom::Array< T, DIM, SPACE >::resize(), axom::Array< T, DIM, SPACE >::size(), and SLIC_INFO.

◆ findTriMeshIntersections()

void axom::quest::findTriMeshIntersections ( mint::UnstructuredMesh< mint::SINGLE_SHAPE > *  surface_mesh,
std::vector< std::pair< int, int >> &  intersections,
std::vector< int > &  degenerateIndices,
int  spatialIndexResolution = 0,
double  intersectionThreshold = 1E-8 
)

Find self-intersections and degenerate triangles in a surface mesh utilizing a Uniform Grid.

Parameters
[in]surface_meshA triangle mesh in three dimensions
[out]intersectionPairs of indices of intersecting mesh triangles
[out]degenerateIndicesindices of degenerate mesh triangles
[in]spatialIndexResolutionThe grid resolution for the index structure (default: 0)
[in]intersectionThresholdTolerance threshold for triangle intersection tests (default: 1E-8)

After running this function over a surface mesh, intersection will be filled with pairs of indices of intersecting triangles and degenerateIndices will be filled with the indices of the degenerate triangles in the mesh. Triangles that share vertex pairs (adjacent triangles in a watertight surface mesh) are not reported as intersecting. Degenerate triangles are not reported as intersecting other triangles.

This function uses a quest::UniformGrid spatial index. Input spatialIndexResolution specifies the bin size for the UniformGrid. The default value of 0 causes this routine to calculate a heuristic bin size based on the cube root of the number of cells in the mesh.

◆ isSurfaceMeshWatertight()

WatertightStatus axom::quest::isSurfaceMeshWatertight ( mint::UnstructuredMesh< mint::SINGLE_SHAPE > *  surface_mesh)

Check a surface mesh for holes using its face relation.

Parameters
[in]surface_meshA surface mesh in three dimensions
Returns
status If the mesh is watertight, is not watertight, or if an error occurred (possibly due to non-manifold mesh).
Note
This method marks the cells on the boundary by creating a new cell-centered field variable, called "boundary", on the given input mesh.
This function computes the mesh's cell-face and face-vertex relations. For large meshes, this can take a long time. The relations are used to check for holes, and remain cached with the mesh after this function finishes.

◆ weldTriMeshVertices()

void axom::quest::weldTriMeshVertices ( mint::UnstructuredMesh< mint::SINGLE_SHAPE > **  surface_mesh,
double  eps 
)

Mesh repair function to weld vertices that are closer than eps.

Parameters
[in,out]surface_meshA pointer to a pointer to a triangle mesh
[in]epsDistance threshold for welding vertices (using the max norm)
Precondition
eps must be greater than zero
surface_mesh is a pointer to a pointer to a non-null triangle mesh.
Postcondition
The triangles of surface_mesh are reindexed using the welded vertices and degenerate triangles are removed. The mesh can still contain vertices that are not referenced by any triangles.

This utility function repairs an input triangle mesh (embedded in three dimensional space) by 'welding' vertices that are closer than eps. The vertices are quantized to an integer lattice with spacing eps and vertices that fall into the same cell on this lattice are identified. All identified vertices are given the coordinates of the first such vertex and all incident triangles use the same index for this vertex.

The input mesh can be a "soup of triangles", where the vertices of adjacent triangles have distinct indices. After running this function, vertices that are closer than eps are welded, and their incident triangles use the new vertex indices. Thus, the output is an "indexed triangle mesh".

This function also removes degenerate triangles from the mesh. These are triangles without three distinct vertices after the welding.

Note
This function is destructive. It modifies the input triangle mesh in place.
The distance metric in this function uses the "max" norm (l_inf).

◆ inout_init() [1/2]

int axom::quest::inout_init ( const std::string &  file,
MPI_Comm  comm = MPI_COMM_SELF 
)

Initializes the quest inout query from a mesh file.

Parameters
[in]filePath to an STL file containing the surface mesh
[in]commThe MPI communicator (when running in parallel)
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.
Precondition
inout_initialized() == false
Postcondition
inout_initialized() == true, when rc is QUEST_INOUT_SUCCESS

◆ inout_init() [2/2]

int axom::quest::inout_init ( mint::Mesh *&  mesh,
MPI_Comm  comm = MPI_COMM_SELF 
)

Initialize the inout query using a pre-loaded mesh.

Parameters
[in,out]meshPointer to the input mesh. This pointer will be updated during this invocation
[in]commThe MPI communicator (when running in parallel)
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.
Precondition
inout_initialized() == false
Postcondition
inout_initialized() == true, when rc is QUEST_INOUT_SUCCESS
Warning
The underlying data structure modifies the input mesh (e.g. by welding vertices) and updates the mesh pointer. It is the user's responsibility to update any other pointers to this same mesh.

◆ inout_finalize()

int axom::quest::inout_finalize ( )

Finalizes the inout query.

Postcondition
inout_initialized() == false

◆ inout_initialized()

bool axom::quest::inout_initialized ( )

Predicate to test whether the inout query has been initialized.

Returns
True if the inout query has been initialized, false otherwise.

◆ inout_evaluate() [1/2]

bool axom::quest::inout_evaluate ( double  x,
double  y,
double  z = 0. 
)

Tests if the point (x, y, z) is inside the contained volume.

Parameters
[in]xThe x-coordinate of the query point
[in]yThe y-coordinate of the query point
[in]zThe z-coordinate of the query point
Returns
True if the point is within the contained volume or on the surface mesh, false otherwise.
Precondition
inout_initialized() == true

◆ inout_evaluate() [2/2]

int axom::quest::inout_evaluate ( const double *  x,
const double *  y,
const double *  z,
int  npoints,
int *  res 
)

Tests an array of points for containment.

Upon successful completion, entries in array res will have the value 1 for points that are inside and value 0 otherwise.

Parameters
[in]xArray of x-coordinates for the query points
[in]yArray of y-coordinates for the query points
[in]zArray of z-coordinates for the query points
[in]npointsThe number of points to test
[out]resAn array of results. Each entry has value 1 if the corresponding point is inside or on the mesh and 0 otherwise.
Returns
Return code is QUEST_INOUT_SUCCESS when all preconditions are satisfied and QUEST_INOUT_FAILED otherwise.
Precondition
inout_initialized() == true
When npoints is greater than zero, arrays x, y, z and res are not nullptr and contain sufficient data/space for npoints points.

◆ inout_mesh_min_bounds()

int axom::quest::inout_mesh_min_bounds ( double *  coords)

Returns the lower coordinates of the mesh's bounding box.

Parameters
[in]coordsA buffer for the returned coordinates
Precondition
coords != nullptr and has sufficient storage for the coordinates
inout_initialized() == true
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.

◆ inout_mesh_max_bounds()

int axom::quest::inout_mesh_max_bounds ( double *  coords)

Returns the upper coordinates of the mesh's bounding box.

Parameters
[in]coordsA buffer for the returned coordinates
Precondition
coords != nullptr and has sufficient storage for the coordinates
inout_initialized() == true
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.

◆ inout_mesh_center_of_mass()

int axom::quest::inout_mesh_center_of_mass ( double *  coords)

Returns the center of mass of the mesh.

The function computes a discrete center of mass defined by the average of the mesh coordinates rather than a continuous center of mass defined by the mesh faces.

Parameters
[in]coordsA buffer for the returned coordinates
Precondition
coords != nullptr and has sufficient storage for the coordinates
inout_initialized() == true
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.

◆ inout_get_dimension()

int axom::quest::inout_get_dimension ( )

Gets the spatial dimension of the query.

Returns
Returns the spatial dimension for the inout query
Note
This determines the number of coordinates for the input mesh
The default dimension is 3
See also
inout_set_dimension

◆ inout_set_dimension()

int axom::quest::inout_set_dimension ( int  dim)

Sets the spatial dimension of the query.

Parameters
dimThe dimension for the inout query
Returns
Return code is QUEST_INOUT_SUCCESS if successful, QUEST_INOUT_FAILED otherwise
Precondition
dim == 2 || dim == 3
inout_initialized() == false

◆ inout_set_verbose()

int axom::quest::inout_set_verbose ( bool  verbosity)

Enables/disables verbose logging output.

By default, the logging verbosity is set to false.

Parameters
verbosityTrue for more verbose, false for less verbose
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.
Precondition
inout_initialized() == false

◆ inout_set_vertex_weld_threshold()

int axom::quest::inout_set_vertex_weld_threshold ( double  thresh)

Sets the cutoff distance for welding vertices during initialization.

By default, the welding threshold is 1E-9.

The inout query requires the input surface to be watertight so this parameter should be set with care. A welding threshold that is too high could unnecessarily merge vertices and create topological defects, while a value that is too low risks leaving gaps in meshes with tolerances between vertices. The default value tends to work well in practice.

Parameters
threshCutoff distance for welding vertices
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.
Precondition
inout_initialized() == false
thresh > 0

◆ inout_set_segments_per_knot_span()

int axom::quest::inout_set_segments_per_knot_span ( int  segmentsPerKnotSpan)

Sets the number of samples for each knot span (2D only)

By default, the welding threshold is 25

Span intervals are the segments of each curve. This parameter controls the number of samples to use when linearizing each knot span of the contour splines

Parameters
segmentsPerKnotSpanThe number of segments for each knot span
Returns
Return code is QUEST_INOUT_SUCCESS if successful and QUEST_INOUT_FAILED otherwise.
Precondition
inout_initialized() == false
segmentsPerKnotSpan >= 1

◆ signed_distance_init() [1/2]

int axom::quest::signed_distance_init ( const std::string &  file,
MPI_Comm  comm = MPI_COMM_SELF 
)

Initializes the Signed Distance Query with a surface given in an STL formatted file.

Parameters
[in]filethe name of the file consisting of the surface mesh.
[in]commthe MPI communicator (applicable when MPI is available)
Returns
status zero on success, or a non-zero value if an error occurs
Note
The Signed Distance Query currently only supports reading in meshes defined in STL file format.
Precondition
file.empty() == false
comm != MPI_COMM_NULL (when MPI is available)
signed_distance_initialized() == false
Postcondition
signed_distance_initialized() == true

◆ signed_distance_init() [2/2]

int axom::quest::signed_distance_init ( const mint::Mesh m,
MPI_Comm  comm = MPI_COMM_SELF 
)

Initializes the Signed Distance Query with the given surface mesh.

Parameters
[in]mpointer to the surface mesh object
[in]commthe MPI communicator (applicable whem MPI is available)
Returns
status zero on success, or a non-zero value if an error occurs.
Note
The Signed Distance Query currently only supports 3-D triangular surface meshes.
Precondition
m != nullptr
comm != MPI_COMM_NULL (when MPI is available)
signed_distance_initialized() == false
Postcondition
signed_distance_initialized() == true
See also
mint::Mesh

◆ signed_distance_initialized()

bool axom::quest::signed_distance_initialized ( )

Checks if the Signed Distance Query has been initialized.

Returns
status true if initialized, else, false.

◆ signed_distance_set_dimension()

void axom::quest::signed_distance_set_dimension ( int  dim)

Sets the dimension for the Signed Distance Query.

Parameters
[in]dimthe dimension, e.g., 2 or 3
Warning
The Signed Distance function is currently supported in 3D
Note
Options must be set before initializing the Signed Distance Query.

◆ signed_distance_set_closed_surface()

void axom::quest::signed_distance_set_closed_surface ( bool  status)

Indicates whether the input to the signed distance consists of a water-tight surface mesh, or not.

Parameters
[in]statusflag indicating whether the input is water-tight
Note
By default the input type is assumed to be a water-tight surface mesh.
Options must be set before initializing the Signed Distance Query.
When the input is not a closed surface mesh, the assumption is that the surface mesh divides the computational mesh domain into two regions. Hence, the surface mesh has to span the entire domain of interest, e.g., the computational mesh at which the signed distance field is evaluated, along some plane.
Warning
The sign of the distance from a given query point is determined by a pseudo-normal which is computed at the closest point on the surface mesh. For a non-watertight mesh, the sign of the distance is not defined everywhere. Specifically, the sign is ambiguous for all points for which a normal projection onto the surface does not exist.

◆ signed_distance_set_compute_signs()

void axom::quest::signed_distance_set_compute_signs ( bool  computeSign)

Sets whether the distance query should compute or ignore the sign.

Parameters
[in]computeSignpredicate indicating if sign should be computed
Note
Options must be set before initializing the Signed Distance Query.

◆ signed_distance_set_allocator()

void axom::quest::signed_distance_set_allocator ( int  allocatorID)

Sets the allocator to use for creating internal signed distance query data structures.

Parameters
[in]allocatorIDthe allocator ID to use
Note
Allocator should be compatible with the set execution space.

◆ signed_distance_set_verbose()

void axom::quest::signed_distance_set_verbose ( bool  status)

Enables/Disables verbose output for the Signed Distance Query.

Parameters
[in]statusflag indicating whether to enable/disable verbose output
Note
Options must be set before initializing the Signed Distance Query.
Currently, this is only applicable when the Signed Distance Query initializes the SLIC logging environment.

◆ signed_distance_use_shared_memory()

void axom::quest::signed_distance_use_shared_memory ( bool  status)

Enable/Disable the use of MPI-3 on-node shared memory for storing the surface mesh. By default this option is disabled.

Parameters
[in]statusflag indicating whether to enable/disable shared memory.
Note
This option utilities MPI-3 features

◆ signed_distance_set_execution_space()

void axom::quest::signed_distance_set_execution_space ( SignedDistExec  execSpace)

Set the execution space in which to run signed distance queries. By default this option is set to SIGNED_DIST_EVAL_CPU.

Parameters
[in]exec_spacethe execution space setting, one of the enum values in SignedDistExec.
Note
This function resets the default allocator based on the requested execution space. If a custom allocator ID is desired, set it in a call to axom::setDefaultAllocator after this call.

◆ signed_distance_evaluate() [1/3]

double axom::quest::signed_distance_evaluate ( double  x,
double  y,
double  z = 0.0 
)

Evaluates the signed distance function at the given point.

Parameters
[in]xthe x-coordinate of the point in query
[in]ythe y-coordinate of the point in query
[in]zthe z-coordinate of the point in query
Returns
d the signed distance evaluated at the specified point.

◆ signed_distance_evaluate() [2/3]

double axom::quest::signed_distance_evaluate ( double  x,
double  y,
double  z,
double &  cp_x,
double &  cp_y,
double &  cp_z,
double &  n_x,
double &  n_y,
double &  n_z 
)

Evaluates the signed distance function at the given 3D point.

Parameters
[in]xthe x-coordinate of the point in query
[in]ythe y-coordinate of the point in query
[in]zthe z-coordinate of the point in query
[out]cp_xthe x-coordinate of the computed closest point on surface to query point
[out]cp_ythe y-coordinate of the computed closest point on surface to query point
[out]cp_zthe z-coordinate of the computed closest point on surface to query point
[out]n_xthe x-coordinate of the surface normal at the computed closest point
[out]n_ythe y-coordinate of the surface normal at the computed closest point
[out]n_zthe z-coordinate of the surface normal at the computed closest point
Returns
d the signed distance evaluated at the specified point. The closest surface point to the query point and the normal at that point are returned as OUT parameters

◆ signed_distance_evaluate() [3/3]

void axom::quest::signed_distance_evaluate ( const double *  x,
const double *  y,
const double *  z,
int  npoints,
double *  phi 
)

Evaluates the signed distance function at the given set of points.

Parameters
[in]xarray consisting of the x-coordinates for each query point
[in]yarray consisting of the y-coordinates for each query point
[in]zarray consisting of the z-coordinates for each query point
[in]npointsthe number of query point
[out]phioutput array storing the signed distance of each point
Precondition
x != nullptr
y != nullptr
z != nullptr
phi != nullptr

◆ signed_distance_get_mesh_bounds()

void axom::quest::signed_distance_get_mesh_bounds ( double *  lo,
double *  hi 
)

Computes the bounds of the specified input mesh supplied to the Signed Distance Query.

Parameters
[out]lobuffer to store the lower bound mesh coordinates.
[out]hibuffer to store the upper bound mesh coordinates.
Precondition
lo != nullptr
hi != nullptr
hi & lo must point to a buffer that is at least ndims long.
signed_distance_initialized() == true

◆ signed_distance_finalize()

void axom::quest::signed_distance_finalize ( )

Finalizes the SignedDistance query.

Postcondition
signed_distance_initialized() == false.

Variable Documentation

◆ QUEST_INOUT_SUCCESS

constexpr int axom::quest::QUEST_INOUT_SUCCESS = 0
constexpr

◆ QUEST_INOUT_FAILED

constexpr int axom::quest::QUEST_INOUT_FAILED = -1
constexpr