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