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