Block-Structured AMR Software Framework
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>
9 #include <AMReX_MultiFabUtil.H>
10 #include <AMReX_LayoutData.H>
11 #include <AMReX_Derive.H>
12 #include <AMReX_BCRec.H>
13 #include <AMReX_Amr.H>
15 #include <AMReX_StateDescriptor.H>
16 #include <AMReX_StateData.H>
17 #include <AMReX_VisMF.H>
18 #include <AMReX_RungeKutta.H>
19 #include <AMReX_FillPatcher.H>
20 #ifdef AMREX_USE_EB
21 #include <AMReX_EBSupport.H>
22 #endif
23 
24 #include <memory>
25 #include <map>
26 
27 namespace amrex {
28 
29 class TagBox;
30 class TagBoxArray;
31 
37 class AmrLevel
38 {
39  friend class Amr;
40  friend class FillPatchIterator;
42 
43 public:
52  virtual ~AmrLevel ();
53 
54  AmrLevel (const AmrLevel&) = delete;
55  AmrLevel (AmrLevel&&) = delete;
56  AmrLevel& operator= (const AmrLevel&) = delete;
58 
60  void LevelDirectoryNames (const std::string &dir,
61  std::string &LevelDir,
62  std::string &FullPath) const;
64  virtual void CreateLevelDirectory (const std::string &dir);
69  void SetLevelDirectoryCreated(bool ldc) noexcept { levelDirectoryCreated = ldc; }
78  virtual std::string thePlotFileType () const
79  {
80  static const std::string the_plot_file_type("HyperCLaw-V1.1");
81  return the_plot_file_type;
82  }
86  virtual void writePlotFile (const std::string& dir,
87  std::ostream& os,
89 
91  virtual void writePlotFilePre (const std::string& dir,
92  std::ostream& os);
93 
95  virtual void writePlotFilePost (const std::string& dir,
96  std::ostream& os);
97 
101  virtual void writeSmallPlotFile (const std::string& /*dir*/,
102  std::ostream& /*os*/,
103  VisMF::How /*how*/ = VisMF::NFiles) {}
105  virtual void checkPoint (const std::string& dir,
106  std::ostream& os,
108  bool dump_old = true);
110  virtual void checkPointPre (const std::string& dir,
111  std::ostream& os);
113  virtual void checkPointPost (const std::string& dir,
114  std::ostream& os);
116  virtual void restart (Amr& papa,
117  std::istream& is,
118  bool bReadSpecial = false);
119 
121  virtual void set_state_in_checkpoint (Vector<int>& state_in_checkpoint);
122 
124  static bool isStateVariable (const std::string& name,
125  int& state_indx,
126  int& ncomp);
127 
128  static void FlushFPICache ();
133  virtual void computeInitialDt (int finest_level,
134  int sub_cycle,
135  Vector<int>& n_cycle,
136  const Vector<IntVect>& ref_ratio,
137  Vector<Real>& dt_level,
138  Real stop_time) = 0;
143  virtual void computeNewDt (int finest_level,
144  int sub_cycle,
145  Vector<int>& n_cycle,
146  const Vector<IntVect>& ref_ratio,
147  Vector<Real>& dt_min,
148  Vector<Real>& dt_level,
149  Real stop_time,
150  int post_regrid_flag) = 0;
156  virtual Real advance (Real time,
157  Real dt,
158  int iteration,
159  int ncycle) = 0;
160 
165  virtual void post_timestep (int iteration);
170  virtual void postCoarseTimeStep (Real time);
174  virtual void post_restart () {}
180  virtual void post_regrid (int lbase,
181  int new_finest) = 0;
187  virtual void post_init (Real stop_time) = 0;
191  virtual int okToContinue () { return 1; }
197  virtual int okToRegrid ();
203  virtual void initData () = 0;
205  virtual void setTimeLevel (Real time,
206  Real dt_old,
207  Real dt_new);
209  virtual void allocOldData ();
211  virtual void removeOldData ();
217  virtual void init (AmrLevel &old) = 0;
223  virtual void init () = 0;
225  void reset ();
227  int Level () const noexcept { return level; }
229  const BoxArray& boxArray () const noexcept { return grids; }
230  const BoxArray& getEdgeBoxArray (int dir) const noexcept;
231  const BoxArray& getNodalBoxArray () const noexcept;
232  //
233  const DistributionMapping& DistributionMap () const noexcept { return dmap; }
234  //
235  const FabFactory<FArrayBox>& Factory () const noexcept { return *m_factory; }
236  //
237 #ifdef AMREX_USE_EB
238  const EBFArrayBoxFactory&
239  EBFactory () const noexcept {
240  return static_cast<amrex::EBFArrayBoxFactory const&>(*m_factory);
241  }
242 #endif
243 
245  int numGrids () const noexcept { return static_cast<int>(grids.size()); }
247  int numStates () const noexcept { return static_cast<int>(state.size()); }
249  const Box& Domain () const noexcept { return geom.Domain(); }
251  int nStep () const noexcept { return parent->levelSteps(level); }
253  const Geometry& Geom () const noexcept { return geom; }
254  //
255  const IntVect& fineRatio () const noexcept { return fine_ratio; }
257  Long countCells () const noexcept;
258 
260  const BoxArray& getAreaNotToTag () noexcept;
261  const Box& getAreaToTag () noexcept;
263  void constructAreaNotToTag ();
265  void setAreaNotToTag (BoxArray& ba) noexcept;
266 
267  void resetFillPatcher ();
268 
273  virtual void errorEst (TagBoxArray& tb,
274  int clearval,
275  int tagval,
276  Real time,
277  int n_error_buf = 0,
278  int ngrow = 0) = 0;
280  void FillCoarsePatch (MultiFab& mf,
281  int dcomp,
282  Real time,
283  int state_idx,
284  int scomp,
285  int ncomp,
286  int nghost = 0);
288  virtual void setPhysBoundaryValues (FArrayBox& dest,
289  int state_indx,
290  Real time,
291  int dest_comp,
292  int src_comp,
293  int num_comp);
300  virtual std::unique_ptr<MultiFab> derive (const std::string& name,
301  Real time,
302  int ngrow);
307  virtual void derive (const std::string& name,
308  Real time,
309  MultiFab& mf,
310  int dcomp);
312  StateData& get_state_data (int state_indx) noexcept { return state[state_indx]; }
314  MultiFab& get_old_data (int state_indx) noexcept { return state[state_indx].oldData(); }
316  const MultiFab& get_old_data (int state_indx) const noexcept { return state[state_indx].oldData(); }
318  MultiFab& get_new_data (int state_indx) noexcept { return state[state_indx].newData(); }
320  const MultiFab& get_new_data (int state_indx) const noexcept { return state[state_indx].newData(); }
322  static const DescriptorList& get_desc_lst () noexcept { return desc_lst; }
324  static DeriveList& get_derive_lst () noexcept;
326  int postStepRegrid () const noexcept { return post_step_regrid; }
328  void setPostStepRegrid (int new_val) noexcept { post_step_regrid = new_val; }
329 
332 
334  Vector<int> getBCArray (int State_Type,
335  int gridno,
336  int strt_comp,
337  int ncomp);
339  MultiFab& get_data (int state_indx, Real time) noexcept;
341  virtual void set_preferred_boundary_values (MultiFab& S,
342  int state_index,
343  int scomp,
344  int dcomp,
345  int ncomp,
346  Real time) const;
351  virtual void manual_tags_placement (TagBoxArray& tags,
352  const Vector<IntVect>& bf_lev);
354  virtual void setPlotVariables ();
356  virtual void setSmallPlotVariables ();
362  virtual Real estimateWork();
363 
365  virtual int WorkEstType () { return -1; }
366 
371  TimeLevel which_time (int state_indx, Real time) const noexcept;
372 
374  virtual bool writePlotNow ();
375 
377  virtual bool writeSmallPlotNow ();
378 
379 #ifdef AMREX_PARTICLES
381  virtual void particle_redistribute (int /*lbase*/ = 0, bool /*a_init*/ = false) {;}
382 #endif
383 
395  void FillPatcherFill (amrex::MultiFab& mf, int dcomp, int ncomp, int nghost,
396  amrex::Real time, int state_index, int scomp);
397 
398  static void FillPatch (AmrLevel& amrlevel,
399  MultiFab& leveldata,
400  int boxGrow,
401  Real time,
402  int index,
403  int scomp,
404  int ncomp,
405  int dcomp=0);
406 
407  static void FillPatchAdd (AmrLevel& amrlevel,
408  MultiFab& leveldata,
409  int boxGrow,
410  Real time,
411  int index,
412  int scomp,
413  int ncomp,
414  int dcomp=0);
415 
437  template <typename F, typename P = RungeKutta::PostStageNoOp>
438  void RK (int order, int state_type, Real time, Real dt, int iteration,
439  int ncycle, F&& f, P&& p = RungeKutta::PostStageNoOp());
440 
441 #ifdef AMREX_USE_EB
442  static void SetEBMaxGrowCells (int nbasic, int nvolume, int nfull) noexcept {
443  m_eb_basic_grow_cells = nbasic;
444  m_eb_volume_grow_cells = nvolume;
445  m_eb_full_grow_cells = nfull;
446  }
447  static int m_eb_basic_grow_cells;
448  static int m_eb_volume_grow_cells;
449  static int m_eb_full_grow_cells;
450  static void SetEBSupportLevel (EBSupport ebs) { m_eb_support_level = ebs; }
451  static EBSupport m_eb_support_level;
452 #endif
453 
455  static IntVect ProperBlockingFactor (AmrLevel const& amr_level, int boxGrow,
456  IndexType const& boxType,
457  StateDescriptor const& desc, int SComp);
458 
459 protected:
461  AmrLevel () noexcept {} // NOLINT
462 
463  AmrLevel (Amr& papa,
464  int lev,
465  const Geometry& level_geom,
466  const BoxArray& ba,
467  const DistributionMapping& dm,
468  Real time);
469 
471  void finishConstructor ();
472 
473  //
474  // The Data.
475  //
476  int level{-1}; // AMR level (0 is coarsest).
477  Geometry geom; // Geom at this level.
478  BoxArray grids; // Cell-centered locations of grids.
479  DistributionMapping dmap; // Distribution of grids among processes
480  Amr* parent{nullptr};// Pointer to parent AMR structure.
481  IntVect crse_ratio; // Refinement ratio to coarser level.
482  IntVect fine_ratio; // Refinement ratio to finer level.
483  static DeriveList derive_lst; // List of derived quantities.
484  static DescriptorList desc_lst; // List of state variables.
485  Vector<StateData> state; // Array of state data.
486 
487  BoxArray m_AreaNotToTag; //Area which shouldn't be tagged on this level.
488  Box m_AreaToTag; //Area which is allowed to be tagged on this level.
489 
490  int post_step_regrid{0}; // Whether or not to do a regrid after the timestep.
491 
492  bool levelDirectoryCreated{false}; // for checkpoints and plotfiles
493 
494  std::unique_ptr<FabFactory<FArrayBox> > m_factory;
495 
497 
498 private:
499 
500  template <std::size_t order>
501  void storeRKCoarseData (int state_type, Real time, Real dt,
502  MultiFab const& S_old,
503  Array<MultiFab,order> const& rkk);
504 
505  void FillRKPatch (int state_index, MultiFab& S, Real time,
506  int stage, int iteration, int ncycle);
507 
508  mutable BoxArray edge_grids[AMREX_SPACEDIM]; // face-centered grids
509  mutable BoxArray nodal_grids; // all nodal grids
510 };
511 
512 //
513 // Forward declaration.
514 //
516 
518  :
519  public MFIter
520 {
521  public:
522 
523  friend class AmrLevel;
524 
525  FillPatchIterator (AmrLevel& amrlevel,
526  MultiFab& leveldata);
527 
528  FillPatchIterator (AmrLevel& amrlevel,
529  MultiFab& leveldata,
530  int boxGrow,
531  Real time,
532  int idx,
533  int scomp,
534  int ncomp);
535 
536  void Initialize (int boxGrow,
537  Real time,
538  int idx,
539  int scomp,
540  int ncomp);
541 
542  FArrayBox& operator() () noexcept { return m_fabs[MFIter::index()]; }
543 
544  Box UngrownBox () const noexcept { return MFIter::validbox(); }
545 
546  MultiFab& get_mf() noexcept { return m_fabs; }
547 
548 private:
549 
550  void FillFromLevel0 (Real time, int idx, int scomp, int dcomp, int ncomp);
551  void FillFromTwoLevels (Real time, int idx, int scomp, int dcomp, int ncomp);
552 
553  //
554  // The data.
555  //
558  std::vector< std::pair<int,int> > m_range;
560  int m_ncomp;
561 };
562 
564 {
565 public:
566 
567  friend class FillPatchIterator;
568 
570  MultiFab& leveldata);
571 
573  MultiFab& leveldata,
574  int boxGrow,
575  Real time,
576  int index,
577  int scomp,
578  int ncomp,
579  InterpBase* mapper);
580 
581  void Initialize (int boxGrow,
582  Real time,
583  int idx,
584  int scomp,
585  int ncomp,
586  InterpBase* mapper);
587 
588  void fill (FArrayBox& fab, int dcomp, int idx);
589 
590 private:
591  //
592  // The data.
593  //
597  Vector< Vector<MultiFabId> > m_mfid; // [level][oldnew]
598  Interpolater* m_map = nullptr;
599  std::map<int,Box> m_ba;
600  Real m_time = std::numeric_limits<Real>::lowest();
601  int m_growsize = std::numeric_limits<int>::lowest();
602  int m_index = -1;
603  int m_scomp = -1;
604  int m_ncomp = -1;
605  bool m_FixUpCorners = false;
606 
607  std::map< int,Vector< Vector<Box> > > m_fbox; // [grid][level][validregion]
608  std::map< int,Vector< Vector<Box> > > m_cbox; // [grid][level][fillablesubbox]
609  std::map< int,Vector< Vector< Vector<FillBoxId> > > > m_fbid; // [grid][level][fillablesubbox][oldnew]
610 };
611 
612 template <typename F, typename P>
613 void AmrLevel::RK (int order, int state_type, Real time, Real dt, int iteration,
614  int ncycle, F&& f, P&& p)
615 {
616  BL_PROFILE("AmrLevel::RK()");
617 
618  AMREX_ASSERT(AmrLevel::desc_lst[state_type].nExtra() > 0); // Need ghost cells in StateData
619 
620  MultiFab& S_old = get_old_data(state_type);
621  MultiFab& S_new = get_new_data(state_type);
622  const Real t_old = state[state_type].prevTime();
623  const Real t_new = state[state_type].curTime();
624  AMREX_ALWAYS_ASSERT(amrex::almostEqual(time,t_old,10) && amrex::almostEqual(time+dt,t_new,10));
625 
626  if (order == 2) {
627  RungeKutta::RK2(S_old, S_new, time, dt, std::forward<F>(f),
628  [&] (int /*stage*/, MultiFab& mf, Real t) {
629  FillPatcherFill(mf, 0, mf.nComp(), mf.nGrow(), t,
630  state_type, 0); },
631  std::forward<P>(p));
632  } else if (order == 3) {
633  RungeKutta::RK3(S_old, S_new, time, dt, std::forward<F>(f),
634  [&] (int stage, MultiFab& mf, Real t) {
635  FillRKPatch(state_type, mf, t, stage, iteration, ncycle);
636  },
637  [&] (Array<MultiFab,3> const& rkk) {
638  if (level < parent->finestLevel()) {
639  storeRKCoarseData(state_type, time, dt, S_old, rkk);
640  }
641  },
642  std::forward<P>(p));
643  } else if (order == 4) {
644  RungeKutta::RK4(S_old, S_new, time, dt, std::forward<F>(f),
645  [&] (int stage, MultiFab& mf, Real t) {
646  FillRKPatch(state_type, mf, t, stage, iteration, ncycle);
647  },
648  [&] (Array<MultiFab,4> const& rkk) {
649  if (level < parent->finestLevel()) {
650  storeRKCoarseData(state_type, time, dt, S_old, rkk);
651  }
652  },
653  std::forward<P>(p));
654  } else {
655  amrex::Abort("AmrLevel::RK: order = "+std::to_string(order)+" is not supported");
656  }
657 }
658 
659 template <std::size_t order>
660 void AmrLevel::storeRKCoarseData (int state_type, Real time, Real dt,
661  MultiFab const& S_old,
662  Array<MultiFab,order> const& rkk)
663 {
664  BL_PROFILE("AmrLevel::storeRKCoarseData()");
665  if (level == parent->finestLevel()) { return; }
666 
667  const StateDescriptor& desc = AmrLevel::desc_lst[state_type];
668 
669  auto& fillpatcher = parent->getLevel(level+1).m_fillpatcher[state_type];
670  fillpatcher = std::make_unique<FillPatcher<MultiFab>>
672  parent->Geom(level+1),
674  parent->Geom(level),
675  IntVect(desc.nExtra()), desc.nComp(), desc.interp(0));
676 
677  fillpatcher->storeRKCoarseData(time, dt, S_old, rkk);
678 }
679 
680 
681 }
682 
683 #endif /*_AmrLevel_H_*/
#define BL_PROFILE(a)
Definition: AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition: AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition: AMReX_BLassert.H:50
Virtual base class for managing individual levels. AmrLevel functions both as a container for state d...
Definition: AMReX_AmrLevel.H:38
virtual void post_regrid(int lbase, int new_finest)=0
Operations to be done after regridding This is a pure virtual function and hence MUST be implemented ...
virtual Real estimateWork()
Estimate the amount of work required to advance Just this level based on the number of cells....
Definition: AMReX_AmrLevel.cpp:2053
virtual void restart(Amr &papa, std::istream &is, bool bReadSpecial=false)
Restart from a checkpoint file.
Definition: AMReX_AmrLevel.cpp:386
const MultiFab & get_new_data(int state_indx) const noexcept
State data at new time.
Definition: AMReX_AmrLevel.H:320
void reset()
Reset data to initial time by swapping new and old time data.
Definition: AMReX_AmrLevel.cpp:596
const DistributionMapping & DistributionMap() const noexcept
Definition: AMReX_AmrLevel.H:233
Long countCells() const noexcept
Returns number of cells on level.
Definition: AMReX_AmrLevel.cpp:504
IntVect fine_ratio
Definition: AMReX_AmrLevel.H:482
int level
Definition: AMReX_AmrLevel.H:476
virtual void checkPoint(const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles, bool dump_old=true)
Write current state to checkpoint file.
Definition: AMReX_AmrLevel.cpp:510
int Level() const noexcept
Returns this AmrLevel.
Definition: AMReX_AmrLevel.H:227
virtual void writePlotFilePre(const std::string &dir, std::ostream &os)
Do pre-plotfile work to avoid synchronizations while writing the amr hierarchy.
Definition: AMReX_AmrLevel.cpp:372
const FabFactory< FArrayBox > & Factory() const noexcept
Definition: AMReX_AmrLevel.H:235
virtual int WorkEstType()
Which state data type is for work estimates? -1 means none.
Definition: AMReX_AmrLevel.H:365
virtual bool writePlotNow()
Does the AmrLevel want Amr to write a plotfile now?
Definition: AMReX_AmrLevel.cpp:2059
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 IntVect ProperBlockingFactor(AmrLevel const &amr_level, int boxGrow, IndexType const &boxType, StateDescriptor const &desc, int SComp)
Recommendation of a proper blocking factor.
Definition: AMReX_AmrLevel.cpp:2276
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. This is a pure virtual function and hence MUST be implemented by deriv...
void FillRKPatch(int state_index, MultiFab &S, Real time, int stage, int iteration, int ncycle)
Definition: AMReX_AmrLevel.cpp:2256
Vector< std::unique_ptr< FillPatcher< MultiFab > > > m_fillpatcher
Definition: AMReX_AmrLevel.H:496
virtual bool writeSmallPlotNow()
Does the AmrLevel want Amr to write a small plotfile now?
Definition: AMReX_AmrLevel.cpp:2065
void LevelDirectoryNames(const std::string &dir, std::string &LevelDir, std::string &FullPath) const
Get the level directory names.
Definition: AMReX_AmrLevel.cpp:2223
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. This is a pure virtual function and hence MUST be implemented by derived ...
virtual std::unique_ptr< MultiFab > derive(const std::string &name, Real time, int ngrow)
Returns a MultiFab containing the derived data for this level. The user is responsible for deleting t...
Definition: AMReX_AmrLevel.cpp:1615
virtual void manual_tags_placement(TagBoxArray &tags, const Vector< IntVect > &bf_lev)
Called in grid_places after other tagging routines to modify the list of tagged points....
Definition: AMReX_AmrLevel.cpp:63
virtual void setPhysBoundaryValues(FArrayBox &dest, int state_indx, Real time, int dest_comp, int src_comp, int num_comp)
Function to set physical boundary conditions.
Definition: AMReX_AmrLevel.cpp:647
void SetLevelDirectoryCreated(bool ldc) noexcept
Set if the Level_ directory was created or to clear the value. CreateLevelDirectory sets levelDirecto...
Definition: AMReX_AmrLevel.H:69
void setPostStepRegrid(int new_val) noexcept
Sets a new value for the post-timestep regrid trigger.
Definition: AMReX_AmrLevel.H:328
MultiFab & get_old_data(int state_indx) noexcept
State data at old time.
Definition: AMReX_AmrLevel.H:314
virtual void writePlotFilePost(const std::string &dir, std::ostream &os)
Do post-plotfile work to avoid synchronizations while writing the amr hierarchy.
Definition: AMReX_AmrLevel.cpp:379
const IntVect & fineRatio() const noexcept
Definition: AMReX_AmrLevel.H:255
static const DescriptorList & get_desc_lst() noexcept
Returns list of Descriptors.
Definition: AMReX_AmrLevel.H:322
AmrLevel & operator=(const AmrLevel &)=delete
BoxArray nodal_grids
Definition: AMReX_AmrLevel.H:509
Vector< int > getBCArray(int State_Type, int gridno, int strt_comp, int ncomp)
Boundary condition access function.
Definition: AMReX_AmrLevel.cpp:1869
Box m_AreaToTag
Definition: AMReX_AmrLevel.H:488
static void FillPatch(AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
Definition: AMReX_AmrLevel.cpp:2187
virtual void removeOldData()
Delete old-time data.
Definition: AMReX_AmrLevel.cpp:587
Geometry geom
Definition: AMReX_AmrLevel.H:477
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:2018
virtual void setTimeLevel(Real time, Real dt_old, Real dt_new)
Set the time levels of state data.
Definition: AMReX_AmrLevel.cpp:476
virtual int okToContinue()
Is it ok to continue the calculation?
Definition: AMReX_AmrLevel.H:191
const BoxArray & getNodalBoxArray() const noexcept
Definition: AMReX_AmrLevel.cpp:637
MultiFab & get_data(int state_indx, Real time) noexcept
Get state data at specified index and time.
Definition: AMReX_AmrLevel.cpp:605
BoxArray grids
Definition: AMReX_AmrLevel.H:478
virtual void allocOldData()
Alloc space for old time data.
Definition: AMReX_AmrLevel.cpp:578
virtual void set_state_in_checkpoint(Vector< int > &state_in_checkpoint)
Old checkpoint may have different number of states than the new source code.
Definition: AMReX_AmrLevel.cpp:467
const MultiFab & get_old_data(int state_indx) const noexcept
State data at old time.
Definition: AMReX_AmrLevel.H:316
void UpdateDistributionMaps(DistributionMapping &dmap)
Update the distribution maps in StateData based on the size of the map.
Definition: AMReX_AmrLevel.cpp:1852
virtual void postCoarseTimeStep(Real time)
Contains operations to be done only after a full coarse timestep. The default implementation does not...
Definition: AMReX_AmrLevel.cpp:43
virtual void post_timestep(int iteration)
Contains operations to be done after a timestep. If this function is overridden, don't forget to rese...
Definition: AMReX_AmrLevel.cpp:35
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:483
void constructAreaNotToTag()
Construct the area not to tag.
Definition: AMReX_AmrLevel.cpp:2085
bool levelDirectoryCreated
Definition: AMReX_AmrLevel.H:492
BoxArray edge_grids[AMREX_SPACEDIM]
Definition: AMReX_AmrLevel.H:508
const BoxArray & boxArray() const noexcept
List of grids at this level.
Definition: AMReX_AmrLevel.H:229
static DeriveList & get_derive_lst() noexcept
Returns list of derived variables.
Definition: AMReX_AmrLevel.cpp:57
IntVect crse_ratio
Definition: AMReX_AmrLevel.H:481
virtual void setPlotVariables()
Modify list of variables to be plotted.
Definition: AMReX_AmrLevel.cpp:1897
const Box & getAreaToTag() noexcept
Definition: AMReX_AmrLevel.cpp:2075
virtual void writeSmallPlotFile(const std::string &, std::ostream &, VisMF::How=VisMF::NFiles)
Write small plot file stuff to specified directory.
Definition: AMReX_AmrLevel.H:101
TimeLevel
What time are we at?
Definition: AMReX_AmrLevel.H:45
@ Amr1QtrTime
Definition: AMReX_AmrLevel.H:48
@ AmrOldTime
Definition: AMReX_AmrLevel.H:45
@ AmrNewTime
Definition: AMReX_AmrLevel.H:47
@ AmrOtherTime
Definition: AMReX_AmrLevel.H:50
@ AmrHalfTime
Definition: AMReX_AmrLevel.H:46
@ Amr3QtrTime
Definition: AMReX_AmrLevel.H:49
virtual void initData()=0
Init grid data at problem start-up. This is a pure virtual function and hence MUST be implemented by ...
virtual void checkPointPost(const std::string &dir, std::ostream &os)
Do post-checkpoint work to avoid synchronizations while writing the amr hierarchy.
Definition: AMReX_AmrLevel.cpp:565
std::unique_ptr< FabFactory< FArrayBox > > m_factory
Definition: AMReX_AmrLevel.H:494
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:613
virtual void CreateLevelDirectory(const std::string &dir)
Create the Level_ directory for checkpoint and plot files.
Definition: AMReX_AmrLevel.cpp:2239
virtual Real advance(Real time, Real dt, int iteration, int ncycle)=0
Do an integration step on this level. Returns maximum safe time step. This is a pure virtual function...
StateData & get_state_data(int state_indx) noexcept
State data object.
Definition: AMReX_AmrLevel.H:312
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:1519
int post_step_regrid
Definition: AMReX_AmrLevel.H:490
BoxArray m_AreaNotToTag
Definition: AMReX_AmrLevel.H:487
virtual int okToRegrid()
Should I regrid with this level as base level? This test is only evaluated if regrid_int > 0 and leve...
Definition: AMReX_AmrLevel.cpp:1891
static void FillPatchAdd(AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
Definition: AMReX_AmrLevel.cpp:2205
virtual void checkPointPre(const std::string &dir, std::ostream &os)
Do pre-checkpoint work to avoid synchronizations while writing the amr hierarchy.
Definition: AMReX_AmrLevel.cpp:557
int numStates() const noexcept
Number of states at this level.
Definition: AMReX_AmrLevel.H:247
MultiFab & get_new_data(int state_indx) noexcept
State data at new time.
Definition: AMReX_AmrLevel.H:318
virtual void writePlotFile(const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles)
Write plot file stuff to specified directory.
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:78
DistributionMapping dmap
Definition: AMReX_AmrLevel.H:479
AmrLevel(AmrLevel &&)=delete
void storeRKCoarseData(int state_type, Real time, Real dt, MultiFab const &S_old, Array< MultiFab, order > const &rkk)
Definition: AMReX_AmrLevel.H:660
Vector< StateData > state
Definition: AMReX_AmrLevel.H:485
int nStep() const noexcept
Timestep n at this level.
Definition: AMReX_AmrLevel.H:251
virtual void post_init(Real stop_time)=0
Operations to be done after initialization. This is a pure virtual function and hence MUST be impleme...
static void FlushFPICache()
void setAreaNotToTag(BoxArray &ba) noexcept
Set the area not to tag.
Definition: AMReX_AmrLevel.cpp:2080
static bool isStateVariable(const std::string &name, int &state_indx, int &ncomp)
Is name a state variable?
Definition: AMReX_AmrLevel.cpp:487
static DescriptorList desc_lst
Definition: AMReX_AmrLevel.H:484
const BoxArray & getEdgeBoxArray(int dir) const noexcept
Definition: AMReX_AmrLevel.cpp:626
const Geometry & Geom() const noexcept
Returns the geometry object.
Definition: AMReX_AmrLevel.H:253
AmrLevel(const AmrLevel &)=delete
virtual void set_preferred_boundary_values(MultiFab &S, int state_index, int scomp, int dcomp, int ncomp, Real time) const
Hack to allow override of (non-fine-fine) fillpatched boundary data.
Definition: AMReX_AmrLevel.cpp:48
virtual void setSmallPlotVariables()
Modify list of variables to be plotted.
Definition: AMReX_AmrLevel.cpp:1957
AmrLevel() noexcept
The constructors – for derived classes.
Definition: AMReX_AmrLevel.H:461
int postStepRegrid() const noexcept
Returns whether or not we want a post-timestep regrid.
Definition: AMReX_AmrLevel.H:326
virtual void init()=0
int numGrids() const noexcept
Number of grids at this level.
Definition: AMReX_AmrLevel.H:245
Amr * parent
Definition: AMReX_AmrLevel.H:480
void resetFillPatcher()
Definition: AMReX_AmrLevel.cpp:2122
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 with FillPatcher on level > 0 and AmrLevel::FillPatch on level 0.
Definition: AMReX_AmrLevel.cpp:2130
const Box & Domain() const noexcept
Returns the indices defining physical domain.
Definition: AMReX_AmrLevel.H:249
virtual void post_restart()
Operations to be done after restart.
Definition: AMReX_AmrLevel.H:174
const BoxArray & getAreaNotToTag() noexcept
Get the area not to tag.
Definition: AMReX_AmrLevel.cpp:2070
const Vector< BoxArray > & boxArray() const noexcept
Definition: AMReX_AmrMesh.H:108
const Vector< Geometry > & Geom() const noexcept
Definition: AMReX_AmrMesh.H:106
int finestLevel() const noexcept
Return the finest level.
Definition: AMReX_AmrMesh.H:95
const Vector< DistributionMapping > & DistributionMap() const noexcept
Definition: AMReX_AmrMesh.H:107
Manage hierarchy of levels for time-dependent AMR computations.
Definition: AMReX_Amr.H:35
AmrLevel & getLevel(int lev) noexcept
AmrLevel lev.
Definition: AMReX_Amr.H:189
int levelSteps(int lev) const noexcept
Number of time steps at specified level.
Definition: AMReX_Amr.H:103
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:550
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition: AMReX_BoxArray.H:597
A list of DeriveRecs.
Definition: AMReX_Derive.H:364
Definition: AMReX_StateDescriptor.H:437
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
Definition: AMReX_EBFabFactory.H:22
A Fortran Array of REALs.
Definition: AMReX_FArrayBox.H:229
int nGrow(int direction=0) const noexcept
Return the grow factor that defines the region of definition.
Definition: AMReX_FabArrayBase.H:77
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition: AMReX_FabArrayBase.H:82
Definition: AMReX_AmrLevel.H:564
void Initialize(int boxGrow, Real time, int idx, int scomp, int ncomp, InterpBase *mapper)
Definition: AMReX_AmrLevel.cpp:739
FillPatchIteratorHelper(AmrLevel &amrlevel, MultiFab &leveldata)
Definition: AMReX_AmrLevel.cpp:664
std::map< int, Vector< Vector< Box > > > m_cbox
Definition: AMReX_AmrLevel.H:608
std::map< int, Vector< Vector< Vector< FillBoxId > > > > m_fbid
Definition: AMReX_AmrLevel.H:609
bool m_FixUpCorners
Definition: AMReX_AmrLevel.H:605
std::map< int, Vector< Vector< Box > > > m_fbox
Definition: AMReX_AmrLevel.H:607
int m_scomp
Definition: AMReX_AmrLevel.H:603
AmrLevel * m_amrlevel
Definition: AMReX_AmrLevel.H:594
int m_growsize
Definition: AMReX_AmrLevel.H:601
void fill(FArrayBox &fab, int dcomp, int idx)
Definition: AMReX_AmrLevel.cpp:1266
int m_index
Definition: AMReX_AmrLevel.H:602
MultiFabCopyDescriptor m_mfcd
Definition: AMReX_AmrLevel.H:596
int m_ncomp
Definition: AMReX_AmrLevel.H:604
Real m_time
Definition: AMReX_AmrLevel.H:600
Vector< Vector< MultiFabId > > m_mfid
Definition: AMReX_AmrLevel.H:597
std::map< int, Box > m_ba
Definition: AMReX_AmrLevel.H:599
MultiFab * m_leveldata
Definition: AMReX_AmrLevel.H:595
Interpolater * m_map
Definition: AMReX_AmrLevel.H:598
Definition: AMReX_AmrLevel.H:520
FArrayBox & operator()() noexcept
Definition: AMReX_AmrLevel.H:542
std::vector< std::pair< int, int > > m_range
Definition: AMReX_AmrLevel.H:558
void FillFromLevel0(Real time, int idx, int scomp, int dcomp, int ncomp)
Definition: AMReX_AmrLevel.cpp:1100
MultiFab & get_mf() noexcept
Definition: AMReX_AmrLevel.H:546
MultiFab * m_leveldata
Definition: AMReX_AmrLevel.H:557
int m_ncomp
Definition: AMReX_AmrLevel.H:560
Box UngrownBox() const noexcept
Definition: AMReX_AmrLevel.H:544
void FillFromTwoLevels(Real time, int idx, int scomp, int dcomp, int ncomp)
Definition: AMReX_AmrLevel.cpp:1119
AmrLevel * m_amrlevel
Definition: AMReX_AmrLevel.H:556
MultiFab m_fabs
Definition: AMReX_AmrLevel.H:559
FillPatchIterator(AmrLevel &amrlevel, MultiFab &leveldata)
Definition: AMReX_AmrLevel.cpp:672
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition: AMReX_Geometry.H:210
Definition: AMReX_InterpBase.H:26
Virtual base class for interpolaters.
Definition: AMReX_Interpolater.H:22
Definition: AMReX_MFIter.H:57
Box validbox() const noexcept
Return the valid Box in which the current tile resides.
Definition: AMReX_MFIter.H:132
int index() const noexcept
The index into the underlying BoxArray of the current FAB.
Definition: AMReX_MFIter.H:144
void Initialize()
Definition: AMReX_MFIter.cpp:274
Definition: AMReX_MFCopyDescriptor.H:46
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
Current and previous level-time data.
Definition: AMReX_StateData.H:33
Attributes of StateData.
Definition: AMReX_StateDescriptor.H:33
InterpBase * interp() const noexcept
Returns the interpolater.
Definition: AMReX_StateDescriptor.cpp:251
int nComp() const noexcept
Returns number of components.
Definition: AMReX_StateDescriptor.cpp:239
int nExtra() const noexcept
Returns the grow factor.
Definition: AMReX_StateDescriptor.cpp:245
An array of TagBoxes.
Definition: AMReX_TagBox.H:150
How
How we write out FabArray<FArrayBox>s. These are deprecated, we always use NFiles....
Definition: AMReX_VisMF.H:41
@ NFiles
Definition: AMReX_VisMF.H:41
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:158
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:246
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:196
static int f(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition: AMReX_SundialsIntegrator.H:44
static constexpr int P
Definition: AMReX_OpenBC.H:14
Definition: AMReX_Amr.cpp:49
EBSupport
Definition: AMReX_EBSupport.H:7
IntVectND< AMREX_SPACEDIM > IntVect
Definition: AMReX_BaseFwd.H:30
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE std::enable_if_t< std::is_floating_point_v< T >, bool > almostEqual(T x, T y, int ulp=2)
Definition: AMReX_Algorithm.H:93
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:225
IndexTypeND< AMREX_SPACEDIM > IndexType
Definition: AMReX_BaseFwd.H:33
std::array< T, N > Array
Definition: AMReX_Array.H:24