Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState > Class Template Reference

Manages source and target sub-meshes for part-by-part remap. It detects boundaries mismatch and provides the necessary fixup for partially filled and empty cells values. More...

#include <parts.h>

Public Member Functions

 PartPair ()=default
 Construct a default source-target mesh parts pair. More...
 
 PartPair (SourceMesh const &source_mesh, SourceState &source_state, TargetMesh const &target_mesh, TargetState &target_state, std::vector< int > const &source_entities, std::vector< int > const &target_entities, Wonton::Executor_type const *executor)
 Construct a source-target mesh parts pair. More...
 
 ~PartPair ()=default
 Delete a source-target mesh parts pair. More...
 
bool has_mismatch () const
 Do source and target meshes have a boundary mismatch? More...
 
bool mismatch_tested () const
 Is mismatch already tested? More...
 
const SourcePartsource () const
 Retrieve a pointer to source mesh part. More...
 
const TargetParttarget () const
 Retrieve a pointer to target mesh part. More...
 
double compute_intersect_volumes (Portage::vector< entity_weights_t > const &source_weights)
 Compute source and target parts intersection volume. More...
 
bool check_mismatch (Portage::vector< entity_weights_t > const &source_weights)
 Check and fix source/target boundaries mismatch. More...
 
bool fix_mismatch (std::string src_var_name, std::string trg_var_name, double global_lower_bound=-infinity_, double global_upper_bound=infinity_, double conservation_tol=tolerance_, int maxiter=5, Partial_fixup_type partial_fixup_type=SHIFTED_CONSERVATIVE, Empty_fixup_type empty_fixup_type=EXTRAPOLATE) const
 Repair the remapped field to account for boundary mismatch. More...
 
bool fix_mismatch_meshvar (std::string const &src_var_name, std::string const &trg_var_name, double global_lower_bound, double global_upper_bound, double conservation_tol, int maxiter, Partial_fixup_type partial_fixup_type, Empty_fixup_type empty_fixup_type) const
 Repair a remapped mesh field to account for boundary mismatch. More...
 

Detailed Description

template<int D, class SourceMesh, class SourceState, class TargetMesh = SourceMesh, class TargetState = SourceState>
class Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >

Manages source and target sub-meshes for part-by-part remap. It detects boundaries mismatch and provides the necessary fixup for partially filled and empty cells values.

Template Parameters
Dthe problem dimension.
SourceMeshthe source mesh wrapper to use
SourceStatethe source state wrapper to use
TargetMeshthe target mesh wrapper to use
TargetStatethe target state wrapper to use

Constructor & Destructor Documentation

◆ PartPair() [1/2]

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::PartPair ( )
default

Construct a default source-target mesh parts pair.

◆ PartPair() [2/2]

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::PartPair ( SourceMesh const &  source_mesh,
SourceState &  source_state,
TargetMesh const &  target_mesh,
TargetState &  target_state,
std::vector< int > const &  source_entities,
std::vector< int > const &  target_entities,
Wonton::Executor_type const *  executor 
)
inline

Construct a source-target mesh parts pair.

Parameters
source_meshthe source mesh
target_meshthe target mesh
source_statethe source mesh data
target_statethe target mesh data
source_entitiesthe list of source entities to remap
target_entitiesthe list of target entities to remap
executorthe MPI executor to use

◆ ~PartPair()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::~PartPair ( )
default

Delete a source-target mesh parts pair.

Member Function Documentation

◆ check_mismatch()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
bool Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::check_mismatch ( Portage::vector< entity_weights_t > const &  source_weights)
inline

Check and fix source/target boundaries mismatch.

T is the target mesh, S is the source mesh T_i is the i'th cell of the target mesh and |*| signifies the extent/volume of an entity

  • if sum_i(|T_i intersect S|) NE |S|, some source cells are not completely covered by the target mesh.
  • if sum_i(|T_i intersect S|) NE |T|, some target cells are not completely covered by the source mesh
  • if there is a mismatch, adjust values for fields to account so that integral quantities are conserved and a constant field is readjusted to a different constant

WARNING: 'source_ents_and_weights' is a GLOBAL list defined on the entire target mesh.

Parameters
source...source entities ID and weights for each target entity.
Returns
true if a mismatch has been identified, false otherwise.

◆ compute_intersect_volumes()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
double Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::compute_intersect_volumes ( Portage::vector< entity_weights_t > const &  source_weights)
inline

Compute source and target parts intersection volume.

Parameters
source_weightscandidate source cells and their intersection moments.
Returns
the total intersection volume.

◆ fix_mismatch()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
bool Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::fix_mismatch ( std::string  src_var_name,
std::string  trg_var_name,
double  global_lower_bound = -infinity_,
double  global_upper_bound = infinity_,
double  conservation_tol = tolerance_,
int  maxiter = 5,
Partial_fixup_type  partial_fixup_type = SHIFTED_CONSERVATIVE,
Empty_fixup_type  empty_fixup_type = EXTRAPOLATE 
) const
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
conservation_tolconservation tolerance treshold
maxitermax number of iterations
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)

◆ fix_mismatch_meshvar()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
bool Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::fix_mismatch_meshvar ( std::string const &  src_var_name,
std::string const &  trg_var_name,
double  global_lower_bound,
double  global_upper_bound,
double  conservation_tol,
int  maxiter,
Partial_fixup_type  partial_fixup_type,
Empty_fixup_type  empty_fixup_type 
) const
inline

Repair a remapped mesh field to account for boundary mismatch.

Parameters
src_var_namefield variable on source mesh
trg_var_namefield variable on source mesh
global_lower_boundlower limit on variable value
global_upper_boundupper limit on variable value
conservation_tolconservation tolerance treshold
maxitermax number of iterations
partial_fixup_typetype of fixup in case of partial mismatch
empty_fixup_typetype of fixup in empty target entities
Returns
true if correctly fixed, false otherwise.

◆ has_mismatch()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
bool Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::has_mismatch ( ) const
inline

Do source and target meshes have a boundary mismatch?

Returns
true if so, false otherwise.

◆ mismatch_tested()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
bool Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::mismatch_tested ( ) const
inline

Is mismatch already tested?

Returns
true if so, false otherwise.

◆ source()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
const SourcePart& Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::source ( ) const
inline

Retrieve a pointer to source mesh part.

Returns
a pointer to source mesh part

◆ target()

template<int D, class SourceMesh , class SourceState , class TargetMesh = SourceMesh, class TargetState = SourceState>
const TargetPart& Portage::PartPair< D, SourceMesh, SourceState, TargetMesh, TargetState >::target ( ) const
inline

Retrieve a pointer to target mesh part.

Returns
a pointer to target mesh part

The documentation for this class was generated from the following file:
  • /home/portage/portage/portage/driver/parts.h