Alamo
Operator.cpp
Go to the documentation of this file.
1
2#include <AMReX_MLNodeLinOp.H>
3#include <AMReX_MLCellLinOp.H>
4#include <AMReX_MLNodeLap_K.H>
5#include <AMReX_MultiFabUtil.H>
6#include "Util/Color.H"
7#include "Set/Set.H"
8#include "Operator.H"
9
10using namespace amrex;
11namespace Operator {
12
13// constexpr amrex::IntVect AMREX_D_DECL(Operator<Grid::Node>::dx,Operator<Grid::Node>::dy,Operator<Grid::Node>::dz);
15
17{
18 BL_PROFILE(Color::FG::Yellow + "Operator::Diagonal()" + Color::Reset);
19 if (!recompute && m_diagonal_computed) return;
20 m_diagonal_computed = true;
21
22 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
23 {
24 for (int mglev = 0; mglev < m_num_mg_levels[amrlev]; ++mglev)
25 {
26 Diagonal(amrlev, mglev, *m_diag[amrlev][mglev]);
27 }
28 }
29}
30
31void Operator<Grid::Node>::Diagonal(int amrlev, int mglev, amrex::MultiFab& diag)
32{
33 BL_PROFILE("Operator::Diagonal()");
34 //Util::Message(INFO);
35
36 int ncomp = diag.nComp();
37 int nghost = 0;
38
39 int sep = 2;
40 int num = AMREX_D_TERM(sep, *sep, *sep);
41 int cntr = 0;
42
43 amrex::MultiFab x(m_diag[amrlev][mglev]->boxArray(), m_diag[amrlev][mglev]->DistributionMap(), ncomp, nghost);
44 amrex::MultiFab Ax(m_diag[amrlev][mglev]->boxArray(), m_diag[amrlev][mglev]->DistributionMap(), ncomp, nghost);
45
46 for (MFIter mfi(x, false); mfi.isValid(); ++mfi)
47 {
48 const Box& bx = mfi.validbox();
49 amrex::FArrayBox& diagfab = diag[mfi];
50 amrex::FArrayBox& xfab = x[mfi];
51 amrex::FArrayBox& Axfab = Ax[mfi];
52
53 diagfab.setVal(0.0);
54
55 for (int i = 0; i < num; i++)
56 {
57 for (int n = 0; n < ncomp; n++)
58 {
59 xfab.setVal(0.0);
60 Axfab.setVal(0.0);
61
62 //BL_PROFILE_VAR("Operator::Part1", part1);
63 AMREX_D_TERM(for (int m1 = bx.loVect()[0]; m1 <= bx.hiVect()[0]; m1++),
64 for (int m2 = bx.loVect()[1]; m2 <= bx.hiVect()[1]; m2++),
65 for (int m3 = bx.loVect()[2]; m3 <= bx.hiVect()[2]; m3++))
66 {
67 amrex::IntVect m(AMREX_D_DECL(m1, m2, m3));
68
69 if (m1 % sep == i / sep && m2 % sep == i % sep) xfab(m, n) = 1.0;
70 else xfab(m, n) = 0.0;
71 }
72 //BL_PROFILE_VAR_STOP(part1);
73
74 BL_PROFILE_VAR("Operator::Part2", part2);
75 Util::Message(INFO, "Calling fapply...", cntr++);
76 Fapply(amrlev, mglev, Ax, x);
77 BL_PROFILE_VAR_STOP(part2);
78
79 //BL_PROFILE_VAR("Operator::Part3", part3);
80 Axfab.mult(xfab, n, n, 1);
81 diagfab.plus(Axfab, n, n, 1);
82 //BL_PROFILE_VAR_STOP(part3);
83 }
84 }
85 }
86}
87
88void Operator<Grid::Node>::Fsmooth(int amrlev, int mglev, amrex::MultiFab& x, const amrex::MultiFab& b) const
89{
90 BL_PROFILE("Operator::Fsmooth()");
91
92 amrex::Box domain(m_geom[amrlev][mglev].growPeriodicDomain(1));
93 domain.convert(amrex::IntVect::TheNodeVector());
94
95 int ncomp = b.nComp();
96 int nghost = 2; //b.nGrow();
97
98
99 amrex::MultiFab Ax(x.boxArray(), x.DistributionMap(), ncomp, nghost);
100 amrex::MultiFab Dx(x.boxArray(), x.DistributionMap(), ncomp, nghost);
101 amrex::MultiFab Rx(x.boxArray(), x.DistributionMap(), ncomp, nghost);
102
103 if (!m_diagonal_computed) Util::Abort(INFO, "Operator::Diagonal() must be called before using Fsmooth");
104
105 // This is a JACOBI iteration, not Gauss-Seidel.
106 // So we need to do twice the number of iterations to get the same behavior as GS.
107 for (int ctr = 0; ctr < 2; ctr++)
108 {
109 Fapply(amrlev, mglev, Ax, x); // find Ax
110
111 amrex::MultiFab::Copy(Dx, x, 0, 0, ncomp, nghost); // Dx = x
112 amrex::MultiFab::Multiply(Dx, *m_diag[amrlev][mglev], 0, 0, ncomp, nghost); // Dx *= diag (Dx = x*diag)
113
114 amrex::MultiFab::Copy(Rx, Ax, 0, 0, ncomp, nghost); // Rx = Ax
115 amrex::MultiFab::Subtract(Rx, Dx, 0, 0, ncomp, nghost); // Rx -= Dx (Rx = Ax - Dx)
116
117 for (MFIter mfi(x, false); mfi.isValid(); ++mfi)
118 {
119 Box bx = mfi.grownnodaltilebox();
120
121 amrex::Array4<amrex::Real> const& xfab = x.array(mfi);
122 amrex::Array4<const amrex::Real> const& bfab = b.array(mfi);
123 amrex::Array4<const amrex::Real> const& Rxfab = Rx.array(mfi);
124 amrex::Array4<const amrex::Real> const& diagfab = (*m_diag[amrlev][mglev]).array(mfi);
125
126
127 for (int n = 0; n < ncomp; n++)
128 {
129 amrex::ParallelFor(bx, [&] AMREX_GPU_DEVICE(int i, int j, int k)
130 {
131
132 // Skip ghost cells outside problem domain
133 if (!domain.contains(i,j,k))
134 {
135 //continue;
136 }
137 else if ( !bx.strictly_contains(i,j,k))
138 {
139 xfab(i, j, k, n) = 0.0;
140 //continue;
141 }
142 else
143 {
144 xfab(i,j,k,n) = (1. - m_omega) * xfab(i,j,k, n) + m_omega * (bfab(i,j,k, n) - Rxfab(i,j,k, n)) / diagfab(i,j,k,n);
145 }
146 });
147 }
148 }
149 amrex::Geometry geom = m_geom[amrlev][mglev];
150 x.setMultiGhost(true);
151 x.FillBoundary(geom.periodicity());
152 nodalSync(amrlev, mglev, x);
153 }
154}
155
156void Operator<Grid::Node>::normalize(int amrlev, int mglev, MultiFab& a_x) const
157{
158 if (!m_diagonal_computed)
159 Util::Abort(INFO, "Operator::Diagonal() must be called before using normalize");
160
161 a_x.divide(*m_diag[amrlev][mglev],0,getNComp(),2);
162
163 a_x.setMultiGhost(true);
164 a_x.FillBoundaryAndSync(Geom(amrlev,mglev).periodicity());
165}
166
167Operator<Grid::Node>::Operator(const Vector<Geometry>& a_geom,
168 const Vector<BoxArray>& a_grids,
169 const Vector<DistributionMapping>& a_dmap,
170 const LPInfo& a_info,
171 const Vector<FabFactory<FArrayBox> const*>& a_factory)
172{
173 BL_PROFILE("Operator::Operator()");
175
176 if (!(a_grids[0].ixType() == amrex::IndexType::TheNodeType()))
177 Util::Abort(INFO, "Operator must be defined using CELL CENTERED boxarrays.");
178
179 define(a_geom, a_grids, a_dmap, a_info, a_factory);
180}
181
184
185void Operator<Grid::Node>::define(const Vector<Geometry>& a_geom,
186 const Vector<BoxArray>& a_grids,
187 const Vector<DistributionMapping>& a_dmap,
188 const LPInfo& a_info,
189 const Vector<FabFactory<FArrayBox> const*>& a_factory)
190{
191 BL_PROFILE("Operator::~Operator()");
192
193 // Make sure we're not trying to parallelize in vain.
194 if (amrex::ParallelDescriptor::NProcs() > a_grids[0].size())
195 {
196 Util::Warning(INFO, "There are more processors than there are boxes in the amrlev=0 boxarray!!\n",
197 "(NProcs = ", amrex::ParallelDescriptor::NProcs(), ", a_grids[0].size() = ", a_grids[0].size(), ")\n",
198 "You should decrease max_grid_size or you will not get proper scaling!");
199 }
200
201 // This makes sure grids are node-centered;
202 Vector<BoxArray> cc_grids = a_grids;
203 for (auto& ba : cc_grids) {
204 ba.enclosedCells();
205 }
206
207 MLNodeLinOp::define(a_geom, a_grids, a_dmap, a_info, a_factory);
208
209 int nghost = 2;
210 // Resize the multifab containing the operator diagonal
211 m_diag.resize(m_num_amr_levels);
212 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
213 {
214 m_diag[amrlev].resize(m_num_mg_levels[amrlev]);
215
216 for (int mglev = 0; mglev < m_num_mg_levels[amrlev]; ++mglev)
217 {
218 m_diag[amrlev][mglev].reset(new MultiFab(amrex::convert(m_grids[amrlev][mglev], amrex::IntVect::TheNodeVector()),
219 m_dmap[amrlev][mglev], getNComp(), nghost));
220 }
221 }
222
223 // We need to instantiate the m_lobc objects.
224 // WE DO NOT USE THEM - our BCs are implemented differently.
225 // But they need to be the right size or the code will segfault.
226 m_lobc.resize(getNComp(), { {AMREX_D_DECL(BCType::bogus,BCType::bogus,BCType::bogus)} });
227 m_hibc.resize(getNComp(), { {AMREX_D_DECL(BCType::bogus,BCType::bogus,BCType::bogus)} });
228}
229
230void Operator<Grid::Node>::fixUpResidualMask(int amrlev, iMultiFab& resmsk)
231{
232 BL_PROFILE("Operator::fixUpResidualMask()");
233
234 if (!m_masks_built) buildMasks();
235
236 const iMultiFab& cfmask = *m_nd_fine_mask[amrlev];
237
238#ifdef _OPENMP
239#pragma omp parallel if (Gpu::notInLaunchRegion())
240#endif
241 for (MFIter mfi(resmsk, TilingIfNotGPU()); mfi.isValid(); ++mfi)
242 {
243 const Box& bx = mfi.tilebox();
244 Array4<int> const& rmsk = resmsk.array(mfi);
245 Array4<int const> const& fmsk = cfmask.const_array(mfi);
246 AMREX_HOST_DEVICE_PARALLEL_FOR_3D(bx, i, j, k,
247 {
248 if (fmsk(i,j,k) == amrex::nodelap_detail::crse_fine_node) rmsk(i,j,k) = 1;
249 });
250 }
251}
252
254{
255 BL_PROFILE("Operator::prepareForSolve()");
256 MLNodeLinOp::prepareForSolve();
257 buildMasks();
258 averageDownCoeffs();
259 Diagonal(true);
260}
261
262void Operator<Grid::Node>::restriction(int amrlev, int cmglev, MultiFab& crse, MultiFab& fine) const
263{
264 BL_PROFILE("Operator::restriction()");
265
266 applyBC(amrlev, cmglev - 1, fine, BCMode::Homogeneous, StateMode::Solution);
267
268 amrex::Box cdomain = m_geom[amrlev][cmglev].growPeriodicDomain(1);
269 cdomain = cdomain.convert(amrex::IntVect::TheNodeVector());
270
271 bool need_parallel_copy = !amrex::isMFIterSafe(crse, fine);
272 MultiFab cfine;
273 if (need_parallel_copy) {
274 const BoxArray& ba = amrex::coarsen(fine.boxArray(), 2);
275 cfine.define(ba, fine.DistributionMap(), fine.nComp(), fine.nGrow());
276 }
277
278 MultiFab* pcrse = (need_parallel_copy) ? &cfine : &crse;
279
280 for (MFIter mfi(*pcrse, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
281 {
282 const Box& bx = mfi.grownnodaltilebox(-1,1) & cdomain;
283
284 amrex::Array4<const amrex::Real> const& fdata = fine.array(mfi);
285 amrex::Array4<amrex::Real> const& cdata = pcrse->array(mfi);
286
287 const Dim3 lo = amrex::lbound(bx), hi = amrex::ubound(bx);
288
289
290 for (int n = 0; n < crse.nComp(); n++)
291 {
292 // I,J,K == coarse coordinates
293 // i,j,k == fine coordinates
294 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int I, int J, int K) {
295 int i = 2 * I, j = 2 * J, k = 2 * K;
296
297 if ((I == lo.x || I == hi.x) &&
298 (J == lo.y || J == hi.y) &&
299 (K == lo.z || K == hi.z)) // Corner
300 {
301 cdata(I, J, K, n) = fdata(i, j, k, n);
302 }
303 else if ((J == lo.y || J == hi.y) &&
304 (K == lo.z || K == hi.z)) // X edge
305 {
306 cdata(I, J, K, n) = 0.25 * fdata(i - 1, j, k, n) + 0.5 * fdata(i, j, k, n) + 0.25 * fdata(i + 1, j, k, n);
307 }
308 else if ((K == lo.z || K == hi.z) &&
309 (I == lo.x || I == hi.x)) // Y edge
310 {
311 cdata(I, J, K, n) = 0.25 * fdata(i, j - 1, k, n) + 0.5 * fdata(i, j, k, n) + 0.25 * fdata(i, j + 1, k, n);
312 }
313 else if ((I == lo.x || I == hi.x) &&
314 (J == lo.y || J == hi.y)) // Z edge
315 {
316 cdata(I, J, K, n) = 0.25 * fdata(i, j, k - 1, n) + 0.5 * fdata(i, j, k, n) + 0.25 * fdata(i, j, k + 1, n);
317 }
318 else if (I == lo.x || I == hi.x) // X face
319 {
320 cdata(I, J, K, n) =
321 (+fdata(i, j - 1, k - 1, n) + 2.0 * fdata(i, j, k - 1, n) + fdata(i, j + 1, k - 1, n)
322 + 2.0 * fdata(i, j - 1, k, n) + 4.0 * fdata(i, j, k, n) + 2.0 * fdata(i, j + 1, k, n)
323 + fdata(i, j - 1, k + 1, n) + 2.0 * fdata(i, j, k + 1, n) + fdata(i, j + 1, k + 1, n)) / 16.0;
324 }
325 else if (J == lo.y || J == hi.y) // Y face
326 {
327 cdata(I, J, K, n) =
328 (+fdata(i - 1, j, k - 1, n) + 2.0 * fdata(i - 1, j, k, n) + fdata(i - 1, j, k + 1, n)
329 + 2.0 * fdata(i, j, k - 1, n) + 4.0 * fdata(i, j, k, n) + 2.0 * fdata(i, j, k + 1, n)
330 + fdata(i + 1, j, k - 1, n) + 2.0 * fdata(i + 1, j, k, n) + fdata(i + 1, j, k + 1, n)) / 16.0;
331 }
332 else if (K == lo.z || K == hi.z) // Z face
333 {
334 cdata(I, J, K, n) =
335 (+fdata(i - 1, j - 1, k, n) + 2.0 * fdata(i, j - 1, k, n) + fdata(i + 1, j - 1, k, n)
336 + 2.0 * fdata(i - 1, j, k, n) + 4.0 * fdata(i, j, k, n) + 2.0 * fdata(i + 1, j, k, n)
337 + fdata(i - 1, j + 1, k, n) + 2.0 * fdata(i, j + 1, k, n) + fdata(i + 1, j + 1, k, n)) / 16.0;
338 }
339 else // Interior
340 cdata(I, J, K, n) =
341 (fdata(i - 1, j - 1, k - 1, n) + fdata(i - 1, j - 1, k + 1, n) + fdata(i - 1, j + 1, k - 1, n) + fdata(i - 1, j + 1, k + 1, n) +
342 fdata(i + 1, j - 1, k - 1, n) + fdata(i + 1, j - 1, k + 1, n) + fdata(i + 1, j + 1, k - 1, n) + fdata(i + 1, j + 1, k + 1, n)) / 64.0
343 +
344 (fdata(i, j - 1, k - 1, n) + fdata(i, j - 1, k + 1, n) + fdata(i, j + 1, k - 1, n) + fdata(i, j + 1, k + 1, n) +
345 fdata(i - 1, j, k - 1, n) + fdata(i + 1, j, k - 1, n) + fdata(i - 1, j, k + 1, n) + fdata(i + 1, j, k + 1, n) +
346 fdata(i - 1, j - 1, k, n) + fdata(i - 1, j + 1, k, n) + fdata(i + 1, j - 1, k, n) + fdata(i + 1, j + 1, k, n)) / 32.0
347 +
348 (fdata(i - 1, j, k, n) + fdata(i, j - 1, k, n) + fdata(i, j, k - 1, n) +
349 fdata(i + 1, j, k, n) + fdata(i, j + 1, k, n) + fdata(i, j, k + 1, n)) / 16.0
350 +
351 fdata(i, j, k, n) / 8.0;
352 });
353 }
354 }
355
356 if (need_parallel_copy) {
357 crse.ParallelCopy(cfine);
358 }
359
360 crse.setMultiGhost(true);
361 crse.FillBoundary(Geom(amrlev,cmglev).periodicity());
362 nodalSync(amrlev, cmglev, crse);
363}
364
365void Operator<Grid::Node>::interpolation(int amrlev, int fmglev, MultiFab& fine, const MultiFab& crse) const
366{
367 BL_PROFILE("Operator::interpolation()");
368 amrex::Box fdomain = m_geom[amrlev][fmglev].growPeriodicDomain(2);
369 fdomain.convert(amrex::IntVect::TheNodeVector());
370
371 bool need_parallel_copy = !amrex::isMFIterSafe(crse, fine);
372 MultiFab cfine;
373 const MultiFab* cmf = &crse;
374 if (need_parallel_copy) {
375 const BoxArray& ba = amrex::coarsen(fine.boxArray(), 2);
376 cfine.define(ba, fine.DistributionMap(), crse.nComp(), crse.nGrow());
377 cfine.ParallelCopy(crse);
378 cmf = &cfine;
379 }
380
381 for (MFIter mfi(fine, false); mfi.isValid(); ++mfi)
382 {
383 Box fine_bx = mfi.validbox() & fdomain;
384
385 const Box& course_bx = amrex::coarsen(fine_bx, 2);
386 const Box& tmpbx = amrex::refine(course_bx, 2);
387 FArrayBox tmpfab;
388 tmpfab.resize(tmpbx, fine.nComp());
389 tmpfab.setVal(0.0);
390 const amrex::FArrayBox& crsefab = (*cmf)[mfi];
391
392 amrex::Array4<const amrex::Real> const& cdata = crsefab.const_array();
393 amrex::Array4<amrex::Real> const& fdata = tmpfab.array();
394
395 for (int n = 0; n < crse.nComp(); n++)
396 {
397 // I,J,K == coarse coordinates
398 // i,j,k == fine coordinates
399 amrex::ParallelFor(fine_bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
400
401 int I = i / 2, J = j / 2, K = k / 2;
402
403 if (i % 2 == 0 && j % 2 == 0 && k % 2 == 0) // Coincident
404 fdata(i, j, k, n) = cdata(I, J, K, n);
405 else if (j % 2 == 0 && k % 2 == 0) // X Edge
406 fdata(i, j, k, n) = 0.5 * (cdata(I, J, K, n) + cdata(I + 1, J, K, n));
407 else if (k % 2 == 0 && i % 2 == 0) // Y Edge
408 fdata(i, j, k, n) = 0.5 * (cdata(I, J, K, n) + cdata(I, J + 1, K, n));
409 else if (i % 2 == 0 && j % 2 == 0) // Z Edge
410 fdata(i, j, k, n) = 0.5 * (cdata(I, J, K, n) + cdata(I, J, K + 1, n));
411 else if (i % 2 == 0) // X Face
412 fdata(i, j, k, n) = 0.25 * (cdata(I, J, K, n) + cdata(I, J + 1, K, n) +
413 cdata(I, J, K + 1, n) + cdata(I, J + 1, K + 1, n));
414 else if (j % 2 == 0) // Y Face
415 fdata(i, j, k, n) = 0.25 * (cdata(I, J, K, n) + cdata(I, J, K + 1, n) +
416 cdata(I + 1, J, K, n) + cdata(I + 1, J, K + 1, n));
417 else if (k % 2 == 0) // Z Face
418 fdata(i, j, k, n) = 0.25 * (cdata(I, J, K, n) + cdata(I + 1, J, K, n) +
419 cdata(I, J + 1, K, n) + cdata(I + 1, J + 1, K, n));
420 else // Center
421 fdata(i, j, k, n) = 0.125 * (cdata(I, J, K, n) +
422 cdata(I + 1, J, K, n) + cdata(I, J + 1, K, n) + cdata(I, J, K + 1, n) +
423 cdata(I, J + 1, K + 1, n) + cdata(I + 1, J, K + 1, n) + cdata(I + 1, J + 1, K, n) +
424 cdata(I + 1, J + 1, K + 1, n));
425
426 });
427 }
428 fine[mfi].plus(tmpfab, fine_bx, fine_bx, 0, 0, fine.nComp());
429 }
430
431 fine.setMultiGhost(true);
432 fine.FillBoundary(Geom(amrlev,fmglev).periodicity());
433 nodalSync(amrlev, fmglev, fine);
434}
435
436void Operator<Grid::Node>::averageDownSolutionRHS(int camrlev, MultiFab& crse_sol, MultiFab& /*crse_rhs*/,
437 const MultiFab& fine_sol, const MultiFab& /*fine_rhs*/)
438{
439 BL_PROFILE("Operator::averageDownSolutionRHS()");
440 const auto& amrrr = AMRRefRatio(camrlev);
441 amrex::average_down(fine_sol, crse_sol, 0, crse_sol.nComp(), amrrr);
442
443 if (isSingular(0))
444 {
445 Util::Abort(INFO, "Singular operators not supported!");
446 }
447
448}
449
450void Operator<Grid::Node>::realFillBoundary(MultiFab& phi, const Geometry& geom)
451{
452 Util::RealFillBoundary(phi, geom);
453}
454
455void Operator<Grid::Node>::applyBC(int amrlev, int mglev, MultiFab& phi, BCMode/* bc_mode*/,
456 amrex::MLLinOp::StateMode /**/, bool skip_fillboundary) const
457{
458 BL_PROFILE("Operator::applyBC()");
459
460 const Geometry& geom = m_geom[amrlev][mglev];
461
462 if (!skip_fillboundary) {
463 //phi.FillBoundary(geom.periodicity());
464 //phi.setMultiGhost(true);
465 phi.FillBoundaryAndSync(geom.periodicity());
466 }
467}
468
469const amrex::FArrayBox&
470Operator<Grid::Node>::GetFab(const int num, const int amrlev, const int mglev, const amrex::MFIter& mfi) const
471{
472 BL_PROFILE("Operator::GetFab()");
474 return m_a_coeffs[num][amrlev][mglev][mfi];
475}
476
477void Operator<Grid::Node>::RegisterNewFab(amrex::Vector<amrex::MultiFab>& input)
478{
479 BL_PROFILE("Operator::RegisterNewFab()");
481 /// \todo assertions here
482 m_a_coeffs.resize(m_a_coeffs.size() + 1);
483 m_a_coeffs[m_num_a_fabs].resize(m_num_amr_levels);
484 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
485 {
486 m_a_coeffs[m_num_a_fabs][amrlev].resize(m_num_mg_levels[amrlev]);
487 for (int mglev = 0; mglev < m_num_mg_levels[amrlev]; ++mglev)
488 m_a_coeffs[m_num_a_fabs][amrlev][mglev].define(m_grids[amrlev][mglev],
489 m_dmap[amrlev][mglev],
490 input[amrlev].nComp(),
491 input[amrlev].nGrow());
492
493 amrex::MultiFab::Copy(m_a_coeffs[m_num_a_fabs][amrlev][0],
494 input[amrlev], 0, 0,
495 input[amrlev].nComp(),
496 input[amrlev].nGrow());
497 }
498 m_num_a_fabs++;
499}
500
501
502void Operator<Grid::Node>::RegisterNewFab(amrex::Vector<std::unique_ptr<amrex::MultiFab> >& input)
503{
504 BL_PROFILE("Operator::RegisterNewFab()");
506 /// \todo assertions here
507 m_a_coeffs.resize(m_a_coeffs.size() + 1);
508 m_a_coeffs[m_num_a_fabs].resize(m_num_amr_levels);
509 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
510 {
511 m_a_coeffs[m_num_a_fabs][amrlev].resize(m_num_mg_levels[amrlev]);
512 for (int mglev = 0; mglev < m_num_mg_levels[amrlev]; ++mglev)
513 m_a_coeffs[m_num_a_fabs][amrlev][mglev].define(m_grids[amrlev][mglev],
514 m_dmap[amrlev][mglev],
515 input[amrlev]->nComp(),
516 input[amrlev]->nGrow());
517
518 amrex::MultiFab::Copy(m_a_coeffs[m_num_a_fabs][amrlev][0],
519 *input[amrlev], 0, 0,
520 input[amrlev]->nComp(),
521 input[amrlev]->nGrow());
522 }
523 m_num_a_fabs++;
524}
525
526void Operator<Grid::Node>::reflux(int crse_amrlev,
527 MultiFab& res, const MultiFab& /*crse_sol*/, const MultiFab& /*crse_rhs*/,
528 MultiFab& fine_res, MultiFab& /*fine_sol*/, const MultiFab& /*fine_rhs*/) const
529{
530 BL_PROFILE("Operator::Elastic::reflux()");
531
532 int ncomp = AMREX_SPACEDIM;
533
534 amrex::Box cdomain(m_geom[crse_amrlev][0].growPeriodicDomain(2));
535 cdomain.convert(amrex::IntVect::TheNodeVector());
536
537 const Geometry& cgeom = m_geom[crse_amrlev][0];
538
539 const BoxArray& fba = fine_res.boxArray();
540 const DistributionMapping& fdm = fine_res.DistributionMap();
541
542 MultiFab fine_res_for_coarse(amrex::coarsen(fba, 2), fdm, ncomp, 2);
543 fine_res_for_coarse.ParallelCopy(res, 0, 0, ncomp, 0, 0, cgeom.periodicity());
544
545 applyBC(crse_amrlev + 1, 0, fine_res, BCMode::Inhomogeneous, StateMode::Solution);
546
547 /// \todo Replace with Enum
548 // const int coarse_coarse_node = 0;
549 const int coarse_fine_node = 1;
550 const int fine_fine_node = 2;
551
552 amrex::iMultiFab nodemask(amrex::coarsen(fba, 2), fdm, 1, 2);
553 nodemask.ParallelCopy(*m_nd_fine_mask[crse_amrlev], 0, 0, 1, 0, 0, cgeom.periodicity());
554
555 for (MFIter mfi(fine_res_for_coarse, false); mfi.isValid(); ++mfi)
556 {
557 const Box& bx = mfi.grownnodaltilebox(-1,1) & cdomain;
558
559 amrex::Array4<const int> const& nmask = nodemask.array(mfi);
560 //amrex::Array4<const int> const& cmask = cellmask.array(mfi);
561
562 amrex::Array4<amrex::Real> const& cdata = fine_res_for_coarse.array(mfi);
563 amrex::Array4<const amrex::Real> const& fdata = fine_res.array(mfi);
564
565 const Dim3 lo = amrex::lbound(bx), hi = amrex::ubound(bx);
566
567 for (int n = 0; n < fine_res.nComp(); n++)
568 {
569 // I,J,K == coarse coordinates
570 // i,j,k == fine coordinates
571 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int I, int J, int K) {
572 int i = I * 2, j = J * 2, k = K * 2;
573
574 if (nmask(I, J, K) == fine_fine_node || nmask(I, J, K) == coarse_fine_node)
575 {
576 if ((I == lo.x || I == hi.x) &&
577 (J == lo.y || J == hi.y) &&
578 (K == lo.z || K == hi.z)) // Corner
579 cdata(I, J, K, n) = fdata(i, j, k, n);
580 else if ((J == lo.y || J == hi.y) &&
581 (K == lo.z || K == hi.z)) // X edge
582 cdata(I, J, K, n) = 0.25 * fdata(i - 1, j, k, n) + 0.5 * fdata(i, j, k, n) + 0.25 * fdata(i + 1, j, k, n);
583 else if ((K == lo.z || K == hi.z) &&
584 (I == lo.x || I == hi.x)) // Y edge
585 cdata(I, J, K, n) = 0.25 * fdata(i, j - 1, k, n) + 0.5 * fdata(i, j, k, n) + 0.25 * fdata(i, j + 1, k, n);
586 else if ((I == lo.x || I == hi.x) &&
587 (J == lo.y || J == hi.y)) // Z edge
588 cdata(I, J, K, n) = 0.25 * fdata(i, j, k - 1, n) + 0.5 * fdata(i, j, k, n) + 0.25 * fdata(i, j, k + 1, n);
589 else if (I == lo.x || I == hi.x) // X face
590 cdata(I, J, K, n) =
591 (+fdata(i, j - 1, k - 1, n) + 2.0 * fdata(i, j, k - 1, n) + fdata(i, j + 1, k - 1, n)
592 + 2.0 * fdata(i, j - 1, k, n) + 4.0 * fdata(i, j, k, n) + 2.0 * fdata(i, j + 1, k, n)
593 + fdata(i, j - 1, k + 1, n) + 2.0 * fdata(i, j, k + 1, n) + fdata(i, j + 1, k + 1, n)) / 16.0;
594 else if (J == lo.y || J == hi.y) // Y face
595 cdata(I, J, K, n) =
596 (+fdata(i - 1, j, k - 1, n) + 2.0 * fdata(i - 1, j, k, n) + fdata(i - 1, j, k + 1, n)
597 + 2.0 * fdata(i, j, k - 1, n) + 4.0 * fdata(i, j, k, n) + 2.0 * fdata(i, j, k + 1, n)
598 + fdata(i + 1, j, k - 1, n) + 2.0 * fdata(i + 1, j, k, n) + fdata(i + 1, j, k + 1, n)) / 16.0;
599 else if (K == lo.z || K == hi.z) // Z face
600 cdata(I, J, K, n) =
601 (+fdata(i - 1, j - 1, k, n) + 2.0 * fdata(i, j - 1, k, n) + fdata(i + 1, j - 1, k, n)
602 + 2.0 * fdata(i - 1, j, k, n) + 4.0 * fdata(i, j, k, n) + 2.0 * fdata(i + 1, j, k, n)
603 + fdata(i - 1, j + 1, k, n) + 2.0 * fdata(i, j + 1, k, n) + fdata(i + 1, j + 1, k, n)) / 16.0;
604 else // Interior
605 cdata(I, J, K, n) =
606 (fdata(i - 1, j - 1, k - 1, n) + fdata(i - 1, j - 1, k + 1, n) + fdata(i - 1, j + 1, k - 1, n) + fdata(i - 1, j + 1, k + 1, n) +
607 fdata(i + 1, j - 1, k - 1, n) + fdata(i + 1, j - 1, k + 1, n) + fdata(i + 1, j + 1, k - 1, n) + fdata(i + 1, j + 1, k + 1, n)) / 64.0
608 +
609 (fdata(i, j - 1, k - 1, n) + fdata(i, j - 1, k + 1, n) + fdata(i, j + 1, k - 1, n) + fdata(i, j + 1, k + 1, n) +
610 fdata(i - 1, j, k - 1, n) + fdata(i + 1, j, k - 1, n) + fdata(i - 1, j, k + 1, n) + fdata(i + 1, j, k + 1, n) +
611 fdata(i - 1, j - 1, k, n) + fdata(i - 1, j + 1, k, n) + fdata(i + 1, j - 1, k, n) + fdata(i + 1, j + 1, k, n)) / 32.0
612 +
613 (fdata(i - 1, j, k, n) + fdata(i, j - 1, k, n) + fdata(i, j, k - 1, n) +
614 fdata(i + 1, j, k, n) + fdata(i, j + 1, k, n) + fdata(i, j, k + 1, n)) / 16.0
615 +
616 fdata(i, j, k, n) / 8.0;
617 }
618
619 });
620 }
621 }
622
623 // Copy the fine residual restricted onto the coarse grid
624 // into the final residual.
625 res.ParallelCopy(fine_res_for_coarse, 0, 0, ncomp, 0, 0, cgeom.periodicity());
626
627 // Sync up ghost nodes
628 res.setMultiGhost(true);
629 res.FillBoundaryAndSync(Geom(crse_amrlev).periodicity());
630 nodalSync(crse_amrlev, 0, res);
631
632 return;
633}
634
635void
636Operator<Grid::Node>::solutionResidual(int amrlev, MultiFab& resid, MultiFab& x, const MultiFab& b,
637 const MultiFab* /*crse_bcdata*/)
638{
639 const int mglev = 0;
640 const int ncomp = b.nComp();
641 apply(amrlev, mglev, resid, x, BCMode::Inhomogeneous, StateMode::Solution);
642 MultiFab::Xpay(resid, -1.0, b, 0, 0, ncomp, 2);
643 resid.setMultiGhost(true);
644 resid.FillBoundaryAndSync(Geom(amrlev).periodicity());
645}
646
647void
648Operator<Grid::Node>::correctionResidual(int amrlev, int mglev, MultiFab& resid, MultiFab& x, const MultiFab& b,
649 BCMode /*bc_mode*/, const MultiFab* /*crse_bcdata*/)
650{
651 resid.setVal(0.0);
652 apply(amrlev, mglev, resid, x, BCMode::Homogeneous, StateMode::Correction);
653 int ncomp = b.nComp();
654 MultiFab::Xpay(resid, -1.0, b, 0, 0, ncomp, resid.nGrow());
655 resid.setMultiGhost(true);
656 resid.FillBoundaryAndSync(Geom(amrlev).periodicity());
657}
658
659
660
661
663{
664 m_ixtype = amrex::IntVect::TheCellVector();
665}
666
667void
668Operator<Grid::Cell>::define(amrex::Vector<amrex::Geometry> a_geom,
669 const amrex::Vector<amrex::BoxArray>& a_grids,
670 const amrex::Vector<amrex::DistributionMapping>& a_dmap,
672 const amrex::LPInfo& a_info,
673 const amrex::Vector<amrex::FabFactory<amrex::FArrayBox> const*>& a_factory)
674{
675 m_bc = &a_bc;
676
677 std::array<int, AMREX_SPACEDIM> is_periodic = m_bc->IsPeriodic();
678
679 MLCellLinOp::define(a_geom, a_grids, a_dmap, a_info, a_factory);
680
681 Util::Warning(INFO, "This section of code has not been tested.");
682 for (int n = 0; n < getNComp(); n++)
683 {
684 m_lobc.push_back({ AMREX_D_DECL(is_periodic[0] ? amrex::LinOpBCType::Periodic : amrex::LinOpBCType::Dirichlet,
685 is_periodic[1] ? amrex::LinOpBCType::Periodic : amrex::LinOpBCType::Dirichlet,
686 is_periodic[2] ? amrex::LinOpBCType::Periodic : amrex::LinOpBCType::Dirichlet) });
687 m_hibc.push_back({ AMREX_D_DECL(is_periodic[0] ? amrex::LinOpBCType::Periodic : amrex::LinOpBCType::Dirichlet,
688 is_periodic[1] ? amrex::LinOpBCType::Periodic : amrex::LinOpBCType::Dirichlet,
689 is_periodic[2] ? amrex::LinOpBCType::Periodic : amrex::LinOpBCType::Dirichlet) });
690 }
691
692 for (int ilev = 0; ilev < a_geom.size(); ++ilev)
693 setLevelBC(ilev, nullptr);
694
695}
696
697
698void
700{
701 MLCellLinOp::prepareForSolve();
702}
703
704Operator<Grid::Cell>::BndryCondLoc::BndryCondLoc(const amrex::BoxArray& ba, const amrex::DistributionMapping& dm)
705 : bcond(ba, dm),
706 bcloc(ba, dm)
707{
708}
709
710void
711Operator<Grid::Cell>::BndryCondLoc::setLOBndryConds(const amrex::Geometry& /*geom*/, const amrex::Real* /*dx*/,
712 const amrex::Array<BCType, AMREX_SPACEDIM>& /*lobc*/,
713 const amrex::Array<BCType, AMREX_SPACEDIM>& /*hibc*/,
714 int /*ratio*/, const amrex::RealVect& /*a_loc*/)
715{
716 Util::Warning(INFO, "This code has not been properlyt tested");
717}
718
719
720void
722{
723 for (int i = 0; i < m_num_a_fabs; i++)
724 {
725 for (int amrlev = m_num_amr_levels - 1; amrlev > 0; --amrlev)
726 {
727 auto& fine_a_coeffs = m_a_coeffs[i][amrlev];
728 averageDownCoeffsSameAmrLevel(fine_a_coeffs);
729 }
730 averageDownCoeffsSameAmrLevel(m_a_coeffs[i][0]);
731 }
732}
733
734void
736{
737 int nmglevs = a.size();
738 for (int mglev = 1; mglev < nmglevs; ++mglev)
739 {
740 amrex::average_down(a[mglev - 1], a[mglev], 0, a[0].nComp(), mg_coarsen_ratio);
741 }
742}
743
744
745
746const amrex::FArrayBox&
747Operator<Grid::Cell>::GetFab(const int num, const int amrlev, const int mglev, const amrex::MFIter& mfi) const
748{
749 return m_a_coeffs[num][amrlev][mglev][mfi];
750}
751
752
753void
754Operator<Grid::Cell>::RegisterNewFab(amrex::Vector<amrex::MultiFab>& input)
755{
756 /// \todo assertions here
757 m_a_coeffs.resize(m_a_coeffs.size() + 1);
758 m_a_coeffs[m_num_a_fabs].resize(m_num_amr_levels);
759 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
760 {
761 m_a_coeffs[m_num_a_fabs][amrlev].resize(m_num_mg_levels[amrlev]);
762 for (int mglev = 0; mglev < m_num_mg_levels[amrlev]; ++mglev)
763 m_a_coeffs[m_num_a_fabs][amrlev][mglev].define(m_grids[amrlev][mglev],
764 m_dmap[amrlev][mglev],
765 input[amrlev].nComp(),
766 input[amrlev].nGrow());
767
768 amrex::MultiFab::Copy(m_a_coeffs[m_num_a_fabs][amrlev][0],
769 input[amrlev], 0, 0,
770 input[amrlev].nComp(),
771 input[amrlev].nGrow());
772 }
773 m_num_a_fabs++;
774}
775
776
777void
778Operator<Grid::Cell>::RegisterNewFab(amrex::Vector<std::unique_ptr<amrex::MultiFab> >& input)
779{
780 /// \todo assertions here
781 m_a_coeffs.resize(m_a_coeffs.size() + 1);
782 m_a_coeffs[m_num_a_fabs].resize(m_num_amr_levels);
783 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
784 {
785 m_a_coeffs[m_num_a_fabs][amrlev].resize(m_num_mg_levels[amrlev]);
786 for (int mglev = 0; mglev < m_num_mg_levels[amrlev]; ++mglev)
787 m_a_coeffs[m_num_a_fabs][amrlev][mglev].define(m_grids[amrlev][mglev],
788 m_dmap[amrlev][mglev],
789 input[amrlev]->nComp(),
790 input[amrlev]->nGrow());
791
792 amrex::MultiFab::Copy(m_a_coeffs[m_num_a_fabs][amrlev][0],
793 *input[amrlev], 0, 0,
794 input[amrlev]->nComp(),
795 input[amrlev]->nGrow());
796 }
797 m_num_a_fabs++;
798}
799
800
801}
#define INFO
Definition Util.H:24
Definition BC.H:43
virtual amrex::Array< int, AMREX_SPACEDIM > IsPeriodic()
Definition BC.H:114
static std::string Yellow
Definition Color.H:22
static std::string Reset
Definition Color.H:8
Documentation for operator namespace.
Definition Diagonal.cpp:14
constexpr amrex::IntVect AMREX_D_DECL(Operator< Grid::Cell >::dx, Operator< Grid::Cell >::dy, Operator< Grid::Cell >::dz)
void Abort(const char *msg)
Definition Util.cpp:268
void Warning(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:202
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:129
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T > > &a_mf, const amrex::Geometry &, const int nghost=2)
Definition Util.H:339