1 #ifndef AMR_MULTI_GRID_LEVEL
2 #define AMR_MULTI_GRID_LEVEL
6 #include <AMReX_IntVect.H>
12 #include <unordered_map>
14 template <
class MatrixType,
class VectorType>
26 typedef amrex::FabArray<basefab_t>
mask_t;
47 typedef std::unordered_map<go_t, scalar_t>
umap_t;
81 const amrex::BoxArray& _grids,
82 const amrex::DistributionMapping& _dmap,
86 const Teuchos::RCP<comm_t>& comm,
87 const Teuchos::RCP<node_t>& node);
182 void buildMap(
const Teuchos::RCP<comm_t>& comm,
183 const Teuchos::RCP<node_t>& node);
187 const amrex::DistributionMapping&
dmap;
193 Teuchos::RCP<matrix_t>
R_p;
194 Teuchos::RCP<matrix_t>
I_p;
200 Teuchos::RCP<matrix_t>
G_p[AMREX_SPACEDIM];
std::unique_ptr< mask_t > refmask
covered (i.e. refined) or not-covered
amrex::Geometry AmrGeometry_t
bool applyBoundary(const AmrIntVect_t &iv, umap_t &map, const scalar_t &value)
std::vector< scalar_t > coefficients_t
Type for matrix entries.
std::unordered_map< go_t, scalar_t > umap_t
Teuchos::RCP< vector_t > rho_p
charge density
amr::AmrField_t AmrField_t
AmrIntVect_t rr_m
refinement
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)
Teuchos::RCP< matrix_t > G_p[AMREX_SPACEDIM]
gradient matrices in x, y, and z to compute electric field
amrex::IntVect AmrIntVect_t
Teuchos::RCP< matrix_t > Awf_p
composite Poisson matrix
const amrex::BoxArray & grids
boxes of this level
scalar_t invdx_m[AMREX_SPACEDIM]
inverse cell size in particle rest frame
KokkosClassic::DefaultNode::DefaultNodeType node_t
const scalar_t * cellSize() const
amr::local_ordinal_t lo_t
bool isValid(const AmrIntVect_t &iv) const
amrex::MultiFab AmrField_t
Teuchos::RCP< matrix_t > I_p
interpolation matrix
boundary_t bc_mp[AMREX_SPACEDIM]
boundary conditions
std::unique_ptr< mask_t > mask
interior, phys boundary, interface, covered
amrex::FabArray< basefab_t > mask_t
Teuchos::RCP< vector_t > phi_p
potential vector
const amrex::DistributionMapping & dmap
AMReX core distribution map.
amr::AmrIntVect_t AmrIntVect_t
void buildMap(const Teuchos::RCP< comm_t > &comm, const Teuchos::RCP< node_t > &node)
Teuchos::RCP< matrix_t > R_p
restriction matrix
scalar_t dx_m[AMREX_SPACEDIM]
cell size in particle rest frame
Teuchos::MpiComm< int > comm_t
Teuchos::RCP< matrix_t > Bfine_p
boundary from fine cells
std::unique_ptr< mask_t > crsemask
std::shared_ptr< AmrBoundary< AmrMultiGridLevel< MatrixType, VectorType > > > boundary_t
Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > dmap_t
std::unique_ptr< AmrField_t > AmrField_u
Teuchos::RCP< matrix_t > UnCovered_p
uncovered cells
bool isBoundary(const AmrIntVect_t &iv) const
Teuchos::RCP< vector_t > residual_p
residual over all cells
amr::global_ordinal_t go_t
Teuchos::RCP< matrix_t > Bcrse_p
boundary from coarse cells
const AmrGeometry_t & geom
geometry of this problem
std::shared_ptr< AmrField_t > AmrField_s
Teuchos::RCP< vector_t > error_p
error over all cells
Teuchos::RCP< matrix_t > Anf_p
no fine Poisson matrix
amrex::BaseFab< int > basefab_t
go_t serialize(const AmrIntVect_t &iv) const
const AmrIntVect_t & refinement() const
Teuchos::RCP< dmap_t > map_p
Tpetra core map.
go_t nr_m[AMREX_SPACEDIM]
number of grid points
const scalar_t * invCellSize() const
std::vector< go_t > indices_t
Type for matrix indices.
amr::AmrGeometry_t AmrGeometry_t