Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_MLABecLaplacian.H
Go to the documentation of this file.
1#ifndef AMREX_ML_ABECLAPLACIAN_H_
2#define AMREX_ML_ABECLAPLACIAN_H_
3#include <AMReX_Config.H>
4
6#include <AMReX_MLABecLap_K.H>
7
8namespace amrex {
9
12template <typename MF>
14 : public MLCellABecLapT<MF>
15{
16public:
17
18 using FAB = typename MF::fab_type;
19 using RT = typename MF::value_type;
20
23
24 MLABecLaplacianT () = default;
26 const Vector<BoxArray>& a_grids,
27 const Vector<DistributionMapping>& a_dmap,
28 const LPInfo& a_info = LPInfo(),
29 const Vector<FabFactory<FAB> const*>& a_factory = {},
30 int a_ncomp = 1);
31
33 const Vector<BoxArray>& a_grids,
34 const Vector<DistributionMapping>& a_dmap,
35 const Vector<iMultiFab const*>& a_overset_mask, // 1: unknown, 0: known
36 const LPInfo& a_info = LPInfo(),
37 const Vector<FabFactory<FAB> const*>& a_factory = {},
38 int a_ncomp = 1);
39
40 ~MLABecLaplacianT () override;
41
46
47 void define (const Vector<Geometry>& a_geom,
48 const Vector<BoxArray>& a_grids,
49 const Vector<DistributionMapping>& a_dmap,
50 const LPInfo& a_info = LPInfo(),
51 const Vector<FabFactory<FAB> const*>& a_factory = {},
52 int a_ncomp = 1);
53
54 void define (const Vector<Geometry>& a_geom,
55 const Vector<BoxArray>& a_grids,
56 const Vector<DistributionMapping>& a_dmap,
57 const Vector<iMultiFab const*>& a_overset_mask,
58 const LPInfo& a_info = LPInfo(),
59 const Vector<FabFactory<FAB> const*>& a_factory = {},
60 int a_ncomp = 1);
61
67 template <typename T1, typename T2>
68 requires (std::is_convertible_v<T1, typename MF::value_type> && std::is_convertible_v<T2, typename MF::value_type>)
69 void setScalars (T1 a, T2 b) noexcept;
70
80 template <FabArrayType AMF>
81 requires (std::is_convertible_v<typename AMF::value_type, typename MF::value_type>)
82 void setACoeffs (int amrlev, const AMF& alpha);
83
93 template <typename T>
94 requires (std::is_convertible_v<T, typename MF::value_type>)
95 void setACoeffs (int amrlev, T alpha);
96
106 template <FabArrayType AMF>
107 requires (std::is_convertible_v<typename AMF::value_type, typename MF::value_type>)
108 void setBCoeffs (int amrlev, const Array<AMF const*,AMREX_SPACEDIM>& beta);
109
119 template <typename T>
120 requires (std::is_convertible_v<T, typename MF::value_type>)
121 void setBCoeffs (int amrlev, T beta);
122
132 template <typename T>
133 requires (std::is_convertible_v<T, typename MF::value_type>)
134 void setBCoeffs (int amrlev, Vector<T> const& beta);
135
136 [[nodiscard]] int getNComp () const override { return m_ncomp; }
137
138 [[nodiscard]] bool needsUpdate () const override {
139 return (m_needs_update || MLCellABecLapT<MF>::needsUpdate());
140 }
141 void update () override;
142
143 void prepareForSolve () override;
144 [[nodiscard]] bool isSingular (int amrlev) const override { return m_is_singular[amrlev]; }
145 [[nodiscard]] bool isBottomSingular () const override { return m_is_singular[0]; }
146 void Fapply (int amrlev, int mglev, MF& out, const MF& in) const override;
147 void Fsmooth (int amrlev, int mglev, MF& sol, const MF& rhs, int redblack) const override;
148 void FFlux (int amrlev, const MFIter& mfi,
149 const Array<FAB*,AMREX_SPACEDIM>& flux,
150 const FAB& sol, Location /* loc */,
151 int face_only=0) const override;
152
153 void normalize (int amrlev, int mglev, MF& mf) const override;
154
155 [[nodiscard]] RT getAScalar () const final { return m_a_scalar; }
156 [[nodiscard]] RT getBScalar () const final { return m_b_scalar; }
157 [[nodiscard]] MF const* getACoeffs (int amrlev, int mglev) const final
158 { return &(m_a_coeffs[amrlev][mglev]); }
159 [[nodiscard]] Array<MF const*,AMREX_SPACEDIM> getBCoeffs (int amrlev, int mglev) const final
160 { return amrex::GetArrOfConstPtrs(m_b_coeffs[amrlev][mglev]); }
161
162 [[nodiscard]] std::unique_ptr<MLLinOpT<MF>> makeNLinOp (int /*grid_size*/) const final;
163
164 [[nodiscard]] bool supportNSolve () const override;
165
166 void copyNSolveSolution (MF& dst, MF const& src) const final;
167
168 void averageDownCoeffsSameAmrLevel (int amrlev, Vector<MF>& a,
170 void averageDownCoeffs ();
171 void averageDownCoeffsToCoarseAmrLevel (int flev);
172
174
176
177 static void FFlux (Box const& box, Real const* dxinv, RT bscalar,
179 Array<FAB*,AMREX_SPACEDIM> const& flux,
180 FAB const& sol, int face_only, int ncomp);
181
182 RT m_a_scalar = std::numeric_limits<RT>::quiet_NaN();
183 RT m_b_scalar = std::numeric_limits<RT>::quiet_NaN();
186
187 bool m_scalars_set = false;
188 bool m_acoef_set = false;
189
190protected:
191
193
194 [[nodiscard]] bool supportRobinBC () const noexcept override { return true; }
195
196private:
197
198 bool m_needs_update = true;
199
200 int m_ncomp = 1;
201
202 void define_ab_coeffs ();
203
204 void update_singular_flags ();
205};
206
207template <typename MF>
209 const Vector<BoxArray>& a_grids,
210 const Vector<DistributionMapping>& a_dmap,
211 const LPInfo& a_info,
212 const Vector<FabFactory<FAB> const*>& a_factory,
213 int a_ncomp)
214{
215 define(a_geom, a_grids, a_dmap, a_info, a_factory, a_ncomp);
216}
217
218template <typename MF>
220 const Vector<BoxArray>& a_grids,
221 const Vector<DistributionMapping>& a_dmap,
222 const Vector<iMultiFab const*>& a_overset_mask,
223 const LPInfo& a_info,
224 const Vector<FabFactory<FAB> const*>& a_factory,
225 int a_ncomp)
226{
227 define(a_geom, a_grids, a_dmap, a_overset_mask, a_info, a_factory, a_ncomp);
228}
229
230template <typename MF> MLABecLaplacianT<MF>::~MLABecLaplacianT () = default;
231
232template <typename MF>
233void
235 const Vector<BoxArray>& a_grids,
236 const Vector<DistributionMapping>& a_dmap,
237 const LPInfo& a_info,
238 const Vector<FabFactory<FAB> const*>& a_factory,
239 int a_ncomp)
240{
241 BL_PROFILE("MLABecLaplacian::define()");
242 this->m_ncomp = a_ncomp;
243 MLCellABecLapT<MF>::define(a_geom, a_grids, a_dmap, a_info, a_factory);
244 define_ab_coeffs();
245}
246
247template <typename MF>
248void
250 const Vector<BoxArray>& a_grids,
251 const Vector<DistributionMapping>& a_dmap,
252 const Vector<iMultiFab const*>& a_overset_mask,
253 const LPInfo& a_info,
254 const Vector<FabFactory<FAB> const*>& a_factory,
255 int a_ncomp)
256{
257 BL_PROFILE("MLABecLaplacian::define(overset)");
258 this->m_ncomp = a_ncomp;
259 MLCellABecLapT<MF>::define(a_geom, a_grids, a_dmap, a_overset_mask, a_info, a_factory);
260 define_ab_coeffs();
261}
262
263template <typename MF>
264void
266{
267 m_a_coeffs.resize(this->m_num_amr_levels);
268 m_b_coeffs.resize(this->m_num_amr_levels);
269 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
270 {
271 m_a_coeffs[amrlev].resize(this->m_num_mg_levels[amrlev]);
272 m_b_coeffs[amrlev].resize(this->m_num_mg_levels[amrlev]);
273 for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
274 {
275 m_a_coeffs[amrlev][mglev].define
276 (this->m_grids[amrlev][mglev], this->m_dmap[amrlev][mglev],
277 1, 0, MFInfo(), *(this->m_factory[amrlev][mglev]));
278 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
279 {
280 const BoxArray& ba = amrex::convert(this->m_grids[amrlev][mglev],
282 m_b_coeffs[amrlev][mglev][idim].define
283 (ba, this->m_dmap[amrlev][mglev], m_ncomp, 0, MFInfo(),
284 *(this->m_factory[amrlev][mglev]));
285 }
286 }
287 }
288}
289
290template <typename MF>
291template <typename T1, typename T2>
292requires (std::is_convertible_v<T1, typename MF::value_type> && std::is_convertible_v<T2, typename MF::value_type>)
293void
295{
296 m_a_scalar = RT(a);
297 m_b_scalar = RT(b);
298 if (m_a_scalar == RT(0.0)) {
299 for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev) {
300 m_a_coeffs[amrlev][0].setVal(RT(0.0));
301 }
302 m_acoef_set = true;
303 }
304 m_scalars_set = true;
305}
306
307template <typename MF>
308template <FabArrayType AMF>
309requires (std::is_convertible_v<typename AMF::value_type, typename MF::value_type>)
310void
311MLABecLaplacianT<MF>::setACoeffs (int amrlev, const AMF& alpha)
312{
313 AMREX_ASSERT_WITH_MESSAGE(alpha.nComp() == 1,
314 "MLABecLaplacian::setACoeffs: alpha is supposed to be single component.");
315 m_a_coeffs[amrlev][0].LocalCopy(alpha, 0, 0, 1, IntVect(0));
316 m_needs_update = true;
317 m_acoef_set = true;
318}
319
320template <typename MF>
321template <typename T>
322requires (std::is_convertible_v<T, typename MF::value_type>)
323void
325{
326 m_a_coeffs[amrlev][0].setVal(RT(alpha));
327 m_needs_update = true;
328 m_acoef_set = true;
329}
330
331
332template <typename MF>
333template <FabArrayType AMF>
334requires (std::is_convertible_v<typename AMF::value_type, typename MF::value_type>)
335void
338{
339 const int ncomp = this->getNComp();
340 AMREX_ASSERT(beta[0]->nComp() == 1 || beta[0]->nComp() == ncomp);
341 if (beta[0]->nComp() == ncomp) {
342 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
343 for (int icomp = 0; icomp < ncomp; ++icomp) {
344 m_b_coeffs[amrlev][0][idim].LocalCopy(*beta[idim], icomp, icomp, 1, IntVect(0));
345 }
346 }
347 } else {
348 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
349 for (int icomp = 0; icomp < ncomp; ++icomp) {
350 m_b_coeffs[amrlev][0][idim].LocalCopy(*beta[idim], 0, icomp, 1, IntVect(0));
351 }
352 }
353 }
354 m_needs_update = true;
355}
356
357template <typename MF>
358template <typename T>
359requires (std::is_convertible_v<T, typename MF::value_type>)
360void
362{
363 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
364 m_b_coeffs[amrlev][0][idim].setVal(RT(beta));
365 }
366 m_needs_update = true;
367}
368
369template <typename MF>
370template <typename T>
371requires (std::is_convertible_v<T, typename MF::value_type>)
372void
374{
375 const int ncomp = this->getNComp();
376 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
377 for (int icomp = 0; icomp < ncomp; ++icomp) {
378 m_b_coeffs[amrlev][0][idim].setVal(RT(beta[icomp]), icomp, 1, 0);
379 }
380 }
381 m_needs_update = true;
382}
383
384template <typename MF>
385void
387{
390 }
391
392#if (AMREX_SPACEDIM != 3)
393 applyMetricTermsCoeffs();
394#endif
395
396 applyRobinBCTermsCoeffs();
397
398 averageDownCoeffs();
399
400 update_singular_flags();
401
402 m_needs_update = false;
403}
404
405template <typename MF>
406void
408{
409 BL_PROFILE("MLABecLaplacian::prepareForSolve()");
410
412
413#if (AMREX_SPACEDIM != 3)
414 applyMetricTermsCoeffs();
415#endif
416
417 applyRobinBCTermsCoeffs();
418
419 averageDownCoeffs();
420
421 update_singular_flags();
422
423 m_needs_update = false;
424}
425
426template <typename MF>
427void
429{
430#if (AMREX_SPACEDIM != 3)
431 for (int alev = 0; alev < this->m_num_amr_levels; ++alev)
432 {
433 const int mglev = 0;
434 this->applyMetricTerm(alev, mglev, m_a_coeffs[alev][mglev]);
435 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
436 {
437 this->applyMetricTerm(alev, mglev, m_b_coeffs[alev][mglev][idim]);
438 }
439 }
440#endif
441}
442
443//
444// Suppose we are solving `alpha u - del (beta grad u) = rhs` (Scalar
445// coefficients can be easily added back in the end) and there is Robin BC
446// `a u + b du/dn = f` at the upper end of the x-direction. The 1D
447// discretization at the last cell i is
448//
449// alpha u_i + (beta_{i-1/2} (du/dx)_{i-1/2} - beta_{i+1/2} (du/dx)_{i+1/2}) / h = rhs_i
450//
451// where h is the cell size. At `i+1/2` (i.e., the boundary), we have
452//
453// a (u_i + u_{i+1})/2 + b (u_{i+1}-u_i)/h = f,
454//
455// according to the Robin BC. This gives
456//
457// u_{i+1} = A + B u_i,
458//
459// where `A = f/(b/h + a/2)` and `B = (b/h - a/2) / (b/h + a/2). We then
460// use `u_i` and `u_{i+1}` to compute `(du/dx)_{i+1/2}`. The discretization
461// at cell i then becomes
462//
463// \tilde{alpha}_i u_i + (beta_{i-1/2} (du/dx)_{i-1/2} - 0) / h = \tilde{rhs}_i
464//
465// This is equivalent to having homogeneous Neumann BC with modified alpha and rhs.
466//
467// \tilde{alpha}_i = alpha_i + (1-B) beta_{i+1/2} / h^2
468// \tilde{rhs}_i = rhs_i + A beta_{i+1/2} / h^2
469//
471namespace detail {
472template <typename LP>
473void applyRobinBCTermsCoeffs (LP& linop)
474{
475 using RT = typename LP::RT;
476
477 const int ncomp = linop.getNComp();
478 bool reset_alpha = false;
479 if (linop.m_a_scalar == RT(0.0)) {
480 linop.m_a_scalar = RT(1.0);
481 reset_alpha = true;
482 }
483 const RT bovera = linop.m_b_scalar/linop.m_a_scalar;
484
485 if (!reset_alpha) {
486 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(linop.m_scalars_set && linop.m_acoef_set,
487 "To reuse solver With Robin BC, one must re-call setScalars (and setACoeffs if the scalar is not zero)");
488 }
489
490 linop.m_scalars_set = false;
491 linop.m_acoef_set = false;
492
493 for (int amrlev = 0; amrlev < linop.NAMRLevels(); ++amrlev) {
494 const int mglev = 0;
495 const Box& domain = linop.Geom(amrlev,mglev).Domain();
496 const RT dxi = static_cast<RT>(linop.Geom(amrlev,mglev).InvCellSize(0));
497 const RT dyi = static_cast<RT>((AMREX_SPACEDIM >= 2) ? linop.Geom(amrlev,mglev).InvCellSize(1) : Real(1.0));
498 const RT dzi = static_cast<RT>((AMREX_SPACEDIM == 3) ? linop.Geom(amrlev,mglev).InvCellSize(2) : Real(1.0));
499
500 if (reset_alpha) {
501 linop.m_a_coeffs[amrlev][mglev].setVal(RT(0.0));
502 }
503
504 MFItInfo mfi_info;
505 if (Gpu::notInLaunchRegion()) { mfi_info.SetDynamic(true); }
506
507#ifdef AMREX_USE_OMP
508#pragma omp parallel if (Gpu::notInLaunchRegion())
509#endif
510 for (MFIter mfi(linop.m_a_coeffs[amrlev][mglev], mfi_info); mfi.isValid(); ++mfi)
511 {
512 const Box& vbx = mfi.validbox();
513 auto const& afab = linop.m_a_coeffs[amrlev][mglev].array(mfi);
514 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
515 auto const& bfab = linop.m_b_coeffs[amrlev][mglev][idim].const_array(mfi);
516 const Box& blo = amrex::adjCellLo(vbx,idim);
517 const Box& bhi = amrex::adjCellHi(vbx,idim);
518 bool outside_domain_lo = !(domain.contains(blo));
519 bool outside_domain_hi = !(domain.contains(bhi));
520 if ((!outside_domain_lo) && (!outside_domain_hi)) { continue; }
521 for (int icomp = 0; icomp < ncomp; ++icomp) {
522 auto const& rbc = (*(linop.m_robin_bcval[amrlev]))[mfi].const_array(icomp*3);
523 if (linop.m_lobc_orig[icomp][idim] == LinOpBCType::Robin && outside_domain_lo)
524 {
525 if (idim == 0) {
526 RT fac = bovera*dxi*dxi;
527 AMREX_HOST_DEVICE_FOR_3D(blo, i, j, k,
528 {
529 RT B = (rbc(i,j,k,1)*dxi - rbc(i,j,k,0)*RT(0.5))
530 / (rbc(i,j,k,1)*dxi + rbc(i,j,k,0)*RT(0.5));
531 afab(i+1,j,k,icomp) += fac*bfab(i+1,j,k,icomp)*(RT(1.0)-B);
532 });
533 } else if (idim == 1) {
534 RT fac = bovera*dyi*dyi;
535 AMREX_HOST_DEVICE_FOR_3D(blo, i, j, k,
536 {
537 RT B = (rbc(i,j,k,1)*dyi - rbc(i,j,k,0)*RT(0.5))
538 / (rbc(i,j,k,1)*dyi + rbc(i,j,k,0)*RT(0.5));
539 afab(i,j+1,k,icomp) += fac*bfab(i,j+1,k,icomp)*(RT(1.0)-B);
540 });
541 } else {
542 RT fac = bovera*dzi*dzi;
543 AMREX_HOST_DEVICE_FOR_3D(blo, i, j, k,
544 {
545 RT B = (rbc(i,j,k,1)*dzi - rbc(i,j,k,0)*RT(0.5))
546 / (rbc(i,j,k,1)*dzi + rbc(i,j,k,0)*RT(0.5));
547 afab(i,j,k+1,icomp) += fac*bfab(i,j,k+1,icomp)*(RT(1.0)-B);
548 });
549 }
550 }
551 if (linop.m_hibc_orig[icomp][idim] == LinOpBCType::Robin && outside_domain_hi)
552 {
553 if (idim == 0) {
554 RT fac = bovera*dxi*dxi;
555 AMREX_HOST_DEVICE_FOR_3D(bhi, i, j, k,
556 {
557 RT B = (rbc(i,j,k,1)*dxi - rbc(i,j,k,0)*RT(0.5))
558 / (rbc(i,j,k,1)*dxi + rbc(i,j,k,0)*RT(0.5));
559 afab(i-1,j,k,icomp) += fac*bfab(i,j,k,icomp)*(RT(1.0)-B);
560 });
561 } else if (idim == 1) {
562 RT fac = bovera*dyi*dyi;
563 AMREX_HOST_DEVICE_FOR_3D(bhi, i, j, k,
564 {
565 RT B = (rbc(i,j,k,1)*dyi - rbc(i,j,k,0)*RT(0.5))
566 / (rbc(i,j,k,1)*dyi + rbc(i,j,k,0)*RT(0.5));
567 afab(i,j-1,k,icomp) += fac*bfab(i,j,k,icomp)*(RT(1.0)-B);
568 });
569 } else {
570 RT fac = bovera*dzi*dzi;
571 AMREX_HOST_DEVICE_FOR_3D(bhi, i, j, k,
572 {
573 RT B = (rbc(i,j,k,1)*dzi - rbc(i,j,k,0)*RT(0.5))
574 / (rbc(i,j,k,1)*dzi + rbc(i,j,k,0)*RT(0.5));
575 afab(i,j,k-1,icomp) += fac*bfab(i,j,k,icomp)*(RT(1.0)-B);
576 });
577 }
578 }
579 }
580 }
581 }
582 }
583}
584} // namespace detail
586
587template <typename MF>
588void
590{
591 if (this->hasRobinBC()) {
592 detail::applyRobinBCTermsCoeffs(*this);
593 }
594}
595
596template <typename MF>
597void
599{
600 BL_PROFILE("MLABecLaplacian::averageDownCoeffs()");
601
602 for (int amrlev = this->m_num_amr_levels-1; amrlev > 0; --amrlev)
603 {
604 auto& fine_a_coeffs = m_a_coeffs[amrlev];
605 auto& fine_b_coeffs = m_b_coeffs[amrlev];
606
607 averageDownCoeffsSameAmrLevel(amrlev, fine_a_coeffs, fine_b_coeffs);
608 averageDownCoeffsToCoarseAmrLevel(amrlev);
609 }
610
611 averageDownCoeffsSameAmrLevel(0, m_a_coeffs[0], m_b_coeffs[0]);
612}
613
614template <typename MF>
615void
618{
619 int nmglevs = a.size();
620 for (int mglev = 1; mglev < nmglevs; ++mglev)
621 {
622 IntVect ratio = (amrlev > 0) ? IntVect(this->mg_coarsen_ratio) : this->mg_coarsen_ratio_vec[mglev-1];
623
624 if (m_a_scalar == 0.0)
625 {
626 a[mglev].setVal(RT(0.0));
627 }
628 else
629 {
630 amrex::average_down(a[mglev-1], a[mglev], 0, 1, ratio);
631 }
632
633 Vector<const MF*> fine {AMREX_D_DECL(&(b[mglev-1][0]),
634 &(b[mglev-1][1]),
635 &(b[mglev-1][2]))};
636 Vector<MF*> crse {AMREX_D_DECL(&(b[mglev][0]),
637 &(b[mglev][1]),
638 &(b[mglev][2]))};
639
641 }
642
643 for (int mglev = 1; mglev < nmglevs; ++mglev)
644 {
645 if (this->m_overset_mask[amrlev][mglev]) {
646 const RT fac = static_cast<RT>(1 << mglev); // 2**mglev
647 const RT osfac = RT(2.0)*fac/(fac+RT(1.0));
648 const int ncomp = this->getNComp();
649#ifdef AMREX_USE_OMP
650#pragma omp parallel if (Gpu::notInLaunchRegion())
651#endif
652 for (MFIter mfi(a[mglev],TilingIfNotGPU()); mfi.isValid(); ++mfi)
653 {
654 AMREX_D_TERM(Box const& xbx = mfi.nodaltilebox(0);,
655 Box const& ybx = mfi.nodaltilebox(1);,
656 Box const& zbx = mfi.nodaltilebox(2));
657 AMREX_D_TERM(auto const& bx = b[mglev][0].array(mfi);,
658 auto const& by = b[mglev][1].array(mfi);,
659 auto const& bz = b[mglev][2].array(mfi));
660 auto const& osm = this->m_overset_mask[amrlev][mglev]->const_array(mfi);
661#if defined(AMREX_USE_CUDA) && defined(_WIN32)
663 (xbx, t_xbx,
664 {
665 overset_rescale_bcoef_x(t_xbx, bx, osm, ncomp, osfac);
666 });
667#if (AMREX_SPACEDIM >= 2)
669 (ybx, t_ybx,
670 {
671 overset_rescale_bcoef_y(t_ybx, by, osm, ncomp, osfac);
672 });
673#endif
674#if (AMREX_SPACEDIM == 3)
676 (zbx, t_zbx,
677 {
678 overset_rescale_bcoef_z(t_zbx, bz, osm, ncomp, osfac);
679 });
680#endif
681#else
683 (xbx, t_xbx,
684 {
685 overset_rescale_bcoef_x(t_xbx, bx, osm, ncomp, osfac);
686 },
687 ybx, t_ybx,
688 {
689 overset_rescale_bcoef_y(t_ybx, by, osm, ncomp, osfac);
690 },
691 zbx, t_zbx,
692 {
693 overset_rescale_bcoef_z(t_zbx, bz, osm, ncomp, osfac);
694 });
695#endif
696 }
697 }
698 }
699}
700
701template <typename MF>
702void
704{
705 auto& fine_a_coeffs = m_a_coeffs[flev ].back();
706 auto& fine_b_coeffs = m_b_coeffs[flev ].back();
707 auto& crse_a_coeffs = m_a_coeffs[flev-1].front();
708 auto& crse_b_coeffs = m_b_coeffs[flev-1].front();
709
710 if (m_a_scalar != 0.0) {
711 // We coarsen from the back of flev to the front of flev-1.
712 // So we use mg_coarsen_ratio.
713 amrex::average_down(fine_a_coeffs, crse_a_coeffs, 0, 1, this->mg_coarsen_ratio);
714 }
715
717 amrex::GetArrOfPtrs(crse_b_coeffs),
718 IntVect(this->mg_coarsen_ratio),
719 this->m_geom[flev-1][0]);
720}
721
722template <typename MF>
723void
725{
726 m_is_singular.clear();
727 m_is_singular.resize(this->m_num_amr_levels, false);
728 auto itlo = std::find(this->m_lobc[0].begin(), this->m_lobc[0].end(), BCType::Dirichlet);
729 auto ithi = std::find(this->m_hibc[0].begin(), this->m_hibc[0].end(), BCType::Dirichlet);
730 if (itlo == this->m_lobc[0].end() && ithi == this->m_hibc[0].end())
731 { // No Dirichlet
732 for (int alev = 0; alev < this->m_num_amr_levels; ++alev)
733 {
734 // For now this assumes that overset regions are treated as Dirichlet bc's
735 if (this->m_domain_covered[alev] && !this->m_overset_mask[alev][0])
736 {
737 if (m_a_scalar == Real(0.0))
738 {
739 m_is_singular[alev] = true;
740 }
741 else
742 {
743 RT asum = m_a_coeffs[alev].back().sum(0,IntVect(0));
744 RT amax = m_a_coeffs[alev].back().norminf(0,1,IntVect(0));
745 m_is_singular[alev] = (std::abs(asum) <= amax * RT(1.e-12));
746 }
747 }
748 }
749 }
750
751 if (!m_is_singular[0] && this->m_needs_coarse_data_for_bc &&
752 this->m_coarse_fine_bc_type == LinOpBCType::Neumann)
753 {
754 AMREX_ASSERT(this->m_overset_mask[0][0] == nullptr);
755
756 bool lev0_a_is_zero = false;
757 if (m_a_scalar == Real(0.0)) {
758 lev0_a_is_zero = true;
759 } else {
760 RT asum = m_a_coeffs[0].back().sum(0,IntVect(0));
761 RT amax = m_a_coeffs[0].back().norminf(0,1,IntVect(0));
762 bool a_is_almost_zero = std::abs(asum) <= amax * RT(1.e-12);
763 if (a_is_almost_zero) { lev0_a_is_zero = true; }
764 }
765
766 if (lev0_a_is_zero) {
767 auto bbox = this->m_grids[0][0].minimalBox();
768 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
769 if (this->m_lobc[0][idim] == LinOpBCType::Dirichlet) {
770 bbox.growLo(idim,1);
771 }
772 if (this->m_hibc[0][idim] == LinOpBCType::Dirichlet) {
773 bbox.growHi(idim,1);
774 }
775 }
776 if (this->m_geom[0][0].Domain().contains(bbox)) {
777 m_is_singular[0] = true;
778 }
779 }
780 }
781}
782
783template <typename MF>
784void
785MLABecLaplacianT<MF>::Fapply (int amrlev, int mglev, MF& out, const MF& in) const
786{
787 BL_PROFILE("MLABecLaplacian::Fapply()");
788
789 const MF& acoef = m_a_coeffs[amrlev][mglev];
790 AMREX_D_TERM(const MF& bxcoef = m_b_coeffs[amrlev][mglev][0];,
791 const MF& bycoef = m_b_coeffs[amrlev][mglev][1];,
792 const MF& bzcoef = m_b_coeffs[amrlev][mglev][2];);
793
795 {AMREX_D_DECL(static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(0)),
796 static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(1)),
797 static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(2)))};
798
799 const RT ascalar = m_a_scalar;
800 const RT bscalar = m_b_scalar;
801
802 const int ncomp = this->getNComp();
803
804#ifdef AMREX_USE_GPU
805 if (Gpu::inLaunchRegion()) {
806 const auto& xma = in.const_arrays();
807 const auto& yma = out.arrays();
808 const auto& ama = acoef.arrays();
809 AMREX_D_TERM(const auto& bxma = bxcoef.const_arrays();,
810 const auto& byma = bycoef.const_arrays();,
811 const auto& bzma = bzcoef.const_arrays(););
812 if (this->m_overset_mask[amrlev][mglev]) {
813 const auto& osmma = this->m_overset_mask[amrlev][mglev]->const_arrays();
814 ParallelFor(out, IntVect(0), ncomp,
815 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
816 {
817 mlabeclap_adotx_os(i,j,k,n, yma[box_no], xma[box_no], ama[box_no],
818 AMREX_D_DECL(bxma[box_no],byma[box_no],bzma[box_no]),
819 osmma[box_no], dxinv, ascalar, bscalar);
820 });
821 } else {
822 ParallelFor(out, IntVect(0), ncomp,
823 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
824 {
825 mlabeclap_adotx(i,j,k,n, yma[box_no], xma[box_no], ama[box_no],
826 AMREX_D_DECL(bxma[box_no],byma[box_no],bzma[box_no]),
827 dxinv, ascalar, bscalar);
828 });
829 }
830 if (!Gpu::inNoSyncRegion()) {
832 }
833 } else
834#endif
835 {
836#ifdef AMREX_USE_OMP
837#pragma omp parallel if (Gpu::notInLaunchRegion())
838#endif
839 for (MFIter mfi(out, TilingIfNotGPU()); mfi.isValid(); ++mfi)
840 {
841 const Box& bx = mfi.tilebox();
842 const auto& xfab = in.array(mfi);
843 const auto& yfab = out.array(mfi);
844 const auto& afab = acoef.array(mfi);
845 AMREX_D_TERM(const auto& bxfab = bxcoef.array(mfi);,
846 const auto& byfab = bycoef.array(mfi);,
847 const auto& bzfab = bzcoef.array(mfi););
848 if (this->m_overset_mask[amrlev][mglev]) {
849 const auto& osm = this->m_overset_mask[amrlev][mglev]->const_array(mfi);
850 AMREX_HOST_DEVICE_PARALLEL_FOR_4D(bx, ncomp, i, j, k, n,
851 {
852 mlabeclap_adotx_os(i,j,k,n, yfab, xfab, afab, AMREX_D_DECL(bxfab,byfab,bzfab),
853 osm, dxinv, ascalar, bscalar);
854 });
855 } else {
856 AMREX_HOST_DEVICE_PARALLEL_FOR_4D(bx, ncomp, i, j, k, n,
857 {
858 mlabeclap_adotx(i,j,k,n, yfab, xfab, afab, AMREX_D_DECL(bxfab,byfab,bzfab),
859 dxinv, ascalar, bscalar);
860 });
861 }
862 }
863 }
864}
865
866template <typename MF>
867void
868MLABecLaplacianT<MF>::Fsmooth (int amrlev, int mglev, MF& sol, const MF& rhs, int redblack) const
869{
870 BL_PROFILE("MLABecLaplacian::Fsmooth()");
871
872 bool regular_coarsening = true;
873 if (amrlev == 0 && mglev > 0) {
874 regular_coarsening = this->mg_coarsen_ratio_vec[mglev-1] == this->mg_coarsen_ratio;
875 }
876
877 MF Ax;
878 if (! this->m_use_gauss_seidel && regular_coarsening) { // jacobi
879 Ax.define(sol.boxArray(), sol.DistributionMap(), sol.nComp(), 0);
880 Fapply(amrlev, mglev, Ax, sol);
881 }
882
883 const MF& acoef = m_a_coeffs[amrlev][mglev];
884 AMREX_ALWAYS_ASSERT(acoef.nGrowVect() == 0);
885 AMREX_D_TERM(const MF& bxcoef = m_b_coeffs[amrlev][mglev][0];,
886 const MF& bycoef = m_b_coeffs[amrlev][mglev][1];,
887 const MF& bzcoef = m_b_coeffs[amrlev][mglev][2];);
888 const auto& undrrelxr = this->m_undrrelxr[amrlev][mglev];
889 const auto& maskvals = this->m_maskvals [amrlev][mglev];
890
891 OrientationIter oitr;
892
893 const auto& f0 = undrrelxr[oitr()]; ++oitr;
894 const auto& f1 = undrrelxr[oitr()]; ++oitr;
895#if (AMREX_SPACEDIM > 1)
896 const auto& f2 = undrrelxr[oitr()]; ++oitr;
897 const auto& f3 = undrrelxr[oitr()]; ++oitr;
898#if (AMREX_SPACEDIM > 2)
899 const auto& f4 = undrrelxr[oitr()]; ++oitr;
900 const auto& f5 = undrrelxr[oitr()]; ++oitr;
901#endif
902#endif
903
904 const MultiMask& mm0 = maskvals[0];
905 const MultiMask& mm1 = maskvals[1];
906#if (AMREX_SPACEDIM > 1)
907 const MultiMask& mm2 = maskvals[2];
908 const MultiMask& mm3 = maskvals[3];
909#if (AMREX_SPACEDIM > 2)
910 const MultiMask& mm4 = maskvals[4];
911 const MultiMask& mm5 = maskvals[5];
912#endif
913#endif
914
915 const int nc = this->getNComp();
916 const Real* h = this->m_geom[amrlev][mglev].CellSize();
917 AMREX_D_TERM(const RT dhx = m_b_scalar/static_cast<RT>(h[0]*h[0]);,
918 const RT dhy = m_b_scalar/static_cast<RT>(h[1]*h[1]);,
919 const RT dhz = m_b_scalar/static_cast<RT>(h[2]*h[2]));
920 const RT alpha = m_a_scalar;
921
922#ifdef AMREX_USE_GPU
924 && (this->m_overset_mask[amrlev][mglev] || regular_coarsening))
925 {
926 const auto& m0ma = mm0.const_arrays();
927 const auto& m1ma = mm1.const_arrays();
928#if (AMREX_SPACEDIM > 1)
929 const auto& m2ma = mm2.const_arrays();
930 const auto& m3ma = mm3.const_arrays();
931#if (AMREX_SPACEDIM > 2)
932 const auto& m4ma = mm4.const_arrays();
933 const auto& m5ma = mm5.const_arrays();
934#endif
935#endif
936
937 const auto& solnma = sol.arrays();
938 const auto& rhsma = rhs.const_arrays();
939 const auto& ama = acoef.const_arrays();
940
941 AMREX_D_TERM(const auto& bxma = bxcoef.const_arrays();,
942 const auto& byma = bycoef.const_arrays();,
943 const auto& bzma = bzcoef.const_arrays(););
944
945 const auto& f0ma = f0.const_arrays();
946 const auto& f1ma = f1.const_arrays();
947#if (AMREX_SPACEDIM > 1)
948 const auto& f2ma = f2.const_arrays();
949 const auto& f3ma = f3.const_arrays();
950#if (AMREX_SPACEDIM > 2)
951 const auto& f4ma = f4.const_arrays();
952 const auto& f5ma = f5.const_arrays();
953#endif
954#endif
955
956 if (this->m_overset_mask[amrlev][mglev]) {
957 const auto& osmma = this->m_overset_mask[amrlev][mglev]->const_arrays();
958 if (this->m_use_gauss_seidel) {
959 ParallelFor(sol, IntVect(0), nc,
960 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
961 {
962 Box vbx(ama[box_no]);
963 abec_gsrb_os(i,j,k,n, solnma[box_no], rhsma[box_no], alpha, ama[box_no],
964 AMREX_D_DECL(dhx, dhy, dhz),
965 AMREX_D_DECL(bxma[box_no],byma[box_no],bzma[box_no]),
966 AMREX_D_DECL(m0ma[box_no],m2ma[box_no],m4ma[box_no]),
967 AMREX_D_DECL(m1ma[box_no],m3ma[box_no],m5ma[box_no]),
968 AMREX_D_DECL(f0ma[box_no],f2ma[box_no],f4ma[box_no]),
969 AMREX_D_DECL(f1ma[box_no],f3ma[box_no],f5ma[box_no]),
970 osmma[box_no], vbx, redblack);
971 });
972 } else {
973 const auto& axma = Ax.const_arrays();
974 ParallelFor(sol, IntVect(0), nc,
975 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
976 {
977 Box vbx(ama[box_no]);
978 abec_jacobi_os(i,j,k,n, solnma[box_no], rhsma[box_no], axma[box_no],
979 alpha, ama[box_no],
980 AMREX_D_DECL(dhx, dhy, dhz),
981 AMREX_D_DECL(bxma[box_no],byma[box_no],bzma[box_no]),
982 AMREX_D_DECL(m0ma[box_no],m2ma[box_no],m4ma[box_no]),
983 AMREX_D_DECL(m1ma[box_no],m3ma[box_no],m5ma[box_no]),
984 AMREX_D_DECL(f0ma[box_no],f2ma[box_no],f4ma[box_no]),
985 AMREX_D_DECL(f1ma[box_no],f3ma[box_no],f5ma[box_no]),
986 osmma[box_no], vbx);
987 });
988 }
989 } else if (regular_coarsening) {
990 if (this->m_use_gauss_seidel) {
991 ParallelFor(sol, IntVect(0), nc,
992 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
993 {
994 Box vbx(ama[box_no]);
995 abec_gsrb(i,j,k,n, solnma[box_no], rhsma[box_no], alpha, ama[box_no],
996 AMREX_D_DECL(dhx, dhy, dhz),
997 AMREX_D_DECL(bxma[box_no],byma[box_no],bzma[box_no]),
998 AMREX_D_DECL(m0ma[box_no],m2ma[box_no],m4ma[box_no]),
999 AMREX_D_DECL(m1ma[box_no],m3ma[box_no],m5ma[box_no]),
1000 AMREX_D_DECL(f0ma[box_no],f2ma[box_no],f4ma[box_no]),
1001 AMREX_D_DECL(f1ma[box_no],f3ma[box_no],f5ma[box_no]),
1002 vbx, redblack);
1003 });
1004 } else {
1005 const auto& axma = Ax.const_arrays();
1006 ParallelFor(sol, IntVect(0), nc,
1007 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1008 {
1009 Box vbx(ama[box_no]);
1010 abec_jacobi(i,j,k,n, solnma[box_no], rhsma[box_no], axma[box_no],
1011 alpha, ama[box_no],
1012 AMREX_D_DECL(dhx, dhy, dhz),
1013 AMREX_D_DECL(bxma[box_no],byma[box_no],bzma[box_no]),
1014 AMREX_D_DECL(m0ma[box_no],m2ma[box_no],m4ma[box_no]),
1015 AMREX_D_DECL(m1ma[box_no],m3ma[box_no],m5ma[box_no]),
1016 AMREX_D_DECL(f0ma[box_no],f2ma[box_no],f4ma[box_no]),
1017 AMREX_D_DECL(f1ma[box_no],f3ma[box_no],f5ma[box_no]),
1018 vbx);
1019 });
1020 }
1021 }
1022 if (Ax.local_size() > 0 || !Gpu::inNoSyncRegion()) {
1024 }
1025 } else
1026#endif
1027 {
1028 MFItInfo mfi_info;
1029 mfi_info.EnableTiling().SetDynamic(true);
1030
1031#ifdef AMREX_USE_OMP
1032#pragma omp parallel if (Gpu::notInLaunchRegion())
1033#endif
1034 for (MFIter mfi(sol,mfi_info); mfi.isValid(); ++mfi)
1035 {
1036 const auto& m0 = mm0.array(mfi);
1037 const auto& m1 = mm1.array(mfi);
1038#if (AMREX_SPACEDIM > 1)
1039 const auto& m2 = mm2.array(mfi);
1040 const auto& m3 = mm3.array(mfi);
1041#if (AMREX_SPACEDIM > 2)
1042 const auto& m4 = mm4.array(mfi);
1043 const auto& m5 = mm5.array(mfi);
1044#endif
1045#endif
1046
1047 const Box& tbx = mfi.tilebox();
1048 const Box& vbx = mfi.validbox();
1049 const auto& solnfab = sol.array(mfi);
1050 const auto& rhsfab = rhs.const_array(mfi);
1051 const auto& afab = acoef.const_array(mfi);
1052
1053 AMREX_D_TERM(const auto& bxfab = bxcoef.const_array(mfi);,
1054 const auto& byfab = bycoef.const_array(mfi);,
1055 const auto& bzfab = bzcoef.const_array(mfi););
1056
1057 const auto& f0fab = f0.const_array(mfi);
1058 const auto& f1fab = f1.const_array(mfi);
1059#if (AMREX_SPACEDIM > 1)
1060 const auto& f2fab = f2.const_array(mfi);
1061 const auto& f3fab = f3.const_array(mfi);
1062#if (AMREX_SPACEDIM > 2)
1063 const auto& f4fab = f4.const_array(mfi);
1064 const auto& f5fab = f5.const_array(mfi);
1065#endif
1066#endif
1067
1068 if (this->m_overset_mask[amrlev][mglev]) {
1069 const auto& osm = this->m_overset_mask[amrlev][mglev]->const_array(mfi);
1070 if (this->m_use_gauss_seidel) {
1071 AMREX_LOOP_4D(tbx, nc, i, j, k, n,
1072 {
1073 abec_gsrb_os(i,j,k,n, solnfab, rhsfab, alpha, afab,
1074 AMREX_D_DECL(dhx, dhy, dhz),
1075 AMREX_D_DECL(bxfab, byfab, bzfab),
1076 AMREX_D_DECL(m0,m2,m4),
1077 AMREX_D_DECL(m1,m3,m5),
1078 AMREX_D_DECL(f0fab,f2fab,f4fab),
1079 AMREX_D_DECL(f1fab,f3fab,f5fab),
1080 osm, vbx, redblack);
1081 });
1082 } else {
1083 const auto& axfab = Ax.const_array(mfi);
1084 AMREX_LOOP_4D(tbx, nc, i, j, k, n,
1085 {
1086 abec_jacobi_os(i,j,k,n, solnfab, rhsfab, axfab,
1087 alpha, afab,
1088 AMREX_D_DECL(dhx, dhy, dhz),
1089 AMREX_D_DECL(bxfab, byfab, bzfab),
1090 AMREX_D_DECL(m0,m2,m4),
1091 AMREX_D_DECL(m1,m3,m5),
1092 AMREX_D_DECL(f0fab,f2fab,f4fab),
1093 AMREX_D_DECL(f1fab,f3fab,f5fab),
1094 osm, vbx);
1095 });
1096 }
1097 } else if (regular_coarsening) {
1098 if (this->m_use_gauss_seidel) {
1099 AMREX_LOOP_4D(tbx, nc, i, j, k, n,
1100 {
1101 abec_gsrb(i,j,k,n, solnfab, rhsfab, alpha, afab,
1102 AMREX_D_DECL(dhx, dhy, dhz),
1103 AMREX_D_DECL(bxfab, byfab, bzfab),
1104 AMREX_D_DECL(m0,m2,m4),
1105 AMREX_D_DECL(m1,m3,m5),
1106 AMREX_D_DECL(f0fab,f2fab,f4fab),
1107 AMREX_D_DECL(f1fab,f3fab,f5fab),
1108 vbx, redblack);
1109 });
1110 } else {
1111 const auto& axfab = Ax.const_array(mfi);
1112 AMREX_LOOP_4D(tbx, nc, i, j, k, n,
1113 {
1114 abec_jacobi(i,j,k,n, solnfab, rhsfab, axfab,
1115 alpha, afab,
1116 AMREX_D_DECL(dhx, dhy, dhz),
1117 AMREX_D_DECL(bxfab, byfab, bzfab),
1118 AMREX_D_DECL(m0,m2,m4),
1119 AMREX_D_DECL(m1,m3,m5),
1120 AMREX_D_DECL(f0fab,f2fab,f4fab),
1121 AMREX_D_DECL(f1fab,f3fab,f5fab),
1122 vbx);
1123 });
1124 }
1125 } else {
1126 // line solve does not with with GPU
1127 abec_gsrb_with_line_solve(tbx, solnfab, rhsfab, alpha, afab,
1128 AMREX_D_DECL(dhx, dhy, dhz),
1129 AMREX_D_DECL(bxfab, byfab, bzfab),
1130 AMREX_D_DECL(m0,m2,m4),
1131 AMREX_D_DECL(m1,m3,m5),
1132 AMREX_D_DECL(f0fab,f2fab,f4fab),
1133 AMREX_D_DECL(f1fab,f3fab,f5fab),
1134 vbx, redblack, nc);
1135 }
1136 }
1137 }
1138}
1139
1140template <typename MF>
1141void
1142MLABecLaplacianT<MF>::FFlux (int amrlev, const MFIter& mfi,
1143 const Array<FAB*,AMREX_SPACEDIM>& flux,
1144 const FAB& sol, Location, int face_only) const
1145{
1146 BL_PROFILE("MLABecLaplacian::FFlux()");
1147
1148 const int mglev = 0;
1149 const Box& box = mfi.tilebox();
1150 const Real* dxinv = this->m_geom[amrlev][mglev].InvCellSize();
1151 const int ncomp = this->getNComp();
1152 FFlux(box, dxinv, m_b_scalar,
1153 Array<FAB const*,AMREX_SPACEDIM>{{AMREX_D_DECL(&(m_b_coeffs[amrlev][mglev][0][mfi]),
1154 &(m_b_coeffs[amrlev][mglev][1][mfi]),
1155 &(m_b_coeffs[amrlev][mglev][2][mfi]))}},
1156 flux, sol, face_only, ncomp);
1157}
1158
1159template <typename MF>
1160void
1161MLABecLaplacianT<MF>::FFlux (Box const& box, Real const* dxinv, RT bscalar,
1163 Array<FAB*,AMREX_SPACEDIM> const& flux,
1164 FAB const& sol, int face_only, int ncomp)
1165{
1166 AMREX_D_TERM(const auto bx = bcoef[0]->const_array();,
1167 const auto by = bcoef[1]->const_array();,
1168 const auto bz = bcoef[2]->const_array(););
1169 AMREX_D_TERM(const auto& fxarr = flux[0]->array();,
1170 const auto& fyarr = flux[1]->array();,
1171 const auto& fzarr = flux[2]->array(););
1172 const auto& solarr = sol.array();
1173
1174 if (face_only)
1175 {
1176 RT fac = bscalar*static_cast<RT>(dxinv[0]);
1177 Box blo = amrex::bdryLo(box, 0);
1178 int blen = box.length(0);
1180 {
1181 mlabeclap_flux_xface(tbox, fxarr, solarr, bx, fac, blen, ncomp);
1182 });
1183#if (AMREX_SPACEDIM >= 2)
1184 fac = bscalar*static_cast<RT>(dxinv[1]);
1185 blo = amrex::bdryLo(box, 1);
1186 blen = box.length(1);
1188 {
1189 mlabeclap_flux_yface(tbox, fyarr, solarr, by, fac, blen, ncomp);
1190 });
1191#endif
1192#if (AMREX_SPACEDIM == 3)
1193 fac = bscalar*static_cast<RT>(dxinv[2]);
1194 blo = amrex::bdryLo(box, 2);
1195 blen = box.length(2);
1197 {
1198 mlabeclap_flux_zface(tbox, fzarr, solarr, bz, fac, blen, ncomp);
1199 });
1200#endif
1201 }
1202 else
1203 {
1204 RT fac = bscalar*static_cast<RT>(dxinv[0]);
1205 Box bflux = amrex::surroundingNodes(box, 0);
1207 {
1208 mlabeclap_flux_x(tbox, fxarr, solarr, bx, fac, ncomp);
1209 });
1210#if (AMREX_SPACEDIM >= 2)
1211 fac = bscalar*static_cast<RT>(dxinv[1]);
1212 bflux = amrex::surroundingNodes(box, 1);
1214 {
1215 mlabeclap_flux_y(tbox, fyarr, solarr, by, fac, ncomp);
1216 });
1217#endif
1218#if (AMREX_SPACEDIM == 3)
1219 fac = bscalar*static_cast<RT>(dxinv[2]);
1220 bflux = amrex::surroundingNodes(box, 2);
1222 {
1223 mlabeclap_flux_z(tbox, fzarr, solarr, bz, fac, ncomp);
1224 });
1225#endif
1226 }
1227}
1228
1229template <typename MF>
1230void
1231MLABecLaplacianT<MF>::normalize (int amrlev, int mglev, MF& mf) const
1232{
1233 BL_PROFILE("MLABecLaplacian::normalize()");
1234
1235 const auto& acoef = m_a_coeffs[amrlev][mglev];
1236 AMREX_D_TERM(const auto& bxcoef = m_b_coeffs[amrlev][mglev][0];,
1237 const auto& bycoef = m_b_coeffs[amrlev][mglev][1];,
1238 const auto& bzcoef = m_b_coeffs[amrlev][mglev][2];);
1239
1240 const GpuArray<RT,AMREX_SPACEDIM> dxinv
1241 {AMREX_D_DECL(static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(0)),
1242 static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(1)),
1243 static_cast<RT>(this->m_geom[amrlev][mglev].InvCellSize(2)))};
1244
1245 const RT ascalar = m_a_scalar;
1246 const RT bscalar = m_b_scalar;
1247
1248 const int ncomp = getNComp();
1249
1250#ifdef AMREX_USE_GPU
1251 if (Gpu::inLaunchRegion()) {
1252 const auto& ma = mf.arrays();
1253 const auto& ama = acoef.const_arrays();
1254 AMREX_D_TERM(const auto& bxma = bxcoef.const_arrays();,
1255 const auto& byma = bycoef.const_arrays();,
1256 const auto& bzma = bzcoef.const_arrays(););
1257 ParallelFor(mf, IntVect(0), ncomp,
1258 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1259 {
1260 mlabeclap_normalize(i,j,k,n, ma[box_no], ama[box_no],
1261 AMREX_D_DECL(bxma[box_no],byma[box_no],bzma[box_no]),
1262 dxinv, ascalar, bscalar);
1263 });
1264 if (!Gpu::inNoSyncRegion()) {
1266 }
1267 } else
1268#endif
1269 {
1270#ifdef AMREX_USE_OMP
1271#pragma omp parallel if (Gpu::notInLaunchRegion())
1272#endif
1273 for (MFIter mfi(mf, TilingIfNotGPU()); mfi.isValid(); ++mfi)
1274 {
1275 const Box& bx = mfi.tilebox();
1276 const auto& fab = mf.array(mfi);
1277 const auto& afab = acoef.array(mfi);
1278 AMREX_D_TERM(const auto& bxfab = bxcoef.array(mfi);,
1279 const auto& byfab = bycoef.array(mfi);,
1280 const auto& bzfab = bzcoef.array(mfi););
1281
1282 AMREX_HOST_DEVICE_PARALLEL_FOR_4D(bx, ncomp, i, j, k, n,
1283 {
1284 mlabeclap_normalize(i,j,k,n, fab, afab, AMREX_D_DECL(bxfab,byfab,bzfab),
1285 dxinv, ascalar, bscalar);
1286 });
1287 }
1288 }
1289}
1290
1291template <typename MF>
1292bool
1294{
1295 bool support = false;
1296 if (this->m_overset_mask[0][0]) {
1297 if (this->m_geom[0].back().Domain().coarsenable(MLLinOp::mg_coarsen_ratio,
1298 this->mg_domain_min_width)
1299 && this->m_grids[0].back().coarsenable(MLLinOp::mg_coarsen_ratio, MLLinOp::mg_box_min_width))
1300 {
1301 support = true;
1302 }
1303 }
1304 return support;
1305}
1306
1307template <typename MF>
1308std::unique_ptr<MLLinOpT<MF>>
1309MLABecLaplacianT<MF>::makeNLinOp (int /*grid_size*/) const
1310{
1311 if (this->m_overset_mask[0][0] == nullptr) { return nullptr; }
1312
1313 const Geometry& geom = this->m_geom[0].back();
1314 const BoxArray& ba = this->m_grids[0].back();
1315 const DistributionMapping& dm = this->m_dmap[0].back();
1316
1317 std::unique_ptr<MLLinOpT<MF>> r
1318 {new MLABecLaplacianT<MF>({geom}, {ba}, {dm}, this->m_lpinfo_arg)};
1319
1320 auto nop = dynamic_cast<MLABecLaplacianT<MF>*>(r.get());
1321 if (!nop) {
1322 return nullptr;
1323 }
1324
1325 nop->m_parent = this;
1326
1327 nop->setMaxOrder(this->maxorder);
1328 nop->setVerbose(this->verbose);
1329
1330 nop->setDomainBC(this->m_lobc, this->m_hibc);
1331
1332 if (this->needsCoarseDataForBC())
1333 {
1334 const Real* dx0 = this->m_geom[0][0].CellSize();
1335 RealVect fac(this->m_coarse_data_crse_ratio);
1336 fac *= Real(0.5);
1337 RealVect cbloc {AMREX_D_DECL(dx0[0]*fac[0], dx0[1]*fac[1], dx0[2]*fac[2])};
1338 nop->setCoarseFineBCLocation(cbloc);
1339 }
1340
1341 nop->setScalars(m_a_scalar, m_b_scalar);
1342
1343 MF const& alpha_bottom = m_a_coeffs[0].back();
1344 iMultiFab const& osm_bottom = *(this->m_overset_mask[0].back());
1345 const int ncomp = alpha_bottom.nComp();
1346 MF alpha(ba, dm, ncomp, 0);
1347
1348 RT a_max = alpha_bottom.norminf(0, ncomp, IntVect(0), true, true);
1349 const int ncomp_b = m_b_coeffs[0].back()[0].nComp();
1350 AMREX_D_TERM(RT bx_max = m_b_coeffs[0].back()[0].norminf(0,ncomp_b,IntVect(0),true,true);,
1351 RT by_max = m_b_coeffs[0].back()[1].norminf(0,ncomp_b,IntVect(0),true,true);,
1352 RT bz_max = m_b_coeffs[0].back()[2].norminf(0,ncomp_b,IntVect(0),true,true));
1353 const GpuArray<RT,AMREX_SPACEDIM> dxinv
1354 {AMREX_D_DECL(static_cast<RT>(geom.InvCellSize(0)),
1355 static_cast<RT>(geom.InvCellSize(1)),
1356 static_cast<RT>(geom.InvCellSize(2)))};
1357 RT huge_alpha = RT(1.e30) *
1358 amrex::max(a_max*std::abs(m_a_scalar),
1359 AMREX_D_DECL(std::abs(m_b_scalar)*bx_max*dxinv[0]*dxinv[0],
1360 std::abs(m_b_scalar)*by_max*dxinv[1]*dxinv[1],
1361 std::abs(m_b_scalar)*bz_max*dxinv[2]*dxinv[2]));
1363
1364#ifdef AMREX_USE_GPU
1365 if (Gpu::inLaunchRegion() && alpha.isFusingCandidate()) {
1366 auto const& ama = alpha.arrays();
1367 auto const& abotma = alpha_bottom.const_arrays();
1368 auto const& mma = osm_bottom.const_arrays();
1369 ParallelFor(alpha, IntVect(0), ncomp,
1370 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1371 {
1372 if (mma[box_no](i,j,k)) {
1373 ama[box_no](i,j,k,n) = abotma[box_no](i,j,k,n);
1374 } else {
1375 ama[box_no](i,j,k,n) = huge_alpha;
1376 }
1377 });
1378 if (alpha.local_size() > 0 || !Gpu::inNoSyncRegion()) {
1380 }
1381 } else
1382#endif
1383 {
1384#ifdef AMREX_USE_OMP
1385#pragma omp parallel if (Gpu::notInLaunchRegion())
1386#endif
1387 for (MFIter mfi(alpha,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
1388 Box const& bx = mfi.tilebox();
1389 auto const& a = alpha.array(mfi);
1390 auto const& abot = alpha_bottom.const_array(mfi);
1391 auto const& m = osm_bottom.const_array(mfi);
1392 AMREX_HOST_DEVICE_PARALLEL_FOR_4D(bx, ncomp, i, j, k, n,
1393 {
1394 if (m(i,j,k)) {
1395 a(i,j,k,n) = abot(i,j,k,n);
1396 } else {
1397 a(i,j,k,n) = huge_alpha;
1398 }
1399 });
1400 }
1401 }
1402
1403 nop->setACoeffs(0, alpha);
1404 nop->setBCoeffs(0, GetArrOfConstPtrs(m_b_coeffs[0].back()));
1405
1406 return r;
1407}
1408
1409template <typename MF>
1410void
1411MLABecLaplacianT<MF>::copyNSolveSolution (MF& dst, MF const& src) const
1412{
1413 if (this->m_overset_mask[0].back() == nullptr) { return; }
1414
1415 const int ncomp = dst.nComp();
1416
1417#ifdef AMREX_USE_GPU
1418 if (Gpu::inLaunchRegion() && dst.isFusingCandidate()) {
1419 auto const& dstma = dst.arrays();
1420 auto const& srcma = src.const_arrays();
1421 auto const& mma = this->m_overset_mask[0].back()->const_arrays();
1422 ParallelFor(dst, IntVect(0), ncomp,
1423 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1424 {
1425 if (mma[box_no](i,j,k)) {
1426 dstma[box_no](i,j,k,n) = srcma[box_no](i,j,k,n);
1427 } else {
1428 dstma[box_no](i,j,k,n) = RT(0.0);
1429 }
1430 });
1431 if (!Gpu::inNoSyncRegion()) {
1433 }
1434 } else
1435#endif
1436 {
1437#ifdef AMREX_USE_OMP
1438#pragma omp parallel if (Gpu::notInLaunchRegion())
1439#endif
1440 for (MFIter mfi(dst,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
1441 Box const& bx = mfi.tilebox();
1442 auto const& dfab = dst.array(mfi);
1443 auto const& sfab = src.const_array(mfi);
1444 auto const& m = this->m_overset_mask[0].back()->const_array(mfi);
1445 AMREX_HOST_DEVICE_PARALLEL_FOR_4D(bx, ncomp, i, j, k, n,
1446 {
1447 if (m(i,j,k)) {
1448 dfab(i,j,k,n) = sfab(i,j,k,n);
1449 } else {
1450 dfab(i,j,k,n) = RT(0.0);
1451 }
1452 });
1453 }
1454 }
1455}
1456
1457extern template class MLABecLaplacianT<MultiFab>;
1458
1461
1462}
1463
1464#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define AMREX_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:37
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_GPU_LAUNCH_HOST_DEVICE_LAMBDA_RANGE(TN, TI, block)
Definition AMReX_GpuLaunchMacrosC.nolint.H:4
#define AMREX_HOST_DEVICE_FOR_3D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:106
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#define AMREX_LAUNCH_HOST_DEVICE_LAMBDA_DIM(...)
Definition AMReX_GpuLaunch.nolint.H:37
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
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
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:155
const Real * InvCellSize() const noexcept
Returns the inverse cellsize for each coordinate direction.
Definition AMReX_CoordSys.H:82
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:83
MultiArray4< typename FabArray< FAB >::value_type const > const_arrays() const noexcept
Definition AMReX_FabArray.H:647
Array4< typename FabArray< FAB >::value_type const > const_array(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:585
Definition AMReX_FabFactory.H:50
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
__host__ static __device__ constexpr IntVectND< dim > TheDimensionVector(int d) noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:790
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
Box tilebox() const noexcept
Return the tile Box at the current index.
Definition AMReX_MFIter.cpp:389
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
Definition AMReX_MLABecLaplacian.H:15
RT getBScalar() const final
Definition AMReX_MLABecLaplacian.H:156
int getNComp() const override
Return number of components.
Definition AMReX_MLABecLaplacian.H:136
void FFlux(int amrlev, const MFIter &mfi, const Array< FAB *, 3 > &flux, const FAB &sol, Location, int face_only=0) const override
Definition AMReX_MLABecLaplacian.H:1142
void setACoeffs(int amrlev, const AMF &alpha)
Definition AMReX_MLABecLaplacian.H:311
bool isSingular(int amrlev) const override
Is it singular on given AMR level?
Definition AMReX_MLABecLaplacian.H:144
void applyRobinBCTermsCoeffs()
Definition AMReX_MLABecLaplacian.H:589
Vector< Vector< Array< MF, 3 > > > m_b_coeffs
Definition AMReX_MLABecLaplacian.H:185
typename MF::fab_type FAB
Definition AMReX_MLABecLaplacian.H:18
void setBCoeffs(int amrlev, const Array< AMF const *, 3 > &beta)
Definition AMReX_MLABecLaplacian.H:336
MLABecLaplacianT< MF > & operator=(const MLABecLaplacianT< MF > &)=delete
bool supportRobinBC() const noexcept override
Definition AMReX_MLABecLaplacian.H:194
RT m_b_scalar
Definition AMReX_MLABecLaplacian.H:183
typename MF::value_type RT
Definition AMReX_MLABecLaplacian.H:19
~MLABecLaplacianT() override
void prepareForSolve() override
Definition AMReX_MLABecLaplacian.H:407
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={}, int a_ncomp=1)
Definition AMReX_MLABecLaplacian.H:234
RT m_a_scalar
Definition AMReX_MLABecLaplacian.H:182
void averageDownCoeffsToCoarseAmrLevel(int flev)
Definition AMReX_MLABecLaplacian.H:703
void averageDownCoeffsSameAmrLevel(int amrlev, Vector< MF > &a, Vector< Array< MF, 3 > > &b)
Definition AMReX_MLABecLaplacian.H:616
bool m_scalars_set
Definition AMReX_MLABecLaplacian.H:187
Vector< int > m_is_singular
Definition AMReX_MLABecLaplacian.H:192
void averageDownCoeffs()
Definition AMReX_MLABecLaplacian.H:598
RT getAScalar() const final
Definition AMReX_MLABecLaplacian.H:155
void normalize(int amrlev, int mglev, MF &mf) const override
Divide mf by the diagonal component of the operator. Used by bicgstab.
Definition AMReX_MLABecLaplacian.H:1231
void Fsmooth(int amrlev, int mglev, MF &sol, const MF &rhs, int redblack) const override
Definition AMReX_MLABecLaplacian.H:868
std::unique_ptr< MLLinOpT< MF > > makeNLinOp(int) const final
Definition AMReX_MLABecLaplacian.H:1309
bool m_acoef_set
Definition AMReX_MLABecLaplacian.H:188
Array< MF const *, 3 > getBCoeffs(int amrlev, int mglev) const final
Definition AMReX_MLABecLaplacian.H:159
MLABecLaplacianT(MLABecLaplacianT< MF > &&)=delete
Vector< Vector< MF > > m_a_coeffs
Definition AMReX_MLABecLaplacian.H:184
void setScalars(T1 a, T2 b) noexcept
Definition AMReX_MLABecLaplacian.H:294
bool isBottomSingular() const override
Is the bottom of MG singular?
Definition AMReX_MLABecLaplacian.H:145
void applyMetricTermsCoeffs()
Definition AMReX_MLABecLaplacian.H:428
void copyNSolveSolution(MF &dst, MF const &src) const final
Definition AMReX_MLABecLaplacian.H:1411
bool supportNSolve() const override
Definition AMReX_MLABecLaplacian.H:1293
typename MLLinOpT< MF >::Location Location
Definition AMReX_MLABecLaplacian.H:22
MF const * getACoeffs(int amrlev, int mglev) const final
Definition AMReX_MLABecLaplacian.H:157
void update() override
Update for reuse.
Definition AMReX_MLABecLaplacian.H:386
bool needsUpdate() const override
Does it need update if it's reused?
Definition AMReX_MLABecLaplacian.H:138
void Fapply(int amrlev, int mglev, MF &out, const MF &in) const override
Definition AMReX_MLABecLaplacian.H:785
MLABecLaplacianT(const MLABecLaplacianT< MF > &)=delete
Definition AMReX_MLCellABecLap.H:14
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_MLCellABecLap.H:95
void prepareForSolve() override
Definition AMReX_MLCellABecLap.H:250
void update() override
Update for reuse.
Definition AMReX_MLCellABecLap.H:243
static constexpr int mg_coarsen_ratio
Definition AMReX_MLLinOp.H:592
static constexpr int mg_box_min_width
Definition AMReX_MLLinOp.H:593
const MLLinOpT< MF > * m_parent
Definition AMReX_MLLinOp.H:608
Definition AMReX_MultiMask.H:23
MultiArray4< int const > const_arrays() const noexcept
Return const multi-array views (alias of arrays()).
Definition AMReX_MultiMask.H:86
Array4< int const > array(const MFIter &mfi) const noexcept
Return an Array4 view (const) for iterator mfi.
Definition AMReX_MultiMask.H:69
An Iterator over the Orientation of Faces of a Box.
Definition AMReX_Orientation.H:135
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:29
Long size() const noexcept
Definition AMReX_Vector.H:54
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:26
void Max(KeyValuePair< K, V > &vi, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:133
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:88
bool notInLaunchRegion() noexcept
Definition AMReX_GpuControl.H:89
bool inNoSyncRegion() noexcept
Definition AMReX_GpuControl.H:148
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
Definition AMReX_Amr.cpp:50
__host__ __device__ 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:1744
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1567
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2852
std::array< T const *, 3 > GetArrOfConstPtrs(const std::array< T, 3 > &a) noexcept
Create an array of const-qualified pointers from an array of objects.
Definition AMReX_Array.H:1047
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:359
__host__ __device__ BoxND< dim > adjCellLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the cell centered BoxND of length len adjacent to b on the low end along the coordinate direct...
Definition AMReX_Box.H:1723
__host__ __device__ BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition AMReX_Box.H:1531
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2018
__host__ __device__ BoxND< dim > bdryLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the edge-centered BoxND (in direction dir) defining the low side of BoxND b.
Definition AMReX_Box.H:1634
void ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:202
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
void average_down_faces(const Vector< const MF * > &fine, const Vector< MF * > &crse, const IntVect &ratio, int ngcrse=0)
Average fine face-based FabArray onto crse face-based FabArray.
Definition AMReX_MultiFabUtil.H:1042
std::array< T *, 3 > GetArrOfPtrs(std::array< T, 3 > &a) noexcept
Create an array of pointers from an array of objects.
Definition AMReX_Array.H:1001
LinOpBCType
Definition AMReX_LO_BCTYPES.H:27
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:45
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
MF::value_type norminf(MF const &mf, int scomp, int ncomp, IntVect const &nghost, bool local=false)
Definition AMReX_FabArrayUtility.H:1961
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2028
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
Definition AMReX_MLLinOp.H:37
Location
Definition AMReX_MLLinOp.H:91
FabArray memory allocation information.
Definition AMReX_FabArray.H:68
Definition AMReX_MFIter.H:20
MFItInfo & SetDynamic(bool f) noexcept
Definition AMReX_MFIter.H:43
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept
Definition AMReX_MFIter.H:31