OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
AmrMultiGridLevel.h
Go to the documentation of this file.
1 #ifndef AMR_MULTI_GRID_LEVEL
2 #define AMR_MULTI_GRID_LEVEL
3 
4 #include <vector>
5 
6 #include <AMReX_IntVect.H>
7 
8 #include "Ippl.h"
9 
10 #include "AmrMultiGridDefs.h"
11 
12 #include <unordered_map>
13 
14 template <class MatrixType, class VectorType>
16 
17 public:
20  typedef std::unique_ptr<AmrField_t> AmrField_u;
21  typedef std::shared_ptr<AmrField_t> AmrField_s;
23  typedef MatrixType matrix_t;
24  typedef VectorType vector_t;
25  typedef amrex::BaseFab<int> basefab_t;
26  typedef amrex::FabArray<basefab_t> mask_t;
27  typedef std::shared_ptr<AmrBoundary<AmrMultiGridLevel<MatrixType,
28  VectorType
29  >
30  >
32 
39 
41  typedef std::vector<go_t> indices_t;
42 
44  typedef std::vector<scalar_t> coefficients_t;
45 
46  // Type with matrix index (column) and coefficient value
47  typedef std::unordered_map<go_t, scalar_t> umap_t;
48 
49  // covered : ghost cells covered by valid cells of this FabArray
50  // (including periodically shifted valid cells)
51  // notcovered: ghost cells not covered by valid cells
52  // (including ghost cells outside periodic boundaries) (BNDRY)
53  // physbnd : boundary cells outside the domain (excluding periodic boundaries)
54  // interior : interior cells (i.e., valid cells)
55  enum Mask {
56  COVERED = -1,
57  INTERIOR = 0,
58  BNDRY = 1,
60  };
61 
62  // NO : not a refined cell
63  // YES : cell got refined
64  enum Refined {
65  YES = 0,
66  NO = 1
67  };
68 
69 public:
80  AmrMultiGridLevel(const Vector_t& meshScaling,
81  const amrex::BoxArray& _grids,
82  const amrex::DistributionMapping& _dmap,
83  const AmrGeometry_t& _geom,
84  const AmrIntVect_t& rr,
85  const boundary_t* bc,
86  const Teuchos::RCP<comm_t>& comm,
87  const Teuchos::RCP<node_t>& node);
88 
90 
95  go_t serialize(const AmrIntVect_t& iv) const;
96 
101  bool isBoundary(const AmrIntVect_t& iv) const;
102 
109  bool applyBoundary(const AmrIntVect_t& iv,
110  umap_t& map,
111  const scalar_t& value);
112 
122  bool applyBoundary(const AmrIntVect_t& iv,
123  const basefab_t& fab,
124  umap_t& map,
125  const scalar_t& value);
126 
134  void applyBoundary(const AmrIntVect_t& iv,
135  const lo_t& dir,
136  umap_t& map,
137  const scalar_t& value);
138 
139  const AmrIntVect_t& refinement() const;
140 
144  const scalar_t* cellSize() const;
145 
150  const scalar_t& cellSize(lo_t dir) const;
151 
155  const scalar_t* invCellSize() const;
156 
161  const scalar_t& invCellSize(lo_t dir) const;
162 
168  bool isValid(const AmrIntVect_t& iv) const;
169 
174  void buildLevelMask();
175 
176 private:
182  void buildMap(const Teuchos::RCP<comm_t>& comm,
183  const Teuchos::RCP<node_t>& node);
184 
185 public:
186  const amrex::BoxArray& grids;
187  const amrex::DistributionMapping& dmap;
189 
190  Teuchos::RCP<dmap_t> map_p;
191 
192  Teuchos::RCP<matrix_t> Anf_p;
193  Teuchos::RCP<matrix_t> R_p;
194  Teuchos::RCP<matrix_t> I_p;
195  Teuchos::RCP<matrix_t> Bcrse_p;
196  Teuchos::RCP<matrix_t> Bfine_p;
197  Teuchos::RCP<matrix_t> Awf_p;
198 
200  Teuchos::RCP<matrix_t> G_p[AMREX_SPACEDIM];
201 
202  Teuchos::RCP<vector_t> rho_p;
203  Teuchos::RCP<vector_t> phi_p;
204  Teuchos::RCP<vector_t> residual_p;
205  Teuchos::RCP<vector_t> error_p;
206  Teuchos::RCP<matrix_t> UnCovered_p;
207 
208  std::unique_ptr<mask_t> mask;
209  std::unique_ptr<mask_t> refmask;
210  std::unique_ptr<mask_t> crsemask;
211 
212 private:
213  go_t nr_m[AMREX_SPACEDIM];
214 
216 
217  boundary_t bc_mp[AMREX_SPACEDIM];
218 
219  scalar_t dx_m[AMREX_SPACEDIM];
220  scalar_t invdx_m[AMREX_SPACEDIM];
221 };
222 
223 
224 #include "AmrMultiGridLevel.hpp"
225 
226 #endif
std::unique_ptr< mask_t > refmask
covered (i.e. refined) or not-covered
amrex::Geometry AmrGeometry_t
Definition: AmrDefs.h:19
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
Definition: TSVMeta.h:24
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
int local_ordinal_t
amrex::IntVect AmrIntVect_t
Definition: AmrDefs.h:28
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
Definition: AmrDefs.h:14
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
long global_ordinal_t
amr::scalar_t scalar_t
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
double scalar_t