OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
AmrMultiGridLevel.h
Go to the documentation of this file.
1//
2// Class AmrMultiGridLevel
3// This class represents a single AMR level, i.e. it stores all matrices
4// and vectors of a level.
5//
6// Copyright (c) 2017 - 2020, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
7// All rights reserved
8//
9// Implemented as part of the PhD thesis
10// "Precise Simulations of Multibunches in High Intensity Cyclotrons"
11//
12// This file is part of OPAL.
13//
14// OPAL is free software: you can redistribute it and/or modify
15// it under the terms of the GNU General Public License as published by
16// the Free Software Foundation, either version 3 of the License, or
17// (at your option) any later version.
18//
19// You should have received a copy of the GNU General Public License
20// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
21//
22#ifndef AMR_MULTI_GRID_LEVEL
23#define AMR_MULTI_GRID_LEVEL
24
25#include <vector>
26
27#include <AMReX_IntVect.H>
28
29#include "Ippl.h"
30
31#include "AmrMultiGridDefs.h"
32
33#include <unordered_map>
34
35template <class MatrixType, class VectorType>
37
38public:
41 typedef std::unique_ptr<AmrField_t> AmrField_u;
42 typedef std::shared_ptr<AmrField_t> AmrField_s;
44 typedef MatrixType matrix_t;
45 typedef VectorType vector_t;
46 typedef amrex::BaseFab<int> basefab_t;
47 typedef amrex::FabArray<basefab_t> mask_t;
48 typedef std::shared_ptr<AmrBoundary<AmrMultiGridLevel<MatrixType,
49 VectorType
50 >
51 >
53
59
61 typedef std::vector<go_t> indices_t;
62
64 typedef std::vector<scalar_t> coefficients_t;
65
66 // Type with matrix index (column) and coefficient value
67 typedef std::unordered_map<go_t, scalar_t> umap_t;
68
69 // covered : ghost cells covered by valid cells of this FabArray
70 // (including periodically shifted valid cells)
71 // notcovered: ghost cells not covered by valid cells
72 // (including ghost cells outside periodic boundaries) (BNDRY)
73 // physbnd : boundary cells outside the domain (excluding periodic boundaries)
74 // interior : interior cells (i.e., valid cells)
75 enum Mask {
76 COVERED = -1,
78 BNDRY = 1,
79 PHYSBNDRY = 2
80 };
81
82 // NO : not a refined cell
83 // YES : cell got refined
84 enum Refined {
85 YES = 0,
86 NO = 1
87 };
88
89public:
99 AmrMultiGridLevel(const Vector_t& meshScaling,
100 const amrex::BoxArray& _grids,
101 const amrex::DistributionMapping& _dmap,
102 const AmrGeometry_t& _geom,
103 const AmrIntVect_t& rr,
104 const boundary_t* bc,
105 const Teuchos::RCP<comm_t>& comm);
106
108
113 go_t serialize(const AmrIntVect_t& iv) const;
114
119 bool isBoundary(const AmrIntVect_t& iv) const;
120
127 bool applyBoundary(const AmrIntVect_t& iv,
128 umap_t& map,
129 const scalar_t& value);
130
140 bool applyBoundary(const AmrIntVect_t& iv,
141 const basefab_t& fab,
142 umap_t& map,
143 const scalar_t& value);
144
152 void applyBoundary(const AmrIntVect_t& iv,
153 const lo_t& dir,
154 umap_t& map,
155 const scalar_t& value);
156
157 const AmrIntVect_t& refinement() const;
158
162 const scalar_t* cellSize() const;
163
168 const scalar_t& cellSize(lo_t dir) const;
169
173 const scalar_t* invCellSize() const;
174
179 const scalar_t& invCellSize(lo_t dir) const;
180
186 bool isValid(const AmrIntVect_t& iv) const;
187
192 void buildLevelMask();
193
194private:
199 void buildMap(const Teuchos::RCP<comm_t>& comm);
200
201public:
202 const amrex::BoxArray& grids;
203 const amrex::DistributionMapping& dmap;
205
206 Teuchos::RCP<dmap_t> map_p;
207
208 Teuchos::RCP<matrix_t> Anf_p;
209 Teuchos::RCP<matrix_t> R_p;
210 Teuchos::RCP<matrix_t> I_p;
211 Teuchos::RCP<matrix_t> Bcrse_p;
212 Teuchos::RCP<matrix_t> Bfine_p;
213 Teuchos::RCP<matrix_t> Awf_p;
214
216 Teuchos::RCP<matrix_t> G_p[AMREX_SPACEDIM];
217
218 Teuchos::RCP<vector_t> rho_p;
219 Teuchos::RCP<vector_t> phi_p;
220 Teuchos::RCP<vector_t> residual_p;
221 Teuchos::RCP<vector_t> error_p;
222 Teuchos::RCP<matrix_t> UnCovered_p;
223
224 std::unique_ptr<mask_t> mask;
225 std::unique_ptr<mask_t> refmask;
226 std::unique_ptr<mask_t> crsemask;
227
228private:
229 go_t nr_m[AMREX_SPACEDIM];
230
232
233 boundary_t bc_mp[AMREX_SPACEDIM];
234
235 scalar_t dx_m[AMREX_SPACEDIM];
236 scalar_t invdx_m[AMREX_SPACEDIM];
237};
238
239
240#include "AmrMultiGridLevel.hpp"
241
242#endif
long global_ordinal_t
amrex::Geometry AmrGeometry_t
Definition: AmrDefs.h:39
Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > dmap_t
amrex::IntVect AmrIntVect_t
Definition: AmrDefs.h:48
int local_ordinal_t
double scalar_t
amrex::MultiFab AmrField_t
Definition: AmrDefs.h:34
Teuchos::MpiComm< int > comm_t
Teuchos::RCP< matrix_t > G_p[AMREX_SPACEDIM]
gradient matrices in x, y, and z to compute electric field
std::unordered_map< go_t, scalar_t > umap_t
std::unique_ptr< AmrField_t > AmrField_u
amrex::FabArray< basefab_t > mask_t
amr::local_ordinal_t lo_t
Teuchos::RCP< vector_t > rho_p
charge density
scalar_t invdx_m[AMREX_SPACEDIM]
inverse cell size in particle rest frame
std::unique_ptr< mask_t > mask
interior, phys boundary, interface, covered
bool isValid(const AmrIntVect_t &iv) const
Teuchos::RCP< matrix_t > UnCovered_p
uncovered cells
bool applyBoundary(const AmrIntVect_t &iv, umap_t &map, const scalar_t &value)
std::shared_ptr< AmrField_t > AmrField_s
amrex::BaseFab< int > basefab_t
Teuchos::RCP< vector_t > residual_p
residual over all cells
AmrIntVect_t rr_m
refinement
Teuchos::RCP< matrix_t > Awf_p
composite Poisson matrix
Teuchos::RCP< vector_t > phi_p
potential vector
std::unique_ptr< mask_t > crsemask
Teuchos::RCP< matrix_t > Bfine_p
boundary from fine cells
Teuchos::RCP< vector_t > error_p
error over all cells
amr::AmrField_t AmrField_t
std::vector< go_t > indices_t
Type for matrix indices.
const amrex::BoxArray & grids
boxes of this level
std::unique_ptr< mask_t > refmask
covered (i.e. refined) or not-covered
amr::AmrIntVect_t AmrIntVect_t
const amrex::DistributionMapping & dmap
AMReX core distribution map.
Teuchos::RCP< matrix_t > Bcrse_p
boundary from coarse cells
bool isBoundary(const AmrIntVect_t &iv) const
Teuchos::RCP< matrix_t > Anf_p
no fine Poisson matrix
amr::AmrGeometry_t AmrGeometry_t
amr::scalar_t scalar_t
const AmrIntVect_t & refinement() const
go_t nr_m[AMREX_SPACEDIM]
number of grid points
const AmrGeometry_t & geom
geometry of this problem
const scalar_t * invCellSize() const
amr::global_ordinal_t go_t
Teuchos::RCP< matrix_t > R_p
restriction matrix
boundary_t bc_mp[AMREX_SPACEDIM]
boundary conditions
std::shared_ptr< AmrBoundary< AmrMultiGridLevel< MatrixType, VectorType > > > boundary_t
std::vector< scalar_t > coefficients_t
Type for matrix entries.
const scalar_t * cellSize() const
Teuchos::RCP< matrix_t > I_p
interpolation matrix
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)
scalar_t dx_m[AMREX_SPACEDIM]
cell size in particle rest frame
go_t serialize(const AmrIntVect_t &iv) const
void buildMap(const Teuchos::RCP< comm_t > &comm)
Teuchos::RCP< dmap_t > map_p
Tpetra core map.