Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_MLLinOp.H
Go to the documentation of this file.
1#ifndef AMREX_ML_LINOP_H_
2#define AMREX_ML_LINOP_H_
3#include <AMReX_Config.H>
4
5#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
6#include <AMReX_Hypre.H>
8#endif
9
10#if defined(AMREX_USE_PETSC) && (AMREX_SPACEDIM > 1)
11#include <AMReX_PETSc.H>
12#endif
13
14#ifdef AMREX_USE_EB
16#include <AMReX_MultiCutFab.H>
17#endif
18
19#include <AMReX_Any.H>
20#include <AMReX_BndryRegister.H>
21#include <AMReX_FabDataType.H>
22#include <AMReX_MLMGBndry.H>
23#include <AMReX_MultiFab.H>
24#include <AMReX_MultiFabUtil.H>
25
26#include <algorithm>
27#include <iterator>
28#include <string>
29
30namespace amrex {
31
32enum class BottomSolver : int {
34};
35
36struct LPInfo
37{
38 bool do_agglomeration = true;
39 bool do_consolidation = true;
40 bool do_semicoarsening = false;
41 int agg_grid_size = -1;
42 int con_grid_size = -1;
43 int con_ratio = 2;
44 int con_strategy = 3;
45 bool has_metric_term = true;
50 bool deterministic = false;
51
52 LPInfo& setAgglomeration (bool x) noexcept { do_agglomeration = x; return *this; }
53 LPInfo& setConsolidation (bool x) noexcept { do_consolidation = x; return *this; }
54 LPInfo& setSemicoarsening (bool x) noexcept { do_semicoarsening = x; return *this; }
55 LPInfo& setAgglomerationGridSize (int x) noexcept { agg_grid_size = x; return *this; }
56 LPInfo& setConsolidationGridSize (int x) noexcept { con_grid_size = x; return *this; }
57 LPInfo& setConsolidationRatio (int x) noexcept { con_ratio = x; return *this; }
58 LPInfo& setConsolidationStrategy (int x) noexcept { con_strategy = x; return *this; }
59 LPInfo& setMetricTerm (bool x) noexcept { has_metric_term = x; return *this; }
60 LPInfo& setMaxCoarseningLevel (int n) noexcept { max_coarsening_level = n; return *this; }
61 LPInfo& setMaxSemicoarseningLevel (int n) noexcept { max_semicoarsening_level = n; return *this; }
62 LPInfo& setSemicoarseningDirection (int n) noexcept { semicoarsening_direction = n; return *this; }
63 LPInfo& setHiddenDirection (int n) noexcept { hidden_direction = n; return *this; }
64 LPInfo& setDeterministic (bool x) noexcept { deterministic = x; return *this; }
65
66 [[nodiscard]] bool hasHiddenDimension () const noexcept {
67 return hidden_direction >=0 && hidden_direction < AMREX_SPACEDIM;
68 }
69
70 static constexpr int getDefaultAgglomerationGridSize () {
71#ifdef AMREX_USE_GPU
72 return 32;
73#else
74 return AMREX_D_PICK(32, 16, 8);
75#endif
76 }
77
78 static constexpr int getDefaultConsolidationGridSize () {
79#ifdef AMREX_USE_GPU
80 return 32;
81#else
82 return AMREX_D_PICK(32, 16, 8);
83#endif
84 }
85};
86
93
94template <typename T> class MLMGT;
95template <typename T> class MLCGSolverT;
96template <typename T> class MLPoissonT;
97template <typename T> class MLABecLaplacianT;
98template <typename T> class GMRESMLMGT;
99
101template <typename MF>
103{
104public:
105
106 template <typename T> friend class MLMGT;
107 template <typename T> friend class MLCGSolverT;
108 template <typename T> friend class MLPoissonT;
109 template <typename T> friend class MLABecLaplacianT;
110 template <typename T> friend class GMRESMLMGT;
111
112 using MFType = MF;
115
120
121 MLLinOpT () = default;
122 virtual ~MLLinOpT () = default;
123
124 MLLinOpT (const MLLinOpT<MF>&) = delete;
125 MLLinOpT (MLLinOpT<MF>&&) = delete;
128
129 void define (const Vector<Geometry>& a_geom,
130 const Vector<BoxArray>& a_grids,
131 const Vector<DistributionMapping>& a_dmap,
132 const LPInfo& a_info,
133 const Vector<FabFactory<FAB> const*>& a_factory,
134 bool eb_limit_coarsening = true);
135
136 [[nodiscard]] virtual std::string name () const { return std::string("Unspecified"); }
137
149 const Array<BCType,AMREX_SPACEDIM>& hibc) noexcept;
150
162
173 const Array<Real,AMREX_SPACEDIM>& hi_bcloc) noexcept;
174
182 [[nodiscard]] bool needsCoarseDataForBC () const noexcept { return m_needs_coarse_data_for_bc; }
183
205 void setCoarseFineBC (const MF* crse, int crse_ratio,
206 LinOpBCType bc_type = LinOpBCType::Dirichlet) noexcept;
207
208 void setCoarseFineBC (const MF* crse, IntVect const& crse_ratio,
209 LinOpBCType bc_type = LinOpBCType::Dirichlet) noexcept;
210
211 template <typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,int> = 0>
212 void setCoarseFineBC (const AMF* crse, int crse_ratio,
213 LinOpBCType bc_type = LinOpBCType::Dirichlet) noexcept;
214
215 template <typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,int> = 0>
216 void setCoarseFineBC (const AMF* crse, IntVect const& crse_ratio,
217 LinOpBCType bc_type = LinOpBCType::Dirichlet) noexcept;
218
237 virtual void setLevelBC (int /*amrlev*/, const MF* /*levelbcdata*/,
238 const MF* /*robinbc_a*/ = nullptr,
239 const MF* /*robinbc_b*/ = nullptr,
240 const MF* /*robinbc_f*/ = nullptr) = 0;
241
242 template <typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF> &&
243 IsMultiFabLike_v<AMF>, int> = 0>
244 void setLevelBC (int amrlev, const AMF* levelbcdata,
245 const AMF* robinbc_a = nullptr,
246 const AMF* robinbc_b = nullptr,
247 const AMF* robinbc_f = nullptr);
248
250 void setVerbose (int v) noexcept { verbose = v; }
251
253 void setMaxOrder (int o) noexcept { maxorder = o; }
255 [[nodiscard]] int getMaxOrder () const noexcept { return maxorder; }
256
262 [[nodiscard]] bool getEnforceSingularSolvable () const noexcept { return enforceSingularSolvable; }
263
264 [[nodiscard]] virtual BottomSolver getDefaultBottomSolver () const { return BottomSolver::bicgstab; }
265
267 [[nodiscard]] virtual int getNComp () const { return 1; }
268
269 [[nodiscard]] virtual int getNGrow (int /*a_lev*/ = 0, int /*mg_lev*/ = 0) const { return 0; }
270
272 [[nodiscard]] virtual bool needsUpdate () const { return false; }
274 virtual void update () {}
275
285 virtual void restriction (int amrlev, int cmglev, MF& crse, MF& fine) const = 0;
286
295 virtual void interpolation (int amrlev, int fmglev, MF& fine, const MF& crse) const = 0;
296
305 virtual void interpAssign (int amrlev, int fmglev, MF& fine, MF& crse) const
306 {
307 amrex::ignore_unused(amrlev, fmglev, fine, crse);
308 amrex::Abort("MLLinOpT::interpAssign: Must be implemented for FMG cycle");
309 }
310
319 virtual void interpolationAmr (int famrlev, MF& fine, const MF& crse,
320 IntVect const& nghost) const
321 {
322 amrex::ignore_unused(famrlev, fine, crse, nghost);
323 amrex::Abort("MLLinOpT::interpolationAmr: Must be implemented for composite solves across multiple AMR levels");
324 }
325
335 virtual void averageDownSolutionRHS (int camrlev, MF& crse_sol, MF& crse_rhs,
336 const MF& fine_sol, const MF& fine_rhs)
337 {
338 amrex::ignore_unused(camrlev, crse_sol, crse_rhs, fine_sol, fine_rhs);
339 amrex::Abort("MLLinOpT::averageDownSolutionRHS: Must be implemented for composite solves across multiple AMR levels");
340 }
341
353 virtual void apply (int amrlev, int mglev, MF& out, MF& in, BCMode bc_mode,
354 StateMode s_mode, const MLMGBndryT<MF>* bndry=nullptr) const = 0;
355
366 virtual void smooth (int amrlev, int mglev, MF& sol, const MF& rhs,
367 bool skip_fillboundary, int niter) const = 0;
368
370 virtual void normalize (int amrlev, int mglev, MF& mf) const {
371 amrex::ignore_unused(amrlev, mglev, mf);
372 }
373
383 virtual void solutionResidual (int amrlev, MF& resid, MF& x, const MF& b,
384 const MF* crse_bcdata=nullptr) = 0;
385
386 virtual void prepareForFluxes (int /*amrlev*/, const MF* /*crse_bcdata*/ = nullptr) {}
387
399 virtual void correctionResidual (int amrlev, int mglev, MF& resid, MF& x, const MF& b,
400 BCMode bc_mode, const MF* crse_bcdata=nullptr) = 0;
401
413 virtual void reflux (int crse_amrlev,
414 MF& res, const MF& crse_sol, const MF& crse_rhs,
415 MF& fine_res, MF& fine_sol, const MF& fine_rhs) const
416 {
417 amrex::ignore_unused(crse_amrlev, res, crse_sol, crse_rhs, fine_res,
418 fine_sol, fine_rhs);
419 amrex::Abort("MLLinOpT::reflux: Must be implemented for composite solves across multiple AMR levels");
420 }
421
430 virtual void compFlux (int amrlev, const Array<MF*,AMREX_SPACEDIM>& fluxes,
431 MF& sol, Location loc) const
432 {
433 amrex::ignore_unused(amrlev, fluxes, sol, loc);
434 amrex::Abort("AMReX_MLLinOp::compFlux::How did we get here?");
435 }
436
445 virtual void compGrad (int amrlev, const Array<MF*,AMREX_SPACEDIM>& grad,
446 MF& sol, Location loc) const
447 {
448 amrex::ignore_unused(amrlev, grad, sol, loc);
449 amrex::Abort("AMReX_MLLinOp::compGrad::How did we get here?");
450 }
451
453 virtual void applyMetricTerm (int /*amrlev*/, int /*mglev*/, MF& /*rhs*/) const {}
455 virtual void unapplyMetricTerm (int /*amrlev*/, int /*mglev*/, MF& /*rhs*/) const {}
456
458 virtual void unimposeNeumannBC (int /*amrlev*/, MF& /*rhs*/) const {}
459
461 virtual void applyInhomogNeumannTerm (int /*amrlev*/, MF& /*rhs*/) const {}
462
464 virtual void applyOverset (int /*amlev*/, MF& /*rhs*/) const {}
465
467 [[nodiscard]] virtual bool scaleRHS (int /*amrlev*/, MF* /*rhs*/) const {
468 return false;
469 }
470
472 virtual Vector<RT> getSolvabilityOffset (int /*amrlev*/, int /*mglev*/,
473 MF const& /*rhs*/) const { return {}; }
474
476 virtual void fixSolvabilityByOffset (int /*amrlev*/, int /*mglev*/, MF& /*rhs*/,
477 Vector<RT> const& /*offset*/) const {}
478
479 virtual void prepareForSolve () = 0;
480
481 virtual void preparePrecond () {}
482
485 virtual void setDirichletNodesToZero (int /*amrlev*/, int /*mglev*/,
486 MF& /*mf*/) const
487 {
488 amrex::Warning("This function might need to be implemented for GMRES to work with this LinOp.");
489 }
490
492 [[nodiscard]] virtual bool isSingular (int amrlev) const = 0;
494 [[nodiscard]] virtual bool isBottomSingular () const = 0;
495
497 virtual RT xdoty (int amrlev, int mglev, const MF& x, const MF& y, bool local) const = 0;
498
500
501 virtual RT norm2Precond (Vector<MF const*> const& x) const;
502
503 virtual std::unique_ptr<MLLinOpT<MF>> makeNLinOp (int /*grid_size*/) const
504 {
505 amrex::Abort("MLLinOp::makeNLinOp: N-Solve not supported");
506 return nullptr;
507 }
508
509 virtual void getFluxes (const Vector<Array<MF*,AMREX_SPACEDIM> >& /*a_flux*/,
510 const Vector<MF*>& /*a_sol*/,
511 Location /*a_loc*/) const {
512 amrex::Abort("MLLinOp::getFluxes: How did we get here?");
513 }
514 virtual void getFluxes (const Vector<MF*>& /*a_flux*/,
515 const Vector<MF*>& /*a_sol*/) const {
516 amrex::Abort("MLLinOp::getFluxes: How did we get here?");
517 }
518
519#ifdef AMREX_USE_EB
520 virtual void getEBFluxes (const Vector<MF*>& /*a_flux*/,
521 const Vector<MF*>& /*a_sol*/) const {
522 amrex::Abort("MLLinOp::getEBFluxes: How did we get here?");
523 }
524#endif
525
526#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
527 [[nodiscard]] virtual std::unique_ptr<Hypre> makeHypre (Hypre::Interface /*hypre_interface*/) const {
528 amrex::Abort("MLLinOp::makeHypre: How did we get here?");
529 return {nullptr};
530 }
531 [[nodiscard]] virtual std::unique_ptr<HypreNodeLap> makeHypreNodeLap(
532 int /*bottom_verbose*/,
533 const std::string& /* options_namespace */) const
534 {
535 amrex::Abort("MLLinOp::makeHypreNodeLap: How did we get here?");
536 return {nullptr};
537 }
538#endif
539
540#if defined(AMREX_USE_PETSC) && (AMREX_SPACEDIM > 1)
541 [[nodiscard]] virtual std::unique_ptr<PETScABecLap> makePETSc () const {
542 amrex::Abort("MLLinOp::makePETSc: How did we get here?");
543 return {nullptr};
544 }
545#endif
546
547 [[nodiscard]] virtual bool supportNSolve () const { return false; }
548
549 virtual void copyNSolveSolution (MF&, MF const&) const {}
550
551 virtual void postSolve (Vector<MF*> const& /*sol*/) const {}
552
553 [[nodiscard]] virtual RT normInf (int amrlev, MF const& mf, bool local) const = 0;
554
555 virtual void averageDownAndSync (Vector<MF>& sol) const = 0;
556
557 virtual void avgDownResAmr (int clev, MF& cres, MF const& fres) const
558 {
559 amrex::ignore_unused(clev, cres, fres);
560 amrex::Abort("MLLinOpT::avgDownResAmr: Must be implemented for composite solves across multiple AMR levels");
561 }
562
563 // This function is needed for FMG cycle, but not V-cycle.
564 virtual void avgDownResMG (int clev, MF& cres, MF const& fres) const;
565
566 virtual void beginPrecondBC () { m_precond_mode = true; }
567 virtual void endPrecondBC () { m_precond_mode = false; }
568
569 [[nodiscard]] bool isMFIterSafe (int amrlev, int mglev1, int mglev2) const;
570
572 [[nodiscard]] int NAMRLevels () const noexcept { return m_num_amr_levels; }
573
575 [[nodiscard]] int NMGLevels (int amrlev) const noexcept { return m_num_mg_levels[amrlev]; }
576
577 [[nodiscard]] const Geometry& Geom (int amr_lev, int mglev=0) const noexcept { return m_geom[amr_lev][mglev]; }
578
579 // BC
582 // Need to save the original copy because we change the BC type to
583 // Neumann for inhomogeneous Neumann and Robin.
586
587protected:
588
589 static constexpr int mg_coarsen_ratio = 2;
590 static constexpr int mg_box_min_width = 2;
592
594
595 int verbose = 0;
596
597 int maxorder = 3;
598
600
603
605 const MLLinOpT<MF>* m_parent = nullptr;
606
608
609 bool m_do_agglomeration = false;
610 bool m_do_consolidation = false;
611
614
621
625 struct CommContainer {
626 MPI_Comm comm;
627 CommContainer (MPI_Comm m) noexcept : comm(m) {}
628 CommContainer (const CommContainer&) = delete;
629 CommContainer (CommContainer&&) = delete;
630 void operator= (const CommContainer&) = delete;
631 void operator= (CommContainer&&) = delete;
632 ~CommContainer () { // NOLINT(modernize-use-equals-default)
633#ifdef BL_USE_MPI
634 if (comm != MPI_COMM_NULL) { MPI_Comm_free(&comm); }
635#endif
636 }
637 };
639 std::unique_ptr<CommContainer> m_raii_comm;
640
643
648 const MF* m_coarse_data_for_bc = nullptr;
650
651 bool m_precond_mode = false;
652
654 [[nodiscard]] const Vector<int>& AMRRefRatio () const noexcept { return m_amr_ref_ratio; }
655
657 [[nodiscard]] int AMRRefRatio (int amr_lev) const noexcept { return m_amr_ref_ratio[amr_lev]; }
658
660 [[nodiscard]] IntVect AMRRefRatioVect (int amr_lev) const noexcept {
661 IntVect rr(m_amr_ref_ratio[amr_lev]);
662 if (info.hasHiddenDimension()) { rr[info.hidden_direction] = 1; }
663 return rr;
664 }
665
666 [[nodiscard]] FabFactory<FAB> const* Factory (int amr_lev, int mglev=0) const noexcept {
667 return m_factory[amr_lev][mglev].get();
668 }
669
670 [[nodiscard]] GpuArray<BCType,AMREX_SPACEDIM> LoBC (int icomp = 0) const noexcept {
672 m_lobc[icomp][1],
673 m_lobc[icomp][2])}};
674 }
675 [[nodiscard]] GpuArray<BCType,AMREX_SPACEDIM> HiBC (int icomp = 0) const noexcept {
677 m_hibc[icomp][1],
678 m_hibc[icomp][2])}};
679 }
680
681 [[nodiscard]] bool hasBC (BCType bct) const noexcept;
682 [[nodiscard]] bool hasInhomogNeumannBC () const noexcept;
683 [[nodiscard]] bool hasRobinBC () const noexcept;
684
685 [[nodiscard]] virtual bool supportRobinBC () const noexcept { return false; }
686 [[nodiscard]] virtual bool supportInhomogNeumannBC () const noexcept { return false; }
687
688#ifdef BL_USE_MPI
689 [[nodiscard]] bool isBottomActive () const noexcept { return m_bottom_comm != MPI_COMM_NULL; }
690#else
691 [[nodiscard]] bool isBottomActive () const noexcept { return true; }
692#endif
693 [[nodiscard]] MPI_Comm BottomCommunicator () const noexcept { return m_bottom_comm; }
694 [[nodiscard]] MPI_Comm Communicator () const noexcept { return m_default_comm; }
695
696 void setCoarseFineBCLocation (const RealVect& cloc) noexcept { m_coarse_bc_loc = cloc; }
697
698 [[nodiscard]] bool doAgglomeration () const noexcept { return m_do_agglomeration; }
699 [[nodiscard]] bool doConsolidation () const noexcept { return m_do_consolidation; }
700 [[nodiscard]] bool doSemicoarsening () const noexcept { return m_do_semicoarsening; }
701
702 [[nodiscard]] bool isCellCentered () const noexcept { return m_ixtype == 0; }
703
704 [[nodiscard]] virtual IntVect getNGrowVectRestriction () const {
705 return isCellCentered() ? IntVect(0) : IntVect(1);
706 }
707
708 virtual void make (Vector<Vector<MF> >& mf, IntVect const& ng) const;
709
710 [[nodiscard]] virtual MF make (int amrlev, int mglev, IntVect const& ng) const;
711
712 [[nodiscard]] virtual MF makeAlias (MF const& mf) const;
713
714 [[nodiscard]] virtual MF makeCoarseMG (int amrlev, int mglev, IntVect const& ng) const;
715
716 [[nodiscard]] virtual MF makeCoarseAmr (int famrlev, IntVect const& ng) const;
717
718 [[nodiscard]] virtual std::unique_ptr<FabFactory<FAB> > makeFactory (int /*amrlev*/, int /*mglev*/) const {
719 return std::make_unique<DefaultFabFactory<FAB>>();
720 }
721
722 virtual void resizeMultiGrid (int new_size);
723
724 [[nodiscard]] bool hasHiddenDimension () const noexcept { return info.hasHiddenDimension(); }
725 [[nodiscard]] int hiddenDirection () const noexcept { return info.hidden_direction; }
726 [[nodiscard]] Box compactify (Box const& b) const noexcept;
727
728 template <typename T>
729 [[nodiscard]] Array4<T> compactify (Array4<T> const& a) const noexcept
730 {
731 if (info.hidden_direction == 0) {
732 return Array4<T>(a.dataPtr(), {a.begin[1],a.begin[2],0}, {a.end[1],a.end[2],1}, a.nComp());
733 } else if (info.hidden_direction == 1) {
734 return Array4<T>(a.dataPtr(), {a.begin[0],a.begin[2],0}, {a.end[0],a.end[2],1}, a.nComp());
735 } else if (info.hidden_direction == 2) {
736 return Array4<T>(a.dataPtr(), {a.begin[0],a.begin[1],0}, {a.end[0],a.end[1],1}, a.nComp());
737 } else {
738 return a;
739 }
740 }
741
742 template <typename T>
743 [[nodiscard]] T get_d0 (T const& dx, T const& dy, T const&) const noexcept
744 {
745 if (info.hidden_direction == 0) {
746 return dy;
747 } else {
748 return dx;
749 }
750 }
751
752 template <typename T>
753 [[nodiscard]] T get_d1 (T const&, T const& dy, T const& dz) const noexcept
754 {
755 if (info.hidden_direction == 0 || info.hidden_direction == 1) {
756 return dz;
757 } else {
758 return dy;
759 }
760 }
761
762private:
763
764 void defineGrids (const Vector<Geometry>& a_geom,
765 const Vector<BoxArray>& a_grids,
766 const Vector<DistributionMapping>& a_dmap,
767 const Vector<FabFactory<FAB> const*>& a_factory);
768 void defineBC ();
769 static void makeAgglomeratedDMap (const Vector<BoxArray>& ba, Vector<DistributionMapping>& dm);
770 static void makeConsolidatedDMap (const Vector<BoxArray>& ba, Vector<DistributionMapping>& dm,
771 int ratio, int strategy);
772 [[nodiscard]] MPI_Comm makeSubCommunicator (const DistributionMapping& dm);
773
774 virtual void checkPoint (std::string const& /*file_name*/) const {
775 amrex::Abort("MLLinOp:checkPoint: not implemented");
776 }
777
778 Vector<std::unique_ptr<MF>> levelbc_raii;
779 Vector<std::unique_ptr<MF>> robin_a_raii;
780 Vector<std::unique_ptr<MF>> robin_b_raii;
781 Vector<std::unique_ptr<MF>> robin_f_raii;
782};
783
784template <typename MF>
785void
787 const Vector<BoxArray>& a_grids,
788 const Vector<DistributionMapping>& a_dmap,
789 const LPInfo& a_info,
790 const Vector<FabFactory<FAB> const*>& a_factory,
791 [[maybe_unused]] bool eb_limit_coarsening)
792{
793 BL_PROFILE("MLLinOp::define()");
794
795 info = a_info;
796#ifdef AMREX_USE_GPU
798 {
799 if (info.agg_grid_size <= 0) { info.agg_grid_size = AMREX_D_PICK(32, 16, 8); }
800 if (info.con_grid_size <= 0) { info.con_grid_size = AMREX_D_PICK(32, 16, 8); }
801 }
802 else
803#endif
804 {
805 if (info.agg_grid_size <= 0) { info.agg_grid_size = LPInfo::getDefaultAgglomerationGridSize(); }
806 if (info.con_grid_size <= 0) { info.con_grid_size = LPInfo::getDefaultConsolidationGridSize(); }
807 }
808
809#ifdef AMREX_USE_EB
810 if (!a_factory.empty() && eb_limit_coarsening) {
811 const auto *f = dynamic_cast<EBFArrayBoxFactory const*>(a_factory[0]);
812 if (f) {
813 info.max_coarsening_level = std::min(info.max_coarsening_level,
814 f->maxCoarseningLevel());
815 }
816 }
817#endif
818 defineGrids(a_geom, a_grids, a_dmap, a_factory);
819 defineBC();
820}
821
822template <typename MF>
823void
825 const Vector<BoxArray>& a_grids,
826 const Vector<DistributionMapping>& a_dmap,
827 const Vector<FabFactory<FAB> const*>& a_factory)
828{
829 BL_PROFILE("MLLinOp::defineGrids()");
830
831#ifdef AMREX_USE_EB
832 if ( ! a_factory.empty() ) {
833 auto const* ebf = dynamic_cast<EBFArrayBoxFactory const*>(a_factory[0]);
834 if (ebf && !(ebf->isAllRegular())) { // Has non-trivial EB
835 mg_domain_min_width = 4;
836 }
837 }
838#endif
839
840 m_num_amr_levels = 0;
841 for (int amrlev = 0; amrlev < std::ssize(a_geom); amrlev++) {
842 if (!a_grids[amrlev].empty()) {
843 m_num_amr_levels++;
844 }
845 }
846
847 m_amr_ref_ratio.resize(m_num_amr_levels);
848 m_num_mg_levels.resize(m_num_amr_levels);
849
850 m_geom.resize(m_num_amr_levels);
851 m_grids.resize(m_num_amr_levels);
852 m_dmap.resize(m_num_amr_levels);
853 m_factory.resize(m_num_amr_levels);
854
855 m_default_comm = ParallelContext::CommunicatorSub();
856
857 const RealBox& rb = a_geom[0].ProbDomain();
858 const int coord = a_geom[0].Coord();
859 const Array<int,AMREX_SPACEDIM>& is_per = a_geom[0].isPeriodic();
860
861 IntVect mg_coarsen_ratio_v(mg_coarsen_ratio);
862 IntVect mg_box_min_width_v(mg_box_min_width);
863 IntVect mg_domain_min_width_v(mg_domain_min_width);
864 if (hasHiddenDimension()) {
865 AMREX_ASSERT_WITH_MESSAGE(AMREX_SPACEDIM == 3,
866 "Hidden direction only supported for 3d");
867 mg_coarsen_ratio_v[info.hidden_direction] = 1;
868 mg_box_min_width_v[info.hidden_direction] = 0;
869 mg_domain_min_width_v[info.hidden_direction] = 0;
870 }
871
872 // fine amr levels
873 for (int amrlev = m_num_amr_levels-1; amrlev > 0; --amrlev)
874 {
875 m_num_mg_levels[amrlev] = 1;
876 m_geom[amrlev].push_back(a_geom[amrlev]);
877 m_grids[amrlev].push_back(a_grids[amrlev]);
878 m_dmap[amrlev].push_back(a_dmap[amrlev]);
879 if (amrlev < std::ssize(a_factory)) {
880 m_factory[amrlev].emplace_back(a_factory[amrlev]->clone());
881 } else {
882 m_factory[amrlev].push_back(std::make_unique<DefaultFabFactory<FAB>>());
883 }
884
885 IntVect rr = mg_coarsen_ratio_v;
886 const Box& dom = a_geom[amrlev].Domain();
887 for (int i = 0; i < 2; ++i)
888 {
889 if (!dom.coarsenable(rr)) { amrex::Abort("MLLinOp: Uncoarsenable domain"); }
890
891 const Box& cdom = amrex::coarsen(dom,rr);
892 if (cdom == a_geom[amrlev-1].Domain()) { break; }
893
894 ++(m_num_mg_levels[amrlev]);
895
896 m_geom[amrlev].emplace_back(cdom, rb, coord, is_per);
897
898 m_grids[amrlev].push_back(a_grids[amrlev]);
899 AMREX_ASSERT(m_grids[amrlev].back().coarsenable(rr));
900 m_grids[amrlev].back().coarsen(rr);
901
902 m_dmap[amrlev].push_back(a_dmap[amrlev]);
903
904 rr *= mg_coarsen_ratio_v;
905 }
906
907#if (AMREX_SPACEDIM > 1)
908 if (hasHiddenDimension()) {
909 m_amr_ref_ratio[amrlev-1] = rr[(info.hidden_direction+1) % AMREX_SPACEDIM];
910 } else
911#endif
912 {
913 m_amr_ref_ratio[amrlev-1] = rr[0];
914 }
915 }
916
917 // coarsest amr level
918 m_num_mg_levels[0] = 1;
919 m_geom[0].push_back(a_geom[0]);
920 m_grids[0].push_back(a_grids[0]);
921 m_dmap[0].push_back(a_dmap[0]);
922 if (!a_factory.empty()) {
923 m_factory[0].emplace_back(a_factory[0]->clone());
924 } else {
925 m_factory[0].push_back(std::make_unique<DefaultFabFactory<FAB>>());
926 }
927
928 m_domain_covered.resize(m_num_amr_levels, false);
929 auto npts0 = m_grids[0][0].numPts();
930 m_domain_covered[0] = (npts0 == compactify(m_geom[0][0].Domain()).numPts());
931 for (int amrlev = 1; amrlev < m_num_amr_levels; ++amrlev)
932 {
933 if (!m_domain_covered[amrlev-1]) { break; }
934 m_domain_covered[amrlev] = (m_grids[amrlev][0].numPts() ==
935 compactify(m_geom[amrlev][0].Domain()).numPts());
936 }
937
938 Box aggbox;
939 bool aggable = false;
940
941 if (m_grids[0][0].size() > 1 && info.do_agglomeration)
942 {
943 if (m_domain_covered[0])
944 {
945 aggbox = m_geom[0][0].Domain();
946 if (hasHiddenDimension()) {
947 aggbox.makeSlab(hiddenDirection(), m_grids[0][0][0].smallEnd(hiddenDirection()));
948 }
949 aggable = true;
950 }
951 else
952 {
953 aggbox = m_grids[0][0].minimalBox();
954 aggable = (aggbox.numPts() == npts0);
955 }
956 }
957
958 bool agged = false;
959 bool coned = false;
960 int agg_lev = 0, con_lev = 0;
961
962 AMREX_ALWAYS_ASSERT( ! (info.do_semicoarsening && info.hasHiddenDimension())
963 && info.semicoarsening_direction >= -1
964 && info.semicoarsening_direction < AMREX_SPACEDIM );
965
966 if (info.do_agglomeration && aggable)
967 {
968 Box dbx = m_geom[0][0].Domain();
969 Box bbx = aggbox;
970 Real const nbxs = static_cast<Real>(m_grids[0][0].size());
971 Real const threshold_npts = static_cast<Real>(AMREX_D_TERM(info.agg_grid_size,
972 *info.agg_grid_size,
973 *info.agg_grid_size));
974 Vector<Box> domainboxes{dbx};
975 Vector<Box> boundboxes{bbx};
976 Vector<int> agg_flag{false};
977 Vector<IntVect> accum_coarsen_ratio{IntVect(1)};
978 int numsclevs = 0;
979
980 for (int lev = 0; lev < info.max_coarsening_level; ++lev)
981 {
982 IntVect rr_level = mg_coarsen_ratio_v;
983 bool const do_semicoarsening_level = info.do_semicoarsening
984 && numsclevs < info.max_semicoarsening_level;
985 if (do_semicoarsening_level
986 && info.semicoarsening_direction != -1)
987 {
988 rr_level[info.semicoarsening_direction] = 1;
989 }
990 IntVect is_coarsenable;
991 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
992 IntVect rr_dir(1);
993 rr_dir[idim] = rr_level[idim];
994 is_coarsenable[idim] = dbx.coarsenable(rr_dir, mg_domain_min_width_v)
995 && bbx.coarsenable(rr_dir, mg_box_min_width_v);
996 if (!is_coarsenable[idim] && do_semicoarsening_level
997 && info.semicoarsening_direction == -1)
998 {
999 is_coarsenable[idim] = true;
1000 rr_level[idim] = 1;
1001 }
1002 }
1003 if (is_coarsenable != IntVect(1) || rr_level == IntVect(1)) {
1004 break;
1005 }
1006 if (do_semicoarsening_level && info.semicoarsening_direction == -1) {
1007 // make sure there is at most one direction that is not coarsened
1008 int n_ones = AMREX_D_TERM( static_cast<int>(rr_level[0] == 1),
1009 + static_cast<int>(rr_level[1] == 1),
1010 + static_cast<int>(rr_level[2] == 1));
1011 if (n_ones > 1) { break; }
1012 }
1013 if (rr_level != mg_coarsen_ratio_v) {
1014 ++numsclevs;
1015 }
1016
1017 accum_coarsen_ratio.push_back(accum_coarsen_ratio.back()*rr_level);
1018 domainboxes.push_back(dbx.coarsen(rr_level));
1019 boundboxes.push_back(bbx.coarsen(rr_level));
1020 bool to_agg = (bbx.d_numPts() / nbxs) < 0.999*threshold_npts;
1021 agg_flag.push_back(to_agg);
1022 }
1023
1024 for (int lev = 1, nlevs = static_cast<int>(domainboxes.size()); lev < nlevs; ++lev) {
1025 if (!agged && !agg_flag[lev] &&
1026 a_grids[0].coarsenable(accum_coarsen_ratio[lev], mg_box_min_width_v))
1027 {
1028 m_grids[0].push_back(amrex::coarsen(a_grids[0], accum_coarsen_ratio[lev]));
1029 m_dmap[0].push_back(a_dmap[0]);
1030 } else {
1031 IntVect cr = domainboxes[lev-1].length() / domainboxes[lev].length();
1032 if (!m_grids[0].back().coarsenable(cr)) {
1033 break; // average_down would fail if fine boxarray is not coarsenable.
1034 }
1035 m_grids[0].emplace_back(boundboxes[lev]);
1036 IntVect max_grid_size(info.agg_grid_size);
1037 if (info.do_semicoarsening && info.max_semicoarsening_level >= lev
1038 && info.semicoarsening_direction != -1)
1039 {
1040 IntVect blen = amrex::enclosedCells(boundboxes[lev]).size();
1041 AMREX_D_TERM(int mgs_0 = (max_grid_size[0]+blen[0]-1) / blen[0];,
1042 int mgs_1 = (max_grid_size[1]+blen[1]-1) / blen[1];,
1043 int mgs_2 = (max_grid_size[2]+blen[2]-1) / blen[2]);
1044 max_grid_size[info.semicoarsening_direction]
1045 *= AMREX_D_TERM(mgs_0, *mgs_1, *mgs_2);
1046 }
1047 m_grids[0].back().maxSize(max_grid_size);
1048 m_dmap[0].push_back(DistributionMapping());
1049 if (!agged) {
1050 agged = true;
1051 agg_lev = lev;
1052 }
1053 }
1054 m_geom[0].emplace_back(domainboxes[lev],rb,coord,is_per);
1055 }
1056 }
1057 else
1058 {
1059 Long consolidation_threshold = 0;
1060 Real avg_npts = 0.0;
1061 if (info.do_consolidation) {
1062 avg_npts = static_cast<Real>(a_grids[0].d_numPts()) / static_cast<Real>(ParallelContext::NProcsSub());
1063 consolidation_threshold = AMREX_D_TERM(Long(info.con_grid_size),
1064 *info.con_grid_size,
1065 *info.con_grid_size);
1066 }
1067
1068 Box const& dom0 = a_geom[0].Domain();
1069 IntVect rr_vec(1);
1070 int numsclevs = 0;
1071 for (int lev = 0; lev < info.max_coarsening_level; ++lev)
1072 {
1073 IntVect rr_level = mg_coarsen_ratio_v;
1074 bool do_semicoarsening_level = info.do_semicoarsening
1075 && numsclevs < info.max_semicoarsening_level;
1076 if (do_semicoarsening_level
1077 && info.semicoarsening_direction != -1)
1078 {
1079 rr_level[info.semicoarsening_direction] = 1;
1080 }
1081 IntVect is_coarsenable;
1082 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1083 IntVect rr_dir(1);
1084 rr_dir[idim] = rr_vec[idim] * rr_level[idim];
1085 is_coarsenable[idim] = dom0.coarsenable(rr_dir, mg_domain_min_width_v)
1086 && a_grids[0].coarsenable(rr_dir, mg_box_min_width_v);
1087 if (!is_coarsenable[idim] && do_semicoarsening_level
1088 && info.semicoarsening_direction == -1)
1089 {
1090 is_coarsenable[idim] = true;
1091 rr_level[idim] = 1;
1092 }
1093 }
1094 if (is_coarsenable != IntVect(1) || rr_level == IntVect(1)) {
1095 break;
1096 }
1097 if (do_semicoarsening_level && info.semicoarsening_direction == -1) {
1098 // make sure there is at most one direction that is not coarsened
1099 int n_ones = AMREX_D_TERM( static_cast<int>(rr_level[0] == 1),
1100 + static_cast<int>(rr_level[1] == 1),
1101 + static_cast<int>(rr_level[2] == 1));
1102 if (n_ones > 1) { break; }
1103 }
1104 if (rr_level != mg_coarsen_ratio_v) {
1105 ++numsclevs;
1106 }
1107 rr_vec *= rr_level;
1108
1109 m_geom[0].emplace_back(amrex::coarsen(dom0, rr_vec), rb, coord, is_per);
1110 m_grids[0].push_back(amrex::coarsen(a_grids[0], rr_vec));
1111
1112 if (info.do_consolidation)
1113 {
1114 if (avg_npts/static_cast<Real>(AMREX_D_TERM(rr_vec[0], *rr_vec[1], *rr_vec[2]))
1115 < Real(0.999)*static_cast<Real>(consolidation_threshold))
1116 {
1117 coned = true;
1118 con_lev = m_dmap[0].size();
1119 m_dmap[0].push_back(DistributionMapping());
1120 }
1121 else
1122 {
1123 m_dmap[0].push_back(m_dmap[0].back());
1124 }
1125 }
1126 else
1127 {
1128 m_dmap[0].push_back(a_dmap[0]);
1129 }
1130 }
1131 }
1132
1133 m_num_mg_levels[0] = m_grids[0].size();
1134
1135 for (int mglev = 0; mglev < m_num_mg_levels[0] - 1; mglev++){
1136 const Box& fine_domain = m_geom[0][mglev].Domain();
1137 const Box& crse_domain = m_geom[0][mglev+1].Domain();
1138 mg_coarsen_ratio_vec.push_back(fine_domain.length()/crse_domain.length());
1139 }
1140
1141 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev) {
1142 if (AMRRefRatio(amrlev) == 4 && mg_coarsen_ratio_vec.empty()) {
1143 mg_coarsen_ratio_vec.push_back(IntVect(2));
1144 }
1145 }
1146
1147 if (agged)
1148 {
1149 makeAgglomeratedDMap(m_grids[0], m_dmap[0]);
1150 }
1151 else if (coned)
1152 {
1153 makeConsolidatedDMap(m_grids[0], m_dmap[0], info.con_ratio, info.con_strategy);
1154 }
1155
1156 if (agged || coned)
1157 {
1158 m_bottom_comm = makeSubCommunicator(m_dmap[0].back());
1159 }
1160 else
1161 {
1162 m_bottom_comm = m_default_comm;
1163 }
1164
1165 m_do_agglomeration = agged;
1166 m_do_consolidation = coned;
1167
1168 if (verbose > 1) {
1169 if (agged) {
1170 Print() << "MLLinOp::defineGrids(): agglomerated AMR level 0 starting at MG level "
1171 << agg_lev << " of " << m_num_mg_levels[0] << "\n";
1172 } else if (coned) {
1173 Print() << "MLLinOp::defineGrids(): consolidated AMR level 0 starting at MG level "
1174 << con_lev << " of " << m_num_mg_levels[0]
1175 << " (ratio = " << info.con_ratio << ")" << "\n";
1176 } else {
1177 Print() << "MLLinOp::defineGrids(): no agglomeration or consolidation of AMR level 0\n";
1178 }
1179 }
1180
1181 for (int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
1182 {
1183 for (int mglev = 1; mglev < m_num_mg_levels[amrlev]; ++mglev)
1184 {
1185 m_factory[amrlev].emplace_back(makeFactory(amrlev,mglev));
1186 }
1187 }
1188
1189 for (int amrlev = 1; amrlev < m_num_amr_levels; ++amrlev)
1190 {
1191 AMREX_ASSERT_WITH_MESSAGE(m_grids[amrlev][0].coarsenable(AMRRefRatioVect(amrlev-1)),
1192 "MLLinOp: grids not coarsenable between AMR levels");
1193 }
1194}
1195
1196template <typename MF>
1197void
1198MLLinOpT<MF>::defineBC ()
1199{
1200 m_needs_coarse_data_for_bc = !m_domain_covered[0];
1201
1202 levelbc_raii.resize(m_num_amr_levels);
1203 robin_a_raii.resize(m_num_amr_levels);
1204 robin_b_raii.resize(m_num_amr_levels);
1205 robin_f_raii.resize(m_num_amr_levels);
1206}
1207
1208template <typename MF>
1209void
1211 const Array<BCType,AMREX_SPACEDIM>& a_hibc) noexcept
1212{
1213 const int ncomp = getNComp();
1214 setDomainBC(Vector<Array<BCType,AMREX_SPACEDIM> >(ncomp,a_lobc),
1215 Vector<Array<BCType,AMREX_SPACEDIM> >(ncomp,a_hibc));
1216}
1217
1218template <typename MF>
1219void
1221 const Vector<Array<BCType,AMREX_SPACEDIM> >& a_hibc)
1222{
1223 const int ncomp = getNComp();
1224 AMREX_ASSERT_WITH_MESSAGE(ncomp == a_lobc.size() && ncomp == a_hibc.size(),
1225 "MLLinOp::setDomainBC: wrong size");
1226 m_lobc = a_lobc;
1227 m_hibc = a_hibc;
1228 m_lobc_orig = m_lobc;
1229 m_hibc_orig = m_hibc;
1230 for (int icomp = 0; icomp < ncomp; ++icomp) {
1231 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1232 if (m_geom[0][0].isPeriodic(idim)) {
1233 AMREX_ALWAYS_ASSERT(m_lobc[icomp][idim] == BCType::Periodic &&
1234 m_hibc[icomp][idim] == BCType::Periodic);
1235 } else {
1236 AMREX_ALWAYS_ASSERT(m_lobc[icomp][idim] != BCType::Periodic &&
1237 m_hibc[icomp][idim] != BCType::Periodic);
1238 }
1239
1240 if (m_lobc[icomp][idim] == LinOpBCType::inhomogNeumann ||
1241 m_lobc[icomp][idim] == LinOpBCType::Robin)
1242 {
1243 m_lobc[icomp][idim] = LinOpBCType::Neumann;
1244 }
1245
1246 if (m_hibc[icomp][idim] == LinOpBCType::inhomogNeumann ||
1247 m_hibc[icomp][idim] == LinOpBCType::Robin)
1248 {
1249 m_hibc[icomp][idim] = LinOpBCType::Neumann;
1250 }
1251 }
1252 }
1253
1254 if (hasHiddenDimension()) {
1255 const int hd = hiddenDirection();
1256 for (int n = 0; n < ncomp; ++n) {
1257 m_lobc[n][hd] = LinOpBCType::Neumann;
1258 m_hibc[n][hd] = LinOpBCType::Neumann;
1259 }
1260 }
1261
1262 if (hasInhomogNeumannBC() && !supportInhomogNeumannBC()) {
1263 amrex::Abort("Inhomogeneous Neumann BC not supported");
1264 }
1265 if (hasRobinBC() && !supportRobinBC()) {
1266 amrex::Abort("Robin BC not supported");
1267 }
1268}
1269
1270template <typename MF>
1271bool
1272MLLinOpT<MF>::hasBC (BCType bct) const noexcept
1273{
1274 int ncomp = m_lobc_orig.size();
1275 for (int n = 0; n < ncomp; ++n) {
1276 for (int idim = 0; idim <AMREX_SPACEDIM; ++idim) {
1277 if (m_lobc_orig[n][idim] == bct || m_hibc_orig[n][idim] == bct) {
1278 return true;
1279 }
1280 }
1281 }
1282 return false;
1283}
1284
1285template <typename MF>
1286bool
1288{
1289 return hasBC(BCType::inhomogNeumann);
1290}
1291
1292template <typename MF>
1293bool
1295{
1296 return hasBC(BCType::Robin);
1297}
1298
1299template <typename MF>
1300Box
1301MLLinOpT<MF>::compactify (Box const& b) const noexcept
1302{
1303#if (AMREX_SPACEDIM == 3)
1304 if (info.hasHiddenDimension()) {
1305 const auto& lo = b.smallEnd();
1306 const auto& hi = b.bigEnd();
1307 if (info.hidden_direction == 0) {
1308 return Box(IntVect(lo[1],lo[2],0), IntVect(hi[1],hi[2],0), b.ixType());
1309 } else if (info.hidden_direction == 1) {
1310 return Box(IntVect(lo[0],lo[2],0), IntVect(hi[0],hi[2],0), b.ixType());
1311 } else {
1312 return Box(IntVect(lo[0],lo[1],0), IntVect(hi[0],hi[1],0), b.ixType());
1313 }
1314 } else
1315#endif
1316 {
1317 return b;
1318 }
1319}
1320
1321template <typename MF>
1322void
1325{
1326 BL_PROFILE("MLLinOp::makeAgglomeratedDMap");
1327
1328 BL_ASSERT(!dm[0].empty());
1329 for (int i = 1, N=static_cast<int>(ba.size()); i < N; ++i)
1330 {
1331 if (dm[i].empty())
1332 {
1333 const std::vector< std::vector<int> >& sfc = DistributionMapping::makeSFC(ba[i]);
1334
1335 const int nprocs = ParallelContext::NProcsSub();
1336 AMREX_ASSERT(std::ssize(sfc) == nprocs);
1337
1338 Vector<int> pmap(ba[i].size());
1339 for (int iproc = 0; iproc < nprocs; ++iproc) {
1340 int grank = ParallelContext::local_to_global_rank(iproc);
1341 for (int ibox : sfc[iproc]) {
1342 pmap[ibox] = grank;
1343 }
1344 }
1345 dm[i].define(std::move(pmap));
1346 }
1347 }
1348}
1349
1350template <typename MF>
1351void
1352MLLinOpT<MF>::makeConsolidatedDMap (const Vector<BoxArray>& ba,
1353 Vector<DistributionMapping>& dm,
1354 int ratio, int strategy)
1355{
1356 BL_PROFILE("MLLinOp::makeConsolidatedDMap()");
1357
1358 int factor = 1;
1359 BL_ASSERT(!dm[0].empty());
1360 for (int i = 1, N=static_cast<int>(ba.size()); i < N; ++i)
1361 {
1362 if (dm[i].empty())
1363 {
1364 factor *= ratio;
1365
1366 const int nprocs = ParallelContext::NProcsSub();
1367 const auto& pmap_fine = dm[i-1].ProcessorMap();
1368 Vector<int> pmap(pmap_fine.size());
1369 ParallelContext::global_to_local_rank(pmap.data(), pmap_fine.data(), static_cast<int>(pmap.size()));
1370 if (strategy == 1) {
1371 for (auto& x: pmap) {
1372 x /= ratio;
1373 }
1374 } else if (strategy == 2) {
1375 int nprocs_con = static_cast<int>(std::ceil(static_cast<Real>(nprocs)
1376 / static_cast<Real>(factor)));
1377 for (auto& x: pmap) {
1378 auto d = std::div(x,nprocs_con);
1379 x = d.rem;
1380 }
1381 } else if (strategy == 3) {
1382 if (factor == ratio) {
1383 const std::vector< std::vector<int> >& sfc = DistributionMapping::makeSFC(ba[i]);
1384 for (int iproc = 0; iproc < nprocs; ++iproc) {
1385 for (int ibox : sfc[iproc]) {
1386 pmap[ibox] = iproc;
1387 }
1388 }
1389 }
1390 for (auto& x: pmap) {
1391 x /= ratio;
1392 }
1393 }
1394
1396 dm[i].define(std::move(pmap));
1397 } else {
1398 Vector<int> pmap_g(pmap.size());
1399 ParallelContext::local_to_global_rank(pmap_g.data(), pmap.data(), static_cast<int>(pmap.size()));
1400 dm[i].define(std::move(pmap_g));
1401 }
1402 }
1403 }
1404}
1405
1406template <typename MF>
1408MLLinOpT<MF>::makeSubCommunicator (const DistributionMapping& dm)
1409{
1410 BL_PROFILE("MLLinOp::makeSubCommunicator()");
1411
1412#ifdef BL_USE_MPI
1413
1414 Vector<int> newgrp_ranks = dm.ProcessorMap();
1415 std::sort(newgrp_ranks.begin(), newgrp_ranks.end());
1416 auto last = std::unique(newgrp_ranks.begin(), newgrp_ranks.end());
1417 newgrp_ranks.erase(last, newgrp_ranks.end());
1418
1419 MPI_Comm newcomm;
1420 MPI_Group defgrp, newgrp;
1421 MPI_Comm_group(m_default_comm, &defgrp);
1423 MPI_Group_incl(defgrp, static_cast<int>(newgrp_ranks.size()), newgrp_ranks.data(), &newgrp);
1424 } else {
1425 Vector<int> local_newgrp_ranks(newgrp_ranks.size());
1426 ParallelContext::global_to_local_rank(local_newgrp_ranks.data(),
1427 newgrp_ranks.data(), static_cast<int>(newgrp_ranks.size()));
1428 MPI_Group_incl(defgrp, static_cast<int>(local_newgrp_ranks.size()), local_newgrp_ranks.data(), &newgrp);
1429 }
1430
1431 MPI_Comm_create(m_default_comm, newgrp, &newcomm);
1432
1433 m_raii_comm = std::make_unique<CommContainer>(newcomm);
1434
1435 MPI_Group_free(&defgrp);
1436 MPI_Group_free(&newgrp);
1437
1438 return newcomm;
1439#else
1441 return m_default_comm;
1442#endif
1443}
1444
1445template <typename MF>
1446void
1448 const Array<Real,AMREX_SPACEDIM>& hi_bcloc) noexcept
1449{
1450 m_domain_bloc_lo = lo_bcloc;
1451 m_domain_bloc_hi = hi_bcloc;
1452}
1453
1454template <typename MF>
1455void
1456MLLinOpT<MF>::setCoarseFineBC (const MF* crse, int crse_ratio,
1457 LinOpBCType bc_type) noexcept
1458{
1459 setCoarseFineBC(crse, IntVect(crse_ratio), bc_type);
1460}
1461
1462template <typename MF>
1463void
1464MLLinOpT<MF>::setCoarseFineBC (const MF* crse, IntVect const& crse_ratio,
1465 LinOpBCType bc_type) noexcept
1466{
1467 m_coarse_data_for_bc = crse;
1468 m_coarse_data_crse_ratio = crse_ratio;
1469 m_coarse_fine_bc_type = bc_type;
1470}
1471
1472template <typename MF>
1473template <typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,int>>
1474void
1475MLLinOpT<MF>::setCoarseFineBC (const AMF* crse, int crse_ratio,
1476 LinOpBCType bc_type) noexcept
1477{
1478 setCoarseFineBC(crse, IntVect(crse_ratio), bc_type);
1479}
1480
1481template <typename MF>
1482template <typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,int>>
1483void
1484MLLinOpT<MF>::setCoarseFineBC (const AMF* crse, IntVect const& crse_ratio,
1485 LinOpBCType bc_type) noexcept
1486{
1487 m_coarse_data_for_bc_raii = MF(crse->boxArray(), crse->DistributionMap(),
1488 crse->nComp(), crse->nGrowVect());
1489 m_coarse_data_for_bc_raii.LocalCopy(*crse, 0, 0, crse->nComp(),
1490 crse->nGrowVect());
1491 m_coarse_data_for_bc = &m_coarse_data_for_bc_raii;
1492 m_coarse_data_crse_ratio = crse_ratio;
1493 m_coarse_fine_bc_type = bc_type;
1494}
1495
1496template <typename MF>
1497void
1499{
1500 mf.clear();
1501 mf.resize(m_num_amr_levels);
1502 for (int alev = 0; alev < m_num_amr_levels; ++alev) {
1503 mf[alev].resize(m_num_mg_levels[alev]);
1504 for (int mlev = 0; mlev < m_num_mg_levels[alev]; ++mlev) {
1505 mf[alev][mlev] = make(alev, mlev, ng);
1506 }
1507 }
1508}
1509
1510template <typename MF>
1511MF
1512MLLinOpT<MF>::make (int amrlev, int mglev, IntVect const& ng) const
1513{
1514 if constexpr (IsMultiFabLike_v<MF>) {
1515 return MF(amrex::convert(m_grids[amrlev][mglev], m_ixtype),
1516 m_dmap[amrlev][mglev], getNComp(), ng, MFInfo(),
1517 *m_factory[amrlev][mglev]);
1518 } else {
1519 amrex::ignore_unused(amrlev, mglev, ng);
1520 amrex::Abort("MLLinOpT::make: how did we get here?");
1521 return {};
1522 }
1523}
1524
1525template <typename MF>
1526MF
1527MLLinOpT<MF>::makeAlias (MF const& mf) const
1528{
1529 if constexpr (IsMultiFabLike_v<MF>) {
1530 return MF(mf, amrex::make_alias, 0, mf.nComp());
1531 } else {
1533 amrex::Abort("MLLinOpT::makeAlias: how did we get here?");
1534 return {};
1535 }
1536}
1537
1538template <typename MF>
1539MF
1540MLLinOpT<MF>::makeCoarseMG (int amrlev, int mglev, IntVect const& ng) const
1541{
1542 if constexpr (IsMultiFabLike_v<MF>) {
1543 BoxArray cba = m_grids[amrlev][mglev];
1544 IntVect ratio = (amrlev > 0) ? IntVect(2) : mg_coarsen_ratio_vec[mglev];
1545 cba.coarsen(ratio);
1546 cba.convert(m_ixtype);
1547 return MF(cba, m_dmap[amrlev][mglev], getNComp(), ng);
1548 } else {
1549 amrex::ignore_unused(amrlev, mglev, ng);
1550 amrex::Abort("MLLinOpT::makeCoarseMG: how did we get here?");
1551 return {};
1552 }
1553}
1554
1555template <typename MF>
1556MF
1557MLLinOpT<MF>::makeCoarseAmr (int famrlev, IntVect const& ng) const
1558{
1559 if constexpr (IsMultiFabLike_v<MF>) {
1560 BoxArray cba = m_grids[famrlev][0];
1561 IntVect ratio(AMRRefRatioVect(famrlev-1));
1562 cba.coarsen(ratio);
1563 cba.convert(m_ixtype);
1564 return MF(cba, m_dmap[famrlev][0], getNComp(), ng);
1565 } else {
1566 amrex::ignore_unused(famrlev, ng);
1567 amrex::Abort("MLLinOpT::makeCoarseAmr: how did we get here?");
1568 return {};
1569 }
1570}
1571
1572template <typename MF>
1573void
1575{
1576 if (new_size <= 0 || new_size >= m_num_mg_levels[0]) { return; }
1577
1578 m_num_mg_levels[0] = new_size;
1579
1580 m_geom[0].resize(new_size);
1581 m_grids[0].resize(new_size);
1582 m_dmap[0].resize(new_size);
1583 m_factory[0].resize(new_size);
1584
1585 if (m_bottom_comm != m_default_comm) {
1586 m_bottom_comm = makeSubCommunicator(m_dmap[0].back());
1587 }
1588}
1589
1590template <typename MF>
1591void
1592MLLinOpT<MF>::avgDownResMG (int clev, MF& cres, MF const& fres) const
1593{
1594 amrex::ignore_unused(clev, cres, fres);
1595 if constexpr (amrex::IsFabArray<MF>::value) {
1596 const int ncomp = this->getNComp();
1597#ifdef AMREX_USE_EB
1598 if (!fres.isAllRegular()) {
1599 if constexpr (std::is_same<MF,MultiFab>()) {
1600 amrex::EB_average_down(fres, cres, 0, ncomp,
1601 mg_coarsen_ratio_vec[clev-1]);
1602 } else {
1603 amrex::Abort("EB_average_down only works with MultiFab");
1604 }
1605 } else
1606#endif
1607 {
1608 amrex::average_down(fres, cres, 0, ncomp, mg_coarsen_ratio_vec[clev-1]);
1609 }
1610 } else {
1611 amrex::Abort("For non-FabArray, MLLinOpT<MF>::avgDownResMG should be overridden.");
1612 }
1613}
1614
1615template <typename MF>
1616bool
1617MLLinOpT<MF>::isMFIterSafe (int amrlev, int mglev1, int mglev2) const
1618{
1619 return m_dmap[amrlev][mglev1] == m_dmap[amrlev][mglev2]
1620 && BoxArray::SameRefs(m_grids[amrlev][mglev1], m_grids[amrlev][mglev2]);
1621}
1622
1623template <typename MF>
1624template <typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF> &&
1625 IsMultiFabLike_v<AMF>, int> >
1626void
1627MLLinOpT<MF>::setLevelBC (int amrlev, const AMF* levelbcdata,
1628 const AMF* robinbc_a, const AMF* robinbc_b,
1629 const AMF* robinbc_f)
1630{
1631 const int ncomp = this->getNComp();
1632 if (levelbcdata) {
1633 levelbc_raii[amrlev] = std::make_unique<MF>(levelbcdata->boxArray(),
1634 levelbcdata->DistributionMap(),
1635 ncomp, levelbcdata->nGrowVect());
1636 levelbc_raii[amrlev]->LocalCopy(*levelbcdata, 0, 0, ncomp,
1637 levelbcdata->nGrowVect());
1638 } else {
1639 levelbc_raii[amrlev].reset();
1640 }
1641
1642 if (robinbc_a) {
1643 robin_a_raii[amrlev] = std::make_unique<MF>(robinbc_a->boxArray(),
1644 robinbc_a->DistributionMap(),
1645 ncomp, robinbc_a->nGrowVect());
1646 robin_a_raii[amrlev]->LocalCopy(*robinbc_a, 0, 0, ncomp,
1647 robinbc_a->nGrowVect());
1648 } else {
1649 robin_a_raii[amrlev].reset();
1650 }
1651
1652 if (robinbc_b) {
1653 robin_b_raii[amrlev] = std::make_unique<MF>(robinbc_b->boxArray(),
1654 robinbc_b->DistributionMap(),
1655 ncomp, robinbc_b->nGrowVect());
1656 robin_b_raii[amrlev]->LocalCopy(*robinbc_b, 0, 0, ncomp,
1657 robinbc_b->nGrowVect());
1658 } else {
1659 robin_b_raii[amrlev].reset();
1660 }
1661
1662 if (robinbc_f) {
1663 robin_f_raii[amrlev] = std::make_unique<MF>(robinbc_f->boxArray(),
1664 robinbc_f->DistributionMap(),
1665 ncomp, robinbc_f->nGrowVect());
1666 robin_f_raii[amrlev]->LocalCopy(*robinbc_f, 0, 0, ncomp,
1667 robinbc_f->nGrowVect());
1668 } else {
1669 robin_f_raii[amrlev].reset();
1670 }
1671
1672 this->setLevelBC(amrlev, levelbc_raii[amrlev].get(), robin_a_raii[amrlev].get(),
1673 robin_b_raii[amrlev].get(), robin_f_raii[amrlev].get());
1674}
1675
1676template <typename MF>
1677auto
1679{
1680 AMREX_ALWAYS_ASSERT(NAMRLevels() == 1);
1681 return xdoty(0,0,*x[0],*y[0],false);
1682}
1683
1684template <typename MF>
1685auto
1687{
1688 AMREX_ALWAYS_ASSERT(NAMRLevels() == 1);
1689 auto r = xdoty(0,0,*x[0],*x[0],false);
1690 return std::sqrt(r);
1691}
1692
1693extern template class MLLinOpT<MultiFab>;
1694
1697
1698}
1699
1700#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#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
Infrastructure for storing per-face boundary data in FabSets.
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_PICK(a, b, c)
Definition AMReX_SPACE.H:173
#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
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:837
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:672
BoxArray & convert(IndexType typ)
Apply Box::convert(IndexType) to each Box in the BoxArray.
Definition AMReX_BoxArray.cpp:813
__host__ __device__ BoxND & makeSlab(int direction, int slab_index) noexcept
Definition AMReX_Box.H:826
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:111
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
static DistributionMapping makeSFC(const MultiFab &weight, bool sort=true)
Definition AMReX_DistributionMapping.cpp:1766
Definition AMReX_EBFabFactory.H:32
Definition AMReX_FabFactory.H:50
Solve using GMRES with multigrid as preconditioner.
Definition AMReX_GMRES_MLMG.H:28
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
Interface
HYPRE interface modes supported.
Definition AMReX_Hypre.H:37
__host__ static __device__ constexpr std::size_t size() noexcept
Definition AMReX_IntVect.H:824
Definition AMReX_MLABecLaplacian.H:15
Definition AMReX_MLCGSolver.H:12
Definition AMReX_MLLinOp.H:103
const MF * m_coarse_data_for_bc
Definition AMReX_MLLinOp.H:648
Vector< Vector< std::unique_ptr< FabFactory< FAB > > > > m_factory
Definition AMReX_MLLinOp.H:619
virtual void avgDownResMG(int clev, MF &cres, MF const &fres) const
Definition AMReX_MLLinOp.H:1592
int NAMRLevels() const noexcept
Return the number of AMR levels.
Definition AMReX_MLLinOp.H:572
bool m_do_consolidation
Definition AMReX_MLLinOp.H:610
bool isCellCentered() const noexcept
Definition AMReX_MLLinOp.H:702
IntVect m_ixtype
Definition AMReX_MLLinOp.H:607
void setVerbose(int v) noexcept
Set verbosity.
Definition AMReX_MLLinOp.H:250
bool isMFIterSafe(int amrlev, int mglev1, int mglev2) const
Definition AMReX_MLLinOp.H:1617
RealVect m_coarse_bc_loc
Definition AMReX_MLLinOp.H:647
virtual bool needsUpdate() const
Does it need update if it's reused?
Definition AMReX_MLLinOp.H:272
virtual void interpolation(int amrlev, int fmglev, MF &fine, const MF &crse) const =0
Add interpolated coarse MG level data to fine MG level data.
virtual void setLevelBC(int, const MF *, const MF *=nullptr, const MF *=nullptr, const MF *=nullptr)=0
Set boundary conditions for given level. For cell-centered solves only.
virtual MF make(int amrlev, int mglev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1512
FabFactory< FAB > const * Factory(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:666
void setDomainBC(const Vector< Array< BCType, 3 > > &lobc, const Vector< Array< BCType, 3 > > &hibc)
Boundary of the whole domain.
Definition AMReX_MLLinOp.H:1220
Array< Real, 3 > m_domain_bloc_hi
Definition AMReX_MLLinOp.H:642
T get_d0(T const &dx, T const &dy, T const &) const noexcept
Definition AMReX_MLLinOp.H:743
MPI_Comm BottomCommunicator() const noexcept
Definition AMReX_MLLinOp.H:693
void setEnforceSingularSolvable(bool o) noexcept
Definition AMReX_MLLinOp.H:259
MPI_Comm Communicator() const noexcept
Definition AMReX_MLLinOp.H:694
int mg_domain_min_width
Definition AMReX_MLLinOp.H:591
void setMaxOrder(int o) noexcept
Set order of interpolation at coarse/fine boundary.
Definition AMReX_MLLinOp.H:253
virtual void compGrad(int amrlev, const Array< MF *, 3 > &grad, MF &sol, Location loc) const
Compute gradients of the solution.
Definition AMReX_MLLinOp.H:445
virtual void interpAssign(int amrlev, int fmglev, MF &fine, MF &crse) const
Overwrite fine MG level data with interpolated coarse data.
Definition AMReX_MLLinOp.H:305
virtual std::string name() const
Definition AMReX_MLLinOp.H:136
void setCoarseFineBC(const AMF *crse, IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Definition AMReX_MLLinOp.H:1484
GpuArray< BCType, 3 > LoBC(int icomp=0) const noexcept
Definition AMReX_MLLinOp.H:670
bool doAgglomeration() const noexcept
Definition AMReX_MLLinOp.H:698
MF m_coarse_data_for_bc_raii
Definition AMReX_MLLinOp.H:649
virtual void setDirichletNodesToZero(int, int, MF &) const
Definition AMReX_MLLinOp.H:485
MLLinOpT< MF > & operator=(const MLLinOpT< MF > &)=delete
std::unique_ptr< CommContainer > m_raii_comm
Definition AMReX_MLLinOp.H:639
bool m_do_semicoarsening
Definition AMReX_MLLinOp.H:612
bool hasRobinBC() const noexcept
Definition AMReX_MLLinOp.H:1294
Vector< Array< BCType, 3 > > m_hibc
Definition AMReX_MLLinOp.H:581
virtual void resizeMultiGrid(int new_size)
Definition AMReX_MLLinOp.H:1574
Vector< Vector< BoxArray > > m_grids
Definition AMReX_MLLinOp.H:617
virtual MF makeCoarseAmr(int famrlev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1557
bool m_do_agglomeration
Definition AMReX_MLLinOp.H:609
void setCoarseFineBC(const AMF *crse, int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Definition AMReX_MLLinOp.H:1475
virtual MF makeAlias(MF const &mf) const
Definition AMReX_MLLinOp.H:1527
Array4< T > compactify(Array4< T > const &a) const noexcept
Definition AMReX_MLLinOp.H:729
static constexpr int mg_coarsen_ratio
Definition AMReX_MLLinOp.H:589
virtual void copyNSolveSolution(MF &, MF const &) const
Definition AMReX_MLLinOp.H:549
virtual void solutionResidual(int amrlev, MF &resid, MF &x, const MF &b, const MF *crse_bcdata=nullptr)=0
Compute residual for solution.
int getMaxOrder() const noexcept
Get order of interpolation at coarse/fine boundary.
Definition AMReX_MLLinOp.H:255
virtual int getNComp() const
Return number of components.
Definition AMReX_MLLinOp.H:267
virtual void getFluxes(const Vector< MF * > &, const Vector< MF * > &) const
Definition AMReX_MLLinOp.H:514
void setCoarseFineBCLocation(const RealVect &cloc) noexcept
Definition AMReX_MLLinOp.H:696
Vector< int > m_amr_ref_ratio
Definition AMReX_MLLinOp.H:602
MPI_Comm m_default_comm
Definition AMReX_MLLinOp.H:622
bool isBottomActive() const noexcept
Definition AMReX_MLLinOp.H:691
virtual Vector< RT > getSolvabilityOffset(int, int, MF const &) const
get offset for fixing solvability
Definition AMReX_MLLinOp.H:472
virtual void getFluxes(const Vector< Array< MF *, 3 > > &, const Vector< MF * > &, Location) const
Definition AMReX_MLLinOp.H:509
virtual BottomSolver getDefaultBottomSolver() const
Definition AMReX_MLLinOp.H:264
typename FabDataType< MF >::fab_type FAB
Definition AMReX_MLLinOp.H:113
virtual RT normInf(int amrlev, MF const &mf, bool local) const =0
Vector< int > m_num_mg_levels
Definition AMReX_MLLinOp.H:604
bool hasBC(BCType bct) const noexcept
Definition AMReX_MLLinOp.H:1272
Vector< Vector< DistributionMapping > > m_dmap
Definition AMReX_MLLinOp.H:618
int verbose
Definition AMReX_MLLinOp.H:595
IntVect m_coarse_data_crse_ratio
Definition AMReX_MLLinOp.H:646
virtual void correctionResidual(int amrlev, int mglev, MF &resid, MF &x, const MF &b, BCMode bc_mode, const MF *crse_bcdata=nullptr)=0
Compute residual for the residual-correction form, resid = b - L(x)
const Vector< int > & AMRRefRatio() const noexcept
Return AMR refinement ratios.
Definition AMReX_MLLinOp.H:654
MLLinOpT(MLLinOpT< MF > &&)=delete
Vector< Array< BCType, 3 > > m_hibc_orig
Definition AMReX_MLLinOp.H:585
void setCoarseFineBC(const MF *crse, int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Set coarse/fine boundary conditions. For cell-centered solves only.
Definition AMReX_MLLinOp.H:1456
virtual void apply(int amrlev, int mglev, MF &out, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr) const =0
Apply the linear operator, out = L(in)
bool needsCoarseDataForBC() const noexcept
Needs coarse data for bc?
Definition AMReX_MLLinOp.H:182
typename FabDataType< MF >::value_type RT
Definition AMReX_MLLinOp.H:114
virtual void update()
Update for reuse.
Definition AMReX_MLLinOp.H:274
Vector< Array< BCType, 3 > > m_lobc_orig
Definition AMReX_MLLinOp.H:584
bool m_precond_mode
Definition AMReX_MLLinOp.H:651
virtual std::unique_ptr< FabFactory< FAB > > makeFactory(int, int) const
Definition AMReX_MLLinOp.H:718
virtual void applyMetricTerm(int, int, MF &) const
apply metric terms if there are any
Definition AMReX_MLLinOp.H:453
virtual void unapplyMetricTerm(int, int, MF &) const
unapply metric terms if there are any
Definition AMReX_MLLinOp.H:455
bool hasHiddenDimension() const noexcept
Definition AMReX_MLLinOp.H:724
virtual bool isBottomSingular() const =0
Is the bottom of MG singular?
virtual void reflux(int crse_amrlev, MF &res, const MF &crse_sol, const MF &crse_rhs, MF &fine_res, MF &fine_sol, const MF &fine_rhs) const
Reflux at AMR coarse/fine boundary.
Definition AMReX_MLLinOp.H:413
virtual IntVect getNGrowVectRestriction() const
Definition AMReX_MLLinOp.H:704
virtual void compFlux(int amrlev, const Array< MF *, 3 > &fluxes, MF &sol, Location loc) const
Compute fluxes.
Definition AMReX_MLLinOp.H:430
static constexpr int mg_box_min_width
Definition AMReX_MLLinOp.H:590
virtual RT dotProductPrecond(Vector< MF const * > const &x, Vector< MF const * > const &y) const
Definition AMReX_MLLinOp.H:1678
virtual MF makeCoarseMG(int amrlev, int mglev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1540
virtual void fixSolvabilityByOffset(int, int, MF &, Vector< RT > const &) const
fix solvability by subtracting offset from RHS
Definition AMReX_MLLinOp.H:476
const Geometry & Geom(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:577
virtual void endPrecondBC()
Definition AMReX_MLLinOp.H:567
int hiddenDirection() const noexcept
Definition AMReX_MLLinOp.H:725
void setDomainBCLoc(const Array< Real, 3 > &lo_bcloc, const Array< Real, 3 > &hi_bcloc) noexcept
Set location of domain boundaries.
Definition AMReX_MLLinOp.H:1447
Vector< Array< BCType, 3 > > m_lobc
Definition AMReX_MLLinOp.H:580
Vector< int > m_domain_covered
Definition AMReX_MLLinOp.H:620
const MLLinOpT< MF > * m_parent
Definition AMReX_MLLinOp.H:605
bool doSemicoarsening() const noexcept
Definition AMReX_MLLinOp.H:700
virtual bool supportNSolve() const
Definition AMReX_MLLinOp.H:547
virtual bool supportRobinBC() const noexcept
Definition AMReX_MLLinOp.H:685
virtual void normalize(int amrlev, int mglev, MF &mf) const
Divide mf by the diagonal component of the operator. Used by bicgstab.
Definition AMReX_MLLinOp.H:370
virtual void applyInhomogNeumannTerm(int, MF &) const
Extra terms introduced when we treat inhomogeneous Nuemann BC as homogeneous.
Definition AMReX_MLLinOp.H:461
virtual void avgDownResAmr(int clev, MF &cres, MF const &fres) const
Definition AMReX_MLLinOp.H:557
Vector< Vector< Geometry > > m_geom
first Vector is for amr level and second is mg level
Definition AMReX_MLLinOp.H:616
virtual RT norm2Precond(Vector< MF const * > const &x) const
Definition AMReX_MLLinOp.H:1686
MLLinOpT(const MLLinOpT< MF > &)=delete
virtual bool scaleRHS(int, MF *) const
scale RHS to fix solvability
Definition AMReX_MLLinOp.H:467
virtual void averageDownAndSync(Vector< MF > &sol) const =0
void setCoarseFineBC(const MF *crse, IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Definition AMReX_MLLinOp.H:1464
MLLinOpT()=default
Box compactify(Box const &b) const noexcept
Definition AMReX_MLLinOp.H:1301
bool m_needs_coarse_data_for_bc
Definition AMReX_MLLinOp.H:644
GpuArray< BCType, 3 > HiBC(int icomp=0) const noexcept
Definition AMReX_MLLinOp.H:675
int maxorder
Definition AMReX_MLLinOp.H:597
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:786
Vector< IntVect > mg_coarsen_ratio_vec
Definition AMReX_MLLinOp.H:613
MF MFType
Definition AMReX_MLLinOp.H:112
virtual void preparePrecond()
Definition AMReX_MLLinOp.H:481
virtual ~MLLinOpT()=default
virtual void averageDownSolutionRHS(int camrlev, MF &crse_sol, MF &crse_rhs, const MF &fine_sol, const MF &fine_rhs)
Average-down data from fine AMR level to coarse AMR level.
Definition AMReX_MLLinOp.H:335
LPInfo info
Definition AMReX_MLLinOp.H:593
int NMGLevels(int amrlev) const noexcept
Return the number of MG levels at given AMR level.
Definition AMReX_MLLinOp.H:575
T get_d1(T const &, T const &dy, T const &dz) const noexcept
Definition AMReX_MLLinOp.H:753
bool enforceSingularSolvable
Definition AMReX_MLLinOp.H:599
void setLevelBC(int amrlev, const AMF *levelbcdata, const AMF *robinbc_a=nullptr, const AMF *robinbc_b=nullptr, const AMF *robinbc_f=nullptr)
Definition AMReX_MLLinOp.H:1627
virtual RT xdoty(int amrlev, int mglev, const MF &x, const MF &y, bool local) const =0
x dot y, used by the bottom solver
virtual void interpolationAmr(int famrlev, MF &fine, const MF &crse, IntVect const &nghost) const
Interpolation between AMR levels.
Definition AMReX_MLLinOp.H:319
bool doConsolidation() const noexcept
Definition AMReX_MLLinOp.H:699
virtual std::unique_ptr< MLLinOpT< MF > > makeNLinOp(int) const
Definition AMReX_MLLinOp.H:503
virtual bool supportInhomogNeumannBC() const noexcept
Definition AMReX_MLLinOp.H:686
virtual void prepareForFluxes(int, const MF *=nullptr)
Definition AMReX_MLLinOp.H:386
LinOpBCType m_coarse_fine_bc_type
Definition AMReX_MLLinOp.H:645
Array< Real, 3 > m_domain_bloc_lo
Definition AMReX_MLLinOp.H:641
virtual bool isSingular(int amrlev) const =0
Is it singular on given AMR level?
virtual void postSolve(Vector< MF * > const &) const
Definition AMReX_MLLinOp.H:551
int AMRRefRatio(int amr_lev) const noexcept
Return AMR refinement ratio at given AMR level.
Definition AMReX_MLLinOp.H:657
MPI_Comm m_bottom_comm
Definition AMReX_MLLinOp.H:623
virtual void restriction(int amrlev, int cmglev, MF &crse, MF &fine) const =0
Restriction onto coarse MG level.
virtual void beginPrecondBC()
Definition AMReX_MLLinOp.H:566
virtual void applyOverset(int, MF &) const
for overset solver only
Definition AMReX_MLLinOp.H:464
virtual void getEBFluxes(const Vector< MF * > &, const Vector< MF * > &) const
Definition AMReX_MLLinOp.H:520
virtual void unimposeNeumannBC(int, MF &) const
This is needed for our nodal projection solver.
Definition AMReX_MLLinOp.H:458
bool getEnforceSingularSolvable() const noexcept
Definition AMReX_MLLinOp.H:262
virtual void prepareForSolve()=0
IntVect AMRRefRatioVect(int amr_lev) const noexcept
Return AMR refinement ratio as IntVect (1 in hidden direction)
Definition AMReX_MLLinOp.H:660
virtual void smooth(int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary, int niter) const =0
Smooth.
virtual void make(Vector< Vector< MF > > &mf, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1498
bool hasInhomogNeumannBC() const noexcept
Definition AMReX_MLLinOp.H:1287
void setDomainBC(const Array< BCType, 3 > &lobc, const Array< BCType, 3 > &hibc) noexcept
Boundary of the whole domain.
Definition AMReX_MLLinOp.H:1210
virtual int getNGrow(int=0, int=0) const
Definition AMReX_MLLinOp.H:269
int m_num_amr_levels
Definition AMReX_MLLinOp.H:601
Definition AMReX_MLMGBndry.H:12
Definition AMReX_MLMG.H:17
Definition AMReX_MLPoisson.H:19
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Long size() const noexcept
Definition AMReX_Vector.H:53
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1409
std::array< T, N > Array
Definition AMReX_Array.H:26
bool notInLaunchRegion() noexcept
Definition AMReX_GpuControl.H:89
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
int local_to_global_rank(int rank) noexcept
translate between local rank and global rank
Definition AMReX_ParallelContext.H:95
int global_to_local_rank(int rank) noexcept
Definition AMReX_ParallelContext.H:98
int NProcsSub() noexcept
number of ranks in current frame
Definition AMReX_ParallelContext.H:74
MPI_Comm Communicator() noexcept
Definition AMReX_ParallelDescriptor.H:223
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
int MPI_Group
Definition AMReX_ccse-mpi.H:52
static constexpr int MPI_COMM_NULL
Definition AMReX_ccse-mpi.H:59
Definition AMReX_Amr.cpp:50
@ make_alias
Definition AMReX_MakeType.H:7
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
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
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
LinOpBCType
Definition AMReX_LO_BCTYPES.H:27
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)
Volume-weighted average-down from fine to coarse using EB volume fractions.
Definition AMReX_EBMultiFabUtil.cpp:336
__host__ __device__ BoxND< dim > enclosedCells(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with CELL based coordinates in direction dir that is enclosed by b....
Definition AMReX_Box.H:1586
BottomSolver
Definition AMReX_MLLinOp.H:32
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
std::unique_ptr< Hypre > makeHypre(const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom, MPI_Comm comm_, Hypre::Interface interface, const iMultiFab *overset_mask)
Factory that instantiates the requested HYPRE interface.
Definition AMReX_Hypre.cpp:12
void Warning(const std::string &msg)
Print out warning message to cerr.
Definition AMReX.cpp:247
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
__host__ __device__ constexpr int get(IntVectND< dim > const &iv) noexcept
Get I'th element of IntVectND<dim>
Definition AMReX_IntVect.H:1334
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_FabDataType.H:9
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
Definition AMReX_TypeTraits.H:29
Definition AMReX_MLLinOp.H:37
LPInfo & setConsolidationRatio(int x) noexcept
Definition AMReX_MLLinOp.H:57
int con_strategy
Definition AMReX_MLLinOp.H:44
bool do_semicoarsening
Definition AMReX_MLLinOp.H:40
bool has_metric_term
Definition AMReX_MLLinOp.H:45
LPInfo & setConsolidationGridSize(int x) noexcept
Definition AMReX_MLLinOp.H:56
int max_semicoarsening_level
Definition AMReX_MLLinOp.H:47
bool hasHiddenDimension() const noexcept
Definition AMReX_MLLinOp.H:66
int con_ratio
Definition AMReX_MLLinOp.H:43
bool do_consolidation
Definition AMReX_MLLinOp.H:39
int con_grid_size
Definition AMReX_MLLinOp.H:42
LPInfo & setSemicoarsening(bool x) noexcept
Definition AMReX_MLLinOp.H:54
LPInfo & setHiddenDirection(int n) noexcept
Definition AMReX_MLLinOp.H:63
LPInfo & setConsolidation(bool x) noexcept
Definition AMReX_MLLinOp.H:53
LPInfo & setSemicoarseningDirection(int n) noexcept
Definition AMReX_MLLinOp.H:62
LPInfo & setMaxSemicoarseningLevel(int n) noexcept
Definition AMReX_MLLinOp.H:61
bool do_agglomeration
Definition AMReX_MLLinOp.H:38
LPInfo & setMetricTerm(bool x) noexcept
Definition AMReX_MLLinOp.H:59
bool deterministic
Enable deterministic mode for GPU operations.
Definition AMReX_MLLinOp.H:50
static constexpr int getDefaultConsolidationGridSize()
Definition AMReX_MLLinOp.H:78
LPInfo & setConsolidationStrategy(int x) noexcept
Definition AMReX_MLLinOp.H:58
int max_coarsening_level
Definition AMReX_MLLinOp.H:46
int agg_grid_size
Definition AMReX_MLLinOp.H:41
static constexpr int getDefaultAgglomerationGridSize()
Definition AMReX_MLLinOp.H:70
int hidden_direction
Definition AMReX_MLLinOp.H:49
LPInfo & setAgglomerationGridSize(int x) noexcept
Definition AMReX_MLLinOp.H:55
LPInfo & setMaxCoarseningLevel(int n) noexcept
Definition AMReX_MLLinOp.H:60
LPInfo & setDeterministic(bool x) noexcept
Definition AMReX_MLLinOp.H:64
LPInfo & setAgglomeration(bool x) noexcept
Definition AMReX_MLLinOp.H:52
int semicoarsening_direction
Definition AMReX_MLLinOp.H:48
Definition AMReX_MLLinOp.H:88
StateMode
Definition AMReX_MLLinOp.H:90
BCMode
Definition AMReX_MLLinOp.H:89
Location
Definition AMReX_MLLinOp.H:91
FabArray memory allocation information.
Definition AMReX_FabArray.H:67