OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
Public Types | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
AmrMultiGridLevel< MatrixType, VectorType > Class Template Reference

#include <AmrMultiGridLevel.h>

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_t
AmrField_u
 
typedef std::shared_ptr
< AmrField_t
AmrField_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_t
mask_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::node_t node_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. More...
 
typedef std::vector< scalar_tcoefficients_t
 Type for matrix entries. More...
 
typedef std::unordered_map
< go_t, scalar_t
umap_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, const Teuchos::RCP< node_t > &node)
 
 ~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 More...
 
const amrex::DistributionMapping & dmap
 AMReX core distribution map. More...
 
const AmrGeometry_tgeom
 geometry of this problem More...
 
Teuchos::RCP< dmap_tmap_p
 Tpetra core map. More...
 
Teuchos::RCP< matrix_tAnf_p
 no fine Poisson matrix More...
 
Teuchos::RCP< matrix_tR_p
 restriction matrix More...
 
Teuchos::RCP< matrix_tI_p
 interpolation matrix More...
 
Teuchos::RCP< matrix_tBcrse_p
 boundary from coarse cells More...
 
Teuchos::RCP< matrix_tBfine_p
 boundary from fine cells More...
 
Teuchos::RCP< matrix_tAwf_p
 composite Poisson matrix More...
 
Teuchos::RCP< matrix_tG_p [AMREX_SPACEDIM]
 gradient matrices in x, y, and z to compute electric field More...
 
Teuchos::RCP< vector_trho_p
 charge density More...
 
Teuchos::RCP< vector_tphi_p
 potential vector More...
 
Teuchos::RCP< vector_tresidual_p
 residual over all cells More...
 
Teuchos::RCP< vector_terror_p
 error over all cells More...
 
Teuchos::RCP< matrix_tUnCovered_p
 uncovered cells More...
 
std::unique_ptr< mask_tmask
 interior, phys boundary, interface, covered More...
 
std::unique_ptr< mask_trefmask
 covered (i.e. refined) or not-covered More...
 
std::unique_ptr< mask_tcrsemask
 

Private Member Functions

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

Private Attributes

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

Detailed Description

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

Definition at line 15 of file AmrMultiGridLevel.h.

Member Typedef Documentation

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

Definition at line 21 of file AmrMultiGridLevel.h.

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

Definition at line 18 of file AmrMultiGridLevel.h.

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

Definition at line 20 of file AmrMultiGridLevel.h.

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

Definition at line 19 of file AmrMultiGridLevel.h.

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

Definition at line 22 of file AmrMultiGridLevel.h.

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

Definition at line 25 of file AmrMultiGridLevel.h.

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

Definition at line 31 of file AmrMultiGridLevel.h.

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

Type for matrix entries.

Definition at line 44 of file AmrMultiGridLevel.h.

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

Definition at line 33 of file AmrMultiGridLevel.h.

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

Definition at line 34 of file AmrMultiGridLevel.h.

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

Definition at line 36 of file AmrMultiGridLevel.h.

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

Type for matrix indices.

Definition at line 41 of file AmrMultiGridLevel.h.

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

Definition at line 38 of file AmrMultiGridLevel.h.

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

Definition at line 26 of file AmrMultiGridLevel.h.

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

Definition at line 23 of file AmrMultiGridLevel.h.

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

Definition at line 35 of file AmrMultiGridLevel.h.

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

Definition at line 37 of file AmrMultiGridLevel.h.

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

Definition at line 47 of file AmrMultiGridLevel.h.

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

Definition at line 24 of file AmrMultiGridLevel.h.

Member Enumeration Documentation

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

Definition at line 55 of file AmrMultiGridLevel.h.

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

Definition at line 64 of file AmrMultiGridLevel.h.

Constructor & Destructor Documentation

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,
const Teuchos::RCP< node_t > &  node 
)
template<class MatrixType , class VectorType >
AmrMultiGridLevel< MatrixType, VectorType >::~AmrMultiGridLevel ( )

Definition at line 62 of file AmrMultiGridLevel.hpp.

Member Function Documentation

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 104 of file AmrMultiGridLevel.hpp.

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) Basefab needs to be a mask with AmrMultiGridLevel::Mask::PHYSBNDRY

Definition at line 120 of file AmrMultiGridLevel.hpp.

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 140 of file AmrMultiGridLevel.hpp.

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 150 of file AmrMultiGridLevel.hpp.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

Build Tpetra::Map of this level

Parameters
commMPI communicator
nodeKokkos node type

Definition at line 198 of file AmrMultiGridLevel.hpp.

References serialize().

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

Here is the call graph for this function:

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 167 of file AmrMultiGridLevel.hpp.

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 173 of file AmrMultiGridLevel.hpp.

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 179 of file AmrMultiGridLevel.hpp.

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 185 of file AmrMultiGridLevel.hpp.

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 97 of file AmrMultiGridLevel.hpp.

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 191 of file AmrMultiGridLevel.hpp.

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

Definition at line 161 of file AmrMultiGridLevel.hpp.

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 87 of file AmrMultiGridLevel.hpp.

Member Data Documentation

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

no fine Poisson matrix

Definition at line 192 of file AmrMultiGridLevel.h.

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

composite Poisson matrix

Definition at line 197 of file AmrMultiGridLevel.h.

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

boundary conditions

Definition at line 217 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

boundary from coarse cells

Definition at line 195 of file AmrMultiGridLevel.h.

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

boundary from fine cells

Definition at line 196 of file AmrMultiGridLevel.h.

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

Definition at line 210 of file AmrMultiGridLevel.h.

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

AMReX core distribution map.

Definition at line 187 of file AmrMultiGridLevel.h.

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

cell size in particle rest frame

Definition at line 219 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

error over all cells

Definition at line 205 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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 200 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

geometry of this problem

Definition at line 188 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

boxes of this level

Definition at line 186 of file AmrMultiGridLevel.h.

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

interpolation matrix

Definition at line 194 of file AmrMultiGridLevel.h.

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 220 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

Tpetra core map.

Definition at line 190 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

interior, phys boundary, interface, covered

Definition at line 208 of file AmrMultiGridLevel.h.

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

number of grid points

Definition at line 213 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

potential vector

Definition at line 203 of file AmrMultiGridLevel.h.

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

restriction matrix

Definition at line 193 of file AmrMultiGridLevel.h.

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

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

Definition at line 209 of file AmrMultiGridLevel.h.

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

residual over all cells

Definition at line 204 of file AmrMultiGridLevel.h.

Referenced by AmrMultiGridLevel< MatrixType, VectorType >::AmrMultiGridLevel().

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

charge density

Definition at line 202 of file AmrMultiGridLevel.h.

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

refinement

Definition at line 215 of file AmrMultiGridLevel.h.

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

uncovered cells

Definition at line 206 of file AmrMultiGridLevel.h.


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