Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys > Class Template Reference

CoreDriver - Core driver that remaps fields on a particular Entity_kind (ONWHAT) like CELL or NODE. More...

#include <coredriver.h>

Inheritance diagram for Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >:
Portage::CoreDriverBase< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >

Public Member Functions

 CoreDriver (SourceMesh const &source_mesh, SourceState const &source_state, TargetMesh const &target_mesh, TargetState &target_state, Wonton::Executor_type const *executor=nullptr)
 Constructor for the CORE remap driver. More...
 
 CoreDriver (const CoreDriver &)=delete
 Copy constructor (disabled) More...
 
CoreDriveroperator= (const CoreDriver &)=delete
 Assignment operator (disabled) More...
 
 ~CoreDriver ()=default
 Destructor. More...
 
Entity_kind onwhat ()
 What entity kind is this defined on? More...
 
template<template< int, Entity_kind, class, class > class Search>
Portage::vector< std::vector< int > > search ()
 
template<template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
Portage::vector< std::vector< Portage::Weights_t > > intersect_meshes (Portage::vector< std::vector< int >> const &candidates)
 
void set_num_tols (const double min_absolute_distance, const double min_absolute_volume)
 Set core numerical tolerances. More...
 
void set_num_tols (const NumericTolerances_t &num_tols)
 Set all numerical tolerances. More...
 
template<template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
std::vector< Portage::vector< std::vector< Weights_t > > > intersect_materials (Portage::vector< std::vector< int >> const &candidates)
 
Portage::vector< Vector< D > > compute_source_gradient (std::string const field_name, Limiter_type limiter_type=NOLIMITER, Boundary_Limiter_type boundary_limiter_type=BND_NOLIMITER, int material_id=0, const Part< SourceMesh, SourceState > *source_part=nullptr) const
 Compute the gradient field of the given variable on source mesh. More...
 
template<typename T = double, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
void interpolate_mesh_var (std::string srcvarname, std::string trgvarname, Portage::vector< std::vector< Weights_t >> const &sources_and_weights, Portage::vector< Vector< D >> *gradients=nullptr)
 Interpolate mesh variable. More...
 
template<typename T = double, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate, Entity_kind ONWHAT1 = ONWHAT, typename = typename std::enable_if<ONWHAT1 == CELL>::type>
void interpolate_mesh_var (std::string srcvarname, std::string trgvarname, Portage::vector< std::vector< Weights_t >> const &sources_and_weights, const PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState > *partition, Portage::vector< Vector< D >> *gradients=nullptr)
 Interpolate mesh variable from source part to target part. More...
 
bool check_mismatch (Portage::vector< std::vector< Weights_t >> const &source_weights)
 
bool has_mismatch ()
 
bool fix_mismatch (std::string const &src_var_name, std::string const &trg_var_name, double global_lower_bound=-std::numeric_limits< double >::max(), double global_upper_bound=std::numeric_limits< double >::max(), double conservation_tol=1e2 *std::numeric_limits< double >::epsilon(), int maxiter=5, Partial_fixup_type partial_fixup_type=Partial_fixup_type::SHIFTED_CONSERVATIVE, Empty_fixup_type empty_fixup_type=Empty_fixup_type::EXTRAPOLATE)
 Repair the remapped field to account for boundary mismatch. More...
 
- Public Member Functions inherited from Portage::CoreDriverBase< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >
 CoreDriverBase ()=default
 
virtual ~CoreDriverBase ()=default
 
template<Entity_kind ONWHAT, template< int, Entity_kind, class, class > class Search>
Portage::vector< std::vector< int > > search ()
 search for candidate source entities whose control volumes (cells, dual cells) overlap the control volumes of target cells More...
 
template<Entity_kind ONWHAT, template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
Portage::vector< std::vector< Portage::Weights_t > > intersect_meshes (Portage::vector< std::vector< int >> const &intersection_candidates)
 intersect target entities with candidate source entities More...
 
template<template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
std::vector< Portage::vector< std::vector< Portage::Weights_t > > > intersect_materials (Portage::vector< std::vector< int >> const &intersection_candidates)
 
template<Entity_kind ONWHAT>
Portage::vector< Vector< D > > compute_source_gradient (std::string const field_name, Limiter_type limiter_type=NOLIMITER, Boundary_Limiter_type boundary_limiter_type=BND_NOLIMITER, int material_id=0, const Part< SourceMesh, SourceState > *source_part=nullptr)
 Compute the gradient field of the given variable on source mesh. More...
 
template<typename T = double, Entity_kind ONWHAT, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
void interpolate_mesh_var (std::string srcvarname, std::string trgvarname, Portage::vector< std::vector< Weights_t >> const &sources_and_weights, Portage::vector< Vector< D >> *gradients=nullptr)
 
template<typename T = double, Entity_kind ONWHAT, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
std::enable_if< ONWHAT==CELL, void > interpolate_mesh_var (std::string srcvarname, std::string trgvarname, Portage::vector< std::vector< Weights_t >> const &sources_and_weights, const PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState > *parts_pair, Portage::vector< Vector< D >> *gradients=nullptr)
 
template<typename T = double, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
void interpolate_mat_var (std::string srcvarname, std::string trgvarname, std::vector< Portage::vector< std::vector< Weights_t >>> const &sources_and_weights_by_mat, std::vector< Portage::vector< Vector< D >>> *gradients=nullptr)
 
template<Entity_kind ONWHAT>
bool check_mismatch (Portage::vector< std::vector< Weights_t >> const &source_weights)
 Check if meshes are mismatched (don't cover identical portions of space) More...
 
template<Entity_kind ONWHAT>
bool has_mismatch ()
 Return if meshes are mismatched (check_mismatch must already have been called) More...
 
template<Entity_kind ONWHAT>
bool fix_mismatch (std::string const &src_var_name, std::string const &trg_var_name, double global_lower_bound=-std::numeric_limits< double >::max(), double global_upper_bound=std::numeric_limits< double >::max(), double conservation_tol=1e2 *std::numeric_limits< double >::epsilon(), int maxiter=5, Partial_fixup_type partial_fixup_type=Partial_fixup_type::SHIFTED_CONSERVATIVE, Empty_fixup_type empty_fixup_type=Empty_fixup_type::EXTRAPOLATE)
 Return if meshes are mismatched (check_mismatch must already have been called) More...
 
template<Entity_kind ONWHAT>
void set_num_tols (const double min_absolute_distance, const double min_absolute_volume)
 Set numerical tolerances for small distances and volumes. More...
 
template<Entity_kind ONWHAT>
void set_num_tols (const NumericTolerances_t &num_tols)
 Set numerical tolerances for small volumes, distances, etc. More...
 

Detailed Description

template<int D, Entity_kind ONWHAT, class SourceMesh, class SourceState, class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
class Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >

CoreDriver - Core driver that remaps fields on a particular Entity_kind (ONWHAT) like CELL or NODE.

NOTE: THIS CLASS ASSUMES THAT ALL SOURCE CELLS OVERLAPPING ANY TARGET CELL ARE AVAILABLE ON THIS PROCESSOR. IT DOES NOT HAVE TO FETCH THE DATA FROM ANYWHERE

Template Parameters
ONWHATOn what kind of entity are we doing the remap
SourceMeshA lightweight wrapper to a specific input mesh implementation that provides certain functionality.
SourceStateA lightweight wrapper to a specific input state manager implementation that provides certain functionality.
TargetMeshA lightweight wrapper to a specific output mesh implementation that provides certain functionality.
TargetStateA lightweight wrapper to a specific output state manager implementation that provides certain functionality.
InterfaceReconstructorTypeThe Interface Reconstructor class we will instantiate
Matpoly_SplitterClass to split a polyhedron into two pieces
Matpoly_ClipperClass to clip a polyhedron with a plane and return the piece behind the plane

CoordSys Coordinate system being used for calculations

Constructor & Destructor Documentation

◆ CoreDriver() [1/2]

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::CoreDriver ( SourceMesh const &  source_mesh,
SourceState const &  source_state,
TargetMesh const &  target_mesh,
TargetState &  target_state,
Wonton::Executor_type const *  executor = nullptr 
)
inline

Constructor for the CORE remap driver.

Parameters
[in]sourceMeshA wrapper to the source mesh (may be native or redistributed source).
[in]sourceStateA wrapper for the data that lives on the source mesh
[in]targetMeshA TargetMesh to the target mesh
[in,out]targetStateA TargetState for the data that will be mapped to the target mesh

◆ CoreDriver() [2/2]

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::CoreDriver ( const CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys > &  )
delete

Copy constructor (disabled)

◆ ~CoreDriver()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::~CoreDriver ( )
default

Destructor.

Member Function Documentation

◆ check_mismatch()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
bool Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::check_mismatch ( Portage::vector< std::vector< Weights_t >> const &  source_weights)
inline

Check mismatch between meshes

Parameters
[in]sources_and_weightsIntersection sources and moments (vols, centroids)
Returns
Whether the meshes are mismatched

◆ compute_source_gradient()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
Portage::vector<Vector<D> > Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::compute_source_gradient ( std::string const  field_name,
Limiter_type  limiter_type = NOLIMITER,
Boundary_Limiter_type  boundary_limiter_type = BND_NOLIMITER,
int  material_id = 0,
const Part< SourceMesh, SourceState > *  source_part = nullptr 
) const
inline

Compute the gradient field of the given variable on source mesh.

Parameters
field_namethe variable name.
limiter_typegradient limiter to use on internal regions.
boundary_limiter_typegradient limiter to use on boundary.
source_partthe source mesh part to consider if any.

◆ fix_mismatch()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
bool Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::fix_mismatch ( std::string const &  src_var_name,
std::string const &  trg_var_name,
double  global_lower_bound = -std::numeric_limits<double>::max(),
double  global_upper_bound = std::numeric_limits<double>::max(),
double  conservation_tol = 1e2*std::numeric_limits<double>::epsilon(),
int  maxiter = 5,
Partial_fixup_type  partial_fixup_type = Partial_fixup_type::SHIFTED_CONSERVATIVE,
Empty_fixup_type  empty_fixup_type = Empty_fixup_type::EXTRAPOLATE 
)
inline

Repair the remapped field to account for boundary mismatch.

Parameters
src_var_namefield variable on source mesh
trg_var_namefield variable on target mesh
global_lower_boundlower limit on variable
global_upper_boundupper limit on variable
partial_fixup_typetype of fixup in case of partial mismatch
empty_fixup_typetype of fixup in empty target entities

partial_fixup_type can be one of three types:

CONSTANT - Fields will see no perturbations BUT REMAP WILL BE NON-CONSERVATIVE (constant preserving, not linearity preserving) LOCALLY_CONSERVATIVE - REMAP WILL BE LOCALLY CONSERVATIVE (target cells will preserve the integral quantities received from source mesh overlap) but perturbations will occur in the field (constant fields may not stay constant if there is mismatch) SHIFTED_CONSERVATIVE - REMAP WILL BE CONSERVATIVE and field perturbations will be minimum but field values may be shifted (Constant fields will be shifted to different constant; no guarantees on linearity preservation)

empty_fixup_type can be one of two types:

LEAVE_EMPTY - Leave empty cells as is EXTRAPOLATE - Fill empty cells with extrapolated values FILL - Fill empty cells with specified values (not yet implemented)

◆ has_mismatch()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
bool Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::has_mismatch ( )
inline

Return mismatch between meshes

Returns
Whether the meshes are mismatched

◆ interpolate_mesh_var() [1/2]

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
template<typename T = double, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
void Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::interpolate_mesh_var ( std::string  srcvarname,
std::string  trgvarname,
Portage::vector< std::vector< Weights_t >> const &  sources_and_weights,
Portage::vector< Vector< D >> *  gradients = nullptr 
)
inline

Interpolate mesh variable.

Parameters
[in]srcvarnamesource mesh variable to remap
[in]trgvarnametarget mesh variable to remap
[in]sources_and_weightsweights for mesh-mesh interpolation
[in]gradientsgradients of variable on source mesh (can be nullptr for 1st order remap)

◆ interpolate_mesh_var() [2/2]

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
template<typename T = double, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate, Entity_kind ONWHAT1 = ONWHAT, typename = typename std::enable_if<ONWHAT1 == CELL>::type>
void Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::interpolate_mesh_var ( std::string  srcvarname,
std::string  trgvarname,
Portage::vector< std::vector< Weights_t >> const &  sources_and_weights,
const PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState > *  partition,
Portage::vector< Vector< D >> *  gradients = nullptr 
)
inline

Interpolate mesh variable from source part to target part.

Parameters
[in]srcvarnamesource mesh variable to remap
[in]trgvarnametarget mesh variable to remap
[in]sources_and_weightsweights for mesh-mesh interpolation
[in]lower_boundlower bound of variable value
[in]upper_boundupper bound of variable value
[in]partitionstructure containing source and target part
[in]gradientsgradients of variable on source mesh (can be nullptr for 1st order remap)

Enable only for cells using SFINAE. Here the class, rather than the function is templated on ONWHAT (as opposed to the equivalent method in the base class); so we have to create a dummy template parameter ONWHAT1 and rely on that to use SFINAE with a second dummy template parameter

 Note ****

If you encounter errors about not being able to find an appropriate overload for interpolate_mesh_var in your application code (particularly something like "no type named 'type' in struct std::enable_if<false, void>"), make sure the compiler does see the possiblity of calling this function with Entity_kinds that are not type CELL (restricting the code flow using 'if' statements will not be enough)

◆ intersect_materials()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
template<template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
std::vector<Portage::vector<std::vector<Weights_t> > > Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::intersect_materials ( Portage::vector< std::vector< int >> const &  candidates)
inline

Intersect target mesh cells with source material polygons and return the intersecting entities and moments of intersection for each entity

Parameters
[in]candidatesVector of intersection candidates for each target entity
Returns
Material-wise vector of intersection moments for each target entity

NOTE: WE COULD SEND IN THE MESH-MESH INTERSECTION MOMENTS AND REUSE THEM IF A SOURCE CELL HAS ONLY ONE MATERIAL

◆ intersect_meshes()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
template<template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
Portage::vector<std::vector<Portage::Weights_t> > Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::intersect_meshes ( Portage::vector< std::vector< int >> const &  candidates)
inline

Intersect source and target mesh entities of kind 'ONWHAT' and return the intersecting entities and moments of intersection for each entity

Parameters
candidatesVector of intersection candidates for each target entity
Returns
vector of intersection moments for each target entity

◆ onwhat()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
Entity_kind Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::onwhat ( )
inlinevirtual

◆ operator=()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
CoreDriver& Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::operator= ( const CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys > &  )
delete

Assignment operator (disabled)

◆ search()

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
template<template< int, Entity_kind, class, class > class Search>
Portage::vector<std::vector<int> > Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::search ( )
inline

Find candidates entities of a particular kind that might intersect each target entity of the same kind

Template Parameters
SearchSearch class templated on dimension, Entity_kind and both meshes
Returns
Vector of intersection candidates for each target entity

◆ set_num_tols() [1/2]

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
void Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::set_num_tols ( const double  min_absolute_distance,
const double  min_absolute_volume 
)
inline

Set core numerical tolerances.

◆ set_num_tols() [2/2]

template<int D, Entity_kind ONWHAT, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState, template< class, int, class, class > class InterfaceReconstructorType = DummyInterfaceReconstructor, class Matpoly_Splitter = void, class Matpoly_Clipper = void, class CoordSys = Wonton::DefaultCoordSys>
void Portage::CoreDriver< D, ONWHAT, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::set_num_tols ( const NumericTolerances_t num_tols)
inline

Set all numerical tolerances.


The documentation for this class was generated from the following file: