OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
AmrMultiGridLevel.hpp
Go to the documentation of this file.
1 #define AMR_NO_SCALE false
2 
3 
4 template <class MatrixType, class VectorType>
5 AmrMultiGridLevel<MatrixType,
6  VectorType>::AmrMultiGridLevel(const Vector_t& meshScaling,
7  const amrex::BoxArray& _grids,
8  const amrex::DistributionMapping& _dmap,
9  const AmrGeometry_t& _geom,
10  const AmrIntVect_t& rr,
11  const boundary_t* bc,
12  const Teuchos::RCP<comm_t>& comm,
13  const Teuchos::RCP<node_t>& node)
14  : grids(_grids),
15  dmap(_dmap),
16  geom(_geom),
17  map_p(Teuchos::null),
18  Anf_p(Teuchos::null),
19  R_p(Teuchos::null),
20  I_p(Teuchos::null),
21  Bcrse_p(Teuchos::null),
22  Bfine_p(Teuchos::null),
23  Awf_p(Teuchos::null),
24  rho_p(Teuchos::null),
25  phi_p(Teuchos::null),
26  residual_p(Teuchos::null),
27  error_p(Teuchos::null),
28  UnCovered_p(Teuchos::null),
29  refmask(nullptr),
30  crsemask(nullptr),
31  rr_m(rr)
32 {
33  for (int j = 0; j < AMREX_SPACEDIM; ++j) {
34  G_p[j] = Teuchos::null;
35 
36  nr_m[j] = _geom.Domain().length(j);
37 
38 #if AMR_NO_SCALE
39  // mesh spacing in particle rest frame
40  dx_m[j] = geom.CellSize(j);
41  invdx_m[j] = geom.InvCellSize(j);
42 #else
43  // mesh spacing in particle rest frame
44  dx_m[j] = meshScaling[j] * geom.CellSize(j);
45  invdx_m[j] = meshScaling[j] * geom.InvCellSize(j);
46 #endif
47 
48  bc_mp[j] = bc[j];
49  }
50 
51  this->buildLevelMask();
52 
53  this->buildMap(comm, node);
54 
55 
56  residual_p = Teuchos::rcp( new vector_t(map_p, false) );
57  error_p = Teuchos::rcp( new vector_t(map_p, false) );
58 }
59 
60 
61 template <class MatrixType, class VectorType>
63 {
64  map_p = Teuchos::null;
65 
66  Anf_p = Teuchos::null;
67  R_p = Teuchos::null;
68  I_p = Teuchos::null;
69  Bcrse_p = Teuchos::null;
70  Bfine_p = Teuchos::null;
71  Awf_p = Teuchos::null;
72 
73  for (int j = 0; j < AMREX_SPACEDIM; ++j)
74  G_p[j] = Teuchos::null;
75 
76  UnCovered_p = Teuchos::null;
77 
78  rho_p = Teuchos::null;
79  phi_p = Teuchos::null;
80  residual_p = Teuchos::null;
81  error_p = Teuchos::null;
82 }
83 
84 
85 template <class MatrixType, class VectorType>
88 #if AMREX_SPACEDIM == 3
89  return iv[0] + (iv[1] + nr_m[1] * iv[2]) * nr_m[0];
90 #else
91  return iv[0] + iv[1] * nr_m[0];
92 #endif
93 }
94 
95 
96 template <class MatrixType, class VectorType>
98  // it doesn't matter with which direction we check, since it checks all
99  return bc_mp[0]->isBoundary(iv, &nr_m[0]);
100 }
101 
102 
103 template <class MatrixType, class VectorType>
105  umap_t& map,
106  const scalar_t& value)
107 {
108  bool applied = false;
109  for (int d = 0; d < AMREX_SPACEDIM; ++d) {
110  if ( bc_mp[d]->isBoundary(iv, d, &nr_m[0]) ) {
111  applied = true;
112  bc_mp[d]->apply(iv, d, map, value, this, &nr_m[0]);
113  }
114  }
115  return applied;
116 }
117 
118 
119 template <class MatrixType, class VectorType>
121  const basefab_t& fab,
122  umap_t& map,
123  const scalar_t& value)
124 {
125  if ( fab(iv) != Mask::PHYSBNDRY )
126  return false;
127 
128  bool applied = false;
129  for (int d = 0; d < AMREX_SPACEDIM; ++d) {
130  if ( bc_mp[d]->isBoundary(iv, d, &nr_m[0]) ) {
131  applied = true;
132  bc_mp[d]->apply(iv, d, map, value, this, &nr_m[0]);
133  }
134  }
135  return applied;
136 }
137 
138 
139 template <class MatrixType, class VectorType>
141  const lo_t& dir,
142  umap_t& map,
143  const scalar_t& value)
144 {
145  bc_mp[dir]->apply(iv, dir, map, value, this, &nr_m[0]);
146 }
147 
148 
149 template <class MatrixType, class VectorType>
151  amrex::Periodicity period(AmrIntVect_t(D_DECL(0, 0, 0)));
152  mask.reset(new mask_t(grids, dmap, 1, 1));
153  mask->BuildMask(geom.Domain(), period,
154  Mask::COVERED, Mask::BNDRY,
155  Mask::PHYSBNDRY, Mask::INTERIOR);
156  mask->FillBoundary(period);
157 }
158 
159 
160 template <class MatrixType, class VectorType>
162  return rr_m;
163 }
164 
165 
166 template <class MatrixType, class VectorType>
168  return dx_m;
169 }
170 
171 
172 template <class MatrixType, class VectorType>
174  return dx_m[dir];
175 }
176 
177 
178 template <class MatrixType, class VectorType>
180  return invdx_m;
181 }
182 
183 
184 template <class MatrixType, class VectorType>
186  return invdx_m[dir];
187 }
188 
189 
190 template <class MatrixType, class VectorType>
192  return ( iv.allGT(AmrIntVect_t(D_DECL(-1, -1, -1))) &&
193  iv.allLT(AmrIntVect_t(D_DECL(nr_m[0], nr_m[1], nr_m[2]))) );
194 }
195 
196 
197 template <class MatrixType, class VectorType>
198 void AmrMultiGridLevel<MatrixType, VectorType>::buildMap(const Teuchos::RCP<comm_t>& comm,
199  const Teuchos::RCP<node_t>& node)
200 {
201 
202  go_t localNumElements = 0;
203 
204  Teuchos::Array<go_t> globalindices;
205 
206  for (amrex::MFIter mfi(grids, dmap, true); mfi.isValid(); ++mfi) {
207  const amrex::Box& tbx = mfi.tilebox();
208  const int* lo = tbx.loVect();
209  const int* hi = tbx.hiVect();
210 
211  for (int i = lo[0]; i <= hi[0]; ++i) {
212  for (int j = lo[1]; j <= hi[1]; ++j) {
213 #if AMREX_SPACEDIM == 3
214  for (int k = lo[2]; k <= hi[2]; ++k) {
215 #endif
216  AmrIntVect_t iv(D_DECL(i, j, k));
217 
218  go_t globalidx = serialize(iv);
219 
220  globalindices.push_back(globalidx);
221 
222  ++localNumElements;
223 #if AMREX_SPACEDIM == 3
224  }
225 #endif
226  }
227  }
228  }
229 
230  /*
231  * create map that specifies which processor gets which data
232  */
233 
234  // get smallest global index of this level
235  amrex::Box bx = grids.minimalBox();
236  const int* lo = bx.loVect();
237  AmrIntVect_t lowcorner(D_DECL(lo[0], lo[1], lo[2]));
238 
239  // where to start indexing
240  go_t baseIndex = serialize(lowcorner);
241 
242  // numGlobalElements == N
243  go_t N = grids.numPts();
244 
245  map_p = Teuchos::rcp( new dmap_t(N, globalindices, baseIndex, comm, node) );
246 }
bool applyBoundary(const AmrIntVect_t &iv, umap_t &map, const scalar_t &value)
std::unordered_map< go_t, scalar_t > umap_t
Teuchos::RCP< matrix_t > G_p[AMREX_SPACEDIM]
gradient matrices in x, y, and z to compute electric field
amrex::IntVect AmrIntVect_t
Definition: AmrDefs.h:28
void serialize(Param_t params, std::ostringstream &os)
serializes params using Boost archive
Definition: MPIHelper.cpp:11
scalar_t invdx_m[AMREX_SPACEDIM]
inverse cell size in particle rest frame
const scalar_t * cellSize() const
amr::local_ordinal_t lo_t
bool isValid(const AmrIntVect_t &iv) const
boundary_t bc_mp[AMREX_SPACEDIM]
boundary conditions
amrex::FabArray< basefab_t > mask_t
amr::AmrIntVect_t AmrIntVect_t
void buildMap(const Teuchos::RCP< comm_t > &comm, const Teuchos::RCP< node_t > &node)
scalar_t dx_m[AMREX_SPACEDIM]
cell size in particle rest frame
std::shared_ptr< AmrBoundary< AmrMultiGridLevel< MatrixType, VectorType > > > boundary_t
Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > dmap_t
bool isBoundary(const AmrIntVect_t &iv) const
Teuchos::RCP< vector_t > residual_p
residual over all cells
amr::global_ordinal_t go_t
amr::scalar_t scalar_t
const AmrGeometry_t & geom
geometry of this problem
Teuchos::RCP< vector_t > error_p
error over all cells
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
amr::AmrGeometry_t AmrGeometry_t
double scalar_t