Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_AmrLevel.H
Go to the documentation of this file.
1
2#ifndef AMREX_AmrLevel_H_
3#define AMREX_AmrLevel_H_
4#include <AMReX_Config.H>
5
6#include <AMReX_REAL.H>
7#include <AMReX_Geometry.H>
8#include <AMReX_MultiFab.H>
10#include <AMReX_LayoutData.H>
11#include <AMReX_Derive.H>
12#include <AMReX_BCRec.H>
13#include <AMReX_Amr.H>
16#include <AMReX_StateData.H>
17#include <AMReX_VisMF.H>
18#include <AMReX_RungeKutta.H>
19#include <AMReX_FillPatcher.H>
20
25#ifdef AMREX_USE_EB
26#include <AMReX_EBSupport.H>
27#endif
28
29#include <memory>
30#include <map>
31
32namespace amrex {
33
34class TagBox;
35class TagBoxArray;
36
44{
45 friend class Amr;
46 friend class FillPatchIterator;
48
49public:
58 virtual ~AmrLevel ();
59
61 AmrLevel (const AmrLevel&) = delete;
63 AmrLevel (AmrLevel&&) = delete;
65 AmrLevel& operator= (const AmrLevel&) = delete;
68
76 void LevelDirectoryNames (const std::string &dir,
77 std::string &LevelDir,
78 std::string &FullPath) const;
84 virtual void CreateLevelDirectory (const std::string &dir);
90 void SetLevelDirectoryCreated(bool ldc) noexcept { levelDirectoryCreated = ldc; }
99 virtual std::string thePlotFileType () const
100 {
101 static const std::string the_plot_file_type("HyperCLaw-V1.1");
102 return the_plot_file_type;
103 }
111 virtual void writePlotFile (const std::string& dir,
112 std::ostream& os,
114
121 virtual void writePlotFilePre (const std::string& dir,
122 std::ostream& os);
123
130 virtual void writePlotFilePost (const std::string& dir,
131 std::ostream& os);
132
140 virtual void writeSmallPlotFile (const std::string& dir,
141 std::ostream& os,
143 {
144 (void)dir;
145 (void)os;
146 (void)how;
147 }
156 virtual void checkPoint (const std::string& dir,
157 std::ostream& os,
159 bool dump_old = true);
166 virtual void checkPointPre (const std::string& dir,
167 std::ostream& os);
174 virtual void checkPointPost (const std::string& dir,
175 std::ostream& os);
183 virtual void restart (Amr& papa,
184 std::istream& is,
185 bool bReadSpecial = false);
186
192 virtual void set_state_in_checkpoint (Vector<int>& state_in_checkpoint);
193
201 static bool isStateVariable (const std::string& name,
202 int& state_indx,
203 int& ncomp);
204
206 static void FlushFPICache ();
217 virtual void computeInitialDt (int finest_level,
218 int sub_cycle,
219 Vector<int>& n_cycle,
220 const Vector<IntVect>& ref_ratio,
221 Vector<Real>& dt_level,
222 Real stop_time) = 0;
235 virtual void computeNewDt (int finest_level,
236 int sub_cycle,
237 Vector<int>& n_cycle,
238 const Vector<IntVect>& ref_ratio,
239 Vector<Real>& dt_min,
240 Vector<Real>& dt_level,
241 Real stop_time,
242 int post_regrid_flag) = 0;
252 virtual Real advance (Real time,
253 Real dt,
254 int iteration,
255 int ncycle) = 0;
256
264 virtual void post_timestep (int iteration);
270 virtual void postCoarseTimeStep (Real time);
272 virtual void post_restart () {}
279 virtual void post_regrid (int lbase,
280 int new_finest) = 0;
286 virtual void post_init (Real stop_time) = 0;
290 virtual int okToContinue () { return 1; }
296 virtual int okToRegrid ();
302 virtual void initData () = 0;
310 virtual void setTimeLevel (Real time,
311 Real dt_old,
312 Real dt_new);
314 virtual void allocOldData ();
316 virtual void removeOldData ();
324 virtual void init (AmrLevel &old) = 0;
330 virtual void init () = 0;
332 void reset ();
334 int Level () const noexcept { return level; }
336 const BoxArray& boxArray () const noexcept { return grids; }
338 const BoxArray& getEdgeBoxArray (int dir) const noexcept;
340 const BoxArray& getNodalBoxArray () const noexcept;
341 //
342 const DistributionMapping& DistributionMap () const noexcept { return dmap; }
343 //
344 const FabFactory<FArrayBox>& Factory () const noexcept { return *m_factory; }
345 //
346#ifdef AMREX_USE_EB
347 const EBFArrayBoxFactory&
348 EBFactory () const noexcept {
349 return static_cast<amrex::EBFArrayBoxFactory const&>(*m_factory);
350 }
351#endif
352 //
354 int numGrids () const noexcept { return static_cast<int>(grids.size()); }
356 int numStates () const noexcept { return static_cast<int>(state.size()); }
358 const Box& Domain () const noexcept { return geom.Domain(); }
360 int nStep () const noexcept { return parent->levelSteps(level); }
362 const Geometry& Geom () const noexcept { return geom; }
363 //
364 const IntVect& fineRatio () const noexcept { return fine_ratio; }
366 Long countCells () const noexcept;
367
369 const BoxArray& getAreaNotToTag () noexcept;
371 const Box& getAreaToTag () noexcept;
373 void constructAreaNotToTag ();
375 void setAreaNotToTag (BoxArray& ba) noexcept;
376
378 void resetFillPatcher ();
379
391 virtual void errorEst (TagBoxArray& tb,
392 int clearval,
393 int tagval,
394 Real time,
395 int n_error_buf = 0,
396 int ngrow = 0) = 0;
398
409 void FillCoarsePatch (MultiFab& mf,
410 int dcomp,
411 Real time,
412 int state_idx,
413 int scomp,
414 int ncomp,
415 int nghost = 0);
426 virtual void setPhysBoundaryValues (FArrayBox& dest,
427 int state_indx,
428 Real time,
429 int dest_comp,
430 int src_comp,
431 int num_comp);
442 virtual std::unique_ptr<MultiFab> derive (const std::string& name,
443 Real time,
444 int ngrow);
454 virtual void derive (const std::string& name,
455 Real time,
456 MultiFab& mf,
457 int dcomp);
459 StateData& get_state_data (int state_indx) noexcept { return state[state_indx]; }
461 MultiFab& get_old_data (int state_indx) noexcept { return state[state_indx].oldData(); }
463 const MultiFab& get_old_data (int state_indx) const noexcept { return state[state_indx].oldData(); }
465 MultiFab& get_new_data (int state_indx) noexcept { return state[state_indx].newData(); }
467 const MultiFab& get_new_data (int state_indx) const noexcept { return state[state_indx].newData(); }
469 static const DescriptorList& get_desc_lst () noexcept { return desc_lst; }
471 static DeriveList& get_derive_lst () noexcept;
473 int postStepRegrid () const noexcept { return post_step_regrid; }
479 void setPostStepRegrid (int new_val) noexcept { post_step_regrid = new_val; }
480
487
496 Vector<int> getBCArray (int State_Type,
497 int gridno,
498 int strt_comp,
499 int ncomp);
507 MultiFab& get_data (int state_indx, Real time) noexcept;
519 int state_index,
520 int scomp,
521 int dcomp,
522 int ncomp,
523 Real time) const;
530 virtual void manual_tags_placement (TagBoxArray& tags,
531 const Vector<IntVect>& bf_lev);
533 virtual void setPlotVariables ();
535 virtual void setSmallPlotVariables ();
541 virtual Real estimateWork();
542
548 virtual int WorkEstType () { return -1; }
549
557 TimeLevel which_time (int state_indx, Real time) const noexcept;
558
564 virtual bool checkPointNow ();
565
569 virtual bool writePlotNow ();
570
574 virtual bool writeSmallPlotNow ();
575
576#ifdef AMREX_PARTICLES
583 virtual void particle_redistribute (int lbase = 0, bool a_init = false)
584 {
585 (void)lbase;
586 (void)a_init;
587 }
588#endif
589
601 void FillPatcherFill (amrex::MultiFab& mf, int dcomp, int ncomp, int nghost,
602 amrex::Real time, int state_index, int scomp);
603
616 static void FillPatch (AmrLevel& amrlevel,
617 MultiFab& leveldata,
618 int boxGrow,
619 Real time,
620 int index,
621 int scomp,
622 int ncomp,
623 int dcomp=0);
624
637 static void FillPatchAdd (AmrLevel& amrlevel,
638 MultiFab& leveldata,
639 int boxGrow,
640 Real time,
641 int index,
642 int scomp,
643 int ncomp,
644 int dcomp=0);
645
667 template <typename F, typename P = RungeKutta::PostStageNoOp>
668 void RK (int order, int state_type, Real time, Real dt, int iteration,
669 int ncycle, F&& f, P&& p = RungeKutta::PostStageNoOp());
670
671#ifdef AMREX_USE_EB
673 static void SetEBMaxGrowCells (int nbasic, int nvolume, int nfull) noexcept {
674 m_eb_basic_grow_cells = nbasic;
675 m_eb_volume_grow_cells = nvolume;
676 m_eb_full_grow_cells = nfull;
677 }
682 static void SetEBSupportLevel (EBSupport ebs) { m_eb_support_level = ebs; }
684#endif
685
695 static IntVect ProperBlockingFactor (AmrLevel const& amr_level, int boxGrow,
696 IndexType const& boxType,
697 StateDescriptor const& desc, int SComp);
698
699protected:
701 AmrLevel () noexcept {} // NOLINT
702
703 AmrLevel (Amr& papa,
704 int lev,
705 const Geometry& level_geom,
706 const BoxArray& ba,
707 const DistributionMapping& dm,
708 Real time);
709
711 void finishConstructor ();
712
713 //
714 // The Data.
715 //
716 int level{-1}; // AMR level (0 is coarsest).
717 Geometry geom; // Geom at this level.
718 BoxArray grids; // Cell-centered locations of grids.
719 DistributionMapping dmap; // Distribution of grids among processes
720 Amr* parent{nullptr};// Pointer to parent AMR structure.
721 IntVect crse_ratio; // Refinement ratio to coarser level.
722 IntVect fine_ratio; // Refinement ratio to finer level.
723 static DeriveList derive_lst; // List of derived quantities.
724 static DescriptorList desc_lst; // List of state variables.
725 Vector<StateData> state; // Array of state data.
726
727 BoxArray m_AreaNotToTag; //Area which shouldn't be tagged on this level.
728 Box m_AreaToTag; //Area which is allowed to be tagged on this level.
729
730 int post_step_regrid{0}; // Whether or not to do a regrid after the timestep.
731
732 bool levelDirectoryCreated{false}; // for checkpoints and plotfiles
733
734 std::unique_ptr<FabFactory<FArrayBox> > m_factory;
735
737
738private:
739
740 template <std::size_t order>
741 void storeRKCoarseData (int state_type, Real time, Real dt,
742 MultiFab const& S_old,
743 Array<MultiFab,order> const& rkk);
744
745 void FillRKPatch (int state_index, MultiFab& S, Real time,
746 int stage, int iteration, int ncycle);
747
748 mutable BoxArray edge_grids[AMREX_SPACEDIM]; // face-centered grids
749 mutable BoxArray nodal_grids; // all nodal grids
750};
751
752//
753// Forward declaration.
754//
756
758 :
759 public MFIter
760{
761 public:
762
763 friend class AmrLevel;
764
771 FillPatchIterator (AmrLevel& amrlevel,
772 MultiFab& leveldata);
773
785 FillPatchIterator (AmrLevel& amrlevel,
786 MultiFab& leveldata,
787 int boxGrow,
788 Real time,
789 int idx,
790 int scomp,
791 int ncomp);
792
802 void Initialize (int boxGrow,
803 Real time,
804 int idx,
805 int scomp,
806 int ncomp);
807
809 FArrayBox& operator() () noexcept { return m_fabs[MFIter::index()]; }
810
812 Box UngrownBox () const noexcept { return MFIter::validbox(); }
813
815 MultiFab& get_mf() noexcept { return m_fabs; }
816
817private:
818
819 void FillFromLevel0 (Real time, int idx, int scomp, int dcomp, int ncomp);
820 void FillFromTwoLevels (Real time, int idx, int scomp, int dcomp, int ncomp);
821
822 //
823 // The data.
824 //
825 AmrLevel* m_amrlevel;
826 MultiFab* m_leveldata;
827 std::vector< std::pair<int,int> > m_range;
828 MultiFab m_fabs;
829 int m_ncomp;
830};
831
833{
834public:
835
836 friend class FillPatchIterator;
837
840 MultiFab& leveldata);
841
855 MultiFab& leveldata,
856 int boxGrow,
857 Real time,
858 int index,
859 int scomp,
860 int ncomp,
861 InterpBase* mapper);
862
873 void Initialize (int boxGrow,
874 Real time,
875 int idx,
876 int scomp,
877 int ncomp,
878 InterpBase* mapper);
879
887 void fill (FArrayBox& fab, int dcomp, int idx);
888
889private:
890 //
891 // The data.
892 //
893 AmrLevel* m_amrlevel;
894 MultiFab* m_leveldata;
896 Vector< Vector<MultiFabId> > m_mfid; // [level][oldnew]
897 Interpolater* m_map = nullptr;
898 std::map<int,Box> m_ba;
899 Real m_time = std::numeric_limits<Real>::lowest();
900 int m_growsize = std::numeric_limits<int>::lowest();
901 int m_index = -1;
902 int m_scomp = -1;
903 int m_ncomp = -1;
904 bool m_FixUpCorners = false;
905
906 std::map< int,Vector< Vector<Box> > > m_fbox; // [grid][level][validregion]
907 std::map< int,Vector< Vector<Box> > > m_cbox; // [grid][level][fillablesubbox]
908 std::map< int,Vector< Vector< Vector<FillBoxId> > > > m_fbid; // [grid][level][fillablesubbox][oldnew]
909};
910
911template <typename F, typename P>
912void AmrLevel::RK (int order, int state_type, Real time, Real dt, int iteration,
913 int ncycle, F&& f, P&& p)
914{
915 BL_PROFILE("AmrLevel::RK()");
916
917 AMREX_ASSERT(AmrLevel::desc_lst[state_type].nExtra() > 0); // Need ghost cells in StateData
918
919 MultiFab& S_old = get_old_data(state_type);
920 MultiFab& S_new = get_new_data(state_type);
921
922 if (order == 2) {
923 RungeKutta::RK2(S_old, S_new, time, dt, std::forward<F>(f),
924 [&] (int /*stage*/, MultiFab& mf, Real t) {
925 FillPatcherFill(mf, 0, mf.nComp(), mf.nGrow(), t,
926 state_type, 0); },
927 std::forward<P>(p));
928 } else if (order == 3) {
929 RungeKutta::RK3(S_old, S_new, time, dt, std::forward<F>(f),
930 [&] (int stage, MultiFab& mf, Real t) {
931 FillRKPatch(state_type, mf, t, stage, iteration, ncycle);
932 },
933 [&] (Array<MultiFab,3> const& rkk) {
934 if (level < parent->finestLevel()) {
935 storeRKCoarseData(state_type, time, dt, S_old, rkk);
936 }
937 },
938 std::forward<P>(p));
939 } else if (order == 4) {
940 RungeKutta::RK4(S_old, S_new, time, dt, std::forward<F>(f),
941 [&] (int stage, MultiFab& mf, Real t) {
942 FillRKPatch(state_type, mf, t, stage, iteration, ncycle);
943 },
944 [&] (Array<MultiFab,4> const& rkk) {
945 if (level < parent->finestLevel()) {
946 storeRKCoarseData(state_type, time, dt, S_old, rkk);
947 }
948 },
949 std::forward<P>(p));
950 } else {
951 amrex::Abort("AmrLevel::RK: order = "+std::to_string(order)+" is not supported");
952 }
953}
954
955template <std::size_t order>
956void AmrLevel::storeRKCoarseData (int state_type, Real time, Real dt,
957 MultiFab const& S_old,
958 Array<MultiFab,order> const& rkk)
959{
960 BL_PROFILE("AmrLevel::storeRKCoarseData()");
961 if (level == parent->finestLevel()) { return; }
962
963 const StateDescriptor& desc = AmrLevel::desc_lst[state_type];
964
965 auto& fillpatcher = parent->getLevel(level+1).m_fillpatcher[state_type];
966 fillpatcher = std::make_unique<FillPatcher<MultiFab>>
968 parent->Geom(level+1),
970 parent->Geom(level),
971 IntVect(desc.nExtra()), desc.nComp(), desc.interp(0));
972
973 fillpatcher->storeRKCoarseData(time, dt, S_old, rkk);
974}
975
976
977}
978
979#endif /*_AmrLevel_H_*/
Driver that coordinates AmrLevel instances, time stepping, and output schedules.
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
Infrastructure for registering and computing derived quantities from state data.
Helper class that caches coarse/fine data for repeated FillPatch operations.
Manages the old/new time slices of a state descriptor on a single AMR level.
Virtual base class for managing individual levels. AmrLevel functions both as a container for state d...
Definition AMReX_AmrLevel.H:44
virtual void post_regrid(int lbase, int new_finest)=0
Operations to be performed immediately after regridding completes.
virtual Real estimateWork()
Estimate the amount of work required to advance just this level, typically based on cell counts.
Definition AMReX_AmrLevel.cpp:2057
virtual void restart(Amr &papa, std::istream &is, bool bReadSpecial=false)
Populate this level from checkpoint input stream is.
Definition AMReX_AmrLevel.cpp:386
virtual void writeSmallPlotFile(const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles)
Write a reduced plot file (same parameters as writePlotFile()) to dir.
Definition AMReX_AmrLevel.H:140
void reset()
Reset data to initial time by swapping new and old time slots.
Definition AMReX_AmrLevel.cpp:596
Long countCells() const noexcept
Returns number of valid cells on this level.
Definition AMReX_AmrLevel.cpp:504
IntVect fine_ratio
Definition AMReX_AmrLevel.H:722
int level
Definition AMReX_AmrLevel.H:716
virtual void checkPoint(const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles, bool dump_old=true)
Write current state to a checkpoint directory.
Definition AMReX_AmrLevel.cpp:510
int Level() const noexcept
Returns this AmrLevel index (0 = coarsest).
Definition AMReX_AmrLevel.H:334
virtual void writePlotFilePre(const std::string &dir, std::ostream &os)
Perform pre-plotfile operations (e.g., gather metadata) for directory dir using os.
Definition AMReX_AmrLevel.cpp:372
virtual int WorkEstType()
Select which state descriptor drives work estimates.
Definition AMReX_AmrLevel.H:548
virtual bool writePlotNow()
Return true if this level requests an immediate plotfile (default false).
Definition AMReX_AmrLevel.cpp:2069
virtual void init(AmrLevel &old)=0
Init data on this level from another AmrLevel (during regrid). This is a pure virtual function and he...
static int m_eb_volume_grow_cells
Definition AMReX_AmrLevel.H:679
const Geometry & Geom() const noexcept
Returns the geometry object.
Definition AMReX_AmrLevel.H:362
static IntVect ProperBlockingFactor(AmrLevel const &amr_level, int boxGrow, IndexType const &boxType, StateDescriptor const &desc, int SComp)
Recommendation of a proper blocking factor for filling.
Definition AMReX_AmrLevel.cpp:2286
virtual void computeInitialDt(int finest_level, int sub_cycle, Vector< int > &n_cycle, const Vector< IntVect > &ref_ratio, Vector< Real > &dt_level, Real stop_time)=0
Compute the initial time step size for each level.
Vector< std::unique_ptr< FillPatcher< MultiFab > > > m_fillpatcher
Definition AMReX_AmrLevel.H:736
virtual bool writeSmallPlotNow()
Return true if this level requests an immediate small-plot file (default false).
Definition AMReX_AmrLevel.cpp:2075
void LevelDirectoryNames(const std::string &dir, std::string &LevelDir, std::string &FullPath) const
Derive the directory name for this level inside a plot/checkpoint tree.
Definition AMReX_AmrLevel.cpp:2233
virtual void computeNewDt(int finest_level, int sub_cycle, Vector< int > &n_cycle, const Vector< IntVect > &ref_ratio, Vector< Real > &dt_min, Vector< Real > &dt_level, Real stop_time, int post_regrid_flag)=0
Compute the next time step sizes after an advance.
virtual std::unique_ptr< MultiFab > derive(const std::string &name, Real time, int ngrow)
Returns a MultiFab containing the derived data for this level. If ngrow>0 the MultiFab is built on th...
Definition AMReX_AmrLevel.cpp:1619
const BoxArray & boxArray() const noexcept
List of grids at this level.
Definition AMReX_AmrLevel.H:336
virtual void manual_tags_placement(TagBoxArray &tags, const Vector< IntVect > &bf_lev)
Adjust tags after core tagging routines run (default: no-op).
Definition AMReX_AmrLevel.cpp:63
virtual void setPhysBoundaryValues(FArrayBox &dest, int state_indx, Real time, int dest_comp, int src_comp, int num_comp)
Apply physical boundary conditions to cell-centered data.
Definition AMReX_AmrLevel.cpp:649
const MultiFab & get_old_data(int state_indx) const noexcept
Const access to the MultiFab at the old time slot for state state_indx.
Definition AMReX_AmrLevel.H:463
void SetLevelDirectoryCreated(bool ldc) noexcept
Record whether CreateLevelDirectory() has already created the folder.
Definition AMReX_AmrLevel.H:90
void setPostStepRegrid(int new_val) noexcept
Set the post-timestep regrid trigger (nonzero requests a regrid).
Definition AMReX_AmrLevel.H:479
const IntVect & fineRatio() const noexcept
Definition AMReX_AmrLevel.H:364
virtual void writePlotFilePost(const std::string &dir, std::ostream &os)
Perform post-plotfile cleanup for directory dir once writes complete.
Definition AMReX_AmrLevel.cpp:379
Vector< int > getBCArray(int State_Type, int gridno, int strt_comp, int ncomp)
Return an array describing boundary types for components [strt_comp,strt_comp+ncomp).
Definition AMReX_AmrLevel.cpp:1873
const FabFactory< FArrayBox > & Factory() const noexcept
Definition AMReX_AmrLevel.H:344
Box m_AreaToTag
Definition AMReX_AmrLevel.H:728
static void FillPatch(AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
Convenience wrapper around FillPatchSingleLevel/TwoLevels.
Definition AMReX_AmrLevel.cpp:2197
virtual void removeOldData()
Delete the old-time slots once they are no longer needed.
Definition AMReX_AmrLevel.cpp:587
Geometry geom
Definition AMReX_AmrLevel.H:717
TimeLevel which_time(int state_indx, Real time) const noexcept
Returns one the TimeLevel enums. Asserts that time is between AmrOldTime and AmrNewTime.
Definition AMReX_AmrLevel.cpp:2022
virtual void setTimeLevel(Real time, Real dt_old, Real dt_new)
Set the time metadata stored in the state data arrays (old/new).
Definition AMReX_AmrLevel.cpp:476
virtual int okToContinue()
Return nonzero to continue advancing; zero to halt the run.
Definition AMReX_AmrLevel.H:290
const BoxArray & getNodalBoxArray() const noexcept
Return the nodal grid layout for this level (cached).
Definition AMReX_AmrLevel.cpp:639
StateData & get_state_data(int state_indx) noexcept
State data object accessor for descriptor state_indx.
Definition AMReX_AmrLevel.H:459
MultiFab & get_data(int state_indx, Real time) noexcept
Fetch the MultiFab (old or new) that best matches time for state state_indx.
Definition AMReX_AmrLevel.cpp:605
const EBFArrayBoxFactory & EBFactory() const noexcept
Definition AMReX_AmrLevel.H:348
BoxArray grids
Definition AMReX_AmrLevel.H:718
virtual void allocOldData()
Allocate storage for "old" time slices before an advance.
Definition AMReX_AmrLevel.cpp:578
virtual void set_state_in_checkpoint(Vector< int > &state_in_checkpoint)
Hook to adjust which state descriptors are expected in checkpoint streams.
Definition AMReX_AmrLevel.cpp:467
static int m_eb_full_grow_cells
Definition AMReX_AmrLevel.H:680
void UpdateDistributionMaps(DistributionMapping &dmap)
Update every StateData's DistributionMapping to match dmap (e.g., after rebalancing).
Definition AMReX_AmrLevel.cpp:1856
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_AmrLevel.H:342
virtual void postCoarseTimeStep(Real time)
Called once per coarse step after fine levels finish refluxing.
Definition AMReX_AmrLevel.cpp:43
virtual void post_timestep(int iteration)
Perform level-specific updates after each subcycle.
Definition AMReX_AmrLevel.cpp:35
const Box & Domain() const noexcept
Returns the indices defining physical domain.
Definition AMReX_AmrLevel.H:358
virtual void errorEst(TagBoxArray &tb, int clearval, int tagval, Real time, int n_error_buf=0, int ngrow=0)=0
Error estimation for regridding. This is a pure virtual function and hence MUST be implemented by der...
static DeriveList derive_lst
Definition AMReX_AmrLevel.H:723
void constructAreaNotToTag()
Build internal data structures defining the area not to tag.
Definition AMReX_AmrLevel.cpp:2095
bool levelDirectoryCreated
Definition AMReX_AmrLevel.H:732
static void SetEBSupportLevel(EBSupport ebs)
Select which level of EB data to build (basic, volume, full) via ebs.
Definition AMReX_AmrLevel.H:682
static DeriveList & get_derive_lst() noexcept
Returns list of derived variables.
Definition AMReX_AmrLevel.cpp:57
IntVect crse_ratio
Definition AMReX_AmrLevel.H:721
virtual void setPlotVariables()
Override to customize which variables appear in full plotfiles on this level (default uses global lis...
Definition AMReX_AmrLevel.cpp:1901
const Box & getAreaToTag() noexcept
Get the Box (single region) that is permitted to be tagged.
Definition AMReX_AmrLevel.cpp:2085
virtual bool checkPointNow()
Return true if this level requests an immediate checkpoint (default false).
Definition AMReX_AmrLevel.cpp:2063
TimeLevel
What time are we at?
Definition AMReX_AmrLevel.H:51
@ Amr1QtrTime
Definition AMReX_AmrLevel.H:54
@ AmrOldTime
Definition AMReX_AmrLevel.H:51
@ AmrNewTime
Definition AMReX_AmrLevel.H:53
@ AmrOtherTime
Definition AMReX_AmrLevel.H:56
@ AmrHalfTime
Definition AMReX_AmrLevel.H:52
@ Amr3QtrTime
Definition AMReX_AmrLevel.H:55
virtual void initData()=0
Init grid data at problem start-up. This is a pure virtual function and hence MUST be implemented by ...
AmrLevel & operator=(const AmrLevel &)=delete
Copy assignment is undefined for the same ownership reasons.
virtual void checkPointPost(const std::string &dir, std::ostream &os)
Do post-checkpoint cleanup for dir once writes complete.
Definition AMReX_AmrLevel.cpp:565
std::unique_ptr< FabFactory< FArrayBox > > m_factory
Definition AMReX_AmrLevel.H:734
void finishConstructor()
Common code used by all constructors.
Definition AMReX_AmrLevel.cpp:473
void RK(int order, int state_type, Real time, Real dt, int iteration, int ncycle, F &&f, P &&p=RungeKutta::PostStageNoOp())
Evolve one step with Runge-Kutta (2, 3, or 4)
Definition AMReX_AmrLevel.H:912
virtual void CreateLevelDirectory(const std::string &dir)
Create the Level_* directory inside dir (if needed) before writing plots/checkpoints.
Definition AMReX_AmrLevel.cpp:2249
virtual Real advance(Real time, Real dt, int iteration, int ncycle)=0
Advance solution on this level by dt and return the next stable dt.
static void SetEBMaxGrowCells(int nbasic, int nvolume, int nfull) noexcept
Configure the ghost-cell requirements (basic/volume/full) for EB FillPatch operations.
Definition AMReX_AmrLevel.H:673
void FillCoarsePatch(MultiFab &mf, int dcomp, Real time, int state_idx, int scomp, int ncomp, int nghost=0)
Interpolate from coarse level to the valid area in mf.
Definition AMReX_AmrLevel.cpp:1523
int post_step_regrid
Definition AMReX_AmrLevel.H:730
static EBSupport m_eb_support_level
Definition AMReX_AmrLevel.H:683
BoxArray m_AreaNotToTag
Definition AMReX_AmrLevel.H:727
virtual int okToRegrid()
Return nonzero if refinement starting from this level is allowed this cycle.
Definition AMReX_AmrLevel.cpp:1895
static void FillPatchAdd(AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
Same as FillPatch(), but adds into leveldata instead of overwriting.
Definition AMReX_AmrLevel.cpp:2215
virtual void checkPointPre(const std::string &dir, std::ostream &os)
Do pre-checkpoint work (e.g., directory creation) for dir using os.
Definition AMReX_AmrLevel.cpp:557
int numStates() const noexcept
Number of states at this level.
Definition AMReX_AmrLevel.H:356
virtual void writePlotFile(const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles)
Write plot file data for this level into directory dir.
Definition AMReX_AmrLevel.cpp:134
virtual std::string thePlotFileType() const
A string written as the first item in writePlotFile() at level zero. It is so we can distinguish betw...
Definition AMReX_AmrLevel.H:99
DistributionMapping dmap
Definition AMReX_AmrLevel.H:719
AmrLevel(AmrLevel &&)=delete
Disable moving; derived classes manage non-movable FillPatch caches.
static int m_eb_basic_grow_cells
Definition AMReX_AmrLevel.H:678
Vector< StateData > state
Definition AMReX_AmrLevel.H:725
const MultiFab & get_new_data(int state_indx) const noexcept
Const access to the MultiFab at the new time slot for state state_indx.
Definition AMReX_AmrLevel.H:467
int nStep() const noexcept
Timestep n at this level.
Definition AMReX_AmrLevel.H:360
virtual void post_init(Real stop_time)=0
Operations to run after initialization finishes.
static void FlushFPICache()
Resize/clear cached FillPatchIterator state in preparation for new fill operations.
void setAreaNotToTag(BoxArray &ba) noexcept
Override the area-not-to-tag cache with user-supplied BoxArray ba.
Definition AMReX_AmrLevel.cpp:2090
virtual void particle_redistribute(int lbase=0, bool a_init=false)
Optional particle redistribution hook called by the parent Amr object.
Definition AMReX_AmrLevel.H:583
static bool isStateVariable(const std::string &name, int &state_indx, int &ncomp)
Return true if name matches a registered state variable, also providing descriptor info.
Definition AMReX_AmrLevel.cpp:487
static DescriptorList desc_lst
Definition AMReX_AmrLevel.H:724
const BoxArray & getEdgeBoxArray(int dir) const noexcept
Return the edge-centered grid layout for direction dir (cached).
Definition AMReX_AmrLevel.cpp:626
AmrLevel(const AmrLevel &)=delete
Disable copying because AmrLevel owns state data and geometry caches.
virtual void set_preferred_boundary_values(MultiFab &S, int state_index, int scomp, int dcomp, int ncomp, Real time) const
Hook to override fillpatched boundary values before applying physical BCs (default no-op).
Definition AMReX_AmrLevel.cpp:48
virtual void setSmallPlotVariables()
Override to customize which variables appear in small plotfiles on this level.
Definition AMReX_AmrLevel.cpp:1961
AmrLevel() noexcept
The constructors – for derived classes.
Definition AMReX_AmrLevel.H:701
static const DescriptorList & get_desc_lst() noexcept
Returns list of Descriptors.
Definition AMReX_AmrLevel.H:469
int postStepRegrid() const noexcept
Returns whether a post-timestep regrid is currently requested (nonzero for true).
Definition AMReX_AmrLevel.H:473
virtual void init()=0
int numGrids() const noexcept
Number of grids at this level.
Definition AMReX_AmrLevel.H:354
Amr * parent
Definition AMReX_AmrLevel.H:720
MultiFab & get_new_data(int state_indx) noexcept
Mutable access to the MultiFab at the new time slot for state state_indx.
Definition AMReX_AmrLevel.H:465
MultiFab & get_old_data(int state_indx) noexcept
Mutable access to the MultiFab at the old time slot for state state_indx.
Definition AMReX_AmrLevel.H:461
void resetFillPatcher()
Clear cached FillPatcher objects (call after reflux/averaging).
Definition AMReX_AmrLevel.cpp:2132
virtual ~AmrLevel()
The destructor.
Definition AMReX_AmrLevel.cpp:572
void FillPatcherFill(amrex::MultiFab &mf, int dcomp, int ncomp, int nghost, amrex::Real time, int state_index, int scomp)
Fill ghost cells using cached FillPatcher (level>0) or FillPatch (level 0).
Definition AMReX_AmrLevel.cpp:2140
virtual void post_restart()
Override to perform any fixups immediately after restart (default: no-op).
Definition AMReX_AmrLevel.H:272
const BoxArray & getAreaNotToTag() noexcept
Get the cached BoxArray describing regions that should not be tagged.
Definition AMReX_AmrLevel.cpp:2080
const Vector< DistributionMapping > & DistributionMap() const noexcept
Definition AMReX_AmrMesh.H:150
int finestLevel() const noexcept
Return the finest level.
Definition AMReX_AmrMesh.H:138
const Vector< Geometry > & Geom() const noexcept
Definition AMReX_AmrMesh.H:149
const Vector< BoxArray > & boxArray() const noexcept
Definition AMReX_AmrMesh.H:151
Manage hierarchy of levels for time-dependent AMR computations.
Definition AMReX_Amr.H:41
AmrLevel & getLevel(int lev) noexcept
Access the AmrLevel object stored at index lev.
Definition AMReX_Amr.H:298
int levelSteps(int lev) const noexcept
Number of time steps completed so far on level lev.
Definition AMReX_Amr.H:146
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:615
A list of DeriveRecs.
Definition AMReX_Derive.H:411
Definition AMReX_StateDescriptor.H:423
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
Definition AMReX_EBFabFactory.H:25
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
int nGrow(int direction=0) const noexcept
Return the grow factor that defines the region of definition.
Definition AMReX_FabArrayBase.H:78
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:83
Definition AMReX_FabFactory.H:50
Definition AMReX_AmrLevel.H:833
void Initialize(int boxGrow, Real time, int idx, int scomp, int ncomp, InterpBase *mapper)
Prepare another fill using cached metadata.
Definition AMReX_AmrLevel.cpp:741
void fill(FArrayBox &fab, int dcomp, int idx)
Fill the provided FAB using cached data.
Definition AMReX_AmrLevel.cpp:1268
Definition AMReX_AmrLevel.H:760
FArrayBox & operator()() noexcept
Access the filled FAB for the current tile.
Definition AMReX_AmrLevel.H:809
MultiFab & get_mf() noexcept
Mutable reference to the internal MultiFab holding filled tiles.
Definition AMReX_AmrLevel.H:815
Box UngrownBox() const noexcept
Return the un-grown valid box for the current tile.
Definition AMReX_AmrLevel.H:812
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:211
Definition AMReX_InterpBase.H:34
Virtual base class for interpolaters.
Definition AMReX_Interpolater.H:27
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
Box validbox() const noexcept
Return the valid Box in which the current tile resides.
Definition AMReX_MFIter.H:160
int index() const noexcept
The index into the underlying BoxArray of the current FAB.
Definition AMReX_MFIter.H:172
void Initialize()
Definition AMReX_MFIter.cpp:278
Definition AMReX_MFCopyDescriptor.H:46
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Current and previous level-time data.
Definition AMReX_StateData.H:38
Attributes of StateData.
Definition AMReX_StateDescriptor.H:33
An array of TagBoxes.
Definition AMReX_TagBox.H:151
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
How
How we write out FabArray<FArrayBox>s. These are deprecated, we always use NFiles....
Definition AMReX_VisMF.H:42
@ NFiles
Definition AMReX_VisMF.H:42
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
std::array< T, N > Array
Definition AMReX_Array.H:26
void RK2(MF &Uold, MF &Unew, Real time, Real dt, F const &frhs, FB const &fillbndry, P const &post_stage=PostStageNoOp())
Time stepping with RK2.
Definition AMReX_RungeKutta.H:170
void RK4(MF &Uold, MF &Unew, Real time, Real dt, F const &frhs, FB const &fillbndry, R const &store_crse_data, P const &post_stage=PostStageNoOp())
Time stepping with RK4.
Definition AMReX_RungeKutta.H:258
void RK3(MF &Uold, MF &Unew, Real time, Real dt, F const &frhs, FB const &fillbndry, R const &store_crse_data, P const &post_stage=PostStageNoOp())
Time stepping with RK3.
Definition AMReX_RungeKutta.H:208
Definition AMReX_Amr.cpp:49
EBSupport
Definition AMReX_EBSupport.H:7
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
Definition AMReX_RungeKutta.H:51