Provides an interface to remap variables from one swarm to another.
More...
#include <driver_swarm.h>
|
| SwarmDriver (SourceSwarm &source_swarm, SourceState &source_state, TargetSwarm const &target_swarm, TargetState &target_state, SmoothingLengths const &smoothing_lengths, Weight::Kernel const kernel_type=Weight::B4, Weight::Geometry const support_geom_type=Weight::ELLIPTIC, WeightCenter const center=Gather) |
| Constructor using a unique kernel_type and support geometry for all particles. More...
|
|
| SwarmDriver (SourceSwarm &source_swarm, SourceState &source_state, TargetSwarm const &target_swarm, TargetState &target_state, SmoothingLengths const &smoothing_lengths, Portage::vector< Weight::Kernel > const &kernel_types, Portage::vector< Weight::Geometry > const &geom_types, WeightCenter const center=Gather) |
| Constructor with search extents computed inline. More...
|
|
| SwarmDriver (SourceSwarm &source_swarm, SourceState &source_state, TargetSwarm const &target_swarm, TargetState &target_state, SmoothingLengths const &smoothing_lengths, Portage::vector< Point< dim >> const &source_extents, Portage::vector< Point< dim >> const &target_extents, WeightCenter const center=Gather) |
| Constructor with search extents specified. More...
|
|
| SwarmDriver (const SwarmDriver &)=delete |
| Disabled copy constructor. More...
|
|
SwarmDriver & | operator= (const SwarmDriver &)=delete |
| Disabled assignment operator. More...
|
|
| ~SwarmDriver ()=default |
| Destructor. More...
|
|
void | set_remap_var_names (std::vector< std::string > const &remap_var_names) |
| Specify the names of the variables to be interpolated. More...
|
|
void | set_remap_var_names (std::vector< std::string > const &source_vars, std::vector< std::string > const &target_vars, EstimateType const estimator_type=LocalRegression, basis::Type const basis_type=basis::Unitary, oper::Type const operator_spec=oper::LastOperator, Portage::vector< oper::Domain > const &operator_domains={}, Portage::vector< std::vector< Point< dim >>> const &operator_data={}, std::string part_field="NONE", double part_tolerance=0.0, SmoothingLengths const &part_smoothing={}) |
| Specify the names of the variables to be interpolated. More...
|
|
std::vector< std::string > | source_vars () const |
| Get all source swarm variables names. More...
|
|
std::vector< std::string > | target_vars () const |
| Get all target swarm variables names. More...
|
|
void | run (Wonton::Executor_type const *executor=nullptr, bool report_time=true) |
| Perform the remap. More...
|
|
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm, class SourceState, class TargetSwarm = SourceSwarm, class TargetState = SourceState>
class Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >
Provides an interface to remap variables from one swarm to another.
- Template Parameters
-
Search | the particle search type to use. |
Accumulate | the particle accumulate type to use. |
Estimate | the particle estimate type to use. |
dim | spatial dimension of the problem. |
SourceSwarm | the wrapper type of the input source swarm. |
SourceState | the wrapper type of the input source state. |
TargetSwarm | the wrapper type of the input target swarm. |
TargetState | the wrapper type of the input target state. |
◆ SwarmDriver() [1/4]
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::SwarmDriver |
( |
SourceSwarm & |
source_swarm, |
|
|
SourceState & |
source_state, |
|
|
TargetSwarm const & |
target_swarm, |
|
|
TargetState & |
target_state, |
|
|
SmoothingLengths const & |
smoothing_lengths, |
|
|
Weight::Kernel const |
kernel_type = Weight::B4 , |
|
|
Weight::Geometry const |
support_geom_type = Weight::ELLIPTIC , |
|
|
WeightCenter const |
center = Gather |
|
) |
| |
|
inline |
Constructor using a unique kernel_type and support geometry for all particles.
- Parameters
-
source_swarm | a reference to the source swarm |
source_state | a reference to the source state |
target_swarm | a reference to the target swarm |
target_state | a reference to the target state |
smoothing_lengths | a vector of smoothing lengths for each target particle - the three levels of vectors allow for a polygonal/polyhedral shaped support around each target particle and a vector of smoothing lengths for each facet of the polygonal/polyhedral support. |
kernel_types | types of weight kernels (B4, SQUARE, etc) to use. |
geom_types | geometry of the supports (ELLIPTIC, TENSOR, FACETED) of the weight functions to use. |
center | specify whether gather-form or scatter-form weights are used. |
NOTE: 'smoothing_lengths' must have the size of 'source_swarm' if 'scatter' and that of 'target_swarm' if 'gather'.
◆ SwarmDriver() [2/4]
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::SwarmDriver |
( |
SourceSwarm & |
source_swarm, |
|
|
SourceState & |
source_state, |
|
|
TargetSwarm const & |
target_swarm, |
|
|
TargetState & |
target_state, |
|
|
SmoothingLengths const & |
smoothing_lengths, |
|
|
Portage::vector< Weight::Kernel > const & |
kernel_types, |
|
|
Portage::vector< Weight::Geometry > const & |
geom_types, |
|
|
WeightCenter const |
center = Gather |
|
) |
| |
|
inline |
Constructor with search extents computed inline.
This is used for any weights except faceted.
- Parameters
-
source_swarm | a reference to the source swarm |
source_state | a reference to the source state |
target_swarm | a reference to the target swarm |
target_state | a reference to the target state |
smoothing_lengths | a vector of smoothing lengths for each target particle - the three levels of vectors allow for a polygonal/polyhedral shaped support around each target particle and a vector of smoothing lengths for each facet of the polygonal/polyhedral support. |
kernel_types | types of weight kernels (B4, SQUARE, etc) to use. |
geom_types | geometry of the supports (ELLIPTIC, TENSOR, FACETED) of the weight functions to use. |
center | specify whether gather-form or scatter-form weights are used. |
NOTE: 'smoothing_lengths' must have the size of 'source_swarm' if 'scatter' and that of 'target_swarm' if 'gather'.
◆ SwarmDriver() [3/4]
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::SwarmDriver |
( |
SourceSwarm & |
source_swarm, |
|
|
SourceState & |
source_state, |
|
|
TargetSwarm const & |
target_swarm, |
|
|
TargetState & |
target_state, |
|
|
SmoothingLengths const & |
smoothing_lengths, |
|
|
Portage::vector< Point< dim >> const & |
source_extents, |
|
|
Portage::vector< Point< dim >> const & |
target_extents, |
|
|
WeightCenter const |
center = Gather |
|
) |
| |
|
inline |
Constructor with search extents specified.
This is used for faceted weights, which will be uniformly used throughout the swarm specified by center. Kernel and geometry type are set internally.
- Parameters
-
source_swarm | a reference to the source swarm |
source_state | a reference to the source state |
target_swarm | a reference to the target swarm |
target_state | a reference to the target state |
smoothing_lengths | a vector of smoothing lengths for each target particle - the three levels of vectors allow for a polygonal/polyhedral shaped support around each target particle and a vector of smoothing lengths for each facet of the polygonal/polyhedral support. |
source_extents | extents for source swarm |
target_extents | extents for target swarm |
center | specify whether gather-form or scatter-form weights are used. |
NOTE: 'smoothing_lengths' must have the size of 'source_swarm' if 'scatter' and that of 'target_swarm' if 'gather'.
◆ SwarmDriver() [4/4]
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::SwarmDriver |
( |
const SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState > & |
| ) |
|
|
delete |
Disabled copy constructor.
◆ ~SwarmDriver()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
◆ check_sizes()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
void Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::check_sizes |
( |
WeightCenter const |
weight_center | ) |
|
|
inlineprotected |
Check sizes according to weight center type.
- Parameters
-
weight_center | weight center type. |
◆ check_sizes_and_set_types()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
Check sizes according to weight center, set kernel and geometry.
- Parameters
-
weight_center | weight center type. |
kernel_type | kernel type |
support_geom_type | support geometry type. |
◆ get_swarm_size()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
Get swarm size according to weight center type.
- Returns
- the number of particles of the swarm.
◆ operator=()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
SwarmDriver& Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::operator= |
( |
const SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState > & |
| ) |
|
|
delete |
Disabled assignment operator.
◆ run()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
void Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::run |
( |
Wonton::Executor_type const * |
executor = nullptr , |
|
|
bool |
report_time = true |
|
) |
| |
|
inline |
Perform the remap.
- Parameters
-
executor | the executor type: serial or parallel. |
report_time | specify if saving timings or not. |
◆ set_extents_from_smoothing_lengths()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
void Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::set_extents_from_smoothing_lengths |
( |
| ) |
|
|
inlineprotected |
Setup search boxes around source and target points for non-faceted weights.
◆ set_remap_var_names() [1/2]
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
void Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::set_remap_var_names |
( |
std::vector< std::string > const & |
remap_var_names | ) |
|
|
inline |
Specify the names of the variables to be interpolated.
- Parameters
-
[in] | remap_var_names | A list of variable names of the variables to interpolate from the source swarm to the target swarm. This variable must exist in both swarms' state manager |
◆ set_remap_var_names() [2/2]
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
void Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::set_remap_var_names |
( |
std::vector< std::string > const & |
source_vars, |
|
|
std::vector< std::string > const & |
target_vars, |
|
|
EstimateType const |
estimator_type = LocalRegression , |
|
|
basis::Type const |
basis_type = basis::Unitary , |
|
|
oper::Type const |
operator_spec = oper::LastOperator , |
|
|
Portage::vector< oper::Domain > const & |
operator_domains = {} , |
|
|
Portage::vector< std::vector< Point< dim >>> const & |
operator_data = {} , |
|
|
std::string |
part_field = "NONE" , |
|
|
double |
part_tolerance = 0.0 , |
|
|
SmoothingLengths const & |
part_smoothing = {} |
|
) |
| |
|
inline |
Specify the names of the variables to be interpolated.
- Parameters
-
source_vars | a list of source swarm variables names. |
target_vars | a list of target swarm variables names. |
estimator_type | estimator to be used (KernelDensity, LocalRegression). |
basis_type | order of the basis used (UNITARY, LINEAR, QUADRATIC). |
operator_spec | operator specification. |
operator_domains | operator domains. |
operator_data | operator data. |
part_field | name of the field for part-by-particle. |
part_tolerance | tolerance threshold for part-by-particle. |
part_smoothing | smoothing lengths matrix for part-by-particle. |
◆ source_vars()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
std::vector<std::string> Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::source_vars |
( |
| ) |
const |
|
inline |
Get all source swarm variables names.
- Returns
- a list of source variables names.
◆ target_vars()
template<template< int, class, class > class Search, template< int, class, class > class Accumulate, template< int, class > class Estimate, int dim, class SourceSwarm , class SourceState , class TargetSwarm = SourceSwarm, class TargetState = SourceState>
std::vector<std::string> Portage::Meshfree::SwarmDriver< Search, Accumulate, Estimate, dim, SourceSwarm, SourceState, TargetSwarm, TargetState >::target_vars |
( |
| ) |
const |
|
inline |
Get all target swarm variables names.
- Returns
- a list of target variables names.
The documentation for this class was generated from the following file: