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#ifdef AMREX_USE_EB
21#include <AMReX_EBSupport.H>
22#endif
23
24#include <memory>
25#include <map>
26
27namespace amrex {
28
29class TagBox;
30class TagBoxArray;
31
39{
40 friend class Amr;
41 friend class FillPatchIterator;
43
44public:
53 virtual ~AmrLevel ();
54
55 AmrLevel (const AmrLevel&) = delete;
56 AmrLevel (AmrLevel&&) = delete;
57 AmrLevel& operator= (const AmrLevel&) = delete;
59
61 void LevelDirectoryNames (const std::string &dir,
62 std::string &LevelDir,
63 std::string &FullPath) const;
65 virtual void CreateLevelDirectory (const std::string &dir);
70 void SetLevelDirectoryCreated(bool ldc) noexcept { levelDirectoryCreated = ldc; }
79 virtual std::string thePlotFileType () const
80 {
81 static const std::string the_plot_file_type("HyperCLaw-V1.1");
82 return the_plot_file_type;
83 }
87 virtual void writePlotFile (const std::string& dir,
88 std::ostream& os,
90
92 virtual void writePlotFilePre (const std::string& dir,
93 std::ostream& os);
94
96 virtual void writePlotFilePost (const std::string& dir,
97 std::ostream& os);
98
102 virtual void writeSmallPlotFile (const std::string& /*dir*/,
103 std::ostream& /*os*/,
104 VisMF::How /*how*/ = VisMF::NFiles) {}
106 virtual void checkPoint (const std::string& dir,
107 std::ostream& os,
109 bool dump_old = true);
111 virtual void checkPointPre (const std::string& dir,
112 std::ostream& os);
114 virtual void checkPointPost (const std::string& dir,
115 std::ostream& os);
117 virtual void restart (Amr& papa,
118 std::istream& is,
119 bool bReadSpecial = false);
120
122 virtual void set_state_in_checkpoint (Vector<int>& state_in_checkpoint);
123
125 static bool isStateVariable (const std::string& name,
126 int& state_indx,
127 int& ncomp);
128
129 static void FlushFPICache ();
134 virtual void computeInitialDt (int finest_level,
135 int sub_cycle,
136 Vector<int>& n_cycle,
137 const Vector<IntVect>& ref_ratio,
138 Vector<Real>& dt_level,
139 Real stop_time) = 0;
144 virtual void computeNewDt (int finest_level,
145 int sub_cycle,
146 Vector<int>& n_cycle,
147 const Vector<IntVect>& ref_ratio,
148 Vector<Real>& dt_min,
149 Vector<Real>& dt_level,
150 Real stop_time,
151 int post_regrid_flag) = 0;
157 virtual Real advance (Real time,
158 Real dt,
159 int iteration,
160 int ncycle) = 0;
161
166 virtual void post_timestep (int iteration);
171 virtual void postCoarseTimeStep (Real time);
175 virtual void post_restart () {}
181 virtual void post_regrid (int lbase,
182 int new_finest) = 0;
188 virtual void post_init (Real stop_time) = 0;
192 virtual int okToContinue () { return 1; }
198 virtual int okToRegrid ();
204 virtual void initData () = 0;
206 virtual void setTimeLevel (Real time,
207 Real dt_old,
208 Real dt_new);
210 virtual void allocOldData ();
212 virtual void removeOldData ();
218 virtual void init (AmrLevel &old) = 0;
224 virtual void init () = 0;
226 void reset ();
228 int Level () const noexcept { return level; }
230 const BoxArray& boxArray () const noexcept { return grids; }
231 const BoxArray& getEdgeBoxArray (int dir) const noexcept;
232 const BoxArray& getNodalBoxArray () const noexcept;
233 //
234 const DistributionMapping& DistributionMap () const noexcept { return dmap; }
235 //
236 const FabFactory<FArrayBox>& Factory () const noexcept { return *m_factory; }
237 //
238#ifdef AMREX_USE_EB
239 const EBFArrayBoxFactory&
240 EBFactory () const noexcept {
241 return static_cast<amrex::EBFArrayBoxFactory const&>(*m_factory);
242 }
243#endif
244
246 int numGrids () const noexcept { return static_cast<int>(grids.size()); }
248 int numStates () const noexcept { return static_cast<int>(state.size()); }
250 const Box& Domain () const noexcept { return geom.Domain(); }
252 int nStep () const noexcept { return parent->levelSteps(level); }
254 const Geometry& Geom () const noexcept { return geom; }
255 //
256 const IntVect& fineRatio () const noexcept { return fine_ratio; }
258 Long countCells () const noexcept;
259
261 const BoxArray& getAreaNotToTag () noexcept;
262 const Box& getAreaToTag () noexcept;
264 void constructAreaNotToTag ();
266 void setAreaNotToTag (BoxArray& ba) noexcept;
267
268 void resetFillPatcher ();
269
274 virtual void errorEst (TagBoxArray& tb,
275 int clearval,
276 int tagval,
277 Real time,
278 int n_error_buf = 0,
279 int ngrow = 0) = 0;
281 void FillCoarsePatch (MultiFab& mf,
282 int dcomp,
283 Real time,
284 int state_idx,
285 int scomp,
286 int ncomp,
287 int nghost = 0);
289 virtual void setPhysBoundaryValues (FArrayBox& dest,
290 int state_indx,
291 Real time,
292 int dest_comp,
293 int src_comp,
294 int num_comp);
301 virtual std::unique_ptr<MultiFab> derive (const std::string& name,
302 Real time,
303 int ngrow);
308 virtual void derive (const std::string& name,
309 Real time,
310 MultiFab& mf,
311 int dcomp);
313 StateData& get_state_data (int state_indx) noexcept { return state[state_indx]; }
315 MultiFab& get_old_data (int state_indx) noexcept { return state[state_indx].oldData(); }
317 const MultiFab& get_old_data (int state_indx) const noexcept { return state[state_indx].oldData(); }
319 MultiFab& get_new_data (int state_indx) noexcept { return state[state_indx].newData(); }
321 const MultiFab& get_new_data (int state_indx) const noexcept { return state[state_indx].newData(); }
323 static const DescriptorList& get_desc_lst () noexcept { return desc_lst; }
325 static DeriveList& get_derive_lst () noexcept;
327 int postStepRegrid () const noexcept { return post_step_regrid; }
329 void setPostStepRegrid (int new_val) noexcept { post_step_regrid = new_val; }
330
333
335 Vector<int> getBCArray (int State_Type,
336 int gridno,
337 int strt_comp,
338 int ncomp);
340 MultiFab& get_data (int state_indx, Real time) noexcept;
343 int state_index,
344 int scomp,
345 int dcomp,
346 int ncomp,
347 Real time) const;
352 virtual void manual_tags_placement (TagBoxArray& tags,
353 const Vector<IntVect>& bf_lev);
355 virtual void setPlotVariables ();
357 virtual void setSmallPlotVariables ();
363 virtual Real estimateWork();
364
366 virtual int WorkEstType () { return -1; }
367
372 TimeLevel which_time (int state_indx, Real time) const noexcept;
373
374 virtual bool checkPointNow ();
375
377 virtual bool writePlotNow ();
378
380 virtual bool writeSmallPlotNow ();
381
382#ifdef AMREX_PARTICLES
384 virtual void particle_redistribute (int /*lbase*/ = 0, bool /*a_init*/ = false) {;}
385#endif
386
398 void FillPatcherFill (amrex::MultiFab& mf, int dcomp, int ncomp, int nghost,
399 amrex::Real time, int state_index, int scomp);
400
401 static void FillPatch (AmrLevel& amrlevel,
402 MultiFab& leveldata,
403 int boxGrow,
404 Real time,
405 int index,
406 int scomp,
407 int ncomp,
408 int dcomp=0);
409
410 static void FillPatchAdd (AmrLevel& amrlevel,
411 MultiFab& leveldata,
412 int boxGrow,
413 Real time,
414 int index,
415 int scomp,
416 int ncomp,
417 int dcomp=0);
418
440 template <typename F, typename P = RungeKutta::PostStageNoOp>
441 void RK (int order, int state_type, Real time, Real dt, int iteration,
442 int ncycle, F&& f, P&& p = RungeKutta::PostStageNoOp());
443
444#ifdef AMREX_USE_EB
445 static void SetEBMaxGrowCells (int nbasic, int nvolume, int nfull) noexcept {
446 m_eb_basic_grow_cells = nbasic;
447 m_eb_volume_grow_cells = nvolume;
448 m_eb_full_grow_cells = nfull;
449 }
453 static void SetEBSupportLevel (EBSupport ebs) { m_eb_support_level = ebs; }
455#endif
456
458 static IntVect ProperBlockingFactor (AmrLevel const& amr_level, int boxGrow,
459 IndexType const& boxType,
460 StateDescriptor const& desc, int SComp);
461
462protected:
464 AmrLevel () noexcept {} // NOLINT
465
466 AmrLevel (Amr& papa,
467 int lev,
468 const Geometry& level_geom,
469 const BoxArray& ba,
470 const DistributionMapping& dm,
471 Real time);
472
474 void finishConstructor ();
475
476 //
477 // The Data.
478 //
479 int level{-1}; // AMR level (0 is coarsest).
480 Geometry geom; // Geom at this level.
481 BoxArray grids; // Cell-centered locations of grids.
482 DistributionMapping dmap; // Distribution of grids among processes
483 Amr* parent{nullptr};// Pointer to parent AMR structure.
484 IntVect crse_ratio; // Refinement ratio to coarser level.
485 IntVect fine_ratio; // Refinement ratio to finer level.
486 static DeriveList derive_lst; // List of derived quantities.
487 static DescriptorList desc_lst; // List of state variables.
488 Vector<StateData> state; // Array of state data.
489
490 BoxArray m_AreaNotToTag; //Area which shouldn't be tagged on this level.
491 Box m_AreaToTag; //Area which is allowed to be tagged on this level.
492
493 int post_step_regrid{0}; // Whether or not to do a regrid after the timestep.
494
495 bool levelDirectoryCreated{false}; // for checkpoints and plotfiles
496
497 std::unique_ptr<FabFactory<FArrayBox> > m_factory;
498
500
501private:
502
503 template <std::size_t order>
504 void storeRKCoarseData (int state_type, Real time, Real dt,
505 MultiFab const& S_old,
506 Array<MultiFab,order> const& rkk);
507
508 void FillRKPatch (int state_index, MultiFab& S, Real time,
509 int stage, int iteration, int ncycle);
510
511 mutable BoxArray edge_grids[AMREX_SPACEDIM]; // face-centered grids
512 mutable BoxArray nodal_grids; // all nodal grids
513};
514
515//
516// Forward declaration.
517//
519
521 :
522 public MFIter
523{
524 public:
525
526 friend class AmrLevel;
527
528 FillPatchIterator (AmrLevel& amrlevel,
529 MultiFab& leveldata);
530
531 FillPatchIterator (AmrLevel& amrlevel,
532 MultiFab& leveldata,
533 int boxGrow,
534 Real time,
535 int idx,
536 int scomp,
537 int ncomp);
538
539 void Initialize (int boxGrow,
540 Real time,
541 int idx,
542 int scomp,
543 int ncomp);
544
545 FArrayBox& operator() () noexcept { return m_fabs[MFIter::index()]; }
546
547 Box UngrownBox () const noexcept { return MFIter::validbox(); }
548
549 MultiFab& get_mf() noexcept { return m_fabs; }
550
551private:
552
553 void FillFromLevel0 (Real time, int idx, int scomp, int dcomp, int ncomp);
554 void FillFromTwoLevels (Real time, int idx, int scomp, int dcomp, int ncomp);
555
556 //
557 // The data.
558 //
559 AmrLevel* m_amrlevel;
560 MultiFab* m_leveldata;
561 std::vector< std::pair<int,int> > m_range;
562 MultiFab m_fabs;
563 int m_ncomp;
564};
565
567{
568public:
569
570 friend class FillPatchIterator;
571
573 MultiFab& leveldata);
574
576 MultiFab& leveldata,
577 int boxGrow,
578 Real time,
579 int index,
580 int scomp,
581 int ncomp,
582 InterpBase* mapper);
583
584 void Initialize (int boxGrow,
585 Real time,
586 int idx,
587 int scomp,
588 int ncomp,
589 InterpBase* mapper);
590
591 void fill (FArrayBox& fab, int dcomp, int idx);
592
593private:
594 //
595 // The data.
596 //
597 AmrLevel* m_amrlevel;
598 MultiFab* m_leveldata;
600 Vector< Vector<MultiFabId> > m_mfid; // [level][oldnew]
601 Interpolater* m_map = nullptr;
602 std::map<int,Box> m_ba;
603 Real m_time = std::numeric_limits<Real>::lowest();
604 int m_growsize = std::numeric_limits<int>::lowest();
605 int m_index = -1;
606 int m_scomp = -1;
607 int m_ncomp = -1;
608 bool m_FixUpCorners = false;
609
610 std::map< int,Vector< Vector<Box> > > m_fbox; // [grid][level][validregion]
611 std::map< int,Vector< Vector<Box> > > m_cbox; // [grid][level][fillablesubbox]
612 std::map< int,Vector< Vector< Vector<FillBoxId> > > > m_fbid; // [grid][level][fillablesubbox][oldnew]
613};
614
615template <typename F, typename P>
616void AmrLevel::RK (int order, int state_type, Real time, Real dt, int iteration,
617 int ncycle, F&& f, P&& p)
618{
619 BL_PROFILE("AmrLevel::RK()");
620
621 AMREX_ASSERT(AmrLevel::desc_lst[state_type].nExtra() > 0); // Need ghost cells in StateData
622
623 MultiFab& S_old = get_old_data(state_type);
624 MultiFab& S_new = get_new_data(state_type);
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
659template <std::size_t order>
660void 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
Virtual base class for managing individual levels. AmrLevel functions both as a container for state d...
Definition AMReX_AmrLevel.H:39
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
void reset()
Reset data to initial time by swapping new and old time data.
Definition AMReX_AmrLevel.cpp:596
Long countCells() const noexcept
Returns number of cells on level.
Definition AMReX_AmrLevel.cpp:504
IntVect fine_ratio
Definition AMReX_AmrLevel.H:485
int level
Definition AMReX_AmrLevel.H:479
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:228
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
virtual int WorkEstType()
Which state data type is for work estimates? -1 means none.
Definition AMReX_AmrLevel.H:366
virtual bool writePlotNow()
Does the AmrLevel want Amr to write a plotfile now?
Definition AMReX_AmrLevel.cpp:2065
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:451
const Geometry & Geom() const noexcept
Returns the geometry object.
Definition AMReX_AmrLevel.H:254
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:2282
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...
Vector< std::unique_ptr< FillPatcher< MultiFab > > > m_fillpatcher
Definition AMReX_AmrLevel.H:499
virtual bool writeSmallPlotNow()
Does the AmrLevel want Amr to write a small plotfile now?
Definition AMReX_AmrLevel.cpp:2071
void LevelDirectoryNames(const std::string &dir, std::string &LevelDir, std::string &FullPath) const
Get the level directory names.
Definition AMReX_AmrLevel.cpp:2229
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
const BoxArray & boxArray() const noexcept
List of grids at this level.
Definition AMReX_AmrLevel.H:230
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
const MultiFab & get_old_data(int state_indx) const noexcept
State data at old time.
Definition AMReX_AmrLevel.H:317
void SetLevelDirectoryCreated(bool ldc) noexcept
Set if the Level_ directory was created or to clear the value. CreateLevelDirectory sets levelDirecto...
Definition AMReX_AmrLevel.H:70
void setPostStepRegrid(int new_val) noexcept
Sets a new value for the post-timestep regrid trigger.
Definition AMReX_AmrLevel.H:329
const IntVect & fineRatio() const noexcept
Definition AMReX_AmrLevel.H:256
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
Vector< int > getBCArray(int State_Type, int gridno, int strt_comp, int ncomp)
Boundary condition access function.
Definition AMReX_AmrLevel.cpp:1869
const FabFactory< FArrayBox > & Factory() const noexcept
Definition AMReX_AmrLevel.H:236
Box m_AreaToTag
Definition AMReX_AmrLevel.H:491
static void FillPatch(AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
Definition AMReX_AmrLevel.cpp:2193
virtual void removeOldData()
Delete old-time data.
Definition AMReX_AmrLevel.cpp:587
Geometry geom
Definition AMReX_AmrLevel.H:480
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:192
const BoxArray & getNodalBoxArray() const noexcept
Definition AMReX_AmrLevel.cpp:637
StateData & get_state_data(int state_indx) noexcept
State data object.
Definition AMReX_AmrLevel.H:313
MultiFab & get_data(int state_indx, Real time) noexcept
Get state data at specified index and time.
Definition AMReX_AmrLevel.cpp:605
const EBFArrayBoxFactory & EBFactory() const noexcept
Definition AMReX_AmrLevel.H:240
BoxArray grids
Definition AMReX_AmrLevel.H:481
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
static int m_eb_full_grow_cells
Definition AMReX_AmrLevel.H:452
void UpdateDistributionMaps(DistributionMapping &dmap)
Update the distribution maps in StateData based on the size of the map.
Definition AMReX_AmrLevel.cpp:1852
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_AmrLevel.H:234
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
const Box & Domain() const noexcept
Returns the indices defining physical domain.
Definition AMReX_AmrLevel.H:250
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:486
void constructAreaNotToTag()
Construct the area not to tag.
Definition AMReX_AmrLevel.cpp:2091
bool levelDirectoryCreated
Definition AMReX_AmrLevel.H:495
static void SetEBSupportLevel(EBSupport ebs)
Definition AMReX_AmrLevel.H:453
static DeriveList & get_derive_lst() noexcept
Returns list of derived variables.
Definition AMReX_AmrLevel.cpp:57
IntVect crse_ratio
Definition AMReX_AmrLevel.H:484
virtual void setPlotVariables()
Modify list of variables to be plotted.
Definition AMReX_AmrLevel.cpp:1897
const Box & getAreaToTag() noexcept
Definition AMReX_AmrLevel.cpp:2081
virtual void writeSmallPlotFile(const std::string &, std::ostream &, VisMF::How=VisMF::NFiles)
Write small plot file stuff to specified directory.
Definition AMReX_AmrLevel.H:102
virtual bool checkPointNow()
Definition AMReX_AmrLevel.cpp:2059
TimeLevel
What time are we at?
Definition AMReX_AmrLevel.H:46
@ Amr1QtrTime
Definition AMReX_AmrLevel.H:49
@ AmrOldTime
Definition AMReX_AmrLevel.H:46
@ AmrNewTime
Definition AMReX_AmrLevel.H:48
@ AmrOtherTime
Definition AMReX_AmrLevel.H:51
@ AmrHalfTime
Definition AMReX_AmrLevel.H:47
@ Amr3QtrTime
Definition AMReX_AmrLevel.H:50
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
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:497
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:616
virtual void CreateLevelDirectory(const std::string &dir)
Create the Level_ directory for checkpoint and plot files.
Definition AMReX_AmrLevel.cpp:2245
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...
static void SetEBMaxGrowCells(int nbasic, int nvolume, int nfull) noexcept
Definition AMReX_AmrLevel.H:445
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:493
static EBSupport m_eb_support_level
Definition AMReX_AmrLevel.H:454
BoxArray m_AreaNotToTag
Definition AMReX_AmrLevel.H:490
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:2211
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:248
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:79
DistributionMapping dmap
Definition AMReX_AmrLevel.H:482
AmrLevel(AmrLevel &&)=delete
static int m_eb_basic_grow_cells
Definition AMReX_AmrLevel.H:450
virtual void particle_redistribute(int=0, bool=false)
This function can be called from the parent.
Definition AMReX_AmrLevel.H:384
Vector< StateData > state
Definition AMReX_AmrLevel.H:488
const MultiFab & get_new_data(int state_indx) const noexcept
State data at new time.
Definition AMReX_AmrLevel.H:321
int nStep() const noexcept
Timestep n at this level.
Definition AMReX_AmrLevel.H:252
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:2086
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:487
const BoxArray & getEdgeBoxArray(int dir) const noexcept
Definition AMReX_AmrLevel.cpp:626
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:464
static const DescriptorList & get_desc_lst() noexcept
Returns list of Descriptors.
Definition AMReX_AmrLevel.H:323
int postStepRegrid() const noexcept
Returns whether or not we want a post-timestep regrid.
Definition AMReX_AmrLevel.H:327
virtual void init()=0
int numGrids() const noexcept
Number of grids at this level.
Definition AMReX_AmrLevel.H:246
Amr * parent
Definition AMReX_AmrLevel.H:483
MultiFab & get_new_data(int state_indx) noexcept
State data at new time.
Definition AMReX_AmrLevel.H:319
MultiFab & get_old_data(int state_indx) noexcept
State data at old time.
Definition AMReX_AmrLevel.H:315
void resetFillPatcher()
Definition AMReX_AmrLevel.cpp:2128
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:2136
virtual void post_restart()
Operations to be done after restart.
Definition AMReX_AmrLevel.H:175
const BoxArray & getAreaNotToTag() noexcept
Get the area not to tag.
Definition AMReX_AmrLevel.cpp:2076
const Vector< DistributionMapping > & DistributionMap() const noexcept
Definition AMReX_AmrMesh.H:109
int finestLevel() const noexcept
Return the finest level.
Definition AMReX_AmrMesh.H:97
const Vector< Geometry > & Geom() const noexcept
Definition AMReX_AmrMesh.H:108
const Vector< BoxArray > & boxArray() const noexcept
Definition AMReX_AmrMesh.H:110
Manage hierarchy of levels for time-dependent AMR computations.
Definition AMReX_Amr.H:36
AmrLevel & getLevel(int lev) noexcept
AmrLevel lev.
Definition AMReX_Amr.H:190
int levelSteps(int lev) const noexcept
Number of time steps at specified level.
Definition AMReX_Amr.H:104
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:567
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:614
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: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:567
void Initialize(int boxGrow, Real time, int idx, int scomp, int ncomp, InterpBase *mapper)
Definition AMReX_AmrLevel.cpp:739
void fill(FArrayBox &fab, int dcomp, int idx)
Definition AMReX_AmrLevel.cpp:1266
Definition AMReX_AmrLevel.H:523
FArrayBox & operator()() noexcept
Definition AMReX_AmrLevel.H:545
MultiFab & get_mf() noexcept
Definition AMReX_AmrLevel.H:549
Box UngrownBox() const noexcept
Definition AMReX_AmrLevel.H:547
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:26
Virtual base class for interpolaters.
Definition AMReX_Interpolater.H:22
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:33
Attributes of StateData.
Definition AMReX_StateDescriptor.H:33
An array of TagBoxes.
Definition AMReX_TagBox.H:150
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:160
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:248
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:198
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:230
Definition AMReX_RungeKutta.H:51