OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
AmrMultiGridLevel< MatrixType, VectorType > Class Template Reference

#include <AmrMultiGridLevel.h>

Inheritance diagram for AmrMultiGridLevel< MatrixType, VectorType >:
Collaboration diagram for AmrMultiGridLevel< MatrixType, VectorType >:

Public Types

enum  Mask { COVERED = -1 , INTERIOR = 0 , BNDRY = 1 , PHYSBNDRY = 2 }
enum  Refined { YES = 0 , NO = 1 }
typedef amr::AmrField_t AmrField_t
typedef amr::AmrGeometry_t AmrGeometry_t
typedef std::unique_ptr< AmrField_tAmrField_u
typedef std::shared_ptr< AmrField_tAmrField_s
typedef amr::AmrIntVect_t AmrIntVect_t
typedef MatrixType matrix_t
typedef VectorType vector_t
typedef amrex::BaseFab< int > basefab_t
typedef amrex::FabArray< basefab_tmask_t
typedef std::shared_ptr< AmrBoundary< AmrMultiGridLevel< MatrixType, VectorType > > > boundary_t
typedef amr::comm_t comm_t
typedef amr::dmap_t dmap_t
typedef amr::global_ordinal_t go_t
typedef amr::scalar_t scalar_t
typedef amr::local_ordinal_t lo_t
typedef std::vector< go_tindices_t
 Type for matrix indices.
typedef std::vector< scalar_tcoefficients_t
 Type for matrix entries.
typedef std::unordered_map< go_t, scalar_tumap_t

Public Member Functions

 AmrMultiGridLevel (const Vector_t &meshScaling, const amrex::BoxArray &_grids, const amrex::DistributionMapping &_dmap, const AmrGeometry_t &_geom, const AmrIntVect_t &rr, const boundary_t *bc, const Teuchos::RCP< comm_t > &comm)
 ~AmrMultiGridLevel ()
go_t serialize (const AmrIntVect_t &iv) const
bool isBoundary (const AmrIntVect_t &iv) const
bool applyBoundary (const AmrIntVect_t &iv, umap_t &map, const scalar_t &value)
bool applyBoundary (const AmrIntVect_t &iv, const basefab_t &fab, umap_t &map, const scalar_t &value)
void applyBoundary (const AmrIntVect_t &iv, const lo_t &dir, umap_t &map, const scalar_t &value)
const AmrIntVect_trefinement () const
const scalar_tcellSize () const
const scalar_tcellSize (lo_t dir) const
const scalar_tinvCellSize () const
const scalar_tinvCellSize (lo_t dir) const
bool isValid (const AmrIntVect_t &iv) const
void buildLevelMask ()

Public Attributes

const amrex::BoxArray & grids
 boxes of this level
const amrex::DistributionMapping & dmap
 AMReX core distribution map.
const AmrGeometry_tgeom
 geometry of this problem
Teuchos::RCP< dmap_tmap_p
 Tpetra core map.
Teuchos::RCP< matrix_tAnf_p
 no fine Poisson matrix
Teuchos::RCP< matrix_tR_p
 restriction matrix
Teuchos::RCP< matrix_tI_p
 interpolation matrix
Teuchos::RCP< matrix_tBcrse_p
 boundary from coarse cells
Teuchos::RCP< matrix_tBfine_p
 boundary from fine cells
Teuchos::RCP< matrix_tAwf_p
 composite Poisson matrix
Teuchos::RCP< matrix_tG_p [AMREX_SPACEDIM]
 gradient matrices in x, y, and z to compute electric field
Teuchos::RCP< vector_trho_p
 charge density
Teuchos::RCP< vector_tphi_p
 potential vector
Teuchos::RCP< vector_tresidual_p
 residual over all cells
Teuchos::RCP< vector_terror_p
 error over all cells
Teuchos::RCP< matrix_tUnCovered_p
 uncovered cells
std::unique_ptr< mask_tmask
 interior, phys boundary, interface, covered
std::unique_ptr< mask_trefmask
 covered (i.e. refined) or not-covered
std::unique_ptr< mask_tcrsemask

Private Member Functions

void buildMap (const Teuchos::RCP< comm_t > &comm)

Private Attributes

go_t nr_m [AMREX_SPACEDIM]
 number of grid points
AmrIntVect_t rr_m
 refinement
boundary_t bc_mp [AMREX_SPACEDIM]
 boundary conditions
scalar_t dx_m [AMREX_SPACEDIM]
 cell size in particle rest frame
scalar_t invdx_m [AMREX_SPACEDIM]
 inverse cell size in particle rest frame

Detailed Description

template<class MatrixType, class VectorType>
class AmrMultiGridLevel< MatrixType, VectorType >

Definition at line 36 of file AmrMultiGridLevel.h.

Member Typedef Documentation

◆ AmrField_s

template<class MatrixType, class VectorType>
typedef std::shared_ptr<AmrField_t> AmrMultiGridLevel< MatrixType, VectorType >::AmrField_s

Definition at line 42 of file AmrMultiGridLevel.h.

◆ AmrField_t

template<class MatrixType, class VectorType>
typedef amr::AmrField_t AmrMultiGridLevel< MatrixType, VectorType >::AmrField_t

Definition at line 39 of file AmrMultiGridLevel.h.

◆ AmrField_u

template<class MatrixType, class VectorType>
typedef std::unique_ptr<AmrField_t> AmrMultiGridLevel< MatrixType, VectorType >::AmrField_u

Definition at line 41 of file AmrMultiGridLevel.h.

◆ AmrGeometry_t

template<class MatrixType, class VectorType>
typedef amr::AmrGeometry_t AmrMultiGridLevel< MatrixType, VectorType >::AmrGeometry_t

Definition at line 40 of file AmrMultiGridLevel.h.

◆ AmrIntVect_t

template<class MatrixType, class VectorType>
typedef amr::AmrIntVect_t AmrMultiGridLevel< MatrixType, VectorType >::AmrIntVect_t

Definition at line 43 of file AmrMultiGridLevel.h.

◆ basefab_t

template<class MatrixType, class VectorType>
typedef amrex::BaseFab<int> AmrMultiGridLevel< MatrixType, VectorType >::basefab_t

Definition at line 46 of file AmrMultiGridLevel.h.

◆ boundary_t

template<class MatrixType, class VectorType>
typedef std::shared_ptr<AmrBoundary<AmrMultiGridLevel<MatrixType, VectorType > > > AmrMultiGridLevel< MatrixType, VectorType >::boundary_t

Definition at line 52 of file AmrMultiGridLevel.h.

◆ coefficients_t

template<class MatrixType, class VectorType>
typedef std::vector<scalar_t> AmrMultiGridLevel< MatrixType, VectorType >::coefficients_t

Type for matrix entries.

Definition at line 64 of file AmrMultiGridLevel.h.

◆ comm_t

template<class MatrixType, class VectorType>
typedef amr::comm_t AmrMultiGridLevel< MatrixType, VectorType >::comm_t

Definition at line 54 of file AmrMultiGridLevel.h.

◆ dmap_t

template<class MatrixType, class VectorType>
typedef amr::dmap_t AmrMultiGridLevel< MatrixType, VectorType >::dmap_t

Definition at line 55 of file AmrMultiGridLevel.h.

◆ go_t

template<class MatrixType, class VectorType>
typedef amr::global_ordinal_t AmrMultiGridLevel< MatrixType, VectorType >::go_t

Definition at line 56 of file AmrMultiGridLevel.h.

◆ indices_t

template<class MatrixType, class VectorType>
typedef std::vector<go_t> AmrMultiGridLevel< MatrixType, VectorType >::indices_t

Type for matrix indices.

Definition at line 61 of file AmrMultiGridLevel.h.

◆ lo_t

template<class MatrixType, class VectorType>
typedef amr::local_ordinal_t AmrMultiGridLevel< MatrixType, VectorType >::lo_t

Definition at line 58 of file AmrMultiGridLevel.h.

◆ mask_t

template<class MatrixType, class VectorType>
typedef amrex::FabArray<basefab_t> AmrMultiGridLevel< MatrixType, VectorType >::mask_t

Definition at line 47 of file AmrMultiGridLevel.h.

◆ matrix_t

template<class MatrixType, class VectorType>
typedef MatrixType AmrMultiGridLevel< MatrixType, VectorType >::matrix_t

Definition at line 44 of file AmrMultiGridLevel.h.

◆ scalar_t

template<class MatrixType, class VectorType>
typedef amr::scalar_t AmrMultiGridLevel< MatrixType, VectorType >::scalar_t

Definition at line 57 of file AmrMultiGridLevel.h.

◆ umap_t

template<class MatrixType, class VectorType>
typedef std::unordered_map<go_t, scalar_t> AmrMultiGridLevel< MatrixType, VectorType >::umap_t

Definition at line 67 of file AmrMultiGridLevel.h.

◆ vector_t

template<class MatrixType, class VectorType>
typedef VectorType AmrMultiGridLevel< MatrixType, VectorType >::vector_t

Definition at line 45 of file AmrMultiGridLevel.h.

Member Enumeration Documentation

◆ Mask

template<class MatrixType, class VectorType>
enum AmrMultiGridLevel::Mask
Enumerator
COVERED 
INTERIOR 
BNDRY 
PHYSBNDRY 

Definition at line 75 of file AmrMultiGridLevel.h.

◆ Refined

template<class MatrixType, class VectorType>
enum AmrMultiGridLevel::Refined
Enumerator
YES 
NO 

Definition at line 84 of file AmrMultiGridLevel.h.

Constructor & Destructor Documentation

◆ AmrMultiGridLevel()

template<class MatrixType, class VectorType>
AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel ( const Vector_t & meshScaling,
const amrex::BoxArray & _grids,
const amrex::DistributionMapping & _dmap,
const AmrGeometry_t & _geom,
const AmrIntVect_t & rr,
const boundary_t * bc,
const Teuchos::RCP< comm_t > & comm )
Parameters
meshscaling due to particle rest frame
_gridsof this level
_dmapAMReX core distribution map
_geomof domain
rrrefinement ratio
bcphysical boundaries (x, y, z)
commMPI communicator

Definition at line 27 of file AmrMultiGridLevel.hpp.

References AmrMultiGridLevel(), Anf_p, Awf_p, bc_mp, Bcrse_p, Bfine_p, buildLevelMask(), buildMap(), crsemask, dmap, dx_m, error_p, G_p, geom, grids, I_p, invdx_m, map_p, nr_m, phi_p, R_p, refmask, residual_p, rho_p, rr_m, and UnCovered_p.

Referenced by AmrMultiGridLevel().

Here is the call graph for this function:

◆ ~AmrMultiGridLevel()

template<class MatrixType, class VectorType>
AmrMultiGridLevel< MatrixType, VectorType >::~AmrMultiGridLevel ( )

Definition at line 82 of file AmrMultiGridLevel.hpp.

References Anf_p, Awf_p, Bcrse_p, Bfine_p, error_p, G_p, I_p, map_p, phi_p, R_p, residual_p, rho_p, and UnCovered_p.

Member Function Documentation

◆ applyBoundary() [1/3]

template<class MatrixType, class VectorType>
bool AmrMultiGridLevel< MatrixType, VectorType >::applyBoundary ( const AmrIntVect_t & iv,
const basefab_t & fab,
umap_t & map,
const scalar_t & value )

Slightly faster version of apply().

Parameters
ivis the cell where we want to have the boundary value
fabis the mask
mapwith indices global matrix indices and matrix values
valuematrix entry (coefficients) @precondition Basefab needs to be a mask with AmrMultiGridLevel::Mask::PHYSBNDRY

Definition at line 140 of file AmrMultiGridLevel.hpp.

References bc_mp, isBoundary(), nr_m, and PHYSBNDRY.

Here is the call graph for this function:

◆ applyBoundary() [2/3]

template<class MatrixType, class VectorType>
void AmrMultiGridLevel< MatrixType, VectorType >::applyBoundary ( const AmrIntVect_t & iv,
const lo_t & dir,
umap_t & map,
const scalar_t & value )

Apply boundary in a certain direction.

Parameters
ivis the cell where we want to have the boundary value
dirdirection of physical / mesh boundary
mapwith indices global matrix indices and matrix values
valuematrix entry (coefficients)

Definition at line 160 of file AmrMultiGridLevel.hpp.

References bc_mp, and nr_m.

◆ applyBoundary() [3/3]

template<class MatrixType, class VectorType>
bool AmrMultiGridLevel< MatrixType, VectorType >::applyBoundary ( const AmrIntVect_t & iv,
umap_t & map,
const scalar_t & value )

Checks all directions if physical / mesh boundary.

Parameters
ivis the cell where we want to have the boundary value
mapwith indices global matrix indices and matrix values
valuematrix entry (coefficients)

Definition at line 124 of file AmrMultiGridLevel.hpp.

References bc_mp, isBoundary(), and nr_m.

Here is the call graph for this function:

◆ buildLevelMask()

template<class MatrixType, class VectorType>
void AmrMultiGridLevel< MatrixType, VectorType >::buildLevelMask ( )

Build a mask specifying if a grid point is covered, an interior cell, at physical boundary or at interior boundary

Definition at line 170 of file AmrMultiGridLevel.hpp.

References BNDRY, COVERED, dmap, geom, grids, INTERIOR, mask, and PHYSBNDRY.

Referenced by AmrMultiGridLevel().

◆ buildMap()

template<class MatrixType, class VectorType>
void AmrMultiGridLevel< MatrixType, VectorType >::buildMap ( const Teuchos::RCP< comm_t > & comm)
private

Build Tpetra::Map of this level

Parameters
commMPI communicator

Definition at line 218 of file AmrMultiGridLevel.hpp.

References dmap, grids, map_p, and serialize().

Referenced by AmrMultiGridLevel().

Here is the call graph for this function:

◆ cellSize() [1/2]

template<class MatrixType, class VectorType>
const amr::scalar_t * AmrMultiGridLevel< MatrixType, VectorType >::cellSize ( ) const
Returns
the mesh spacing in particle rest frame

Definition at line 187 of file AmrMultiGridLevel.hpp.

References dx_m.

◆ cellSize() [2/2]

template<class MatrixType, class VectorType>
const amr::scalar_t & AmrMultiGridLevel< MatrixType, VectorType >::cellSize ( lo_t dir) const
Returns
the mesh spacing in particle rest frame for a certain direction

Definition at line 193 of file AmrMultiGridLevel.hpp.

References dx_m.

◆ invCellSize() [1/2]

template<class MatrixType, class VectorType>
const amr::scalar_t * AmrMultiGridLevel< MatrixType, VectorType >::invCellSize ( ) const
Returns
the inverse mesh spacing in particle rest frame

Definition at line 199 of file AmrMultiGridLevel.hpp.

References invdx_m.

◆ invCellSize() [2/2]

template<class MatrixType, class VectorType>
const amr::scalar_t & AmrMultiGridLevel< MatrixType, VectorType >::invCellSize ( lo_t dir) const
Returns
the inverse mesh spacing in particle rest frame for a certain direction

Definition at line 205 of file AmrMultiGridLevel.hpp.

References invdx_m.

◆ isBoundary()

template<class MatrixType, class VectorType>
bool AmrMultiGridLevel< MatrixType, VectorType >::isBoundary ( const AmrIntVect_t & iv) const

Checks if grid point is on the physical / mesh boundary

Parameters
ivgrid point (i, j, k)

Definition at line 117 of file AmrMultiGridLevel.hpp.

References bc_mp, and nr_m.

Referenced by applyBoundary(), and applyBoundary().

◆ isValid()

template<class MatrixType, class VectorType>
bool AmrMultiGridLevel< MatrixType, VectorType >::isValid ( const AmrIntVect_t & iv) const

Check if a cell on that level is valid

Parameters
ivis the cell
Returns
true if the cell is valid

Definition at line 211 of file AmrMultiGridLevel.hpp.

References nr_m.

◆ refinement()

template<class MatrixType, class VectorType>
const amr::AmrIntVect_t & AmrMultiGridLevel< MatrixType, VectorType >::refinement ( ) const

Definition at line 181 of file AmrMultiGridLevel.hpp.

References rr_m.

◆ serialize()

template<class MatrixType, class VectorType>
AmrMultiGridLevel< MatrixType, VectorType >::go_t AmrMultiGridLevel< MatrixType, VectorType >::serialize ( const AmrIntVect_t & iv) const

Map a 2D / 3D grid point to an array index

Parameters
ivgrid point (i, j, k)

Definition at line 107 of file AmrMultiGridLevel.hpp.

References nr_m.

Member Data Documentation

◆ Anf_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::Anf_p

no fine Poisson matrix

Definition at line 208 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ Awf_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::Awf_p

composite Poisson matrix

Definition at line 213 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ bc_mp

template<class MatrixType, class VectorType>
boundary_t AmrMultiGridLevel< MatrixType, VectorType >::bc_mp[AMREX_SPACEDIM]
private

boundary conditions

Definition at line 233 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), applyBoundary(), applyBoundary(), applyBoundary(), and isBoundary().

◆ Bcrse_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::Bcrse_p

boundary from coarse cells

Definition at line 211 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ Bfine_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::Bfine_p

boundary from fine cells

Definition at line 212 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ crsemask

template<class MatrixType, class VectorType>
std::unique_ptr<mask_t> AmrMultiGridLevel< MatrixType, VectorType >::crsemask

Definition at line 226 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel().

◆ dmap

template<class MatrixType, class VectorType>
const amrex::DistributionMapping& AmrMultiGridLevel< MatrixType, VectorType >::dmap

AMReX core distribution map.

Definition at line 203 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), buildLevelMask(), and buildMap().

◆ dx_m

template<class MatrixType, class VectorType>
scalar_t AmrMultiGridLevel< MatrixType, VectorType >::dx_m[AMREX_SPACEDIM]
private

cell size in particle rest frame

Definition at line 235 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), cellSize(), and cellSize().

◆ error_p

template<class MatrixType, class VectorType>
Teuchos::RCP<vector_t> AmrMultiGridLevel< MatrixType, VectorType >::error_p

error over all cells

Definition at line 221 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ G_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::G_p[AMREX_SPACEDIM]

gradient matrices in x, y, and z to compute electric field

Definition at line 216 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ geom

template<class MatrixType, class VectorType>
const AmrGeometry_t& AmrMultiGridLevel< MatrixType, VectorType >::geom

geometry of this problem

Definition at line 204 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and buildLevelMask().

◆ grids

template<class MatrixType, class VectorType>
const amrex::BoxArray& AmrMultiGridLevel< MatrixType, VectorType >::grids

boxes of this level

Definition at line 202 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), buildLevelMask(), and buildMap().

◆ I_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::I_p

interpolation matrix

Definition at line 210 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ invdx_m

template<class MatrixType, class VectorType>
scalar_t AmrMultiGridLevel< MatrixType, VectorType >::invdx_m[AMREX_SPACEDIM]
private

inverse cell size in particle rest frame

Definition at line 236 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), invCellSize(), and invCellSize().

◆ map_p

template<class MatrixType, class VectorType>
Teuchos::RCP<dmap_t> AmrMultiGridLevel< MatrixType, VectorType >::map_p

Tpetra core map.

Definition at line 206 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), buildMap(), and ~AmrMultiGridLevel().

◆ mask

template<class MatrixType, class VectorType>
std::unique_ptr<mask_t> AmrMultiGridLevel< MatrixType, VectorType >::mask

interior, phys boundary, interface, covered

Definition at line 224 of file AmrMultiGridLevel.h.

Referenced by buildLevelMask().

◆ nr_m

template<class MatrixType, class VectorType>
go_t AmrMultiGridLevel< MatrixType, VectorType >::nr_m[AMREX_SPACEDIM]
private

number of grid points

Definition at line 229 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), applyBoundary(), applyBoundary(), applyBoundary(), isBoundary(), isValid(), and serialize().

◆ phi_p

template<class MatrixType, class VectorType>
Teuchos::RCP<vector_t> AmrMultiGridLevel< MatrixType, VectorType >::phi_p

potential vector

Definition at line 219 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ R_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::R_p

restriction matrix

Definition at line 209 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ refmask

template<class MatrixType, class VectorType>
std::unique_ptr<mask_t> AmrMultiGridLevel< MatrixType, VectorType >::refmask

covered (i.e. refined) or not-covered

Definition at line 225 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel().

◆ residual_p

template<class MatrixType, class VectorType>
Teuchos::RCP<vector_t> AmrMultiGridLevel< MatrixType, VectorType >::residual_p

residual over all cells

Definition at line 220 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ rho_p

template<class MatrixType, class VectorType>
Teuchos::RCP<vector_t> AmrMultiGridLevel< MatrixType, VectorType >::rho_p

charge density

Definition at line 218 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().

◆ rr_m

template<class MatrixType, class VectorType>
AmrIntVect_t AmrMultiGridLevel< MatrixType, VectorType >::rr_m
private

refinement

Definition at line 231 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and refinement().

◆ UnCovered_p

template<class MatrixType, class VectorType>
Teuchos::RCP<matrix_t> AmrMultiGridLevel< MatrixType, VectorType >::UnCovered_p

uncovered cells

Definition at line 222 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel(), and ~AmrMultiGridLevel().


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