OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
AmrBoxLib.h
Go to the documentation of this file.
1//
2// Class AmrBoxLib
3// Concrete AMR object. It is based on the AMReX library
4// (cf. https://amrex-codes.github.io/ or https://ccse.lbl.gov/AMReX/).
5// AMReX is the successor of BoxLib. This class represents the interface
6// to AMReX and the AMR framework in OPAL. It implements the functions of
7// the AmrObject class.
8//
9// Copyright (c) 2016 - 2020, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
10// All rights reserved
11//
12// Implemented as part of the PhD thesis
13// "Precise Simulations of Multibunches in High Intensity Cyclotrons"
14//
15// This file is part of OPAL.
16//
17// OPAL is free software: you can redistribute it and/or modify
18// it under the terms of the GNU General Public License as published by
19// the Free Software Foundation, either version 3 of the License, or
20// (at your option) any later version.
21//
22// You should have received a copy of the GNU General Public License
23// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
24//
25#ifndef AMR_BOXLIB_H
26#define AMR_BOXLIB_H
27
28#include "Amr/AmrObject.h"
29
30// AMReX headers
31#include <AMReX_AmrMesh.H>
32#include <AMReX.H>
33
34#include <map>
35#include <vector>
36
37class AmrPartBunch;
38
39class AmrBoxLib: public AmrObject,
40 public amrex::AmrMesh {
41
42public:
56
57 typedef amrex::FArrayBox FArrayBox_t;
58 typedef amrex::Box Box_t;
59 typedef amrex::TagBox TagBox_t;
60 typedef amrex::TagBoxArray TagBoxArray_t;
61 typedef amrex::MFIter MFIter_t;
62
66 enum Strategy {
67 RANK_ZERO = 0, // all grids to processor zero
68 PFC = 1,
69 RANDOM = 2,
71 };
72
73public:
81 AmrBoxLib(const AmrDomain_t& domain,
82 const AmrIntArray_t& nGridPts,
83 int maxLevel,
84 AmrPartBunch* bunch_p);
85
93 static std::unique_ptr<AmrBoxLib> create(const AmrInfo& info,
94 AmrPartBunch* bunch_p);
95
99 void regrid(double time);
100
101 void getGridStatistics(std::map<int, long>& gridPtsPerCore,
102 std::vector<int>& gridsPerLevel) const;
103
107 void initFineLevels();
108
110
111 double getRho(int x, int y, int z);
112
113 void computeSelfFields();
114
115 void computeSelfFields(int bin);
116
117 void computeSelfFields_cycl(double gamma);
118
119 void computeSelfFields_cycl(int bin);
120
121 void updateMesh();
122
127 const Vector_t& getMeshScaling() const;
128
130
131 const int& maxLevel() const;
132 const int& finestLevel() const;
133
137 double getT() const;
138
139 void redistributeGrids(int how);
140
141
142protected:
143 /*
144 * AmrMesh functions
145 */
146
155 void RemakeLevel (int lev, AmrReal_t time,
156 const AmrGrid_t& new_grids, const AmrProcMap_t& new_dmap);
157
165 void MakeNewLevel (int lev, AmrReal_t time,
166 const AmrGrid_t& new_grids, const AmrProcMap_t& new_dmap);
167
172 void ClearLevel(int lev);
173
177 virtual void ErrorEst(int lev, TagBoxArray_t& tags,
178 AmrReal_t time, int ngrow) override;
179
192 void MakeNewLevelFromScratch(int lev, AmrReal_t time,
193 const AmrGrid_t& ba,
194 const AmrProcMap_t& dm);
195
207 void MakeNewLevelFromCoarse(int lev, AmrReal_t time,
208 const AmrGrid_t& ba,
209 const AmrProcMap_t& dm);
210
211
212private:
219 void doRegrid_m(int lbase, double time);
220
226 void preRegrid_m();
227
233 void postRegrid_m(int old_finest);
234
235 double solvePoisson_m();
236
237 /* ATTENTION
238 * The tagging routines assume the particles to be in the
239 * AMR domain, i.e. [-1, 1]^3
240 */
241
251 void tagForChargeDensity_m(int lev, TagBoxArray_t& tags,
252 AmrReal_t time, int ngrow);
253
264 void tagForPotentialStrength_m(int lev, TagBoxArray_t& tags,
265 AmrReal_t time, int ngrow);
266
277 void tagForEfield_m(int lev, TagBoxArray_t& tags,
278 AmrReal_t time, int ngrow);
279
291 void tagForMomenta_m(int lev, TagBoxArray_t& tags,
292 AmrReal_t time, int ngrow);
293
303 void tagForMaxNumParticles_m(int lev, TagBoxArray_t& tags,
304 AmrReal_t time, int ngrow);
305
315 void tagForMinNumParticles_m(int lev, TagBoxArray_t& tags,
316 AmrReal_t time, int ngrow);
317
324 void initBaseLevel_m(const AmrIntArray_t& nGridPts);
325
335 static void initParmParse_m(const AmrInfo& info, AmrLayout_t* layout_p);
336
342 void fillPhysbc_m(AmrField_t& mf, int lev = 0);
343
344// void gradient(int lev) {
345//
346// phi_m[lev]->FillBoundary(geom[lev].periodicity());
347//
348// for (MFIter_t mfi(*phi_m[lev], true); mfi.isValid(); ++mfi) {
349// const Box_t& tilebx = mfi.tilebox();
350// FArrayBox_t& fab = (*phi_m[lev])[mfi];
351// FArrayBox_t& efab = (*efield_m[lev])[mfi];
352//
353// const int* tlo = tilebx.loVect();
354// const int* thi = tilebx.hiVect();
355//
356// for (int i = tlo[0]; i <= thi[0]; ++i) {
357// for (int j = tlo[1]; j <= thi[1]; ++j) {
358// for (int k = tlo[2]; k <= thi[2]; ++k) {
359//
360// amrex::IntVect iv(D_DECL(i,j,k));
361//
362// // x-field
363// amrex::IntVect liv(D_DECL(i-1,j,k));
364// amrex::IntVect riv(D_DECL(i+1,j,k));
365// efab(iv, 0) = 0.5 * (fab(liv) - fab(riv));
366//
367// // y-field
368// liv = amrex::IntVect(D_DECL(i,j-1,k));
369// riv = amrex::IntVect(D_DECL(i,j+1,k));
370// efab(iv, 1) = 0.5 * (fab(liv) - fab(riv));
371//
372// // z-field
373// liv = amrex::IntVect(D_DECL(i,j-1,k));
374// riv = amrex::IntVect(D_DECL(i,j+1,k));
375// efab(iv, 2) = 0.5 * (fab(liv) - fab(riv));
376// }
377// }
378// }
379// }
380// }
381
382
383private:
386
387 // the layout of the bunch
389
392
395
398
401
402 // only used in case of potential and efield tagging
403 std::vector<bool> isFirstTagging_m;
404
406};
407
408#endif
BoxLibLayout< double, 3 > AmrLayout_t
Definition PBunchDefs.h:37
amrex::Vector< AmrVectorField_t > AmrVectorFieldContainer_t
Definition AmrDefs.h:42
amrex::DistributionMapping AmrProcMap_t
Definition AmrDefs.h:38
amrex::Geometry AmrGeometry_t
Definition AmrDefs.h:39
amrex::RealBox AmrDomain_t
Definition AmrDefs.h:46
amrex::Vector< AmrGeometry_t > AmrGeomContainer_t
Definition AmrDefs.h:43
amrex::IntVect AmrIntVect_t
Definition AmrDefs.h:48
amrex::Real AmrReal_t
Definition AmrDefs.h:51
amrex::MultiFab AmrField_t
Definition AmrDefs.h:34
amrex::Vector< std::unique_ptr< AmrField_t > > AmrScalarFieldContainer_t
Definition AmrDefs.h:41
amrex::Vector< int > AmrIntArray_t
Definition AmrDefs.h:47
amrex::BoxArray AmrGrid_t
Definition AmrDefs.h:40
amrex::Vector< AmrGrid_t > AmrGridContainer_t
Definition AmrDefs.h:44
amrex::Vector< AmrProcMap_t > AmrProcMapContainer_t
Definition AmrDefs.h:45
amr::AmrGeometry_t AmrGeometry_t
Definition AmrBoxLib.h:54
amr::AmrProcMap_t AmrProcMap_t
Definition AmrBoxLib.h:53
Vector_t meshScaling_m
in particle rest frame, the longitudinal length enlarged
Definition AmrBoxLib.h:400
const Vector_t & getMeshScaling() const
void redistributeGrids(int how)
void MakeNewLevelFromScratch(int lev, AmrReal_t time, const AmrGrid_t &ba, const AmrProcMap_t &dm)
virtual void ErrorEst(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow) override
const int & maxLevel() const
void preRegrid_m()
AmrBoxLib(const AmrDomain_t &domain, const AmrIntArray_t &nGridPts, int maxLevel, AmrPartBunch *bunch_p)
Definition AmrBoxLib.cpp:47
amrex::Box Box_t
Definition AmrBoxLib.h:58
AmrScalarFieldContainer_t rho_m
charge density on the grid for all levels
Definition AmrBoxLib.h:391
void tagForPotentialStrength_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
VectorPair_t getEExtrema()
void getGridStatistics(std::map< int, long > &gridPtsPerCore, std::vector< int > &gridsPerLevel) const
void tagForMaxNumParticles_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
amrex::MFIter MFIter_t
Definition AmrBoxLib.h:61
void postRegrid_m(int old_finest)
amr::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
Definition AmrBoxLib.h:45
void RemakeLevel(int lev, AmrReal_t time, const AmrGrid_t &new_grids, const AmrProcMap_t &new_dmap)
AmrScalarFieldContainer_t phi_m
scalar potential on the grid for all levels
Definition AmrBoxLib.h:394
amr::AmrField_t AmrField_t
Definition AmrBoxLib.h:43
void computeSelfFields()
amr::AmrIntArray_t AmrIntArray_t
Definition AmrBoxLib.h:50
void updateMesh()
std::vector< bool > isFirstTagging_m
Definition AmrBoxLib.h:403
amr::AmrReal_t AmrReal_t
Definition AmrBoxLib.h:51
void computeSelfFields_cycl(double gamma)
amrex::FArrayBox FArrayBox_t
Definition AmrBoxLib.h:57
double solvePoisson_m()
double getT() const
amr::AmrProcMapContainer_t AmrProcMapContainer_t
Definition AmrBoxLib.h:48
AmrLayout_t * layout_mp
Definition AmrBoxLib.h:388
amr::AmrGrid_t AmrGrid_t
Definition AmrBoxLib.h:52
AmrPartBunch * bunch_mp
bunch used for tagging strategies
Definition AmrBoxLib.h:385
amr::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
Definition AmrBoxLib.h:44
void MakeNewLevel(int lev, AmrReal_t time, const AmrGrid_t &new_grids, const AmrProcMap_t &new_dmap)
void tagForEfield_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
void tagForMinNumParticles_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
amrex::TagBoxArray TagBoxArray_t
Definition AmrBoxLib.h:60
void initFineLevels()
void regrid(double time)
void initBaseLevel_m(const AmrIntArray_t &nGridPts)
void MakeNewLevelFromCoarse(int lev, AmrReal_t time, const AmrGrid_t &ba, const AmrProcMap_t &dm)
amr::AmrGeomContainer_t AmrGeomContainer_t
Definition AmrBoxLib.h:46
const int & finestLevel() const
void doRegrid_m(int lbase, double time)
void tagForMomenta_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
void fillPhysbc_m(AmrField_t &mf, int lev=0)
amrex::TagBox TagBox_t
Definition AmrBoxLib.h:59
static std::unique_ptr< AmrBoxLib > create(const AmrInfo &info, AmrPartBunch *bunch_p)
Definition AmrBoxLib.cpp:74
void tagForChargeDensity_m(int lev, TagBoxArray_t &tags, AmrReal_t time, int ngrow)
amr::AmrDomain_t AmrDomain_t
Definition AmrBoxLib.h:49
AmrVectorFieldContainer_t efield_m
vector field on the grid for all levels
Definition AmrBoxLib.h:397
static void initParmParse_m(const AmrInfo &info, AmrLayout_t *layout_p)
void ClearLevel(int lev)
Vektor< int, 3 > getBaseLevelGridPoints() const
bool isPoissonSolved_m
Definition AmrBoxLib.h:405
amr::AmrGridContainer_t AmrGridContainer_t
Definition AmrBoxLib.h:47
amr::AmrIntVect_t AmrIntVect_t
Definition AmrBoxLib.h:55
double getRho(int x, int y, int z)
std::pair< Vector_t, Vector_t > VectorPair_t
Definition AmrObject.h:35
Vektor< double, 3 > Vector_t