OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
AmrBoxLib.h
Go to the documentation of this file.
1 #ifndef AMR_BOXLIB_H
2 #define AMR_BOXLIB_H
3 
4 #include "Amr/AmrObject.h"
5 
6 class AmrPartBunch;
7 
8 // AMReX headers
9 #include <AMReX_AmrMesh.H>
10 #include <AMReX.H>
11 
18 class AmrBoxLib : public AmrObject,
19  public amrex::AmrMesh
20 {
21 
22 public:
36 
37  typedef amrex::FArrayBox FArrayBox_t;
38  typedef amrex::Box Box_t;
39  typedef amrex::TagBox TagBox_t;
40  typedef amrex::TagBoxArray TagBoxArray_t;
41  typedef amrex::MFIter MFIter_t;
42 
43 
47  enum Strategy {
48  RANK_ZERO = 0, // all grids to processor zero
49  PFC = 1,
50  RANDOM = 2,
52  };
53 
54 public:
55 
63  AmrBoxLib(const AmrDomain_t& domain,
64  const AmrIntArray_t& nGridPts,
65  int maxLevel,
66  AmrPartBunch* bunch_p);
67 
75  static std::unique_ptr<AmrBoxLib> create(const AmrInfo& info,
76  AmrPartBunch* bunch_p);
77 
81  void regrid(double time);
82 
83  void getGridStatistics(std::map<int, long>& gridPtsPerCore,
84  std::vector<int>& gridsPerLevel) const;
85 
89  void initFineLevels();
90 
92 
93  double getRho(int x, int y, int z);
94 
95  void computeSelfFields();
96 
97  void computeSelfFields(int bin);
98 
99  void computeSelfFields_cycl(double gamma);
100 
101  void computeSelfFields_cycl(int bin);
102 
103  void updateMesh();
104 
109  const Vector_t& getMeshScaling() const;
110 
112 
113  const int& maxLevel() const;
114  const int& finestLevel() const;
115 
119  double getT() const;
120 
121  void redistributeGrids(int how);
122 
123 
124 protected:
125  /*
126  * AmrMesh functions
127  */
128 
137  void RemakeLevel (int lev, AmrReal_t time,
138  const AmrGrid_t& new_grids, const AmrProcMap_t& new_dmap);
139 
147  void MakeNewLevel (int lev, AmrReal_t time,
148  const AmrGrid_t& new_grids, const AmrProcMap_t& new_dmap);
149 
154  void ClearLevel(int lev);
155 
159  virtual void ErrorEst(int lev, TagBoxArray_t& tags,
160  AmrReal_t time, int ngrow) override;
161 
162 
175  void MakeNewLevelFromScratch (int lev, AmrReal_t time,
176  const AmrGrid_t& ba,
177  const AmrProcMap_t& dm);
178 
190  void MakeNewLevelFromCoarse (int lev, AmrReal_t time,
191  const AmrGrid_t& ba,
192  const AmrProcMap_t& dm);
193 
194 private:
201  void doRegrid_m(int lbase, double time);
202 
208  void postRegrid_m(int old_finest);
209 
210 
211  /* ATTENTION
212  * The tagging routines assume the particles to be in the
213  * AMR domain, i.e. [-1, 1]^3
214  */
215 
225  void tagForChargeDensity_m(int lev, TagBoxArray_t& tags,
226  AmrReal_t time, int ngrow);
227 
238  void tagForPotentialStrength_m(int lev, TagBoxArray_t& tags,
239  AmrReal_t time, int ngrow);
240 
251  void tagForEfield_m(int lev, TagBoxArray_t& tags,
252  AmrReal_t time, int ngrow);
253 
265  void tagForMomenta_m(int lev, TagBoxArray_t& tags,
266  AmrReal_t time, int ngrow);
267 
277  void tagForMaxNumParticles_m(int lev, TagBoxArray_t& tags,
278  AmrReal_t time, int ngrow);
279 
289  void tagForMinNumParticles_m(int lev, TagBoxArray_t& tags,
290  AmrReal_t time, int ngrow);
291 
298  void initBaseLevel_m(const AmrIntArray_t& nGridPts);
299 
309  static void initParmParse_m(const AmrInfo& info, AmrLayout_t* layout_p);
310 
316  void fillPhysbc_m(AmrField_t& mf, int lev = 0);
317 
318 
319 // void gradient(int lev) {
320 //
321 // phi_m[lev]->FillBoundary(geom[lev].periodicity());
322 //
323 // for (MFIter_t mfi(*phi_m[lev], true); mfi.isValid(); ++mfi) {
324 // const Box_t& tilebx = mfi.tilebox();
325 // FArrayBox_t& fab = (*phi_m[lev])[mfi];
326 // FArrayBox_t& efab = (*efield_m[lev])[mfi];
327 //
328 // const int* tlo = tilebx.loVect();
329 // const int* thi = tilebx.hiVect();
330 //
331 // for (int i = tlo[0]; i <= thi[0]; ++i) {
332 // for (int j = tlo[1]; j <= thi[1]; ++j) {
333 // for (int k = tlo[2]; k <= thi[2]; ++k) {
334 //
335 // amrex::IntVect iv(D_DECL(i,j,k));
336 //
337 // // x-field
338 // amrex::IntVect liv(D_DECL(i-1,j,k));
339 // amrex::IntVect riv(D_DECL(i+1,j,k));
340 // efab(iv, 0) = 0.5 * (fab(liv) - fab(riv));
341 //
342 // // y-field
343 // liv = amrex::IntVect(D_DECL(i,j-1,k));
344 // riv = amrex::IntVect(D_DECL(i,j+1,k));
345 // efab(iv, 1) = 0.5 * (fab(liv) - fab(riv));
346 //
347 // // z-field
348 // liv = amrex::IntVect(D_DECL(i,j-1,k));
349 // riv = amrex::IntVect(D_DECL(i,j+1,k));
350 // efab(iv, 2) = 0.5 * (fab(liv) - fab(riv));
351 // }
352 // }
353 // }
354 // }
355 // }
356 
357 private:
360 
361  // the layout of the bunch
363 
366 
369 
372 
375 
376  // only used in case of potential and efield tagging
377  std::vector<bool> isFirstTagging_m;
378 
380 };
381 
382 #endif
void RemakeLevel(int lev, AmrReal_t time, const AmrGrid_t &new_grids, const AmrProcMap_t &new_dmap)
Definition: AmrBoxLib.cpp:603
amr::AmrGrid_t AmrGrid_t
Definition: AmrBoxLib.h:32
amr::AmrIntVect_t AmrIntVect_t
Definition: AmrBoxLib.h:35
amrex::Geometry AmrGeometry_t
Definition: AmrDefs.h:19
amr::AmrGeometry_t AmrGeometry_t
Definition: AmrBoxLib.h:34
void tagForChargeDensity_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
Definition: AmrBoxLib.cpp:808
void computeSelfFields()
Definition: AmrBoxLib.cpp:199
void initFineLevels()
Definition: AmrBoxLib.cpp:135
static std::unique_ptr< AmrBoxLib > create(const AmrInfo &info, AmrPartBunch *bunch_p)
Definition: AmrBoxLib.cpp:46
amrex::MFIter MFIter_t
Definition: AmrBoxLib.h:41
std::vector< bool > isFirstTagging_m
Definition: AmrBoxLib.h:377
Definition: TSVMeta.h:24
amr::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
Definition: AmrBoxLib.h:25
amrex::Box Box_t
Definition: AmrBoxLib.h:38
AmrScalarFieldContainer_t phi_m
scalar potential on the grid for all levels
Definition: AmrBoxLib.h:368
amr::AmrIntArray_t AmrIntArray_t
Definition: AmrBoxLib.h:30
void doRegrid_m(int lbase, double time)
Definition: AmrBoxLib.cpp:733
Vektor< int, 3 > getBaseLevelGridPoints() const
Definition: AmrBoxLib.cpp:515
amrex::IntVect AmrIntVect_t
Definition: AmrDefs.h:28
void tagForMomenta_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
Definition: AmrBoxLib.cpp:980
void tagForMinNumParticles_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
Definition: AmrBoxLib.cpp:1089
void MakeNewLevel(int lev, AmrReal_t time, const AmrGrid_t &new_grids, const AmrProcMap_t &new_dmap)
Definition: AmrBoxLib.cpp:636
void updateMesh()
Definition: AmrBoxLib.cpp:498
bool info
Info flag.
Definition: Options.cpp:8
amrex::Vector< AmrGeometry_t > AmrGeomContainer_t
Definition: AmrDefs.h:23
void tagForEfield_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
Definition: AmrBoxLib.cpp:915
void getGridStatistics(std::map< int, long > &gridPtsPerCore, std::vector< int > &gridsPerLevel) const
Definition: AmrBoxLib.cpp:108
void regrid(double time)
Definition: AmrBoxLib.cpp:79
amrex::MultiFab AmrField_t
Definition: AmrDefs.h:14
amr::AmrProcMap_t AmrProcMap_t
Definition: AmrBoxLib.h:33
void ClearLevel(int lev)
Definition: AmrBoxLib.cpp:671
amrex::TagBox TagBox_t
Definition: AmrBoxLib.h:39
amrex::DistributionMapping AmrProcMap_t
Definition: AmrDefs.h:18
AmrBoxLib(const AmrDomain_t &domain, const AmrIntArray_t &nGridPts, int maxLevel, AmrPartBunch *bunch_p)
Definition: AmrBoxLib.cpp:19
amrex::Vector< int > AmrIntArray_t
Definition: AmrDefs.h:27
const int & finestLevel() const
Definition: AmrBoxLib.cpp:532
amrex::FArrayBox FArrayBox_t
Definition: AmrBoxLib.h:37
AmrScalarFieldContainer_t rho_m
charge density on the grid for all levels
Definition: AmrBoxLib.h:365
amrex::RealBox AmrDomain_t
Definition: AmrDefs.h:26
amr::AmrField_t AmrField_t
Definition: AmrBoxLib.h:23
amr::AmrGridContainer_t AmrGridContainer_t
Definition: AmrBoxLib.h:27
void MakeNewLevelFromCoarse(int lev, AmrReal_t time, const AmrGrid_t &ba, const AmrProcMap_t &dm)
Definition: AmrBoxLib.cpp:725
amr::AmrReal_t AmrReal_t
Definition: AmrBoxLib.h:31
amr::AmrProcMapContainer_t AmrProcMapContainer_t
Definition: AmrBoxLib.h:28
const int & maxLevel() const
Definition: AmrBoxLib.cpp:527
double getRho(int x, int y, int z)
Definition: AmrBoxLib.cpp:192
amrex::Vector< AmrGrid_t > AmrGridContainer_t
Definition: AmrDefs.h:24
static void initParmParse_m(const AmrInfo &info, AmrLayout_t *layout_p)
Definition: AmrBoxLib.cpp:1167
virtual void ErrorEst(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow) override
Definition: AmrBoxLib.cpp:684
amrex::BoxArray AmrGrid_t
Definition: AmrDefs.h:20
bool isPoissonSolved_m
Definition: AmrBoxLib.h:379
VectorPair_t getEExtrema()
Definition: AmrBoxLib.cpp:169
amr::AmrDomain_t AmrDomain_t
Definition: AmrBoxLib.h:29
Vector_t meshScaling_m
in particle rest frame, the longitudinal length enlarged
Definition: AmrBoxLib.h:374
double getT() const
Definition: AmrBoxLib.cpp:537
void MakeNewLevelFromScratch(int lev, AmrReal_t time, const AmrGrid_t &ba, const AmrProcMap_t &dm)
Definition: AmrBoxLib.cpp:717
void tagForMaxNumParticles_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
Definition: AmrBoxLib.cpp:1039
AmrVectorFieldContainer_t efield_m
vector field on the grid for all levels
Definition: AmrBoxLib.h:371
void postRegrid_m(int old_finest)
Definition: AmrBoxLib.cpp:769
void fillPhysbc_m(AmrField_t &mf, int lev=0)
Definition: AmrBoxLib.cpp:1468
void computeSelfFields_cycl(double gamma)
Definition: AmrBoxLib.cpp:211
AmrLayout_t * layout_mp
Definition: AmrBoxLib.h:362
amrex::TagBoxArray TagBoxArray_t
Definition: AmrBoxLib.h:40
amr::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
Definition: AmrBoxLib.h:24
amrex::Vector< AmrProcMap_t > AmrProcMapContainer_t
Definition: AmrDefs.h:25
std::pair< Vector_t, Vector_t > VectorPair_t
Definition: AmrObject.h:23
AmrPartBunch * bunch_mp
bunch used for tagging strategies
Definition: AmrBoxLib.h:359
const Vector_t & getMeshScaling() const
Definition: AmrBoxLib.cpp:510
amrex::Vector< AmrVectorField_t > AmrVectorFieldContainer_t
Definition: AmrDefs.h:22
amrex::Vector< std::unique_ptr< AmrField_t > > AmrScalarFieldContainer_t
Definition: AmrDefs.h:21
void redistributeGrids(int how)
Definition: AmrBoxLib.cpp:542
void tagForPotentialStrength_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
Definition: AmrBoxLib.cpp:861
amrex::Real AmrReal_t
Definition: AmrDefs.h:31
void initBaseLevel_m(const AmrIntArray_t &nGridPts)
Definition: AmrBoxLib.cpp:1141
amr::AmrField_t AmrField_t
Definition: PBunchDefs.h:59
amr::AmrGeomContainer_t AmrGeomContainer_t
Definition: AmrBoxLib.h:26