Block-Structured AMR Software Framework
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
AMReX_MLCellLinOp.H
Go to the documentation of this file.
1#ifndef AMREX_ML_CELL_LINOP_H_
2#define AMREX_ML_CELL_LINOP_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_MLLinOp.H>
6#include <AMReX_iMultiFab.H>
9#include <AMReX_MLLinOp_K.H>
10#include <AMReX_MLMG_K.H>
11
12#ifndef BL_NO_FORT
13#include <AMReX_MLLinOp_F.H>
14#endif
15
16namespace amrex {
17
18template <typename MF>
19class MLCellLinOpT // NOLINT(cppcoreguidelines-virtual-class-destructor)
20 : public MLLinOpT<MF>
21{
22public:
23
24 using FAB = typename MF::fab_type;
25 using RT = typename MF::value_type;
26
27 using BCType = LinOpBCType;
28 using BCMode = typename MLLinOpT<MF>::BCMode;
31
32 MLCellLinOpT ();
33 ~MLCellLinOpT () override = default;
34
35 MLCellLinOpT (const MLCellLinOpT<MF>&) = delete;
39
40 void define (const Vector<Geometry>& a_geom,
41 const Vector<BoxArray>& a_grids,
42 const Vector<DistributionMapping>& a_dmap,
43 const LPInfo& a_info = LPInfo(),
44 const Vector<FabFactory<FAB> const*>& a_factory = {});
45
46 void setLevelBC (int amrlev, const MF* levelbcdata,
47 const MF* robinbc_a = nullptr,
48 const MF* robinbc_b = nullptr,
49 const MF* robinbc_f = nullptr) final;
50
51 using MLLinOpT<MF>::setLevelBC;
52
53 bool needsUpdate () const override {
55 }
56 void update () override;
57
58 void setGaussSeidel (bool flag) noexcept { m_use_gauss_seidel = flag; }
59
60 virtual bool isCrossStencil () const { return true; }
61 virtual bool isTensorOp () const { return false; }
62
63 void updateSolBC (int amrlev, const MF& crse_bcdata) const;
64 void updateCorBC (int amrlev, const MF& crse_bcdata) const;
65
66 virtual void applyBC (int amrlev, int mglev, MF& in, BCMode bc_mode, StateMode s_mode,
67 const MLMGBndryT<MF>* bndry=nullptr, bool skip_fillboundary=false) const;
68
69 BoxArray makeNGrids (int grid_size) const;
70
71 void restriction (int, int, MF& crse, MF& fine) const override;
72
73 void interpolation (int amrlev, int fmglev, MF& fine, const MF& crse) const override;
74
75 void interpAssign (int amrlev, int fmglev, MF& fine, MF& crse) const override;
76
77 void interpolationAmr (int famrlev, MF& fine, const MF& crse,
78 IntVect const& nghost) const override;
79
80 void averageDownSolutionRHS (int camrlev, MF& crse_sol, MF& crse_rhs,
81 const MF& fine_sol, const MF& fine_rhs) override;
82
83 void apply (int amrlev, int mglev, MF& out, MF& in, BCMode bc_mode,
84 StateMode s_mode, const MLMGBndryT<MF>* bndry=nullptr) const override;
85 void smooth (int amrlev, int mglev, MF& sol, const MF& rhs,
86 bool skip_fillboundary, int niter) const final;
87
88 void solutionResidual (int amrlev, MF& resid, MF& x, const MF& b,
89 const MF* crse_bcdata=nullptr) override;
90
91 void prepareForFluxes (int amrlev, const MF* crse_bcdata = nullptr) override;
92
93 void correctionResidual (int amrlev, int mglev, MF& resid, MF& x, const MF& b,
94 BCMode bc_mode, const MF* crse_bcdata=nullptr) final;
95
96 // The assumption is crse_sol's boundary has been filled, but not fine_sol.
97 void reflux (int crse_amrlev,
98 MF& res, const MF& crse_sol, const MF&,
99 MF&, MF& fine_sol, const MF&) const final;
100 void compFlux (int amrlev, const Array<MF*,AMREX_SPACEDIM>& fluxes,
101 MF& sol, Location loc) const override;
102 void compGrad (int amrlev, const Array<MF*,AMREX_SPACEDIM>& grad,
103 MF& sol, Location loc) const override;
104
105 void applyMetricTerm (int amrlev, int mglev, MF& rhs) const final;
106 void unapplyMetricTerm (int amrlev, int mglev, MF& rhs) const final;
107 Vector<RT> getSolvabilityOffset (int amrlev, int mglev,
108 MF const& rhs) const override;
109 void fixSolvabilityByOffset (int amrlev, int mglev, MF& rhs,
110 Vector<RT> const& offset) const override;
111
112 void prepareForSolve () override;
113
114 RT xdoty (int amrlev, int mglev, const MF& x, const MF& y, bool local) const final;
115
116 RT dotProductPrecond (Vector<MF const*> const& x, Vector<MF const*> const& y) const final;
117
118 RT norm2Precond (Vector<MF const*> const& x) const final;
119
120 virtual void Fapply (int amrlev, int mglev, MF& out, const MF& in) const = 0;
121 virtual void Fsmooth (int amrlev, int mglev, MF& sol, const MF& rhs, int redblack) const = 0;
122 virtual void FFlux (int amrlev, const MFIter& mfi,
123 const Array<FAB*,AMREX_SPACEDIM>& flux,
124 const FAB& sol, Location loc, int face_only=0) const = 0;
125
126 virtual void addInhomogNeumannFlux (int /*amrlev*/,
127 const Array<MF*,AMREX_SPACEDIM>& /*grad*/,
128 MF const& /*sol*/,
129 bool /*mult_bcoef*/) const {}
130
131 RT normInf (int amrlev, MF const& mf, bool local) const override;
132
133 void averageDownAndSync (Vector<MF>& sol) const override;
134
135 void avgDownResAmr (int clev, MF& cres, MF const& fres) const override;
136
137 void beginPrecondBC () override;
138 void endPrecondBC () override;
139
144
146
148
149protected:
150
151 bool m_has_metric_term = false;
152
155
158
160
161 // In case of agglomeration, coarse MG grids on amr level 0 are
162 // not simply coarsened from fine MG grids. So we need to build
163 // bcond and bcloc for each MG level.
167 {
168 public:
169 BndryCondLoc (const BoxArray& ba, const DistributionMapping& dm, int ncomp);
170
171 void setLOBndryConds (const Geometry& geom, const Real* dx,
174 IntVect const& ratio, const RealVect& interior_bloc,
175 const Array<Real,AMREX_SPACEDIM>& domain_bloc_lo,
176 const Array<Real,AMREX_SPACEDIM>& domain_bloc_hi,
177 LinOpBCType crse_fine_bc_type);
178
179 const Vector<BCTuple>& bndryConds (const MFIter& mfi) const noexcept {
180 return bcond[mfi];
181 }
182 const Vector<RealTuple>& bndryLocs (const MFIter& mfi) const noexcept {
183 return bcloc[mfi];
184 }
185 const BCTuple& bndryConds (const MFIter& mfi, int icomp) const noexcept {
186 return bcond[mfi][icomp];
187 }
188 const RealTuple& bndryLocs (const MFIter& mfi, int icomp) const noexcept {
189 return bcloc[mfi][icomp];
190 }
191 GpuArray<BCTL,2*AMREX_SPACEDIM> const* getBCTLPtr (const MFIter& mfi) const noexcept {
192 return bctl[mfi];
193 }
194 private:
200 };
202
203 // used to save interpolation coefficients of the first interior cells
205
206 // boundary cell flags for covered, not_covered, outside_domain
208
210
212
213 bool m_use_gauss_seidel = true; // use red-black Gauss-Seidel by default
214
215private:
216
217 void defineAuxData ();
218 void defineBC ();
219
220 void computeVolInv () const;
221 mutable Vector<Vector<RT> > m_volinv; // used by solvability fix
222
224};
225
226template <typename T>
245
246template <typename T>
264
265#ifdef AMREX_USE_EB
266template <typename T>
267struct MLMGPSEBTag {
268 Array4<T> flo;
269 Array4<T> fhi;
270 Array4<T const> ap;
271 Array4<int const> mlo;
272 Array4<int const> mhi;
273 T bcllo;
274 T bclhi;
275 Box bx;
276 BoundCond bctlo;
277 BoundCond bcthi;
278 int blen;
279 int comp;
280 int dir;
281
283 Box const& box() const noexcept { return bx; }
284};
285#endif
286
287template <typename MF>
289 const DistributionMapping& dm,
290 int ncomp)
291 : bcond(ba, dm),
292 bcloc(ba, dm),
293 bctl(ba, dm),
294 bctl_dv(bctl.local_size()*ncomp),
295 m_ncomp(ncomp)
296{
297 auto* dp = bctl_dv.data();
298 for (MFIter mfi(bcloc); mfi.isValid(); ++mfi) {
299 bcond[mfi].resize(ncomp);
300 bcloc[mfi].resize(ncomp);
301 bctl[mfi] = dp;
302 dp += ncomp;
303 }
304}
305
306template <typename MF>
307void
309setLOBndryConds (const Geometry& geom, const Real* dx,
312 IntVect const& ratio, const RealVect& interior_bloc,
313 const Array<Real,AMREX_SPACEDIM>& domain_bloc_lo,
314 const Array<Real,AMREX_SPACEDIM>& domain_bloc_hi,
315 LinOpBCType crse_fine_bc_type)
316{
317 const Box& domain = geom.Domain();
318
319#ifdef AMREX_USE_OMP
320#pragma omp parallel
321#endif
322 for (MFIter mfi(bcloc); mfi.isValid(); ++mfi)
323 {
324 const Box& bx = mfi.validbox();
325 for (int icomp = 0; icomp < m_ncomp; ++icomp) {
326 RealTuple & bloc = bcloc[mfi][icomp];
327 BCTuple & bctag = bcond[mfi][icomp];
328 MLMGBndryT<MF>::setBoxBC(bloc, bctag, bx, domain,
329 lobc[icomp], hibc[icomp],
330 dx, ratio, interior_bloc,
331 domain_bloc_lo, domain_bloc_hi,
332 geom.isPeriodicArray(),
333 crse_fine_bc_type);
334 }
335 }
336
338 hv.reserve(bctl_dv.size());
339 for (MFIter mfi(bctl); mfi.isValid(); ++mfi)
340 {
341 for (int icomp = 0; icomp < m_ncomp; ++icomp) {
343 for (int m = 0; m < 2*AMREX_SPACEDIM; ++m) {
344 tmp[m].type = bcond[mfi][icomp][m];
345 tmp[m].location = bcloc[mfi][icomp][m];
346 }
347 hv.push_back(std::move(tmp));
348 }
349 }
350 Gpu::copyAsync(Gpu::hostToDevice, hv.begin(), hv.end(), bctl_dv.begin());
352}
353
354template <typename MF>
359
360template <typename MF>
361void
363 const Vector<BoxArray>& a_grids,
364 const Vector<DistributionMapping>& a_dmap,
365 const LPInfo& a_info,
366 const Vector<FabFactory<FAB> const*>& a_factory)
367{
368 MLLinOpT<MF>::define(a_geom, a_grids, a_dmap, a_info, a_factory);
370 defineBC();
371}
372
373template <typename MF>
374void
376{
377 BL_PROFILE("MLCellLinOp::defineAuxData()");
378
379 m_undrrelxr.resize(this->m_num_amr_levels);
380 m_maskvals.resize(this->m_num_amr_levels);
381 m_fluxreg.resize(this->m_num_amr_levels-1);
382 m_norm_fine_mask.resize(this->m_num_amr_levels-1);
383
384 const int ncomp = this->getNComp();
385
386 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
387 {
388 m_undrrelxr[amrlev].resize(this->m_num_mg_levels[amrlev]);
389 for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
390 {
391 m_undrrelxr[amrlev][mglev].define(this->m_grids[amrlev][mglev],
392 this->m_dmap[amrlev][mglev],
393 1, 0, 0, ncomp);
394 }
395 }
396
397 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
398 {
399 m_maskvals[amrlev].resize(this->m_num_mg_levels[amrlev]);
400 for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
401 {
402 for (OrientationIter oitr; oitr; ++oitr)
403 {
404 const Orientation face = oitr();
405 const int ngrow = 1;
406 const int extent = this->isCrossStencil() ? 0 : 1; // extend to corners
407 m_maskvals[amrlev][mglev][face].define(this->m_grids[amrlev][mglev],
408 this->m_dmap[amrlev][mglev],
409 this->m_geom[amrlev][mglev],
410 face, 0, ngrow, extent, 1, true);
411 }
412 }
413 }
414
415 for (int amrlev = 0; amrlev < this->m_num_amr_levels-1; ++amrlev)
416 {
417 const IntVect ratio{this->m_amr_ref_ratio[amrlev]};
418 m_fluxreg[amrlev].define(this->m_grids[amrlev+1][0],
419 this->m_grids[amrlev][0],
420 this->m_dmap[amrlev+1][0],
421 this->m_dmap[amrlev][0],
422 this->m_geom[amrlev+1][0],
423 this->m_geom[amrlev][0],
424 ratio, amrlev+1, ncomp);
425 m_norm_fine_mask[amrlev] = std::make_unique<iMultiFab>
426 (makeFineMask(this->m_grids[amrlev][0], this->m_dmap[amrlev][0],
427 this->m_grids[amrlev+1][0],
428 ratio, 1, 0));
429 }
430
431#if (AMREX_SPACEDIM != 3)
432 m_has_metric_term = !this->m_geom[0][0].IsCartesian() && this->info.has_metric_term;
433#endif
434}
435
436template <typename MF>
437void
439{
440 BL_PROFILE("MLCellLinOp::defineBC()");
441
442 const int ncomp = this->getNComp();
443
444 m_bndry_sol.resize(this->m_num_amr_levels);
445 m_crse_sol_br.resize(this->m_num_amr_levels);
446
447 m_bndry_cor.resize(this->m_num_amr_levels);
448 m_crse_cor_br.resize(this->m_num_amr_levels);
449
450 m_robin_bcval.resize(this->m_num_amr_levels);
451
452 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
453 {
454 m_bndry_sol[amrlev] = std::make_unique<MLMGBndryT<MF>>(this->m_grids[amrlev][0],
455 this->m_dmap[amrlev][0],
456 ncomp,
457 this->m_geom[amrlev][0]);
458 }
459
460 for (int amrlev = 1; amrlev < this->m_num_amr_levels; ++amrlev)
461 {
462 const int in_rad = 0;
463 const int out_rad = 1;
464 const int extent_rad = 2;
465 const int crse_ratio = this->m_amr_ref_ratio[amrlev-1];
466 BoxArray cba = this->m_grids[amrlev][0];
467 cba.coarsen(crse_ratio);
468 m_crse_sol_br[amrlev] = std::make_unique<BndryRegisterT<MF>>
469 (cba, this->m_dmap[amrlev][0], in_rad, out_rad, extent_rad, ncomp);
470 }
471
472 for (int amrlev = 1; amrlev < this->m_num_amr_levels; ++amrlev)
473 {
474 const int in_rad = 0;
475 const int out_rad = 1;
476 const int extent_rad = 2;
477 const int crse_ratio = this->m_amr_ref_ratio[amrlev-1];
478 BoxArray cba = this->m_grids[amrlev][0];
479 cba.coarsen(crse_ratio);
480 m_crse_cor_br[amrlev] = std::make_unique<BndryRegisterT<MF>>
481 (cba, this->m_dmap[amrlev][0], in_rad, out_rad, extent_rad, ncomp);
482 m_crse_cor_br[amrlev]->setVal(RT(0.0));
483 }
484
485 // This has be to done after m_crse_cor_br is defined.
486 for (int amrlev = 1; amrlev < this->m_num_amr_levels; ++amrlev)
487 {
488 m_bndry_cor[amrlev] = std::make_unique<MLMGBndryT<MF>>
489 (this->m_grids[amrlev][0], this->m_dmap[amrlev][0], ncomp, this->m_geom[amrlev][0]);
490 MF bc_data(this->m_grids[amrlev][0], this->m_dmap[amrlev][0], ncomp, 1);
491 bc_data.setVal(0.0);
492
493 m_bndry_cor[amrlev]->setBndryValues(*m_crse_cor_br[amrlev], 0, bc_data, 0, 0, ncomp,
494 IntVect(this->m_amr_ref_ratio[amrlev-1]),
497
499 (ncomp,Array<LinOpBCType,AMREX_SPACEDIM>{{AMREX_D_DECL(BCType::Dirichlet,
500 BCType::Dirichlet,
501 BCType::Dirichlet)}});
502 m_bndry_cor[amrlev]->setLOBndryConds(bclohi, bclohi, this->m_amr_ref_ratio[amrlev-1], RealVect{});
503 }
504
505 m_bcondloc.resize(this->m_num_amr_levels);
506 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
507 {
508 m_bcondloc[amrlev].resize(this->m_num_mg_levels[amrlev]);
509 for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
510 {
511 m_bcondloc[amrlev][mglev] = std::make_unique<BndryCondLoc>(this->m_grids[amrlev][mglev],
512 this->m_dmap[amrlev][mglev],
513 ncomp);
514 }
515 }
516}
517
518template <typename MF>
519void
520MLCellLinOpT<MF>::setLevelBC (int amrlev, const MF* a_levelbcdata, const MF* robinbc_a,
521 const MF* robinbc_b, const MF* robinbc_f)
522{
523 BL_PROFILE("MLCellLinOp::setLevelBC()");
524
525 AMREX_ALWAYS_ASSERT(amrlev >= 0 && amrlev < this->m_num_amr_levels);
526
527 const int ncomp = this->getNComp();
528
529 MF zero;
530 IntVect ng(1);
531 if (this->hasHiddenDimension()) { ng[this->hiddenDirection()] = 0; }
532 if (a_levelbcdata == nullptr) {
533 zero.define(this->m_grids[amrlev][0], this->m_dmap[amrlev][0], ncomp, ng);
534 zero.setVal(RT(0.0));
535 } else {
536 AMREX_ALWAYS_ASSERT(a_levelbcdata->nGrowVect().allGE(ng));
537 }
538 const MF& bcdata = (a_levelbcdata == nullptr) ? zero : *a_levelbcdata;
539
540 IntVect br_ref_ratio(-1);
541
542 if (amrlev == 0)
543 {
544 if (this->needsCoarseDataForBC())
545 {
546 // AMREX_ALWAYS_ASSERT(!this->hasHiddenDimension());
547 if (this->hasHiddenDimension()) {
548 int hidden_dir = this->hiddenDirection();
549 AMREX_ALWAYS_ASSERT(this->m_coarse_data_crse_ratio[hidden_dir] == 1);
550 }
551 br_ref_ratio = this->m_coarse_data_crse_ratio.allGT(0) ? this->m_coarse_data_crse_ratio : IntVect(2);
552 if (m_crse_sol_br[amrlev] == nullptr && br_ref_ratio.allGT(0))
553 {
554 const int in_rad = 0;
555 const int out_rad = 1;
556 const int extent_rad = 2;
557 const IntVect crse_ratio = br_ref_ratio;
558 BoxArray cba = this->m_grids[amrlev][0];
559 cba.coarsen(crse_ratio);
560 m_crse_sol_br[amrlev] = std::make_unique<BndryRegisterT<MF>>
561 (cba, this->m_dmap[amrlev][0], in_rad, out_rad, extent_rad, ncomp);
562 }
563 if (this->m_coarse_data_for_bc != nullptr) {
565 const Box& cbx = amrex::coarsen(this->m_geom[0][0].Domain(), this->m_coarse_data_crse_ratio);
566 m_crse_sol_br[amrlev]->copyFrom(*this->m_coarse_data_for_bc, 0, 0, 0, ncomp,
567 this->m_geom[0][0].periodicity(cbx));
568 } else {
569 m_crse_sol_br[amrlev]->setVal(RT(0.0));
570 }
571 m_bndry_sol[amrlev]->setBndryValues(*m_crse_sol_br[amrlev], 0,
572 bcdata, 0, 0, ncomp, br_ref_ratio,
575 br_ref_ratio = this->m_coarse_data_crse_ratio;
576 }
577 else
578 {
579 m_bndry_sol[amrlev]->setPhysBndryValues(bcdata,0,0,ncomp);
580 br_ref_ratio = IntVect(1);
581 }
582 }
583 else
584 {
585 m_bndry_sol[amrlev]->setPhysBndryValues(bcdata,0,0,ncomp);
586 br_ref_ratio = IntVect(this->m_amr_ref_ratio[amrlev-1]);
587 }
588
589 auto crse_fine_bc_type = (amrlev == 0) ? this->m_coarse_fine_bc_type : LinOpBCType::Dirichlet;
590 m_bndry_sol[amrlev]->setLOBndryConds(this->m_lobc, this->m_hibc, br_ref_ratio,
591 this->m_coarse_bc_loc, crse_fine_bc_type);
592
593 const Real* dx = this->m_geom[amrlev][0].CellSize();
594 for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
595 {
596 m_bcondloc[amrlev][mglev]->setLOBndryConds(this->m_geom[amrlev][mglev], dx,
597 this->m_lobc, this->m_hibc,
598 br_ref_ratio, this->m_coarse_bc_loc,
600 crse_fine_bc_type);
601 }
602
603 if (this->hasRobinBC()) {
604 AMREX_ASSERT(robinbc_a != nullptr && robinbc_b != nullptr && robinbc_f != nullptr);
605 m_robin_bcval[amrlev] = std::make_unique<MF>(this->m_grids[amrlev][0],
606 this->m_dmap[amrlev][0],
607 ncomp*3, 1);
608 const Box& domain = this->m_geom[amrlev][0].Domain();
609 MFItInfo mfi_info;
610 if (Gpu::notInLaunchRegion()) { mfi_info.SetDynamic(true); }
611#ifdef AMREX_USE_OMP
612#pragma omp parallel if (Gpu::notInLaunchRegion())
613#endif
614 for (MFIter mfi(*m_robin_bcval[amrlev], mfi_info); mfi.isValid(); ++mfi) {
615 Box const& vbx = mfi.validbox();
616 Array4<RT const> const& ra = robinbc_a->const_array(mfi);
617 Array4<RT const> const& rb = robinbc_b->const_array(mfi);
618 Array4<RT const> const& rf = robinbc_f->const_array(mfi);
619 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
620 const Box& blo = amrex::adjCellLo(vbx, idim);
621 const Box& bhi = amrex::adjCellHi(vbx, idim);
622 bool outside_domain_lo = !(domain.contains(blo));
623 bool outside_domain_hi = !(domain.contains(bhi));
624 if ((!outside_domain_lo) && (!outside_domain_hi)) { continue; }
625 for (int icomp = 0; icomp < ncomp; ++icomp) {
626 Array4<RT> const& rbc = (*m_robin_bcval[amrlev])[mfi].array(icomp*3);
627 if (this->m_lobc_orig[icomp][idim] == LinOpBCType::Robin && outside_domain_lo)
628 {
630 {
631 rbc(i,j,k,0) = ra(i,j,k,icomp);
632 rbc(i,j,k,1) = rb(i,j,k,icomp);
633 rbc(i,j,k,2) = rf(i,j,k,icomp);
634 });
635 }
636 if (this->m_hibc_orig[icomp][idim] == LinOpBCType::Robin && outside_domain_hi)
637 {
639 {
640 rbc(i,j,k,0) = ra(i,j,k,icomp);
641 rbc(i,j,k,1) = rb(i,j,k,icomp);
642 rbc(i,j,k,2) = rf(i,j,k,icomp);
643 });
644 }
645 }
646 }
647 }
648 }
649}
650
651template <typename MF>
652void
657
658template <typename MF>
659void
660MLCellLinOpT<MF>::updateSolBC (int amrlev, const MF& crse_bcdata) const
661{
662 BL_PROFILE("MLCellLinOp::updateSolBC()");
663
664 AMREX_ALWAYS_ASSERT(amrlev > 0);
665 const int ncomp = this->getNComp();
666 m_crse_sol_br[amrlev]->copyFrom(crse_bcdata, 0, 0, 0, ncomp,
667 this->m_geom[amrlev-1][0].periodicity());
668 m_bndry_sol[amrlev]->updateBndryValues(*m_crse_sol_br[amrlev], 0, 0, ncomp,
669 IntVect(this->m_amr_ref_ratio[amrlev-1]),
672}
673
674template <typename MF>
675void
676MLCellLinOpT<MF>::updateCorBC (int amrlev, const MF& crse_bcdata) const
677{
678 BL_PROFILE("MLCellLinOp::updateCorBC()");
679 AMREX_ALWAYS_ASSERT(amrlev > 0);
680 const int ncomp = this->getNComp();
681 m_crse_cor_br[amrlev]->copyFrom(crse_bcdata, 0, 0, 0, ncomp,
682 this->m_geom[amrlev-1][0].periodicity());
683 m_bndry_cor[amrlev]->updateBndryValues(*m_crse_cor_br[amrlev], 0, 0, ncomp,
684 IntVect(this->m_amr_ref_ratio[amrlev-1]),
687}
688
689template <typename MF>
690void
691MLCellLinOpT<MF>::applyBC (int amrlev, int mglev, MF& in, BCMode bc_mode, StateMode,
692 const MLMGBndryT<MF>* bndry, bool skip_fillboundary) const
693{
694 BL_PROFILE("MLCellLinOp::applyBC()");
695 // No coarsened boundary values, cannot apply inhomog at mglev>0.
696 BL_ASSERT(mglev == 0 || bc_mode == BCMode::Homogeneous);
697 BL_ASSERT(bndry != nullptr || bc_mode == BCMode::Homogeneous);
698
699 const int ncomp = this->getNComp();
700 const int cross = isCrossStencil();
701 const int tensorop = isTensorOp();
702 if (!skip_fillboundary) {
703 in.FillBoundary(0, ncomp, this->m_geom[amrlev][mglev].periodicity(), cross);
704 }
705
706 int flagbc = bc_mode == BCMode::Inhomogeneous;
707 const int imaxorder = this->maxorder;
708
709 const Real* dxinv = this->m_geom[amrlev][mglev].InvCellSize();
710 const RT dxi = static_cast<RT>(dxinv[0]);
711 const RT dyi = (AMREX_SPACEDIM >= 2) ? static_cast<RT>(dxinv[1]) : RT(1.0);
712 const RT dzi = (AMREX_SPACEDIM == 3) ? static_cast<RT>(dxinv[2]) : RT(1.0);
713
714 const auto& maskvals = m_maskvals[amrlev][mglev];
715 const auto& bcondloc = *m_bcondloc[amrlev][mglev];
716
717 FAB foofab(Box::TheUnitBox(),ncomp);
718 const auto& foo = foofab.const_array();
719
720 MFItInfo mfi_info;
721 if (Gpu::notInLaunchRegion()) { mfi_info.SetDynamic(true); }
722
724 "non-cross stencil not support for gpu");
725
726 const int hidden_direction = this->hiddenDirection();
727
728#ifdef AMREX_USE_GPU
729 if ((cross || tensorop) && Gpu::inLaunchRegion())
730 {
732 tags.reserve(in.local_size()*AMREX_SPACEDIM*ncomp);
733
734 for (MFIter mfi(in); mfi.isValid(); ++mfi) {
735 const Box& vbx = mfi.validbox();
736 const auto& iofab = in.array(mfi);
737 const auto & bdlv = bcondloc.bndryLocs(mfi);
738 const auto & bdcv = bcondloc.bndryConds(mfi);
739
740 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
741 if (idim != hidden_direction) {
742 const Orientation olo(idim,Orientation::low);
743 const Orientation ohi(idim,Orientation::high);
744 const auto& bvlo = (bndry != nullptr) ?
745 bndry->bndryValues(olo).const_array(mfi) : foo;
746 const auto& bvhi = (bndry != nullptr) ?
747 bndry->bndryValues(ohi).const_array(mfi) : foo;
748 for (int icomp = 0; icomp < ncomp; ++icomp) {
749 tags.emplace_back(MLMGABCTag<RT>{iofab, bvlo, bvhi,
750 maskvals[olo].const_array(mfi),
751 maskvals[ohi].const_array(mfi),
752 bdlv[icomp][olo], bdlv[icomp][ohi],
753 amrex::adjCell(vbx,olo),
754 bdcv[icomp][olo], bdcv[icomp][ohi],
755 vbx.length(idim), icomp, idim});
756 }
757 }
758 }
759 }
760
761 ParallelFor(tags,
762 [=] AMREX_GPU_DEVICE (int i, int j, int k, MLMGABCTag<RT> const& tag) noexcept
763 {
764 if (tag.dir == 0)
765 {
766 mllinop_apply_bc_x(0, i, j, k, tag.blen, tag.fab,
767 tag.mask_lo, tag.bctype_lo, tag.bcloc_lo, tag.bcval_lo,
768 imaxorder, dxi, flagbc, tag.comp);
769 mllinop_apply_bc_x(1, i+tag.blen+1, j, k, tag.blen, tag.fab,
770 tag.mask_hi, tag.bctype_hi, tag.bcloc_hi, tag.bcval_hi,
771 imaxorder, dxi, flagbc, tag.comp);
772 }
773#if (AMREX_SPACEDIM > 1)
774 else
775#if (AMREX_SPACEDIM > 2)
776 if (tag.dir == 1)
777#endif
778 {
779 mllinop_apply_bc_y(0, i, j, k, tag.blen, tag.fab,
780 tag.mask_lo, tag.bctype_lo, tag.bcloc_lo, tag.bcval_lo,
781 imaxorder, dyi, flagbc, tag.comp);
782 mllinop_apply_bc_y(1, i, j+tag.blen+1, k, tag.blen, tag.fab,
783 tag.mask_hi, tag.bctype_hi, tag.bcloc_hi, tag.bcval_hi,
784 imaxorder, dyi, flagbc, tag.comp);
785 }
786#if (AMREX_SPACEDIM > 2)
787 else {
788 mllinop_apply_bc_z(0, i, j, k, tag.blen, tag.fab,
789 tag.mask_lo, tag.bctype_lo, tag.bcloc_lo, tag.bcval_lo,
790 imaxorder, dzi, flagbc, tag.comp);
791 mllinop_apply_bc_z(1, i, j, k+tag.blen+1, tag.blen, tag.fab,
792 tag.mask_hi, tag.bctype_hi, tag.bcloc_hi, tag.bcval_hi,
793 imaxorder, dzi, flagbc, tag.comp);
794 }
795#endif
796#endif
797 });
798 } else
799#endif
800 if (cross || tensorop)
801 {
802#ifdef AMREX_USE_OMP
803#pragma omp parallel if (Gpu::notInLaunchRegion())
804#endif
805 for (MFIter mfi(in, mfi_info); mfi.isValid(); ++mfi)
806 {
807 const Box& vbx = mfi.validbox();
808 const auto& iofab = in.array(mfi);
809
810 const auto & bdlv = bcondloc.bndryLocs(mfi);
811 const auto & bdcv = bcondloc.bndryConds(mfi);
812
813 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
814 {
815 if (hidden_direction == idim) { continue; }
816 const Orientation olo(idim,Orientation::low);
817 const Orientation ohi(idim,Orientation::high);
818 const Box blo = amrex::adjCellLo(vbx, idim);
819 const Box bhi = amrex::adjCellHi(vbx, idim);
820 const int blen = vbx.length(idim);
821 const auto& mlo = maskvals[olo].array(mfi);
822 const auto& mhi = maskvals[ohi].array(mfi);
823 const auto& bvlo = (bndry != nullptr) ? bndry->bndryValues(olo).const_array(mfi) : foo;
824 const auto& bvhi = (bndry != nullptr) ? bndry->bndryValues(ohi).const_array(mfi) : foo;
825 for (int icomp = 0; icomp < ncomp; ++icomp) {
826 const BoundCond bctlo = bdcv[icomp][olo];
827 const BoundCond bcthi = bdcv[icomp][ohi];
828 const RT bcllo = bdlv[icomp][olo];
829 const RT bclhi = bdlv[icomp][ohi];
830 if (idim == 0) {
831 mllinop_apply_bc_x(0, blo, blen, iofab, mlo,
832 bctlo, bcllo, bvlo,
833 imaxorder, dxi, flagbc, icomp);
834 mllinop_apply_bc_x(1, bhi, blen, iofab, mhi,
835 bcthi, bclhi, bvhi,
836 imaxorder, dxi, flagbc, icomp);
837 } else if (idim == 1) {
838 mllinop_apply_bc_y(0, blo, blen, iofab, mlo,
839 bctlo, bcllo, bvlo,
840 imaxorder, dyi, flagbc, icomp);
841 mllinop_apply_bc_y(1, bhi, blen, iofab, mhi,
842 bcthi, bclhi, bvhi,
843 imaxorder, dyi, flagbc, icomp);
844 } else {
845 mllinop_apply_bc_z(0, blo, blen, iofab, mlo,
846 bctlo, bcllo, bvlo,
847 imaxorder, dzi, flagbc, icomp);
848 mllinop_apply_bc_z(1, bhi, blen, iofab, mhi,
849 bcthi, bclhi, bvhi,
850 imaxorder, dzi, flagbc, icomp);
851 }
852 }
853 }
854 }
855 }
856 else
857 {
858#ifdef BL_NO_FORT
859 amrex::Abort("amrex_mllinop_apply_bc not available when BL_NO_FORT=TRUE");
860#else
861 if constexpr (std::is_same_v<Real,RT>) {
862#ifdef AMREX_USE_OMP
863#pragma omp parallel
864#endif
865 for (MFIter mfi(in, mfi_info); mfi.isValid(); ++mfi)
866 {
867 const Box& vbx = mfi.validbox();
868
869 const auto & bdlv = bcondloc.bndryLocs(mfi);
870 const auto & bdcv = bcondloc.bndryConds(mfi);
871
872 const RealTuple & bdl = bdlv[0];
873 const BCTuple & bdc = bdcv[0];
874
875 for (OrientationIter oitr; oitr; ++oitr)
876 {
877 const Orientation ori = oitr();
878
879 int cdr = ori;
880 RT bcl = bdl[ori];
881 int bct = bdc[ori];
882
883 const auto& fsfab = (bndry != nullptr) ? bndry->bndryValues(ori)[mfi] : foofab;
884
885 const Mask& m = maskvals[ori][mfi];
886
887 amrex_mllinop_apply_bc(BL_TO_FORTRAN_BOX(vbx),
888 BL_TO_FORTRAN_ANYD(in[mfi]),
889 BL_TO_FORTRAN_ANYD(m),
890 cdr, bct, bcl,
891 BL_TO_FORTRAN_ANYD(fsfab),
892 imaxorder, dxinv, flagbc, ncomp, cross);
893 }
894 }
895 } else {
896 amrex::Abort("Not supported");
897 }
898#endif
899 }
900}
901
902template <typename MF>
904MLCellLinOpT<MF>::makeNGrids (int grid_size) const
905{
906 const Box& dombx = this->m_geom[0].back().Domain();
907
908 const BoxArray& old_ba = this->m_grids[0].back();
909 const int N = old_ba.size();
910 Vector<Box> bv;
911 bv.reserve(N);
912 for (int i = 0; i < N; ++i)
913 {
914 Box b = old_ba[i];
915 b.coarsen(grid_size);
916 b.refine(grid_size);
917 IntVect sz = b.size();
918 const IntVect nblks {AMREX_D_DECL(sz[0]/grid_size, sz[1]/grid_size, sz[2]/grid_size)};
919
920 IntVect big = b.smallEnd() + grid_size - 1;
921 b.setBig(big);
922
923#if (AMREX_SPACEDIM == 3)
924 for (int kk = 0; kk < nblks[2]; ++kk) {
925#endif
926#if (AMREX_SPACEDIM >= 2)
927 for (int jj = 0; jj < nblks[1]; ++jj) {
928#endif
929 for (int ii = 0; ii < nblks[0]; ++ii)
930 {
931 IntVect shft{AMREX_D_DECL(ii*grid_size,jj*grid_size,kk*grid_size)};
932 Box bb = amrex::shift(b,shft);
933 bb &= dombx;
934 bv.push_back(bb);
935 }
936#if (AMREX_SPACEDIM >= 2)
937 }
938#endif
939#if (AMREX_SPACEDIM == 3)
940 }
941#endif
942 }
943
944 std::sort(bv.begin(), bv.end());
945 bv.erase(std::unique(bv.begin(), bv.end()), bv.end());
946
947 BoxList bl(std::move(bv));
948
949 return BoxArray{std::move(bl)};
950}
951
952template <typename MF>
953void
954MLCellLinOpT<MF>::restriction (int amrlev, int cmglev, MF& crse, MF& fine) const
955{
956 const int ncomp = this->getNComp();
957 IntVect ratio = (amrlev > 0) ? IntVect(2) : this->mg_coarsen_ratio_vec[cmglev-1];
958 amrex::average_down(fine, crse, 0, ncomp, ratio);
959}
960
961template <typename MF>
962void
963MLCellLinOpT<MF>::interpolation (int amrlev, int fmglev, MF& fine, const MF& crse) const
964{
965 const int ncomp = this->getNComp();
966
967 Dim3 ratio3 = {2,2,2};
968 IntVect ratio = (amrlev > 0) ? IntVect(2) : this->mg_coarsen_ratio_vec[fmglev];
969 AMREX_D_TERM(ratio3.x = ratio[0];,
970 ratio3.y = ratio[1];,
971 ratio3.z = ratio[2];);
972
973#ifdef AMREX_USE_GPU
974 if (Gpu::inLaunchRegion() && fine.isFusingCandidate()) {
975 auto const& finema = fine.arrays();
976 auto const& crsema = crse.const_arrays();
977 ParallelFor(fine, IntVect(0), ncomp,
978 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
979 {
980 int ic = amrex::coarsen(i,ratio3.x);
981 int jc = amrex::coarsen(j,ratio3.y);
982 int kc = amrex::coarsen(k,ratio3.z);
983 finema[box_no](i,j,k,n) += crsema[box_no](ic,jc,kc,n);
984 });
986 } else
987#endif
988 {
989#ifdef AMREX_USE_OMP
990#pragma omp parallel if (Gpu::notInLaunchRegion())
991#endif
992 for (MFIter mfi(fine,TilingIfNotGPU()); mfi.isValid(); ++mfi)
993 {
994 const Box& bx = mfi.tilebox();
995 Array4<RT const> const& cfab = crse.const_array(mfi);
996 Array4<RT> const& ffab = fine.array(mfi);
997 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
998 {
999 int ic = amrex::coarsen(i,ratio3.x);
1000 int jc = amrex::coarsen(j,ratio3.y);
1001 int kc = amrex::coarsen(k,ratio3.z);
1002 ffab(i,j,k,n) += cfab(ic,jc,kc,n);
1003 });
1004 }
1005 }
1006}
1007
1008template <typename MF>
1009void
1010MLCellLinOpT<MF>::interpAssign (int amrlev, int fmglev, MF& fine, MF& crse) const
1011{
1012 const int ncomp = this->getNComp();
1013
1014 const Geometry& crse_geom = this->Geom(amrlev,fmglev+1);
1015 const IntVect refratio = (amrlev > 0) ? IntVect(2) : this->mg_coarsen_ratio_vec[fmglev];
1016 const IntVect ng = crse.nGrowVect();
1017
1018 MF cfine;
1019 const MF* cmf;
1020
1022 {
1023 crse.FillBoundary(crse_geom.periodicity());
1024 cmf = &crse;
1025 }
1026 else
1027 {
1028 BoxArray cba = fine.boxArray();
1029 cba.coarsen(refratio);
1030 cfine.define(cba, fine.DistributionMap(), ncomp, ng);
1031 cfine.setVal(RT(0.0));
1032 cfine.ParallelCopy(crse, 0, 0, ncomp, IntVect(0), ng, crse_geom.periodicity());
1033 cmf = & cfine;
1034 }
1035
1036 bool isEB = fine.hasEBFabFactory();
1037 ignore_unused(isEB);
1038
1039#ifdef AMREX_USE_EB
1040 const auto *factory = dynamic_cast<EBFArrayBoxFactory const*>(&(fine.Factory()));
1041 const FabArray<EBCellFlagFab>* flags = (factory) ? &(factory->getMultiEBCellFlagFab()) : nullptr;
1042#endif
1043
1044 MFItInfo mfi_info;
1045 if (Gpu::notInLaunchRegion()) { mfi_info.EnableTiling().SetDynamic(true); }
1046#ifdef AMREX_USE_OMP
1047#pragma omp parallel if (Gpu::notInLaunchRegion())
1048#endif
1049 for (MFIter mfi(fine, mfi_info); mfi.isValid(); ++mfi)
1050 {
1051 const Box& bx = mfi.tilebox();
1052 const auto& ff = fine.array(mfi);
1053 const auto& cc = cmf->array(mfi);
1054#ifdef AMREX_USE_EB
1055 bool call_lincc;
1056 if (isEB)
1057 {
1058 const auto& flag = (*flags)[mfi];
1059 if (flag.getType(amrex::grow(bx,1)) == FabType::regular) {
1060 call_lincc = true;
1061 } else {
1062 Array4<EBCellFlag const> const& flg = flag.const_array();
1064 {
1065 mlmg_eb_cc_interp_r<2>(tbx, ff, cc, flg, ncomp);
1066 });
1067
1068 call_lincc = false;
1069 }
1070 }
1071 else
1072 {
1073 call_lincc = true;
1074 }
1075#else
1076 const bool call_lincc = true;
1077#endif
1078 if (call_lincc)
1079 {
1080#if (AMREX_SPACEDIM == 3)
1081 if (this->hasHiddenDimension()) {
1082 Box const& bx_2d = this->compactify(bx);
1083 auto const& ff_2d = this->compactify(ff);
1084 auto const& cc_2d = this->compactify(cc);
1086 {
1087 TwoD::mlmg_lin_cc_interp_r2(tbx, ff_2d, cc_2d, ncomp);
1088 });
1089 } else
1090#endif
1091 {
1093 {
1094 mlmg_lin_cc_interp_r2(tbx, ff, cc, ncomp);
1095 });
1096 }
1097 }
1098 }
1099}
1100
1101template <typename MF>
1102void
1103MLCellLinOpT<MF>::interpolationAmr (int famrlev, MF& fine, const MF& crse,
1104 IntVect const& /*nghost*/) const
1105{
1106 const int ncomp = this->getNComp();
1107 const int refratio = this->AMRRefRatio(famrlev-1);
1108
1109#ifdef AMREX_USE_EB
1110 const auto *factory = dynamic_cast<EBFArrayBoxFactory const*>(this->Factory(famrlev));
1111 const FabArray<EBCellFlagFab>* flags = (factory) ? &(factory->getMultiEBCellFlagFab()) : nullptr;
1112#endif
1113
1114 MFItInfo mfi_info;
1115 if (Gpu::notInLaunchRegion()) { mfi_info.EnableTiling().SetDynamic(true); }
1116#ifdef AMREX_USE_OMP
1117#pragma omp parallel if (Gpu::notInLaunchRegion())
1118#endif
1119 for (MFIter mfi(fine, mfi_info); mfi.isValid(); ++mfi)
1120 {
1121 const Box& bx = mfi.tilebox();
1122 auto const& ff = fine.array(mfi);
1123 auto const& cc = crse.const_array(mfi);
1124#ifdef AMREX_USE_EB
1125 bool call_lincc;
1126 if (factory)
1127 {
1128 const auto& flag = (*flags)[mfi];
1129 if (flag.getType(amrex::grow(bx,1)) == FabType::regular) {
1130 call_lincc = true;
1131 } else {
1132 Array4<EBCellFlag const> const& flg = flag.const_array();
1133 switch(refratio) {
1134 case 2:
1135 {
1137 {
1138 mlmg_eb_cc_interp_r<2>(tbx, ff, cc, flg, ncomp);
1139 });
1140 break;
1141 }
1142 case 4:
1143 {
1145 {
1146 mlmg_eb_cc_interp_r<4>(tbx, ff, cc, flg, ncomp);
1147 });
1148 break;
1149 }
1150 default:
1151 amrex::Abort("mlmg_eb_cc_interp: only refratio 2 and 4 are supported");
1152 }
1153
1154 call_lincc = false;
1155 }
1156 }
1157 else
1158 {
1159 call_lincc = true;
1160 }
1161#else
1162 const bool call_lincc = true;
1163#endif
1164 if (call_lincc)
1165 {
1166 switch(refratio) {
1167 case 2:
1168 {
1170 {
1171 mlmg_lin_cc_interp_r2(tbx, ff, cc, ncomp);
1172 });
1173 break;
1174 }
1175 case 4:
1176 {
1178 {
1179 mlmg_lin_cc_interp_r4(tbx, ff, cc, ncomp);
1180 });
1181 break;
1182 }
1183 default:
1184 amrex::Abort("mlmg_lin_cc_interp: only refratio 2 and 4 are supported");
1185 }
1186 }
1187 }
1188}
1189
1190template <typename MF>
1191void
1192MLCellLinOpT<MF>::averageDownSolutionRHS (int camrlev, MF& crse_sol, MF& crse_rhs,
1193 const MF& fine_sol, const MF& fine_rhs)
1194{
1195 const auto amrrr = this->AMRRefRatio(camrlev);
1196 const int ncomp = this->getNComp();
1197 amrex::average_down(fine_sol, crse_sol, 0, ncomp, amrrr);
1198 amrex::average_down(fine_rhs, crse_rhs, 0, ncomp, amrrr);
1199}
1200
1201template <typename MF>
1202void
1203MLCellLinOpT<MF>::apply (int amrlev, int mglev, MF& out, MF& in, BCMode bc_mode,
1204 StateMode s_mode, const MLMGBndryT<MF>* bndry) const
1205{
1206 BL_PROFILE("MLCellLinOp::apply()");
1207 applyBC(amrlev, mglev, in, bc_mode, s_mode, bndry);
1208 Fapply(amrlev, mglev, out, in);
1209}
1210
1211template <typename MF>
1212void
1213MLCellLinOpT<MF>::smooth (int amrlev, int mglev, MF& sol, const MF& rhs,
1214 bool skip_fillboundary, int niter) const
1215{
1216 BL_PROFILE("MLCellLinOp::smooth()");
1217 for (int i = 0; i < niter; ++i) {
1218 for (int redblack = 0; redblack < 2; ++redblack)
1219 {
1220 applyBC(amrlev, mglev, sol, BCMode::Homogeneous, StateMode::Solution,
1221 nullptr, skip_fillboundary);
1222 Fsmooth(amrlev, mglev, sol, rhs, redblack);
1223 skip_fillboundary = false;
1224 }
1225 }
1226}
1227
1228template <typename MF>
1229void
1230MLCellLinOpT<MF>::solutionResidual (int amrlev, MF& resid, MF& x, const MF& b,
1231 const MF* crse_bcdata)
1232{
1233 BL_PROFILE("MLCellLinOp::solutionResidual()");
1234 const int ncomp = this->getNComp();
1235 if (crse_bcdata != nullptr) {
1236 updateSolBC(amrlev, *crse_bcdata);
1237 }
1238 const int mglev = 0;
1239 apply(amrlev, mglev, resid, x, BCMode::Inhomogeneous, StateMode::Solution,
1240 m_bndry_sol[amrlev].get());
1241
1242 AMREX_ASSERT(resid.nComp() == b.nComp());
1243 MF::Xpay(resid, RT(-1.0), b, 0, 0, ncomp, IntVect(0));
1244}
1245
1246template <typename MF>
1247void
1248MLCellLinOpT<MF>::prepareForFluxes (int amrlev, const MF* crse_bcdata)
1249{
1250 if (crse_bcdata != nullptr) {
1251 updateSolBC(amrlev, *crse_bcdata);
1252 }
1253}
1254
1255template <typename MF>
1256void
1257MLCellLinOpT<MF>::correctionResidual (int amrlev, int mglev, MF& resid, MF& x, const MF& b,
1258 BCMode bc_mode, const MF* crse_bcdata)
1259{
1260 BL_PROFILE("MLCellLinOp::correctionResidual()");
1261 const int ncomp = this->getNComp();
1262 if (bc_mode == BCMode::Inhomogeneous)
1263 {
1264 if (crse_bcdata)
1265 {
1266 AMREX_ASSERT(mglev == 0 && amrlev > 0);
1267 updateCorBC(amrlev, *crse_bcdata);
1268 }
1269 apply(amrlev, mglev, resid, x, BCMode::Inhomogeneous, StateMode::Correction,
1270 m_bndry_cor[amrlev].get());
1271 }
1272 else
1273 {
1274 AMREX_ASSERT(crse_bcdata == nullptr);
1275 apply(amrlev, mglev, resid, x, BCMode::Homogeneous, StateMode::Correction, nullptr);
1276 }
1277
1278 MF::Xpay(resid, Real(-1.0), b, 0, 0, ncomp, IntVect(0));
1279}
1280
1281template <typename MF>
1282void
1283MLCellLinOpT<MF>::reflux (int crse_amrlev, MF& res, const MF& crse_sol, const MF&,
1284 MF&, MF& fine_sol, const MF&) const
1285{
1286 BL_PROFILE("MLCellLinOp::reflux()");
1287
1288 auto& fluxreg = m_fluxreg[crse_amrlev];
1289 fluxreg.reset();
1290
1291 const int ncomp = this->getNComp();
1292
1293 const int fine_amrlev = crse_amrlev+1;
1294
1295 Real dt = Real(1.0);
1296 const Real* crse_dx = this->m_geom[crse_amrlev][0].CellSize();
1297 const Real* fine_dx = this->m_geom[fine_amrlev][0].CellSize();
1298
1299 const int mglev = 0;
1300 applyBC(fine_amrlev, mglev, fine_sol, BCMode::Inhomogeneous, StateMode::Solution,
1301 m_bndry_sol[fine_amrlev].get());
1302
1303 MFItInfo mfi_info;
1304 if (Gpu::notInLaunchRegion()) { mfi_info.EnableTiling().SetDynamic(true); }
1305
1306#ifdef AMREX_USE_OMP
1307#pragma omp parallel if (Gpu::notInLaunchRegion())
1308#endif
1309 {
1311 Array<FAB*,AMREX_SPACEDIM> pflux {{ AMREX_D_DECL(flux.data(), flux.data()+1, flux.data()+2) }};
1312 Array<FAB const*,AMREX_SPACEDIM> cpflux {{ AMREX_D_DECL(flux.data(), flux.data()+1, flux.data()+2) }};
1313
1314 for (MFIter mfi(crse_sol, mfi_info); mfi.isValid(); ++mfi)
1315 {
1316 if (fluxreg.CrseHasWork(mfi))
1317 {
1318 const Box& tbx = mfi.tilebox();
1319 AMREX_D_TERM(flux[0].resize(amrex::surroundingNodes(tbx,0),ncomp,The_Async_Arena());,
1320 flux[1].resize(amrex::surroundingNodes(tbx,1),ncomp,The_Async_Arena());,
1321 flux[2].resize(amrex::surroundingNodes(tbx,2),ncomp,The_Async_Arena()););
1322 FFlux(crse_amrlev, mfi, pflux, crse_sol[mfi], Location::FaceCentroid);
1323 fluxreg.CrseAdd(mfi, cpflux, crse_dx, dt, RunOn::Gpu);
1324 }
1325 }
1326
1327#ifdef AMREX_USE_OMP
1328#pragma omp barrier
1329#endif
1330
1331 for (MFIter mfi(fine_sol, mfi_info); mfi.isValid(); ++mfi)
1332 {
1333 if (fluxreg.FineHasWork(mfi))
1334 {
1335 const Box& tbx = mfi.tilebox();
1336 const int face_only = true;
1337 AMREX_D_TERM(flux[0].resize(amrex::surroundingNodes(tbx,0),ncomp,The_Async_Arena());,
1338 flux[1].resize(amrex::surroundingNodes(tbx,1),ncomp,The_Async_Arena());,
1339 flux[2].resize(amrex::surroundingNodes(tbx,2),ncomp,The_Async_Arena()););
1340 FFlux(fine_amrlev, mfi, pflux, fine_sol[mfi], Location::FaceCentroid, face_only);
1341 fluxreg.FineAdd(mfi, cpflux, fine_dx, dt, RunOn::Gpu);
1342 }
1343 }
1344 }
1345
1346 fluxreg.Reflux(res);
1347}
1348
1349template <typename MF>
1350void
1352 MF& sol, Location loc) const
1353{
1354 BL_PROFILE("MLCellLinOp::compFlux()");
1355
1356 const int mglev = 0;
1357 const int ncomp = this->getNComp();
1358 applyBC(amrlev, mglev, sol, BCMode::Inhomogeneous, StateMode::Solution,
1359 m_bndry_sol[amrlev].get());
1360
1361 MFItInfo mfi_info;
1362 if (Gpu::notInLaunchRegion()) { mfi_info.EnableTiling().SetDynamic(true); }
1363
1364#ifdef AMREX_USE_OMP
1365#pragma omp parallel if (Gpu::notInLaunchRegion())
1366#endif
1367 {
1369 Array<FAB*,AMREX_SPACEDIM> pflux {{ AMREX_D_DECL(flux.data(), flux.data()+1, flux.data()+2) }};
1370 for (MFIter mfi(sol, mfi_info); mfi.isValid(); ++mfi)
1371 {
1372 const Box& tbx = mfi.tilebox();
1373 AMREX_D_TERM(flux[0].resize(amrex::surroundingNodes(tbx,0),ncomp,The_Async_Arena());,
1374 flux[1].resize(amrex::surroundingNodes(tbx,1),ncomp,The_Async_Arena());,
1375 flux[2].resize(amrex::surroundingNodes(tbx,2),ncomp,The_Async_Arena()););
1376 FFlux(amrlev, mfi, pflux, sol[mfi], loc);
1377 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1378 const Box& nbx = mfi.nodaltilebox(idim);
1379 auto const& dst = fluxes[idim]->array(mfi);
1380 auto const& src = pflux[idim]->const_array();
1381 AMREX_HOST_DEVICE_PARALLEL_FOR_4D (nbx, ncomp, i, j, k, n,
1382 {
1383 dst(i,j,k,n) = src(i,j,k,n);
1384 });
1385 }
1386 }
1387 }
1388}
1389
1390template <typename MF>
1391void
1393 MF& sol, Location /*loc*/) const
1394{
1395 BL_PROFILE("MLCellLinOp::compGrad()");
1396
1397 if (sol.nComp() > 1) {
1398 amrex::Abort("MLCellLinOp::compGrad called, but only works for single-component solves");
1399 }
1400
1401 const int mglev = 0;
1402 applyBC(amrlev, mglev, sol, BCMode::Inhomogeneous, StateMode::Solution,
1403 m_bndry_sol[amrlev].get());
1404
1405 const int ncomp = this->getNComp();
1406
1407 AMREX_D_TERM(const RT dxi = static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(0));,
1408 const RT dyi = static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(1));,
1409 const RT dzi = static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(2)););
1410#ifdef AMREX_USE_OMP
1411#pragma omp parallel if (Gpu::notInLaunchRegion())
1412#endif
1413 for (MFIter mfi(sol, TilingIfNotGPU()); mfi.isValid(); ++mfi)
1414 {
1415 AMREX_D_TERM(const Box& xbx = mfi.nodaltilebox(0);,
1416 const Box& ybx = mfi.nodaltilebox(1);,
1417 const Box& zbx = mfi.nodaltilebox(2););
1418 const auto& s = sol.array(mfi);
1419 AMREX_D_TERM(const auto& gx = grad[0]->array(mfi);,
1420 const auto& gy = grad[1]->array(mfi);,
1421 const auto& gz = grad[2]->array(mfi););
1422
1423 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( xbx, ncomp, i, j, k, n,
1424 {
1425 gx(i,j,k,n) = dxi*(s(i,j,k,n) - s(i-1,j,k,n));
1426 });
1427#if (AMREX_SPACEDIM >= 2)
1428 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( ybx, ncomp, i, j, k, n,
1429 {
1430 gy(i,j,k,n) = dyi*(s(i,j,k,n) - s(i,j-1,k,n));
1431 });
1432#endif
1433#if (AMREX_SPACEDIM == 3)
1434 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( zbx, ncomp, i, j, k, n,
1435 {
1436 gz(i,j,k,n) = dzi*(s(i,j,k,n) - s(i,j,k-1,n));
1437 });
1438#endif
1439 }
1440
1441 addInhomogNeumannFlux(amrlev, grad, sol, false);
1442}
1443
1444template <typename MF>
1445void
1446MLCellLinOpT<MF>::applyMetricTerm (int amrlev, int mglev, MF& rhs) const
1447{
1448 amrex::ignore_unused(amrlev,mglev,rhs);
1449#if (AMREX_SPACEDIM != 3)
1450 if (!m_has_metric_term) { return; }
1451
1452 const int ncomp = rhs.nComp();
1453
1454 bool cc = rhs.ixType().cellCentered(0);
1455
1456 const Geometry& geom = this->m_geom[amrlev][mglev];
1457 const RT dx = static_cast<RT>(geom.CellSize(0));
1458 const RT probxlo = static_cast<RT>(geom.ProbLo(0));
1459
1460#ifdef AMREX_USE_OMP
1461#pragma omp parallel if (Gpu::notInLaunchRegion())
1462#endif
1463 for (MFIter mfi(rhs,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1464 {
1465 const Box& tbx = mfi.tilebox();
1466 auto const& rhsarr = rhs.array(mfi);
1467#if (AMREX_SPACEDIM == 1)
1468 if (cc) {
1469 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1470 {
1471 RT rc = probxlo + (i+RT(0.5))*dx;
1472 rhsarr(i,j,k,n) *= rc*rc;
1473 });
1474 } else {
1475 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1476 {
1477 RT re = probxlo + i*dx;
1478 rhsarr(i,j,k,n) *= re*re;
1479 });
1480 }
1481#elif (AMREX_SPACEDIM == 2)
1482 if (cc) {
1483 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1484 {
1485 RT rc = probxlo + (i+RT(0.5))*dx;
1486 rhsarr(i,j,k,n) *= rc;
1487 });
1488 } else {
1489 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1490 {
1491 RT re = probxlo + i*dx;
1492 rhsarr(i,j,k,n) *= re;
1493 });
1494 }
1495#endif
1496 }
1497#endif
1498}
1499
1500template <typename MF>
1501void
1502MLCellLinOpT<MF>::unapplyMetricTerm (int amrlev, int mglev, MF& rhs) const
1503{
1504 amrex::ignore_unused(amrlev,mglev,rhs);
1505#if (AMREX_SPACEDIM != 3)
1506 if (!m_has_metric_term) { return; }
1507
1508 const int ncomp = rhs.nComp();
1509
1510 bool cc = rhs.ixType().cellCentered(0);
1511
1512 const Geometry& geom = this->m_geom[amrlev][mglev];
1513 const RT dx = static_cast<RT>(geom.CellSize(0));
1514 const RT probxlo = static_cast<RT>(geom.ProbLo(0));
1515
1516#ifdef AMREX_USE_OMP
1517#pragma omp parallel if (Gpu::notInLaunchRegion())
1518#endif
1519 for (MFIter mfi(rhs,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1520 {
1521 const Box& tbx = mfi.tilebox();
1522 auto const& rhsarr = rhs.array(mfi);
1523#if (AMREX_SPACEDIM == 1)
1524 if (cc) {
1525 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1526 {
1527 RT rcinv = RT(1.0)/(probxlo + (i+RT(0.5))*dx);
1528 rhsarr(i,j,k,n) *= rcinv*rcinv;
1529 });
1530 } else {
1531 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1532 {
1533 RT re = probxlo + i*dx;
1534 RT reinv = (re==RT(0.0)) ? RT(0.0) : RT(1.)/re;
1535 rhsarr(i,j,k,n) *= reinv*reinv;
1536 });
1537 }
1538#elif (AMREX_SPACEDIM == 2)
1539 if (cc) {
1540 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1541 {
1542 RT rcinv = RT(1.0)/(probxlo + (i+RT(0.5))*dx);
1543 rhsarr(i,j,k,n) *= rcinv;
1544 });
1545 } else {
1546 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
1547 {
1548 RT re = probxlo + i*dx;
1549 RT reinv = (re==RT(0.0)) ? RT(0.0) : RT(1.)/re;
1550 rhsarr(i,j,k,n) *= reinv;
1551 });
1552 }
1553#endif
1554 }
1555#endif
1556}
1557
1558template <typename MF>
1559auto
1560MLCellLinOpT<MF>::getSolvabilityOffset (int amrlev, int mglev, MF const& rhs) const
1561 -> Vector<RT>
1562{
1563 computeVolInv();
1564
1565 const int ncomp = this->getNComp();
1566 Vector<RT> offset(ncomp);
1567
1568#ifdef AMREX_USE_EB
1569 const auto *factory = dynamic_cast<EBFArrayBoxFactory const*>(this->Factory(amrlev,mglev));
1570 if (factory && !factory->isAllRegular())
1571 {
1572 if constexpr (std::is_same<MF,MultiFab>()) {
1573 const MultiFab& vfrac = factory->getVolFrac();
1574 for (int c = 0; c < ncomp; ++c) {
1575 offset[c] = amrex::Dot(rhs, c, vfrac, 0, 1, IntVect(0), true)
1576 * m_volinv[amrlev][mglev];
1577 }
1578 } else {
1579 amrex::Abort("TODO: MLMG with EB only works with MultiFab");
1580 }
1581 }
1582 else
1583#endif
1584 {
1585 for (int c = 0; c < ncomp; ++c) {
1586 offset[c] = rhs.sum(c,IntVect(0),true) * m_volinv[amrlev][mglev];
1587 }
1588 }
1589
1591
1592 return offset;
1593}
1594
1595template <typename MF>
1596void
1597MLCellLinOpT<MF>::fixSolvabilityByOffset (int /*amrlev*/, int /*mglev*/, MF& rhs,
1598 Vector<RT> const& offset) const
1599{
1600 const int ncomp = this->getNComp();
1601 for (int c = 0; c < ncomp; ++c) {
1602 rhs.plus(-offset[c], c, 1);
1603 }
1604#ifdef AMREX_USE_EB
1605 if (!rhs.isAllRegular()) {
1606 if constexpr (std::is_same<MF,MultiFab>()) {
1607 amrex::EB_set_covered(rhs, 0, ncomp, 0, 0.0_rt);
1608 } else {
1609 amrex::Abort("amrex::EB_set_covered only works with MultiFab");
1610 }
1611 }
1612#endif
1613}
1614
1615template <typename MF>
1616void
1618{
1619 BL_PROFILE("MLCellLinOp::prepareForSolve()");
1620
1621 const int imaxorder = this->maxorder;
1622 const int ncomp = this->getNComp();
1623 const int hidden_direction = this->hiddenDirection();
1624 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
1625 {
1626 for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
1627 {
1628 const auto& bcondloc = *m_bcondloc[amrlev][mglev];
1629 const auto& maskvals = m_maskvals[amrlev][mglev];
1630
1631 const RT dxi = static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(0));
1632 const RT dyi = static_cast<RT>((AMREX_SPACEDIM >= 2) ? this->m_geom[amrlev][mglev].InvCellSize(1) : Real(1.0));
1633 const RT dzi = static_cast<RT>((AMREX_SPACEDIM == 3) ? this->m_geom[amrlev][mglev].InvCellSize(2) : Real(1.0));
1634
1635 auto& undrrelxr = this->m_undrrelxr[amrlev][mglev];
1636 MF foo(this->m_grids[amrlev][mglev], this->m_dmap[amrlev][mglev], ncomp, 0, MFInfo().SetAlloc(false));
1637
1638#ifdef AMREX_USE_EB
1639 const auto *factory = dynamic_cast<EBFArrayBoxFactory const*>(this->m_factory[amrlev][mglev].get());
1640 const FabArray<EBCellFlagFab>* flags =
1641 (factory) ? &(factory->getMultiEBCellFlagFab()) : nullptr;
1642 auto area = (factory) ? factory->getAreaFrac()
1643 : Array<const MultiCutFab*,AMREX_SPACEDIM>{AMREX_D_DECL(nullptr,nullptr,nullptr)};
1645#endif
1646
1647#ifdef AMREX_USE_GPU
1648 if (Gpu::inLaunchRegion()) {
1649#ifdef AMREX_USE_EB
1650 if (factory && !factory->isAllRegular()) {
1651#if defined(AMREX_USE_CUDA) && defined(_WIN32)
1652 if (!std::is_same<MF,MultiFab>()) {
1653#else
1654 if constexpr (!std::is_same<MF,MultiFab>()) {
1655#endif
1656 amrex::Abort("MLCellLinOp with EB only works with MultiFab");
1657 } else {
1659 tags.reserve(foo.local_size()*AMREX_SPACEDIM*ncomp);
1660
1661 for (MFIter mfi(foo); mfi.isValid(); ++mfi)
1662 {
1663 const Box& vbx = mfi.validbox();
1664
1665 const auto & bdlv = bcondloc.bndryLocs(mfi);
1666 const auto & bdcv = bcondloc.bndryConds(mfi);
1667
1668 auto fabtyp = (flags) ? (*flags)[mfi].getType(vbx) : FabType::regular;
1669
1670 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
1671 {
1672 if (idim != hidden_direction && fabtyp != FabType::covered) {
1673 const Orientation olo(idim,Orientation::low);
1674 const Orientation ohi(idim,Orientation::high);
1675 auto const& ap = (fabtyp == FabType::singlevalued)
1676 ? area[idim]->const_array(mfi) : Array4<Real const>{};
1677 for (int icomp = 0; icomp < ncomp; ++icomp) {
1678 tags.emplace_back(MLMGPSEBTag<RT>{undrrelxr[olo].array(mfi),
1679 undrrelxr[ohi].array(mfi),
1680 ap,
1681 maskvals[olo].const_array(mfi),
1682 maskvals[ohi].const_array(mfi),
1683 bdlv[icomp][olo], bdlv[icomp][ohi],
1684 amrex::adjCell(vbx,olo),
1685 bdcv[icomp][olo], bdcv[icomp][ohi],
1686 vbx.length(idim), icomp, idim});
1687 }
1688 }
1689 }
1690 }
1691
1692 ParallelFor(tags,
1693 [=] AMREX_GPU_DEVICE (int i, int j, int k, MLMGPSEBTag<RT> const& tag) noexcept
1694 {
1695 if (tag.ap) {
1696 if (tag.dir == 0)
1697 {
1698 mllinop_comp_interp_coef0_x_eb
1699 (0, i , j, k, tag.blen, tag.flo, tag.mlo, tag.ap,
1700 tag.bctlo, tag.bcllo, imaxorder, dxi, tag.comp);
1701 mllinop_comp_interp_coef0_x_eb
1702 (1, i+tag.blen+1, j, k, tag.blen, tag.fhi, tag.mhi, tag.ap,
1703 tag.bcthi, tag.bclhi, imaxorder, dxi, tag.comp);
1704 }
1705#if (AMREX_SPACEDIM > 1)
1706 else
1707#if (AMREX_SPACEDIM > 2)
1708 if (tag.dir == 1)
1709#endif
1710 {
1711 mllinop_comp_interp_coef0_y_eb
1712 (0, i, j , k, tag.blen, tag.flo, tag.mlo, tag.ap,
1713 tag.bctlo, tag.bcllo, imaxorder, dyi, tag.comp);
1714 mllinop_comp_interp_coef0_y_eb
1715 (1, i, j+tag.blen+1, k, tag.blen, tag.fhi, tag.mhi, tag.ap,
1716 tag.bcthi, tag.bclhi, imaxorder, dyi, tag.comp);
1717 }
1718#if (AMREX_SPACEDIM > 2)
1719 else {
1720 mllinop_comp_interp_coef0_z_eb
1721 (0, i, j, k , tag.blen, tag.flo, tag.mlo, tag.ap,
1722 tag.bctlo, tag.bcllo, imaxorder, dzi, tag.comp);
1723 mllinop_comp_interp_coef0_z_eb
1724 (1, i, j, k+tag.blen+1, tag.blen, tag.fhi, tag.mhi, tag.ap,
1725 tag.bcthi, tag.bclhi, imaxorder, dzi, tag.comp);
1726 }
1727#endif
1728#endif
1729 } else {
1730 if (tag.dir == 0)
1731 {
1733 (0, i , j, k, tag.blen, tag.flo, tag.mlo,
1734 tag.bctlo, tag.bcllo, imaxorder, dxi, tag.comp);
1736 (1, i+tag.blen+1, j, k, tag.blen, tag.fhi, tag.mhi,
1737 tag.bcthi, tag.bclhi, imaxorder, dxi, tag.comp);
1738 }
1739#if (AMREX_SPACEDIM > 1)
1740 else
1741#if (AMREX_SPACEDIM > 2)
1742 if (tag.dir == 1)
1743#endif
1744 {
1746 (0, i, j , k, tag.blen, tag.flo, tag.mlo,
1747 tag.bctlo, tag.bcllo, imaxorder, dyi, tag.comp);
1749 (1, i, j+tag.blen+1, k, tag.blen, tag.fhi, tag.mhi,
1750 tag.bcthi, tag.bclhi, imaxorder, dyi, tag.comp);
1751 }
1752#if (AMREX_SPACEDIM > 2)
1753 else {
1755 (0, i, j, k , tag.blen, tag.flo, tag.mlo,
1756 tag.bctlo, tag.bcllo, imaxorder, dzi, tag.comp);
1758 (1, i, j, k+tag.blen+1, tag.blen, tag.fhi, tag.mhi,
1759 tag.bcthi, tag.bclhi, imaxorder, dzi, tag.comp);
1760 }
1761#endif
1762#endif
1763 }
1764 });
1765 }
1766 } else
1767#endif
1768 {
1770 tags.reserve(foo.local_size()*AMREX_SPACEDIM*ncomp);
1771
1772 for (MFIter mfi(foo); mfi.isValid(); ++mfi)
1773 {
1774 const Box& vbx = mfi.validbox();
1775
1776 const auto & bdlv = bcondloc.bndryLocs(mfi);
1777 const auto & bdcv = bcondloc.bndryConds(mfi);
1778
1779 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
1780 {
1781 if (idim != hidden_direction) {
1782 const Orientation olo(idim,Orientation::low);
1783 const Orientation ohi(idim,Orientation::high);
1784 for (int icomp = 0; icomp < ncomp; ++icomp) {
1785 tags.emplace_back(MLMGPSTag<RT>{undrrelxr[olo].array(mfi),
1786 undrrelxr[ohi].array(mfi),
1787 maskvals[olo].const_array(mfi),
1788 maskvals[ohi].const_array(mfi),
1789 bdlv[icomp][olo], bdlv[icomp][ohi],
1790 amrex::adjCell(vbx,olo),
1791 bdcv[icomp][olo], bdcv[icomp][ohi],
1792 vbx.length(idim), icomp, idim});
1793 }
1794 }
1795 }
1796 }
1797
1798 ParallelFor(tags,
1799 [=] AMREX_GPU_DEVICE (int i, int j, int k, MLMGPSTag<RT> const& tag) noexcept
1800 {
1801 if (tag.dir == 0)
1802 {
1803 mllinop_comp_interp_coef0_x
1804 (0, i , j, k, tag.blen, tag.flo, tag.mlo,
1805 tag.bctlo, tag.bcllo, imaxorder, dxi, tag.comp);
1806 mllinop_comp_interp_coef0_x
1807 (1, i+tag.blen+1, j, k, tag.blen, tag.fhi, tag.mhi,
1808 tag.bcthi, tag.bclhi, imaxorder, dxi, tag.comp);
1809 }
1810#if (AMREX_SPACEDIM > 1)
1811 else
1812#if (AMREX_SPACEDIM > 2)
1813 if (tag.dir == 1)
1814#endif
1815 {
1816 mllinop_comp_interp_coef0_y
1817 (0, i, j , k, tag.blen, tag.flo, tag.mlo,
1818 tag.bctlo, tag.bcllo, imaxorder, dyi, tag.comp);
1819 mllinop_comp_interp_coef0_y
1820 (1, i, j+tag.blen+1, k, tag.blen, tag.fhi, tag.mhi,
1821 tag.bcthi, tag.bclhi, imaxorder, dyi, tag.comp);
1822 }
1823#if (AMREX_SPACEDIM > 2)
1824 else {
1826 (0, i, j, k , tag.blen, tag.flo, tag.mlo,
1827 tag.bctlo, tag.bcllo, imaxorder, dzi, tag.comp);
1829 (1, i, j, k+tag.blen+1, tag.blen, tag.fhi, tag.mhi,
1830 tag.bcthi, tag.bclhi, imaxorder, dzi, tag.comp);
1831 }
1832#endif
1833#endif
1834 });
1835 }
1836 } else
1837#endif
1838 {
1839#ifdef AMREX_USE_OMP
1840#pragma omp parallel
1841#endif
1842 for (MFIter mfi(foo, MFItInfo{}.SetDynamic(true)); mfi.isValid(); ++mfi)
1843 {
1844 const Box& vbx = mfi.validbox();
1845
1846 const auto & bdlv = bcondloc.bndryLocs(mfi);
1847 const auto & bdcv = bcondloc.bndryConds(mfi);
1848
1849#ifdef AMREX_USE_EB
1850 auto fabtyp = (flags) ? (*flags)[mfi].getType(vbx) : FabType::regular;
1851#endif
1852 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
1853 {
1854 if (idim == hidden_direction) { continue; }
1855 const Orientation olo(idim,Orientation::low);
1856 const Orientation ohi(idim,Orientation::high);
1857 const Box blo = amrex::adjCellLo(vbx, idim);
1858 const Box bhi = amrex::adjCellHi(vbx, idim);
1859 const int blen = vbx.length(idim);
1860 const auto& mlo = maskvals[olo].array(mfi);
1861 const auto& mhi = maskvals[ohi].array(mfi);
1862 const auto& flo = undrrelxr[olo].array(mfi);
1863 const auto& fhi = undrrelxr[ohi].array(mfi);
1864 for (int icomp = 0; icomp < ncomp; ++icomp) {
1865 const BoundCond bctlo = bdcv[icomp][olo];
1866 const BoundCond bcthi = bdcv[icomp][ohi];
1867 const auto bcllo = bdlv[icomp][olo];
1868 const auto bclhi = bdlv[icomp][ohi];
1869#ifdef AMREX_USE_EB
1870 if (fabtyp == FabType::singlevalued) {
1871 if constexpr (!std::is_same<MF,MultiFab>()) {
1872 amrex::Abort("MLCellLinOp with EB only works with MultiFab");
1873 } else {
1874 auto const& ap = area[idim]->const_array(mfi);
1875 if (idim == 0) {
1876 mllinop_comp_interp_coef0_x_eb
1877 (0, blo, blen, flo, mlo, ap, bctlo, bcllo,
1878 imaxorder, dxi, icomp);
1879 mllinop_comp_interp_coef0_x_eb
1880 (1, bhi, blen, fhi, mhi, ap, bcthi, bclhi,
1881 imaxorder, dxi, icomp);
1882 } else if (idim == 1) {
1883 mllinop_comp_interp_coef0_y_eb
1884 (0, blo, blen, flo, mlo, ap, bctlo, bcllo,
1885 imaxorder, dyi, icomp);
1886 mllinop_comp_interp_coef0_y_eb
1887 (1, bhi, blen, fhi, mhi, ap, bcthi, bclhi,
1888 imaxorder, dyi, icomp);
1889 } else {
1890 mllinop_comp_interp_coef0_z_eb
1891 (0, blo, blen, flo, mlo, ap, bctlo, bcllo,
1892 imaxorder, dzi, icomp);
1893 mllinop_comp_interp_coef0_z_eb
1894 (1, bhi, blen, fhi, mhi, ap, bcthi, bclhi,
1895 imaxorder, dzi, icomp);
1896 }
1897 }
1898 } else if (fabtyp == FabType::regular)
1899#endif
1900 {
1901 if (idim == 0) {
1903 (0, blo, blen, flo, mlo, bctlo, bcllo,
1904 imaxorder, dxi, icomp);
1906 (1, bhi, blen, fhi, mhi, bcthi, bclhi,
1907 imaxorder, dxi, icomp);
1908 } else if (idim == 1) {
1910 (0, blo, blen, flo, mlo, bctlo, bcllo,
1911 imaxorder, dyi, icomp);
1913 (1, bhi, blen, fhi, mhi, bcthi, bclhi,
1914 imaxorder, dyi, icomp);
1915 } else {
1917 (0, blo, blen, flo, mlo, bctlo, bcllo,
1918 imaxorder, dzi, icomp);
1920 (1, bhi, blen, fhi, mhi, bcthi, bclhi,
1921 imaxorder, dzi, icomp);
1922 }
1923 }
1924 }
1925 }
1926 }
1927 }
1928 }
1929 }
1930}
1931
1932template <typename MF>
1933auto
1934MLCellLinOpT<MF>::xdoty (int /*amrlev*/, int /*mglev*/, const MF& x, const MF& y, bool local) const
1935 -> RT
1936{
1937 const int ncomp = this->getNComp();
1938 const IntVect nghost(0);
1939 RT result = amrex::Dot(x,0,y,0,ncomp,nghost,true);
1940 if (!local) {
1942 }
1943 return result;
1944}
1945
1946template <typename MF>
1947auto
1949{
1950 const int ncomp = this->getNComp();
1951 const IntVect nghost(0);
1952 RT result = 0;
1953 for (int ilev = 0; ilev < this->NAMRLevels()-1; ++ilev) {
1954 result += amrex::Dot(*m_norm_fine_mask[ilev], *x[ilev], 0, *y[ilev], 0, ncomp, nghost, true);
1955 }
1956 result += amrex::Dot(*x[this->NAMRLevels()-1], 0,
1957 *y[this->NAMRLevels()-1], 0, ncomp, nghost, true);
1959 return result;
1960}
1961
1962template <typename MF>
1963auto
1965{
1966 const int ncomp = this->getNComp();
1967 const IntVect nghost(0);
1968 RT result = 0;
1969 for (int ilev = 0; ilev < this->NAMRLevels()-1; ++ilev) {
1970 result += amrex::Dot(*m_norm_fine_mask[ilev], *x[ilev], 0, ncomp, nghost, true);
1971 }
1972 result += amrex::Dot(*x[this->NAMRLevels()-1], 0, ncomp, nghost, true);
1974 return std::sqrt(result);
1975}
1976
1977template <typename MF>
1978void
1980{
1981 if (!m_volinv.empty()) { return; }
1982
1983 m_volinv.resize(this->m_num_amr_levels);
1984 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev) {
1985 m_volinv[amrlev].resize(this->NMGLevels(amrlev));
1986 }
1987
1988 // We don't need to compute for every level
1989
1990 auto f = [&] (int amrlev, int mglev) {
1991#ifdef AMREX_USE_EB
1992 const auto *factory = dynamic_cast<EBFArrayBoxFactory const*>(this->Factory(amrlev,mglev));
1993 if (factory && !factory->isAllRegular())
1994 {
1995 if constexpr (std::is_same<MF,MultiFab>()) {
1996 const auto& vfrac = factory->getVolFrac();
1997 m_volinv[amrlev][mglev] = vfrac.sum(0,true);
1998 } else {
1999 amrex::Abort("MLCellLinOp with EB only works with MultiFab");
2000 }
2001 }
2002 else
2003#endif
2004 {
2005 auto const npts = (this->m_coarse_fine_bc_type == LinOpBCType::Dirichlet)
2006 ? this->compactify(this->Geom(amrlev,mglev).Domain()).d_numPts()
2007 : this->m_grids[amrlev][mglev].d_numPts();
2008 AMREX_ASSERT(npts > 0.);
2009 m_volinv[amrlev][mglev] = RT(1.0 / npts);
2010 }
2011 };
2012
2013 // amrlev = 0, mglev = 0
2014 f(0,0);
2015
2016 int mgbottom = this->NMGLevels(0)-1;
2017 f(0,mgbottom);
2018
2019#ifdef AMREX_USE_EB
2020 RT temp1, temp2;
2021 const auto *factory = dynamic_cast<EBFArrayBoxFactory const*>(this->Factory(0,0));
2022 if (factory && !factory->isAllRegular())
2023 {
2024 ParallelAllReduce::Sum<RT>({m_volinv[0][0], m_volinv[0][mgbottom]},
2026 temp1 = RT(1.0)/m_volinv[0][0];
2027 temp2 = RT(1.0)/m_volinv[0][mgbottom];
2028 }
2029 else
2030 {
2031 temp1 = m_volinv[0][0];
2032 temp2 = m_volinv[0][mgbottom];
2033 }
2034 m_volinv[0][0] = temp1;
2035 m_volinv[0][mgbottom] = temp2;
2036#endif
2037}
2038
2039template <typename MF>
2040auto
2041MLCellLinOpT<MF>::normInf (int amrlev, MF const& mf, bool local) const -> RT
2042{
2043 const int ncomp = this->getNComp();
2044 const int finest_level = this->NAMRLevels() - 1;
2045 RT norm = RT(0.0);
2046#ifdef AMREX_USE_EB
2047 const auto *factory = dynamic_cast<EBFArrayBoxFactory const*>(this->Factory(amrlev));
2048 if (factory && !factory->isAllRegular()) {
2049#if defined(AMREX_USE_CUDA) && defined(_WIN32)
2050 if (!std::is_same<MF,MultiFab>()) {
2051#else
2052 if constexpr (!std::is_same<MF,MultiFab>()) {
2053#endif
2054 amrex::Abort("MLCellLinOpT with EB only works with MultiFab");
2055 } else {
2056 const MultiFab& vfrac = factory->getVolFrac();
2057 if (amrlev == finest_level) {
2058#ifdef AMREX_USE_GPU
2059 if (Gpu::inLaunchRegion()) {
2060 auto const& ma = mf.const_arrays();
2061 auto const& vfrac_ma = vfrac.const_arrays();
2063 mf, IntVect(0), ncomp,
2064 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n)
2066 {
2067 return std::abs(ma[box_no](i,j,k,n)
2068 *vfrac_ma[box_no](i,j,k));
2069 });
2070 } else
2071#endif
2072 {
2073#ifdef AMREX_USE_OMP
2074#pragma omp parallel reduction(max:norm)
2075#endif
2076 for (MFIter mfi(mf,true); mfi.isValid(); ++mfi) {
2077 Box const& bx = mfi.tilebox();
2078 auto const& fab = mf.const_array(mfi);
2079 auto const& v = vfrac.const_array(mfi);
2080 AMREX_LOOP_4D(bx, ncomp, i, j, k, n,
2081 {
2082 norm = std::max(norm, std::abs(fab(i,j,k,n)*v(i,j,k)));
2083 });
2084 }
2085 }
2086 } else {
2087#ifdef AMREX_USE_GPU
2088 if (Gpu::inLaunchRegion()) {
2089 auto const& ma = mf.const_arrays();
2090 auto const& mask_ma = m_norm_fine_mask[amrlev]->const_arrays();
2091 auto const& vfrac_ma = vfrac.const_arrays();
2093 mf, IntVect(0), ncomp,
2094 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n)
2096 {
2097 if (mask_ma[box_no](i,j,k)) {
2098 return std::abs(ma[box_no](i,j,k,n)
2099 *vfrac_ma[box_no](i,j,k));
2100 } else {
2101 return Real(0.0);
2102 }
2103 });
2104 } else
2105#endif
2106 {
2107#ifdef AMREX_USE_OMP
2108#pragma omp parallel reduction(max:norm)
2109#endif
2110 for (MFIter mfi(mf,true); mfi.isValid(); ++mfi) {
2111 Box const& bx = mfi.tilebox();
2112 auto const& fab = mf.const_array(mfi);
2113 auto const& mask = m_norm_fine_mask[amrlev]->const_array(mfi);
2114 auto const& v = vfrac.const_array(mfi);
2115 AMREX_LOOP_4D(bx, ncomp, i, j, k, n,
2116 {
2117 if (mask(i,j,k)) {
2118 norm = std::max(norm, std::abs(fab(i,j,k,n)*v(i,j,k)));
2119 }
2120 });
2121 }
2122 }
2123 }
2124 }
2125 } else
2126#endif
2127 {
2128 if (amrlev == finest_level) {
2129 norm = mf.norminf(0, ncomp, IntVect(0), true);
2130 } else {
2131 norm = mf.norminf(*m_norm_fine_mask[amrlev], 0, ncomp, IntVect(0), true);
2132 }
2133 }
2134
2136 return norm;
2137}
2138
2139template <typename MF>
2140void
2142{
2143 int ncomp = this->getNComp();
2144 for (int falev = this->NAMRLevels()-1; falev > 0; --falev)
2145 {
2146#ifdef AMREX_USE_EB
2147 if (!sol[falev].isAllRegular()) {
2148 if constexpr (std::is_same<MF,MultiFab>()) {
2149 amrex::EB_average_down(sol[falev], sol[falev-1], 0, ncomp, this->AMRRefRatio(falev-1));
2150 } else {
2151 amrex::Abort("EB_average_down only works with MultiFab");
2152 }
2153 } else
2154#endif
2155 {
2156 amrex::average_down(sol[falev], sol[falev-1], 0, ncomp, this->AMRRefRatio(falev-1));
2157 }
2158 }
2159}
2160
2161template <typename MF>
2162void
2163MLCellLinOpT<MF>::avgDownResAmr (int clev, MF& cres, MF const& fres) const
2164{
2165#ifdef AMREX_USE_EB
2166 if (!fres.isAllRegular()) {
2167 if constexpr (std::is_same<MF,MultiFab>()) {
2168 amrex::EB_average_down(fres, cres, 0, this->getNComp(),
2169 this->AMRRefRatio(clev));
2170 } else {
2171 amrex::Abort("EB_average_down only works with MultiFab");
2172 }
2173 } else
2174#endif
2175 {
2176 amrex::average_down(fres, cres, 0, this->getNComp(),
2177 this->AMRRefRatio(clev));
2178 }
2179}
2180
2181template <typename MF>
2182void
2184{
2185 this->m_precond_mode = true;
2186
2187 if (m_bndry_sol_zero.empty()) {
2188 m_bndry_sol_zero.resize(m_bndry_sol.size());
2189 const int ncomp = this->getNComp();
2190 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev) {
2191 m_bndry_sol_zero[amrlev] = std::make_unique<MLMGBndryT<MF>>
2192 (this->m_grids[amrlev][0],
2193 this->m_dmap[amrlev][0],
2194 ncomp,
2195 this->m_geom[amrlev][0]);
2196 }
2197 std::swap(m_bndry_sol, m_bndry_sol_zero);
2198 MF const* coarse_data_for_bc_save = this->m_coarse_data_for_bc;
2199 this->m_coarse_data_for_bc = nullptr;
2200 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev) {
2201 this->setLevelBC(amrlev, nullptr);
2202 }
2203 this->m_coarse_data_for_bc = coarse_data_for_bc_save;
2204 } else {
2205 std::swap(m_bndry_sol, m_bndry_sol_zero);
2206 }
2207}
2208
2209template <typename MF>
2210void
2212{
2213 this->m_precond_mode = false;
2214 std::swap(m_bndry_sol, m_bndry_sol_zero);
2215}
2216
2217extern template class MLCellLinOpT<MultiFab>;
2218
2220
2221}
2222
2223#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_HOST_DEVICE_PARALLEL_FOR_3D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:110
#define AMREX_GPU_LAUNCH_HOST_DEVICE_LAMBDA_RANGE(TN, TI, block)
Definition AMReX_GpuLaunchMacrosC.nolint.H:4
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1089
Box cbx
Definition AMReX_HypreMLABecLap.cpp:1091
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< int const > mask
Definition AMReX_InterpFaceRegister.cpp:93
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
#define AMREX_LOOP_4D(bx, ncomp, i, j, k, n, block)
Definition AMReX_Loop.nolint.H:16
void amrex_mllinop_apply_bc(const int *lo, const int *hi, amrex_real *phi, const int *philo, const int *phihi, const int *mask, const int *mlo, const int *mhi, int cdir, int bct, amrex_real bcl, const amrex_real *bcval, const int *blo, const int *bhi, int maxorder, const amrex_real *dxinv, int inhomog, int nc, int cross)
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:129
if(!(yy_init))
Definition amrex_iparser.lex.nolint.H:935
const FabSetT< MF > & bndryValues(Orientation face) const noexcept
Return FabSet on given face.
Definition AMReX_BndryData.H:77
Maintain an identifier for boundary condition types.
Definition AMReX_BoundCond.H:20
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:550
void define(const Box &bx)
Initialize the BoxArray from a single box. It is an error if the BoxArray has already been initialize...
BoxArray & refine(int refinement_ratio)
Refine each Box in the BoxArray to the specified ratio.
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:597
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
static AMREX_GPU_HOST_DEVICE BoxND TheUnitBox() noexcept
This static member function returns a constant reference to an object of type BoxND representing the ...
Definition AMReX_Box.H:737
AMREX_GPU_HOST_DEVICE double d_numPts() const noexcept
Returns the number of points contained in the BoxND. This is intended for use only in diagnostic mess...
Definition AMReX_Box.H:366
AMREX_GPU_HOST_DEVICE IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:146
AMREX_GPU_HOST_DEVICE bool contains(const IntVectND< dim > &p) const noexcept
Returns true if argument is contained within BoxND.
Definition AMReX_Box.H:204
const Real * CellSize() const noexcept
Returns the cellsize for each coordinate direction.
Definition AMReX_CoordSys.H:71
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
Definition AMReX_EBFabFactory.H:24
const MultiFab & getVolFrac() const noexcept
Definition AMReX_EBFabFactory.H:55
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:344
Array4< typename FabArray< FAB >::value_type const > const_array(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:584
MultiArray4< typename FabArray< FAB >::value_type const > const_arrays() const noexcept
Definition AMReX_FabArray.H:646
Definition AMReX_FabFactory.H:50
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:210
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:355
GpuArray< int, AMREX_SPACEDIM > isPeriodicArray() const noexcept
Definition AMReX_Geometry.H:347
const Real * ProbLo() const noexcept
Returns the lo end of the problem domain in each dimension.
Definition AMReX_Geometry.H:178
Definition AMReX_Tuple.H:93
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool allGT(const IntVectND< dim > &rhs) const noexcept
Returns true if this is greater than argument for all components. NOTE: This is NOT a strict weak ord...
Definition AMReX_IntVect.H:416
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr std::size_t size() noexcept
Definition AMReX_IntVect.H:723
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE IntVectND & setVal(int i, int val) noexcept
Set i'th coordinate of IntVectND to val.
Definition AMReX_IntVect.H:272
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr IntVectND< dim > TheCellVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:709
An InterpBndryData object adds to a BndryData object the ability to manipulate and set the data store...
Definition AMReX_InterpBndryData.H:40
a one-thingy-per-box distributed object
Definition AMReX_LayoutData.H:13
Definition AMReX_MFIter.H:57
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:141
Definition AMReX_MLCellLinOp.H:167
const Vector< BCTuple > & bndryConds(const MFIter &mfi) const noexcept
Definition AMReX_MLCellLinOp.H:179
GpuArray< BCTL, 2 *AMREX_SPACEDIM > const * getBCTLPtr(const MFIter &mfi) const noexcept
Definition AMReX_MLCellLinOp.H:191
const BCTuple & bndryConds(const MFIter &mfi, int icomp) const noexcept
Definition AMReX_MLCellLinOp.H:185
Gpu::DeviceVector< GpuArray< BCTL, 2 *AMREX_SPACEDIM > > bctl_dv
Definition AMReX_MLCellLinOp.H:198
BndryCondLoc(const BoxArray &ba, const DistributionMapping &dm, int ncomp)
Definition AMReX_MLCellLinOp.H:288
LayoutData< GpuArray< BCTL, 2 *AMREX_SPACEDIM > * > bctl
Definition AMReX_MLCellLinOp.H:197
const Vector< RealTuple > & bndryLocs(const MFIter &mfi) const noexcept
Definition AMReX_MLCellLinOp.H:182
int m_ncomp
Definition AMReX_MLCellLinOp.H:199
void setLOBndryConds(const Geometry &geom, const Real *dx, const Vector< Array< BCType, AMREX_SPACEDIM > > &lobc, const Vector< Array< BCType, AMREX_SPACEDIM > > &hibc, IntVect const &ratio, const RealVect &interior_bloc, const Array< Real, AMREX_SPACEDIM > &domain_bloc_lo, const Array< Real, AMREX_SPACEDIM > &domain_bloc_hi, LinOpBCType crse_fine_bc_type)
Definition AMReX_MLCellLinOp.H:309
LayoutData< Vector< RealTuple > > bcloc
Definition AMReX_MLCellLinOp.H:196
LayoutData< Vector< BCTuple > > bcond
Definition AMReX_MLCellLinOp.H:195
const RealTuple & bndryLocs(const MFIter &mfi, int icomp) const noexcept
Definition AMReX_MLCellLinOp.H:188
Definition AMReX_MLCellLinOp.H:21
virtual void Fsmooth(int amrlev, int mglev, MF &sol, const MF &rhs, int redblack) const =0
Vector< RT > getSolvabilityOffset(int amrlev, int mglev, MF const &rhs) const override
get offset for fixing solvability
Definition AMReX_MLCellLinOp.H:1560
typename MF::fab_type FAB
Definition AMReX_MLCellLinOp.H:24
void averageDownSolutionRHS(int camrlev, MF &crse_sol, MF &crse_rhs, const MF &fine_sol, const MF &fine_rhs) override
Average-down data from fine AMR level to coarse AMR level.
Definition AMReX_MLCellLinOp.H:1192
void averageDownAndSync(Vector< MF > &sol) const override
Definition AMReX_MLCellLinOp.H:2141
BoxArray makeNGrids(int grid_size) const
Definition AMReX_MLCellLinOp.H:904
Vector< YAFluxRegisterT< MF > > m_fluxreg
Definition AMReX_MLCellLinOp.H:211
virtual void applyBC(int amrlev, int mglev, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr, bool skip_fillboundary=false) const
Definition AMReX_MLCellLinOp.H:691
void updateSolBC(int amrlev, const MF &crse_bcdata) const
Definition AMReX_MLCellLinOp.H:660
Vector< std::unique_ptr< MLMGBndryT< MF > > > m_bndry_sol
Definition AMReX_MLCellLinOp.H:153
MLCellLinOpT< MF > & operator=(const MLCellLinOpT< MF > &)=delete
void compGrad(int amrlev, const Array< MF *, AMREX_SPACEDIM > &grad, MF &sol, Location loc) const override
Compute gradients of the solution.
Definition AMReX_MLCellLinOp.H:1392
void avgDownResAmr(int clev, MF &cres, MF const &fres) const override
Definition AMReX_MLCellLinOp.H:2163
void reflux(int crse_amrlev, MF &res, const MF &crse_sol, const MF &, MF &, MF &fine_sol, const MF &) const final
Reflux at AMR coarse/fine boundary.
Definition AMReX_MLCellLinOp.H:1283
RT dotProductPrecond(Vector< MF const * > const &x, Vector< MF const * > const &y) const final
Definition AMReX_MLCellLinOp.H:1948
virtual void Fapply(int amrlev, int mglev, MF &out, const MF &in) const =0
void defineBC()
Definition AMReX_MLCellLinOp.H:438
void endPrecondBC() override
Definition AMReX_MLCellLinOp.H:2211
void update() override
Update for reuse.
Definition AMReX_MLCellLinOp.H:653
typename MF::value_type RT
Definition AMReX_MLCellLinOp.H:25
LinOpBCType BCType
Definition AMReX_MLCellLinOp.H:27
void correctionResidual(int amrlev, int mglev, MF &resid, MF &x, const MF &b, BCMode bc_mode, const MF *crse_bcdata=nullptr) final
Compute residual for the residual-correction form, resid = b - L(x)
Definition AMReX_MLCellLinOp.H:1257
MLCellLinOpT(const MLCellLinOpT< MF > &)=delete
typename MLLinOpT< MF >::BCMode BCMode
Definition AMReX_MLCellLinOp.H:28
void smooth(int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary, int niter) const final
Smooth.
Definition AMReX_MLCellLinOp.H:1213
RT normInf(int amrlev, MF const &mf, bool local) const override
Definition AMReX_MLCellLinOp.H:2041
Array< RT, 2 *BL_SPACEDIM > RealTuple
Definition AMReX_MLCellLinOp.H:164
virtual bool isCrossStencil() const
Definition AMReX_MLCellLinOp.H:60
void updateCorBC(int amrlev, const MF &crse_bcdata) const
Definition AMReX_MLCellLinOp.H:676
Array< BoundCond, 2 *BL_SPACEDIM > BCTuple
Definition AMReX_MLCellLinOp.H:165
void interpolation(int amrlev, int fmglev, MF &fine, const MF &crse) const override
Add interpolated coarse MG level data to fine MG level data.
Definition AMReX_MLCellLinOp.H:963
virtual void addInhomogNeumannFlux(int, const Array< MF *, AMREX_SPACEDIM > &, MF const &, bool) const
Definition AMReX_MLCellLinOp.H:126
RT xdoty(int amrlev, int mglev, const MF &x, const MF &y, bool local) const final
x dot y, used by the bottom solver
Definition AMReX_MLCellLinOp.H:1934
void prepareForSolve() override
Definition AMReX_MLCellLinOp.H:1617
void unapplyMetricTerm(int amrlev, int mglev, MF &rhs) const final
unapply metric terms if there are any
Definition AMReX_MLCellLinOp.H:1502
void apply(int amrlev, int mglev, MF &out, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr) const override
Apply the linear operator, out = L(in)
Definition AMReX_MLCellLinOp.H:1203
void applyMetricTerm(int amrlev, int mglev, MF &rhs) const final
apply metric terms if there are any
Definition AMReX_MLCellLinOp.H:1446
Vector< std::unique_ptr< BndryRegisterT< MF > > > m_crse_cor_br
Definition AMReX_MLCellLinOp.H:157
RT norm2Precond(Vector< MF const * > const &x) const final
Definition AMReX_MLCellLinOp.H:1964
void beginPrecondBC() override
Definition AMReX_MLCellLinOp.H:2183
Vector< Vector< std::unique_ptr< BndryCondLoc > > > m_bcondloc
Definition AMReX_MLCellLinOp.H:201
void setGaussSeidel(bool flag) noexcept
Definition AMReX_MLCellLinOp.H:58
void defineAuxData()
Definition AMReX_MLCellLinOp.H:375
Vector< std::unique_ptr< MF > > m_robin_bcval
Definition AMReX_MLCellLinOp.H:145
Vector< Vector< Array< MultiMask, 2 *AMREX_SPACEDIM > > > m_maskvals
Definition AMReX_MLCellLinOp.H:207
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info=LPInfo(), const Vector< FabFactory< FAB > const * > &a_factory={})
Definition AMReX_MLCellLinOp.H:362
~MLCellLinOpT() override=default
virtual bool isTensorOp() const
Definition AMReX_MLCellLinOp.H:61
void setInterpBndryHalfWidth(int w)
Definition AMReX_MLCellLinOp.H:147
void restriction(int, int, MF &crse, MF &fine) const override
Restriction onto coarse MG level.
Definition AMReX_MLCellLinOp.H:954
Vector< Vector< BndryRegisterT< MF > > > m_undrrelxr
Definition AMReX_MLCellLinOp.H:204
int m_interpbndry_halfwidth
Definition AMReX_MLCellLinOp.H:223
MLCellLinOpT()
Definition AMReX_MLCellLinOp.H:355
MLCellLinOpT(MLCellLinOpT< MF > &&)=delete
void interpAssign(int amrlev, int fmglev, MF &fine, MF &crse) const override
Overwrite fine MG level data with interpolated coarse data.
Definition AMReX_MLCellLinOp.H:1010
void fixSolvabilityByOffset(int amrlev, int mglev, MF &rhs, Vector< RT > const &offset) const override
Definition AMReX_MLCellLinOp.H:1597
Vector< std::unique_ptr< MLMGBndryT< MF > > > m_bndry_cor
Definition AMReX_MLCellLinOp.H:156
typename MLLinOpT< MF >::Location Location
Definition AMReX_MLCellLinOp.H:30
void prepareForFluxes(int amrlev, const MF *crse_bcdata=nullptr) override
Definition AMReX_MLCellLinOp.H:1248
void solutionResidual(int amrlev, MF &resid, MF &x, const MF &b, const MF *crse_bcdata=nullptr) override
Compute residual for solution.
Definition AMReX_MLCellLinOp.H:1230
Vector< Vector< RT > > m_volinv
Definition AMReX_MLCellLinOp.H:221
bool m_has_metric_term
Definition AMReX_MLCellLinOp.H:151
bool m_use_gauss_seidel
Definition AMReX_MLCellLinOp.H:213
Vector< std::unique_ptr< MLMGBndryT< MF > > > m_bndry_sol_zero
Definition AMReX_MLCellLinOp.H:159
Vector< std::unique_ptr< BndryRegisterT< MF > > > m_crse_sol_br
Definition AMReX_MLCellLinOp.H:154
void interpolationAmr(int famrlev, MF &fine, const MF &crse, IntVect const &nghost) const override
Interpolation between AMR levels.
Definition AMReX_MLCellLinOp.H:1103
void compFlux(int amrlev, const Array< MF *, AMREX_SPACEDIM > &fluxes, MF &sol, Location loc) const override
Compute fluxes.
Definition AMReX_MLCellLinOp.H:1351
void computeVolInv() const
Definition AMReX_MLCellLinOp.H:1979
Vector< std::unique_ptr< iMultiFab > > m_norm_fine_mask
Definition AMReX_MLCellLinOp.H:209
virtual void FFlux(int amrlev, const MFIter &mfi, const Array< FAB *, AMREX_SPACEDIM > &flux, const FAB &sol, Location loc, int face_only=0) const =0
typename MLLinOpT< MF >::StateMode StateMode
Definition AMReX_MLCellLinOp.H:29
void setLevelBC(int amrlev, const MF *levelbcdata, const MF *robinbc_a=nullptr, const MF *robinbc_b=nullptr, const MF *robinbc_f=nullptr) final
Set boundary conditions for given level. For cell-centered solves only.
Definition AMReX_MLCellLinOp.H:520
bool needsUpdate() const override
Does it need update if it's reused?
Definition AMReX_MLCellLinOp.H:53
Definition AMReX_MLLinOp.H:98
const MF * m_coarse_data_for_bc
Definition AMReX_MLLinOp.H:638
Vector< Vector< std::unique_ptr< FabFactory< FAB > > > > m_factory
Definition AMReX_MLLinOp.H:611
int NAMRLevels() const noexcept
Return the number of AMR levels.
Definition AMReX_MLLinOp.H:564
IntVect m_ixtype
Definition AMReX_MLLinOp.H:599
Array< Real, AMREX_SPACEDIM > m_domain_bloc_hi
Definition AMReX_MLLinOp.H:632
RealVect m_coarse_bc_loc
Definition AMReX_MLLinOp.H:637
virtual bool needsUpdate() const
Does it need update if it's reused?
Definition AMReX_MLLinOp.H:266
FabFactory< FAB > const * Factory(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:649
Vector< Array< BCType, AMREX_SPACEDIM > > m_hibc_orig
Definition AMReX_MLLinOp.H:577
bool hasRobinBC() const noexcept
Definition AMReX_MLLinOp.H:1277
Array< Real, AMREX_SPACEDIM > m_domain_bloc_lo
Definition AMReX_MLLinOp.H:631
Vector< Vector< BoxArray > > m_grids
Definition AMReX_MLLinOp.H:609
virtual int getNComp() const
Return number of components.
Definition AMReX_MLLinOp.H:261
Vector< int > m_amr_ref_ratio
Definition AMReX_MLLinOp.H:594
Vector< int > m_num_mg_levels
Definition AMReX_MLLinOp.H:596
Vector< Vector< DistributionMapping > > m_dmap
Definition AMReX_MLLinOp.H:610
IntVect m_coarse_data_crse_ratio
Definition AMReX_MLLinOp.H:636
const Vector< int > & AMRRefRatio() const noexcept
Return AMR refinement ratios.
Definition AMReX_MLLinOp.H:644
Vector< Array< BCType, AMREX_SPACEDIM > > m_lobc_orig
Definition AMReX_MLLinOp.H:576
bool needsCoarseDataForBC() const noexcept
Needs coarse data for bc?
Definition AMReX_MLLinOp.H:177
virtual void update()
Update for reuse.
Definition AMReX_MLLinOp.H:268
bool m_precond_mode
Definition AMReX_MLLinOp.H:641
bool hasHiddenDimension() const noexcept
Definition AMReX_MLLinOp.H:707
const Geometry & Geom(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:569
Vector< Array< BCType, AMREX_SPACEDIM > > m_lobc
Definition AMReX_MLLinOp.H:572
int hiddenDirection() const noexcept
Definition AMReX_MLLinOp.H:708
Vector< Vector< Geometry > > m_geom
first Vector is for amr level and second is mg level
Definition AMReX_MLLinOp.H:608
Box compactify(Box const &b) const noexcept
Definition AMReX_MLLinOp.H:1284
int maxorder
Definition AMReX_MLLinOp.H:589
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info, const Vector< FabFactory< FAB > const * > &a_factory, bool eb_limit_coarsening=true)
Definition AMReX_MLLinOp.H:769
Vector< IntVect > mg_coarsen_ratio_vec
Definition AMReX_MLLinOp.H:605
LPInfo info
Definition AMReX_MLLinOp.H:585
Vector< Array< BCType, AMREX_SPACEDIM > > m_hibc
Definition AMReX_MLLinOp.H:573
int NMGLevels(int amrlev) const noexcept
Return the number of MG levels at given AMR level.
Definition AMReX_MLLinOp.H:567
LinOpBCType m_coarse_fine_bc_type
Definition AMReX_MLLinOp.H:635
int m_num_amr_levels
Definition AMReX_MLLinOp.H:593
Definition AMReX_MLMGBndry.H:12
static void setBoxBC(RealTuple &bloc, BCTuple &bctag, const Box &bx, const Box &domain, const Array< LinOpBCType, AMREX_SPACEDIM > &lo, const Array< LinOpBCType, AMREX_SPACEDIM > &hi, const Real *dx, IntVect const &ratio, const RealVect &interior_bloc, const Array< Real, AMREX_SPACEDIM > &domain_bloc_lo, const Array< Real, AMREX_SPACEDIM > &domain_bloc_hi, const GpuArray< int, AMREX_SPACEDIM > &is_periodic, LinOpBCType a_crse_fine_bc_type)
Definition AMReX_MLMGBndry.H:107
Definition AMReX_Mask.H:28
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
Real sum(int comp=0, bool local=false) const
Returns the sum of component "comp" over the MultiFab – no ghost cells are included.
Definition AMReX_MultiFab.cpp:1213
An Iterator over the Orientation of Faces of a Box.
Definition AMReX_Orientation.H:135
Encapsulation of the Orientation of the Faces of a Box.
Definition AMReX_Orientation.H:29
@ low
Definition AMReX_Orientation.H:34
@ high
Definition AMReX_Orientation.H:34
Definition AMReX_PODVector.H:262
void reserve(size_type a_capacity)
Definition AMReX_PODVector.H:663
iterator begin() noexcept
Definition AMReX_PODVector.H:617
iterator end() noexcept
Definition AMReX_PODVector.H:621
void push_back(const T &a_value)
Definition AMReX_PODVector.H:572
A Real vector in SpaceDim-dimensional space.
Definition AMReX_RealVect.H:32
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:27
void copyAsync(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
A host-to-device copy routine. Note this is just a wrapper around memcpy, so it assumes contiguous st...
Definition AMReX_GpuContainers.H:233
static constexpr HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:98
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:237
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:86
bool notInLaunchRegion() noexcept
Definition AMReX_GpuControl.H:87
void Sum(T &v, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:204
void Max(KeyValuePair< K, V > &vi, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:126
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlmg_lin_cc_interp_r2(Box const &bx, Array4< T > const &ff, Array4< T const > const &cc, int nc) noexcept
Definition AMReX_MLMG_2D_K.H:9
Definition AMReX_Amr.cpp:49
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:191
void average_down(const MultiFab &S_fine, MultiFab &S_crse, const Geometry &fgeom, const Geometry &cgeom, int scomp, int ncomp, int rr)
Definition AMReX_MultiFabUtil.cpp:314
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlmg_lin_cc_interp_r2(Box const &bx, Array4< T > const &ff, Array4< T const > const &cc, int nc) noexcept
Definition AMReX_MLMG_1D_K.H:9
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Returns a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition AMReX_Box.H:1399
ReduceData< Ts... >::Type ParReduce(TypeList< Ops... > operation_list, TypeList< Ts... > type_list, FabArray< FAB > const &fa, IntVect const &nghost, F &&f)
Parallel reduce for MultiFab/FabArray.
Definition AMReX_ParReduce.H:47
void EB_set_covered(MultiFab &mf, Real val)
Definition AMReX_EBMultiFabUtil.cpp:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlmg_lin_cc_interp_r4(Box const &bx, Array4< T > const &ff, Array4< T const > const &cc, int nc) noexcept
Definition AMReX_MLMG_1D_K.H:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > adjCellLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Returns the cell centered BoxND of length len adjacent to b on the low end along the coordinate direc...
Definition AMReX_Box.H:1591
bool isMFIterSafe(const FabArrayBase &x, const FabArrayBase &y)
Definition AMReX_MFIter.H:219
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:1304
void EB_average_down(const MultiFab &S_fine, MultiFab &S_crse, const MultiFab &vol_fine, const MultiFab &vfrac_fine, int scomp, int ncomp, const IntVect &ratio)
Definition AMReX_EBMultiFabUtil.cpp:336
IntVectND< AMREX_SPACEDIM > IntVect
Definition AMReX_BaseFwd.H:30
AMREX_GPU_HOST_DEVICE constexpr GpuTupleElement< I, GpuTuple< Ts... > >::type & get(GpuTuple< Ts... > &tup) noexcept
Definition AMReX_Tuple.H:179
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1211
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > adjCell(const BoxND< dim > &b, Orientation face, int len=1) noexcept
Similar to adjCellLo and adjCellHi; operates on given face.
Definition AMReX_Box.H:1634
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:127
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE T norm(const GpuComplex< T > &a_z) noexcept
Return the norm (magnitude squared) of a complex number.
Definition AMReX_GpuComplex.H:344
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mllinop_apply_bc_x(int side, Box const &box, int blen, Array4< T > const &phi, Array4< int const > const &mask, BoundCond bct, T bcl, Array4< T const > const &bcval, int maxorder, T dxinv, int inhomog, int icomp) noexcept
Definition AMReX_MLLinOp_K.H:14
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mllinop_comp_interp_coef0_x(int side, Box const &box, int blen, Array4< T > const &f, Array4< int const > const &mask, BoundCond bct, T bcl, int maxorder, T dxinv, int icomp) noexcept
Definition AMReX_MLLinOp_K.H:329
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mllinop_apply_bc_z(int side, Box const &box, int blen, Array4< T > const &phi, Array4< int const > const &mask, BoundCond bct, T bcl, Array4< T const > const &bcval, int maxorder, T dzinv, int inhomog, int icomp) noexcept
Definition AMReX_MLLinOp_K.H:224
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
Arena * The_Async_Arena()
Definition AMReX_Arena.cpp:626
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mllinop_comp_interp_coef0_y(int side, Box const &box, int blen, Array4< T > const &f, Array4< int const > const &mask, BoundCond bct, T bcl, int maxorder, T dyinv, int icomp) noexcept
Definition AMReX_MLLinOp_K.H:410
FAB::value_type Dot(FabArray< FAB > const &x, int xcomp, FabArray< FAB > const &y, int ycomp, int ncomp, IntVect const &nghost, bool local=false)
Compute dot products of two FabArrays.
Definition AMReX_FabArrayUtility.H:1621
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > adjCellHi(const BoxND< dim > &b, int dir, int len=1) noexcept
Similar to adjCellLo but builds an adjacent BoxND on the high end.
Definition AMReX_Box.H:1612
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > shift(const BoxND< dim > &b, int dir, int nzones) noexcept
Return a BoxND with indices shifted by nzones in dir direction.
Definition AMReX_Box.H:1372
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mllinop_comp_interp_coef0_z(int side, Box const &box, int blen, Array4< T > const &f, Array4< int const > const &mask, BoundCond bct, T bcl, int maxorder, T dzinv, int icomp) noexcept
Definition AMReX_MLLinOp_K.H:491
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mllinop_apply_bc_y(int side, Box const &box, int blen, Array4< T > const &phi, Array4< int const > const &mask, BoundCond bct, T bcl, Array4< T const > const &bcval, int maxorder, T dyinv, int inhomog, int icomp) noexcept
Definition AMReX_MLLinOp_K.H:119
iMultiFab makeFineMask(const BoxArray &cba, const DistributionMapping &cdm, const BoxArray &fba, const IntVect &ratio, int crse_value, int fine_value)
Definition AMReX_MultiFabUtil.cpp:607
std::array< T, N > Array
Definition AMReX_Array.H:24
Definition AMReX_Array4.H:61
Definition AMReX_Dim3.H:12
int x
Definition AMReX_Dim3.H:12
int z
Definition AMReX_Dim3.H:12
int y
Definition AMReX_Dim3.H:12
Definition AMReX_Array.H:34
Definition AMReX_MLLinOp.H:35
bool has_metric_term
Definition AMReX_MLLinOp.H:43
StateMode
Definition AMReX_MLLinOp.H:86
BCMode
Definition AMReX_MLLinOp.H:85
Location
Definition AMReX_MLLinOp.H:87
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_MFIter.H:20
MFItInfo & SetDynamic(bool f) noexcept
Definition AMReX_MFIter.H:34
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept
Definition AMReX_MFIter.H:29
Definition AMReX_MLCellLinOp.H:140
BoundCond type
Definition AMReX_MLCellLinOp.H:141
RT location
Definition AMReX_MLCellLinOp.H:142
Definition AMReX_MLCellLinOp.H:227
T bcloc_hi
Definition AMReX_MLCellLinOp.H:234
Box bx
Definition AMReX_MLCellLinOp.H:235
Array4< T const > bcval_hi
Definition AMReX_MLCellLinOp.H:230
Array4< int const > mask_hi
Definition AMReX_MLCellLinOp.H:232
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box const & box() const noexcept
Definition AMReX_MLCellLinOp.H:243
int comp
Definition AMReX_MLCellLinOp.H:239
int dir
Definition AMReX_MLCellLinOp.H:240
Array4< T const > bcval_lo
Definition AMReX_MLCellLinOp.H:229
int blen
Definition AMReX_MLCellLinOp.H:238
Array4< T > fab
Definition AMReX_MLCellLinOp.H:228
BoundCond bctype_lo
Definition AMReX_MLCellLinOp.H:236
T bcloc_lo
Definition AMReX_MLCellLinOp.H:233
Array4< int const > mask_lo
Definition AMReX_MLCellLinOp.H:231
BoundCond bctype_hi
Definition AMReX_MLCellLinOp.H:237
Definition AMReX_MLCellLinOp.H:247
int blen
Definition AMReX_MLCellLinOp.H:257
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box const & box() const noexcept
Definition AMReX_MLCellLinOp.H:262
BoundCond bcthi
Definition AMReX_MLCellLinOp.H:256
int comp
Definition AMReX_MLCellLinOp.H:258
Array4< T > flo
Definition AMReX_MLCellLinOp.H:248
Array4< int const > mlo
Definition AMReX_MLCellLinOp.H:250
int dir
Definition AMReX_MLCellLinOp.H:259
Array4< T > fhi
Definition AMReX_MLCellLinOp.H:249
T bcllo
Definition AMReX_MLCellLinOp.H:252
Box bx
Definition AMReX_MLCellLinOp.H:254
BoundCond bctlo
Definition AMReX_MLCellLinOp.H:255
Array4< int const > mhi
Definition AMReX_MLCellLinOp.H:251
T bclhi
Definition AMReX_MLCellLinOp.H:253
Struct for holding types.
Definition AMReX_TypeList.H:12