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