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

UberDriver provides the API to mapping multi-material data from one mesh to another in a general way. More...

#include <uberdriver.h>

Public Types

using SerialDriverType = CoreDriverBase< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >
 
using ParallelDriverType = CoreDriverBase< D, Flat_Mesh_Wrapper<>, Flat_State_Wrapper< Flat_Mesh_Wrapper<> >, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >
 

Public Member Functions

 UberDriver (SourceMesh const &source_mesh, SourceState const &source_state, TargetMesh const &target_mesh, TargetState &target_state, std::vector< std::string > const &source_vars_to_remap, Wonton::Executor_type const *executor=nullptr, std::string *errmsg=nullptr)
 Constructor for the remap driver. More...
 
 UberDriver (SourceMesh const &source_mesh, SourceState const &source_state, TargetMesh const &target_mesh, TargetState &target_state, Wonton::Executor_type const *executor=nullptr, std::string *errmsg=nullptr)
 Constructor for the remap driver. More...
 
 UberDriver (const UberDriver &)=delete
 Copy constructor (disabled) More...
 
UberDriveroperator= (const UberDriver &)=delete
 Assignment operator (disabled) More...
 
 UberDriver (UberDriver &&) noexcept=default
 Enable move semantics. More...
 
 ~UberDriver ()=default
 Destructor. More...
 
bool is_distributed_run (Wonton::Executor_type const *executor=nullptr)
 Is this a distributed (multi-rank) run? More...
 
bool source_needs_redistribution (Wonton::Executor_type const *executor=nullptr)
 
template<template< int, Entity_kind, class, class > class Search, template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
void compute_interpolation_weights ()
 Compute interpolation weights in advance of actual interpolation of variables. More...
 
void set_num_tols (const double min_absolute_distance, const double min_absolute_volume)
 Set numerical tolerances for small distances and volumes in core driver. More...
 
void set_num_tols (const NumericTolerances_t &num_tols)
 set numerical tolerances in core driver More...
 
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 &candidates)
 intersect target control volumes with source control volumes 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 &candidates)
 
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 (std::string srcvarname, T lower_bound, T upper_bound, Limiter_type limiter=DEFAULT_LIMITER, Boundary_Limiter_type bnd_limiter=DEFAULT_BND_LIMITER, Partial_fixup_type partial_fixup_type=DEFAULT_PARTIAL_FIXUP_TYPE, Empty_fixup_type empty_fixup_type=DEFAULT_EMPTY_FIXUP_TYPE, double conservation_tol=DEFAULT_NUMERIC_TOLERANCES< D >.relative_conservation_eps, int max_fixup_iter=DEFAULT_NUMERIC_TOLERANCES< D >.max_num_fixup_iter)
 
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 (std::string srcvarname, std::string trgvarname, T lower_bound, T upper_bound, Limiter_type limiter=DEFAULT_LIMITER, Boundary_Limiter_type bnd_limiter=DEFAULT_BND_LIMITER, Partial_fixup_type partial_fixup_type=DEFAULT_PARTIAL_FIXUP_TYPE, Empty_fixup_type empty_fixup_type=DEFAULT_EMPTY_FIXUP_TYPE, double conservation_tol=DEFAULT_NUMERIC_TOLERANCES< D >.relative_conservation_eps, int max_fixup_iter=DEFAULT_NUMERIC_TOLERANCES< D >.max_num_fixup_iter)
 
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_in, T lower_bound, T upper_bound, Limiter_type limiter, Boundary_Limiter_type bnd_limiter, Partial_fixup_type partial_fixup_type, Empty_fixup_type empty_fixup_type, double conservation_tol, int max_fixup_iter)
 
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_in, T lower_bound, T upper_bound, Limiter_type limiter, Boundary_Limiter_type bnd_limiter, Partial_fixup_type partial_fixup_type, Empty_fixup_type empty_fixup_type, double conservation_tol, int max_fixup_iter)
 

Detailed Description

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >

UberDriver provides the API to mapping multi-material data from one mesh to another in a general way.

Template Parameters
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 target mesh implementation that provides certain functionality.
TargetStateA lightweight wrapper to a specific target state manager implementation that provides certain functionality.
InterfaceReconstructorTypeAn interface reconstruction class that takes the raw interface reconstruction method, the dimension of the problem and the source mesh class as template parameters
Matpoly_SplitterA polygon/polyhedron splitting class (returns both pieces of the polygon)
Matpoly_ClipperA polygon/polyhedron clipping class (returns only the piece below/behind the clipping plane)
CoordinateSystem

Member Typedef Documentation

◆ ParallelDriverType

template<int D, 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>
using Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::ParallelDriverType = CoreDriverBase<D, Flat_Mesh_Wrapper<>, Flat_State_Wrapper<Flat_Mesh_Wrapper<> >, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys>

◆ SerialDriverType

template<int D, 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>
using Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::SerialDriverType = CoreDriverBase<D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys>

Constructor & Destructor Documentation

◆ UberDriver() [1/4]

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::UberDriver ( SourceMesh const &  source_mesh,
SourceState const &  source_state,
TargetMesh const &  target_mesh,
TargetState &  target_state,
std::vector< std::string > const &  source_vars_to_remap,
Wonton::Executor_type const *  executor = nullptr,
std::string *  errmsg = nullptr 
)
inline

Constructor for the remap driver.

Parameters
[in]source_meshA wrapper to the source mesh.
[in]source_stateA wrapper for the data that lives on the source mesh
[in]target_meshA TargetMesh to the target mesh
[in,out]target_stateA TargetState for the data that will be mapped to the target mesh
[in]source_vars_to_remapOptional list of source variables to remap (if not everything will be remapped)
[in]executorpointer to an executor allowing us choose between serial and parallel runs

◆ UberDriver() [2/4]

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::UberDriver ( SourceMesh const &  source_mesh,
SourceState const &  source_state,
TargetMesh const &  target_mesh,
TargetState &  target_state,
Wonton::Executor_type const *  executor = nullptr,
std::string *  errmsg = nullptr 
)
inline

Constructor for the remap driver.

Parameters
[in]sourceMeshA wrapper to the source mesh.
[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

◆ UberDriver() [3/4]

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::UberDriver ( const UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys > &  )
delete

Copy constructor (disabled)

◆ UberDriver() [4/4]

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::UberDriver ( UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys > &&  )
defaultnoexcept

Enable move semantics.

◆ ~UberDriver()

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::~UberDriver ( )
default

Destructor.

Member Function Documentation

◆ compute_interpolation_weights()

template<int D, 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, template< Entity_kind, class, class, class, template< class, int, class, class > class, class, class > class Intersect>
void Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::compute_interpolation_weights ( )
inline

Compute interpolation weights in advance of actual interpolation of variables.

Template Parameters
SearchA search method that takes the dimension, source mesh class and target mesh class as template parameters
IntersectA polyhedron-polyhedron intersection class that takes the source and taget mesh classes as template parameters

◆ interpolate() [1/2]

template<int D, 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, Entity_kind ONWHAT, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
void Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::interpolate ( std::string  srcvarname,
lower_bound,
upper_bound,
Limiter_type  limiter = DEFAULT_LIMITER,
Boundary_Limiter_type  bnd_limiter = DEFAULT_BND_LIMITER,
Partial_fixup_type  partial_fixup_type = DEFAULT_PARTIAL_FIXUP_TYPE,
Empty_fixup_type  empty_fixup_type = DEFAULT_EMPTY_FIXUP_TYPE,
double  conservation_tol = DEFAULT_NUMERIC_TOLERANCES<D>.relative_conservation_eps,
int  max_fixup_iter = DEFAULT_NUMERIC_TOLERANCES<D>.max_num_fixup_iter 
)
inline

Interpolate a mesh variable of type T residing on entity kind ONWHAT using previously computed intersection weights (same variable name on source and target)

Template Parameters
Ttype of variable being remapped (default double) - underlying interpolator must be able to handle this type
ONWHATEntity kind on which variable resides
InterpolateInterpolate functor to do the actual interpolation
Parameters
[in]srcvarnameVariable name on source and target meshes
[in]lower_boundLower bound for variable
[in]upper_boundUpper bound for variable
[in]limiterLimiter to use for second order reconstruction
[in]bnd_limiterBoundary limiter to use for second order reconstruction
[in]partial_fixup_typeMethod to populate fields on partially filled target entities (cells or dual cells)
[in]empty_fixup_typeMethod to populate fields on empty target entities (cells or dual cells)
[in]conservation_tolTolerance to which source and target integral quantities are to be matched
[in]max_fixup_iterMax number of iterations for global repair

See support/portage.h for options on limiter, partial_fixup_type and empty_fixup_type

◆ interpolate() [2/2]

template<int D, 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, Entity_kind ONWHAT, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
void Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::interpolate ( std::string  srcvarname,
std::string  trgvarname,
lower_bound,
upper_bound,
Limiter_type  limiter = DEFAULT_LIMITER,
Boundary_Limiter_type  bnd_limiter = DEFAULT_BND_LIMITER,
Partial_fixup_type  partial_fixup_type = DEFAULT_PARTIAL_FIXUP_TYPE,
Empty_fixup_type  empty_fixup_type = DEFAULT_EMPTY_FIXUP_TYPE,
double  conservation_tol = DEFAULT_NUMERIC_TOLERANCES<D>.relative_conservation_eps,
int  max_fixup_iter = DEFAULT_NUMERIC_TOLERANCES<D>.max_num_fixup_iter 
)
inline

Interpolate a mesh variable of type T residing on entity kind ONWHAT using previously computed intersection weights (different variable name on source and target)

Template Parameters
Ttype of variable being remapped (default double) - underlying interpolator must be able to handle this type
Parameters
[in]srcvarnameVariable name on source mesh
[in]trgvarnameVariable name on target mesh
[in]lower_boundLower bound for variable
[in]upper_boundUpper bound for variable
[in]limiterLimiter to use for second order reconstruction
[in]bnd_limiterBoundary limiter to use for second order reconstruction
[in]partial_fixup_typeMethod to populate fields on partially filled target entities (cells or dual cells)
[in]empty_fixup_typeMethod to populate fields on empty target entities (cells or dual cells)
[in]conservation_tolTolerance to which source and target integral quantities are to be matched
[in]max_fixup_iterMax number of iterations for global repair

See support/portage.h for options on limiter, partial_fixup_type and empty_fixup_type

◆ interpolate_mat_var()

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::interpolate_mat_var ( std::string  srcvarname,
std::string  trgvarname,
std::vector< Portage::vector< std::vector< Weights_t >>> const &  sources_and_weights_by_mat_in,
lower_bound,
upper_bound,
Limiter_type  limiter,
Boundary_Limiter_type  bnd_limiter,
Partial_fixup_type  partial_fixup_type,
Empty_fixup_type  empty_fixup_type,
double  conservation_tol,
int  max_fixup_iter 
)
inline

Interpolate a (multi-)material variable of type T residing on CELLs

Parameters
[in]srcvarnameVariable name on source mesh
[in]trgvarnameVariable name on target mesh
[in]lower_boundLower bound for variable
[in]upper_boundUpper bound for variable
[in]limiterLimiter to use for second order reconstruction
[in]bnd_limiterBoundary limiter to use for second order reconstruction
[in]partial_fixup_typeMethod to populate fields on partially filled target entities (cells or dual cells)
[in]empty_fixup_typeMethod to populate fields on empty target entities (cells or dual cells)
[in]conservation_tolTolerance to which source and target integral quantities are to be matched
[in]max_fixup_iterMax number of iterations for global repair

See support/portage.h for options on limiter, partial_fixup_type and empty_fixup_type

Since this call explicitly takes intersection weights we don't have to check if intersection step is complete

◆ interpolate_mesh_var()

template<int D, 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, Entity_kind ONWHAT, template< int, Entity_kind, class, class, class, class, class, template< class, int, class, class > class, class, class, class > class Interpolate>
void Portage::UberDriver< D, 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_in,
lower_bound,
upper_bound,
Limiter_type  limiter,
Boundary_Limiter_type  bnd_limiter,
Partial_fixup_type  partial_fixup_type,
Empty_fixup_type  empty_fixup_type,
double  conservation_tol,
int  max_fixup_iter 
)
inline

Interpolate a mesh variable of type T residing on entity kind ONWHAT using previously computed intersection weights

Template Parameters
Ttype of variable
ONWHATEntity_kind that field resides on
InterpolateFunctor for doing the interpolate from mesh to mesh
Parameters
[in]srcvarnameVariable name on source mesh
[in]trgvarnameVariable name on target mesh
[in]lower_boundLower bound for variable
[in]upper_boundUpper bound for variable
[in]limiterLimiter to use for second order reconstruction
[in]bnd_limiterBoundary limiter to use for second order reconstruction
[in]partial_fixup_typeMethod to populate fields on partially filled target entities (cells or dual cells)
[in]empty_fixup_typeMethod to populate fields on empty target entities (cells or dual cells)
[in]conservation_tolTolerance to which source and target integral quantities are to be matched
[in]max_fixup_iterMax number of iterations for global repair

See support/portage.h for options on limiter, partial_fixup_type and empty_fixup_type

Since this call explicitly takes intersection weights we don't have to check if intersection step is complete

◆ intersect_materials()

template<int D, 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<Portage::Weights_t> > > Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::intersect_materials ( Portage::vector< std::vector< int >> const &  candidates)
inline

◆ intersect_meshes()

template<int D, 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<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> > Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::intersect_meshes ( Portage::vector< std::vector< int >> const &  candidates)
inline

intersect target control volumes with source control volumes

Template Parameters
Entity_kindwhat kind of entities are we intersecting
Intersectintersect functor
Parameters
[in]candidatesIntersection candidates for each target cell
Returns
vector of weights for each target cell

◆ is_distributed_run()

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::is_distributed_run ( Wonton::Executor_type const *  executor = nullptr)
inline

Is this a distributed (multi-rank) run?

◆ operator=()

template<int D, 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>
UberDriver& Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::operator= ( const UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys > &  )
delete

Assignment operator (disabled)

◆ search()

template<int D, 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<Entity_kind ONWHAT, template< int, Entity_kind, class, class > class Search>
Portage::vector<std::vector<int> > Portage::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::search ( )
inline

search for candidate source entities whose control volumes (cells, dual cells) overlap the control volumes of target cells

Template Parameters
Entity_kindwhat kind of entity are we searching on/for
Searchsearch functor
Returns
vector of candidate cells for each target cell

◆ set_num_tols() [1/2]

template<int D, 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::UberDriver< D, 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 numerical tolerances for small distances and volumes in core driver.

Template Parameters
Entity_kindwhat kind of entity are we setting for
Parameters
min_absolute_distanceselected minimal distance
min_absolute_volumeselected minimal volume

◆ set_num_tols() [2/2]

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::set_num_tols ( const NumericTolerances_t num_tols)
inline

set numerical tolerances in core driver

Template Parameters
Entity_kindwhat kind of entity are we setting for
Parameters
num_tolsstruct of selected numerical tolerances

◆ source_needs_redistribution()

template<int D, 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::UberDriver< D, SourceMesh, SourceState, TargetMesh, TargetState, InterfaceReconstructorType, Matpoly_Splitter, Matpoly_Clipper, CoordSys >::source_needs_redistribution ( Wonton::Executor_type const *  executor = nullptr)
inline

Does the source mesh need redistribution due to geometric mismatch of partitons (different from mismatch of overall domain geometry)


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