WarpX
Loading...
Searching...
No Matches
ablastr::fields::MultiFabRegister Struct Reference

#include <MultiFabRegister.H>

Public Member Functions

 MultiFabRegister ()=default
 
 MultiFabRegister (MultiFabRegister const &)=delete
 
 MultiFabRegister (MultiFabRegister &&)=delete
 
MultiFabRegisteroperator= (MultiFabRegister const &)=delete
 
MultiFabRegisteroperator= (MultiFabRegister &&)=delete
 
 ~MultiFabRegister ()=default
 
template<typename T>
amrex::MultiFaballoc_init (T name, int level, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, int ncomp, amrex::IntVect const &ngrow, std::optional< amrex::Real const > initial_value=std::nullopt, bool remake=true, bool redistribute_on_remake=true)
 
template<typename T>
amrex::MultiFaballoc_init (T name, Direction dir, int level, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, int ncomp, amrex::IntVect const &ngrow, std::optional< amrex::Real const > initial_value=std::nullopt, bool remake=true, bool redistribute_on_remake=true)
 
template<typename N, typename A>
amrex::MultiFabalias_init (N new_name, A alias_name, int level, std::optional< amrex::Real const > initial_value=std::nullopt)
 
template<typename N, typename A>
amrex::MultiFabalias_init (N new_name, A alias_name, Direction dir, int level, std::optional< amrex::Real const > initial_value=std::nullopt)
 
template<typename T>
bool has (T name, int level) const
 
template<typename T>
bool has (T name, Direction dir, int level) const
 
template<typename T>
bool has_vector (T name, int level) const
 
template<typename T>
amrex::MultiFabget (T name, int level)
 
template<typename T>
amrex::MultiFabget (T name, Direction dir, int level)
 
template<typename T>
amrex::MultiFab const * get (T name, int level) const
 
template<typename T>
amrex::MultiFab const * get (T name, Direction dir, int level) const
 
template<typename T>
MultiLevelScalarField get_mr_levels (T name, int finest_level, bool skip_level_0=false)
 
template<typename T>
ConstMultiLevelScalarField get_mr_levels (T name, int finest_level, bool skip_level_0=false) const
 
template<typename T>
VectorField get_alldirs (T name, int level)
 
template<typename T>
ConstVectorField get_alldirs (T name, int level) const
 
template<typename T>
MultiLevelVectorField get_mr_levels_alldirs (T name, int finest_level, bool skip_level_0=false)
 
template<typename T>
ConstMultiLevelVectorField get_mr_levels_alldirs (T name, int finest_level, bool skip_level_0=false) const
 
std::vector< std::string > list () const
 
template<typename T>
void erase (T name, int level)
 
template<typename T>
void erase (T name, Direction dir, int level)
 
void clear_level (int level)
 
void remake_level (int other_level, amrex::DistributionMapping const &new_dm)
 
std::string mf_name (std::string name, int level) const
 
std::string mf_name (std::string name, Direction dir, int level) const
 
bool internal_has (std::string const &internal_name)
 
amrex::MultiFabinternal_get (std::string const &internal_name)
 

Static Public Attributes

static std::vector< Directionm_all_dirs
 

Private Member Functions

amrex::MultiFab const * internal_get (std::string const &internal_name) const
 
amrex::MultiFabinternal_alloc_init (std::string const &name, int level, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, int ncomp, amrex::IntVect const &ngrow, std::optional< amrex::Real const > initial_value=std::nullopt, bool remake=true, bool redistribute_on_remake=true)
 
amrex::MultiFabinternal_alloc_init (std::string const &name, Direction dir, int level, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, int ncomp, amrex::IntVect const &ngrow, std::optional< amrex::Real const > initial_value=std::nullopt, bool remake=true, bool redistribute_on_remake=true)
 
amrex::MultiFabinternal_alias_init (std::string const &new_name, std::string const &alias_name, int level, std::optional< amrex::Real const > initial_value=std::nullopt)
 
amrex::MultiFabinternal_alias_init (std::string const &new_name, std::string const &alias_name, Direction dir, int level, std::optional< amrex::Real const > initial_value=std::nullopt)
 
bool internal_has (std::string const &name, int level) const
 
bool internal_has (std::string const &name, Direction dir, int level) const
 
bool internal_has_vector (std::string const &name, int level) const
 
amrex::MultiFabinternal_get (std::string const &name, int level)
 
amrex::MultiFab const * internal_get (std::string const &name, int level) const
 
amrex::MultiFabinternal_get (std::string const &name, Direction dir, int level)
 
amrex::MultiFab const * internal_get (std::string const &name, Direction dir, int level) const
 
MultiLevelScalarField internal_get_mr_levels (std::string const &name, int finest_level, bool skip_level_0)
 
ConstMultiLevelScalarField internal_get_mr_levels (std::string const &name, int finest_level, bool skip_level_0) const
 
VectorField internal_get_alldirs (std::string const &name, int level)
 
ConstVectorField internal_get_alldirs (std::string const &name, int level) const
 
MultiLevelVectorField internal_get_mr_levels_alldirs (std::string const &name, int finest_level, bool skip_level_0)
 
ConstMultiLevelVectorField internal_get_mr_levels_alldirs (std::string const &name, int finest_level, bool skip_level_0) const
 
void internal_erase (std::string const &name, int level)
 
void internal_erase (std::string const &name, Direction dir, int level)
 

Private Attributes

std::map< std::string, MultiFabOwnerm_mf_register
 

Detailed Description

This is a register of fields aka amrex::MultiFabs.

This is owned by a simulation instance. All used fields should be registered here. Internally, this contains

See also
MultiFabOwner values.

Constructor & Destructor Documentation

◆ MultiFabRegister() [1/3]

ablastr::fields::MultiFabRegister::MultiFabRegister ( )
default

◆ MultiFabRegister() [2/3]

ablastr::fields::MultiFabRegister::MultiFabRegister ( MultiFabRegister const & )
delete

◆ MultiFabRegister() [3/3]

ablastr::fields::MultiFabRegister::MultiFabRegister ( MultiFabRegister && )
delete

◆ ~MultiFabRegister()

ablastr::fields::MultiFabRegister::~MultiFabRegister ( )
default

Member Function Documentation

◆ alias_init() [1/2]

template<typename N, typename A>
amrex::MultiFab * ablastr::fields::MultiFabRegister::alias_init ( N new_name,
A alias_name,
Direction dir,
int level,
std::optional< amrex::Real const > initial_value = std::nullopt )
inline

Create an alias of a MultiFab (field)

Registers a new name for an existing MultiFab (field) and optionally assigning a value.

Parameters
new_namenew name
alias_nameowner name to alias
dirthe field component for vector fields ("direction" of the unit vector) both in the alias and aliased
levelthe MR level to represent
initial_valuethe optional value to assign
Returns
the newly aliased MultiFab

◆ alias_init() [2/2]

template<typename N, typename A>
amrex::MultiFab * ablastr::fields::MultiFabRegister::alias_init ( N new_name,
A alias_name,
int level,
std::optional< amrex::Real const > initial_value = std::nullopt )
inline

Create an alias of a MultiFab (field)

Registers a new name for an existing MultiFab (field) and optionally assigning a value.

Parameters
new_namenew name
alias_nameowner name to alias
levelthe MR level to represent
initial_valuethe optional value to assign
Returns
the newly aliased MultiFab

◆ alloc_init() [1/2]

template<typename T>
amrex::MultiFab * ablastr::fields::MultiFabRegister::alloc_init ( T name,
Direction dir,
int level,
amrex::BoxArray const & ba,
amrex::DistributionMapping const & dm,
int ncomp,
amrex::IntVect const & ngrow,
std::optional< amrex::Real const > initial_value = std::nullopt,
bool remake = true,
bool redistribute_on_remake = true )
inline

Allocate and optionally initialize a MultiFab (field)

This registers a new MultiFab under a unique name, allocates it and optionally assigns it an initial value.

Parameters
namea unique name for this field
dirthe field component for vector fields ("direction" of the unit vector)
levelthe MR level to represent
bathe list of boxes to cover the field
dmthe distribution mapping for load balancing with MPI
ncompthe number of components of the field (all with the same staggering)
ngrowthe number of guard (ghost, halo) cells
initial_valuethe optional initial value
remakefollow the default domain decomposition of the simulation
redistribute_on_remakeredistribute on
See also
amrex::AmrCore::RemakeLevel
Returns
pointer to newly allocated MultiFab

◆ alloc_init() [2/2]

template<typename T>
amrex::MultiFab * ablastr::fields::MultiFabRegister::alloc_init ( T name,
int level,
amrex::BoxArray const & ba,
amrex::DistributionMapping const & dm,
int ncomp,
amrex::IntVect const & ngrow,
std::optional< amrex::Real const > initial_value = std::nullopt,
bool remake = true,
bool redistribute_on_remake = true )
inline

Allocate and optionally initialize a MultiFab (field)

This registers a new MultiFab under a unique name, allocates it and optionally assigns it an initial value.

Parameters
namea unique name for this field
levelthe MR level to represent
bathe list of boxes to cover the field
dmthe distribution mapping for load balancing with MPI
ncompthe number of components of the field (all with the same staggering)
ngrowthe number of guard (ghost, halo) cells
initial_valuethe optional initial value
remakefollow the default domain decomposition of the simulation
redistribute_on_remakeredistribute on
See also
amrex::AmrCore::RemakeLevel
Returns
pointer to newly allocated MultiFab

◆ clear_level()

void ablastr::fields::MultiFabRegister::clear_level ( int level)

Erase all MultiFabs on a specific MR level.

Calls

See also
erase for all MultiFabs on a specific level.
Parameters
levelthe MR level to erase all MultiFabs from

◆ erase() [1/2]

template<typename T>
void ablastr::fields::MultiFabRegister::erase ( T name,
Direction dir,
int level )
inline

Deallocate and remove a vector field component.

Parameters
namethe name of the field
dirthe field component for vector fields ("direction" of the unit vector)
levelthe MR level

◆ erase() [2/2]

template<typename T>
void ablastr::fields::MultiFabRegister::erase ( T name,
int level )
inline

Deallocate and remove a scalar field.

Parameters
namethe name of the field
levelthe MR level

◆ get() [1/4]

template<typename T>
amrex::MultiFab * ablastr::fields::MultiFabRegister::get ( T name,
Direction dir,
int level )
inlinenodiscard

Return a MultiFab that is part of a vector/tensor field.

This throws a runtime error if the requested field is not present.

Parameters
namethe name of the field
dirthe field component for vector fields ("direction" of the unit vector)
levelthe MR level
Returns
a non-owning pointer to the MultiFab (field)

◆ get() [2/4]

template<typename T>
amrex::MultiFab const * ablastr::fields::MultiFabRegister::get ( T name,
Direction dir,
int level ) const
inlinenodiscard

Return a MultiFab that is part of a vector/tensor field.

This throws a runtime error if the requested field is not present.

Parameters
namethe name of the field
dirthe field component for vector fields ("direction" of the unit vector)
levelthe MR level
Returns
a non-owning pointer to the MultiFab (field)

◆ get() [3/4]

template<typename T>
amrex::MultiFab * ablastr::fields::MultiFabRegister::get ( T name,
int level )
inlinenodiscard

Return a scalar MultiFab (field).

This throws a runtime error if the requested field is not present.

Parameters
namethe name of the field
levelthe MR level
Returns
a non-owning pointer to the MultiFab (field)

◆ get() [4/4]

template<typename T>
amrex::MultiFab const * ablastr::fields::MultiFabRegister::get ( T name,
int level ) const
inlinenodiscard

Return a scalar MultiFab (field).

This throws a runtime error if the requested field is not present.

Parameters
namethe name of the field
levelthe MR level
Returns
a non-owning pointer to the MultiFab (field)

◆ get_alldirs() [1/2]

template<typename T>
VectorField ablastr::fields::MultiFabRegister::get_alldirs ( T name,
int level )
inlinenodiscard

title

Same as get above, but returns all levels for a name.

Parameters
namethe name of the field
levelthe MR level
Returns
non-owning pointers to all components of a vector field

◆ get_alldirs() [2/2]

template<typename T>
ConstVectorField ablastr::fields::MultiFabRegister::get_alldirs ( T name,
int level ) const
inlinenodiscard

◆ get_mr_levels() [1/2]

template<typename T>
MultiLevelScalarField ablastr::fields::MultiFabRegister::get_mr_levels ( T name,
int finest_level,
bool skip_level_0 = false )
inlinenodiscard

Return the MultiFab of a scalar field on all MR levels.

This throws a runtime error if the requested field is not present.

Parameters
namethe name of the field
finest_levelthe highest MR level to return
skip_level_0return a nullptr for level 0
Returns
non-owning pointers to the MultiFab (field) on all levels

◆ get_mr_levels() [2/2]

template<typename T>
ConstMultiLevelScalarField ablastr::fields::MultiFabRegister::get_mr_levels ( T name,
int finest_level,
bool skip_level_0 = false ) const
inlinenodiscard

◆ get_mr_levels_alldirs() [1/2]

template<typename T>
MultiLevelVectorField ablastr::fields::MultiFabRegister::get_mr_levels_alldirs ( T name,
int finest_level,
bool skip_level_0 = false )
inlinenodiscard

Return a vector field on all MR levels.

Out loop: MR levels. Inner loop: directions (components).

Parameters
namethe name of the field
finest_levelthe highest MR level to return
skip_level_0return a nullptr for level 0
Returns
non-owning pointers to all components of a vector field on all MR levels

◆ get_mr_levels_alldirs() [2/2]

template<typename T>
ConstMultiLevelVectorField ablastr::fields::MultiFabRegister::get_mr_levels_alldirs ( T name,
int finest_level,
bool skip_level_0 = false ) const
inlinenodiscard

◆ has() [1/2]

template<typename T>
bool ablastr::fields::MultiFabRegister::has ( T name,
Direction dir,
int level ) const
inlinenodiscard

Check if a MultiFab that is part of a vector/tensor field is registered.

Parameters
namethe name to check if registered
dirthe field component for vector fields ("direction" of the unit vector)
levelthe MR level to check
Returns
true if contained, otherwise false

◆ has() [2/2]

template<typename T>
bool ablastr::fields::MultiFabRegister::has ( T name,
int level ) const
inlinenodiscard

Check if a scalar MultiFab (field) is registered.

Parameters
namethe name to check if registered
levelthe MR level to check
Returns
true if contained, otherwise false

◆ has_vector()

template<typename T>
bool ablastr::fields::MultiFabRegister::has_vector ( T name,
int level ) const
inlinenodiscard

Check if a MultiFab vector field is registered.

Parameters
namethe name to check if registered
levelthe MR level to check
Returns
true if contained, otherwise false

◆ internal_alias_init() [1/2]

amrex::MultiFab * ablastr::fields::MultiFabRegister::internal_alias_init ( std::string const & new_name,
std::string const & alias_name,
Direction dir,
int level,
std::optional< amrex::Real const > initial_value = std::nullopt )
private

◆ internal_alias_init() [2/2]

amrex::MultiFab * ablastr::fields::MultiFabRegister::internal_alias_init ( std::string const & new_name,
std::string const & alias_name,
int level,
std::optional< amrex::Real const > initial_value = std::nullopt )
private

◆ internal_alloc_init() [1/2]

amrex::MultiFab * ablastr::fields::MultiFabRegister::internal_alloc_init ( std::string const & name,
Direction dir,
int level,
amrex::BoxArray const & ba,
amrex::DistributionMapping const & dm,
int ncomp,
amrex::IntVect const & ngrow,
std::optional< amrex::Real const > initial_value = std::nullopt,
bool remake = true,
bool redistribute_on_remake = true )
private

◆ internal_alloc_init() [2/2]

amrex::MultiFab * ablastr::fields::MultiFabRegister::internal_alloc_init ( std::string const & name,
int level,
amrex::BoxArray const & ba,
amrex::DistributionMapping const & dm,
int ncomp,
amrex::IntVect const & ngrow,
std::optional< amrex::Real const > initial_value = std::nullopt,
bool remake = true,
bool redistribute_on_remake = true )
private

◆ internal_erase() [1/2]

void ablastr::fields::MultiFabRegister::internal_erase ( std::string const & name,
Direction dir,
int level )
private

◆ internal_erase() [2/2]

void ablastr::fields::MultiFabRegister::internal_erase ( std::string const & name,
int level )
private

◆ internal_get() [1/6]

amrex::MultiFab * ablastr::fields::MultiFabRegister::internal_get ( std::string const & internal_name)
nodiscard

◆ internal_get() [2/6]

amrex::MultiFab const * ablastr::fields::MultiFabRegister::internal_get ( std::string const & internal_name) const
nodiscardprivate

◆ internal_get() [3/6]

amrex::MultiFab * ablastr::fields::MultiFabRegister::internal_get ( std::string const & name,
Direction dir,
int level )
nodiscardprivate

◆ internal_get() [4/6]

amrex::MultiFab const * ablastr::fields::MultiFabRegister::internal_get ( std::string const & name,
Direction dir,
int level ) const
nodiscardprivate

◆ internal_get() [5/6]

amrex::MultiFab * ablastr::fields::MultiFabRegister::internal_get ( std::string const & name,
int level )
nodiscardprivate

◆ internal_get() [6/6]

amrex::MultiFab const * ablastr::fields::MultiFabRegister::internal_get ( std::string const & name,
int level ) const
nodiscardprivate

◆ internal_get_alldirs() [1/2]

VectorField ablastr::fields::MultiFabRegister::internal_get_alldirs ( std::string const & name,
int level )
nodiscardprivate

◆ internal_get_alldirs() [2/2]

ConstVectorField ablastr::fields::MultiFabRegister::internal_get_alldirs ( std::string const & name,
int level ) const
nodiscardprivate

◆ internal_get_mr_levels() [1/2]

MultiLevelScalarField ablastr::fields::MultiFabRegister::internal_get_mr_levels ( std::string const & name,
int finest_level,
bool skip_level_0 )
nodiscardprivate

◆ internal_get_mr_levels() [2/2]

ConstMultiLevelScalarField ablastr::fields::MultiFabRegister::internal_get_mr_levels ( std::string const & name,
int finest_level,
bool skip_level_0 ) const
nodiscardprivate

◆ internal_get_mr_levels_alldirs() [1/2]

MultiLevelVectorField ablastr::fields::MultiFabRegister::internal_get_mr_levels_alldirs ( std::string const & name,
int finest_level,
bool skip_level_0 )
nodiscardprivate

◆ internal_get_mr_levels_alldirs() [2/2]

ConstMultiLevelVectorField ablastr::fields::MultiFabRegister::internal_get_mr_levels_alldirs ( std::string const & name,
int finest_level,
bool skip_level_0 ) const
nodiscardprivate

◆ internal_has() [1/3]

bool ablastr::fields::MultiFabRegister::internal_has ( std::string const & internal_name)
nodiscard

Temporary test function for legacy Python bindings

◆ internal_has() [2/3]

bool ablastr::fields::MultiFabRegister::internal_has ( std::string const & name,
Direction dir,
int level ) const
nodiscardprivate

◆ internal_has() [3/3]

bool ablastr::fields::MultiFabRegister::internal_has ( std::string const & name,
int level ) const
nodiscardprivate

◆ internal_has_vector()

bool ablastr::fields::MultiFabRegister::internal_has_vector ( std::string const & name,
int level ) const
nodiscardprivate

◆ list()

std::vector< std::string > ablastr::fields::MultiFabRegister::list ( ) const
nodiscard

List the internal names of all registered fields.

Returns
all currently allocated and registered fields

◆ mf_name() [1/2]

std::string ablastr::fields::MultiFabRegister::mf_name ( std::string name,
Direction dir,
int level ) const
nodiscard

Create the register name of vector field, component direction and MR level

Parameters
namethe name of the field
dirthe field component for vector fields ("direction" of the unit vector)
levelthe MR level
Returns
internal name of the field in the register

◆ mf_name() [2/2]

std::string ablastr::fields::MultiFabRegister::mf_name ( std::string name,
int level ) const
nodiscard

Create the register name of scalar field and MR level

Parameters
namethe name of the field
levelthe MR level
Returns
internal name of the field in the register

◆ operator=() [1/2]

MultiFabRegister & ablastr::fields::MultiFabRegister::operator= ( MultiFabRegister && )
delete

◆ operator=() [2/2]

MultiFabRegister & ablastr::fields::MultiFabRegister::operator= ( MultiFabRegister const & )
delete

◆ remake_level()

void ablastr::fields::MultiFabRegister::remake_level ( int other_level,
amrex::DistributionMapping const & new_dm )

Remake all (i)MultiFab with a new distribution mapping.

If redistribute is true, we also copy from the old data into the new.

Parameters
other_levelthe MR level to erase all MultiFabs from
new_dmnew distribution mapping

Member Data Documentation

◆ m_all_dirs

std::vector<Direction> ablastr::fields::MultiFabRegister::m_all_dirs
inlinestatic
Initial value:
=
static const Direction z
Definition MultiFabRegister.H:89
static const Direction y
Definition MultiFabRegister.H:86
static const Direction x
Definition MultiFabRegister.H:85

The directions of a vector field as stored in the simulation.

Cartesian: always x,y,z in any ND RZ: r,t(heta),z (r,z over azimuthal modes) RCylinder: r,t(heta),z RSphere: r,t(heta),p(hi)

See also
https://warpx.readthedocs.io/en/latest/developers/dimensionality.html

◆ m_mf_register

std::map< std::string, MultiFabOwner > ablastr::fields::MultiFabRegister::m_mf_register
private

data storage: ownership and lifetime control


The documentation for this struct was generated from the following files: