OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
AmrBoxLib.cpp
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#include "Amr/AmrBoxLib.h"
26
28#include "Amr/AmrYtWriter.h"
29#include "Physics/Physics.h"
30#include "Physics/Units.h"
34#include "Utilities/Options.h"
35#include "Utility/PAssert.h"
36
37#include <AMReX_BCUtil.H>
38#include <AMReX_MultiFabUtil.H>
39#include <AMReX_ParmParse.H> // used in initialize function
40
41#include <algorithm>
42#include <cmath>
43
44extern Inform* gmsg;
45
46
48 const AmrIntArray_t& nGridPts,
49 int maxLevel,
50 AmrPartBunch* bunch_p)
51 : AmrObject()
52 , amrex::AmrMesh(&domain, maxLevel, nGridPts, 0 /* cartesian */)
53 , bunch_mp(bunch_p)
54 , layout_mp(static_cast<AmrLayout_t*>(&bunch_p->getLayout()))
55 , rho_m(maxLevel + 1)
56 , phi_m(maxLevel + 1)
57 , efield_m(maxLevel + 1)
58 , meshScaling_m(Vector_t(1.0, 1.0, 1.0))
59 , isFirstTagging_m(maxLevel + 1, true)
60 , isPoissonSolved_m(false)
61{
62 /*
63 * The layout needs to know how many levels we can make.
64 */
65 layout_mp->resize(maxLevel);
66
67 initBaseLevel_m(nGridPts);
68
69 // set mesh spacing of bunch
70 updateMesh();
71}
72
73
74std::unique_ptr<AmrBoxLib> AmrBoxLib::create(const AmrInfo& info,
75 AmrPartBunch* bunch_p) {
76 /* The bunch is initialized first with a Geometry,
77 * BoxArray and DistributionMapping on
78 * the base level (level = 0). Thus, we take the domain specified there in
79 * order to create the Amr object.
80 */
81 AmrLayout_t* layout_p = static_cast<AmrLayout_t*>(&bunch_p->getLayout());
82 AmrDomain_t domain = layout_p->Geom(0).ProbDomain();
83
84 AmrIntArray_t nGridPts = {
85 info.grid[0],
86 info.grid[1],
87 info.grid[2]
88 };
89
90 int maxlevel = info.maxlevel;
91
92 /*
93 * further attributes are given by the BoxLib's ParmParse class.
94 */
95 initParmParse_m(info, layout_p);
96
97 return std::unique_ptr<AmrBoxLib>(new AmrBoxLib(domain,
98 nGridPts,
99 maxlevel,
100 bunch_p
101 )
102 );
103}
104
105
106void AmrBoxLib::regrid(double time) {
108
109 *gmsg << "* Start regriding:" << endl
110 << "* Old finest level: "
111 << finest_level << endl;
112
113 this->preRegrid_m();
114
115 /* ATTENTION: The bunch might be updated during
116 * the regrid process!
117 * We regrid from base level 0 up to the finest level.
118 */
119 int old_finest = finest_level;
120
121 int lev_top = std::min(finest_level, max_level - 1);
122 for (int i = 0; i <= lev_top; ++i) {
123 this->doRegrid_m(i, time);
124 lev_top = std::min(finest_level, max_level - 1);
125 }
126
127 this->postRegrid_m(old_finest);
128
129 *gmsg << "* New finest level: "
130 << finest_level << endl
131 << "* Finished regriding" << endl;
132
134}
135
136
137void AmrBoxLib::getGridStatistics(std::map<int, long>& gridPtsPerCore,
138 std::vector<int>& gridsPerLevel) const {
139
140 typedef std::vector<int> container_t;
141
142 gridPtsPerCore.clear();
143 gridsPerLevel.clear();
144
145 gridsPerLevel.resize(max_level + 1);
146
147 for (int lev = 0; lev <= finest_level; ++lev) {
148 /* container index: box
149 * container value: cores that owns box
150 */
151 const container_t& pmap = this->dmap[lev].ProcessorMap();
152 const AmrGrid_t& ba = this->grids[lev];
153
154 gridsPerLevel[lev] = pmap.size();
155
156 // iterate over all boxes
157 for (unsigned int i = 0; i < ba.size(); ++i) {
158 gridPtsPerCore[pmap[i]] += ba[i].numPts();
159 }
160 }
161}
162
163
165 if ( bunch_mp->getNumBunch() > 1 && !refined_m ) {
166 *gmsg << "* Initialization of all levels" << endl;
167
168 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
169
170 /* we do an explicit domain mapping of the particles and then
171 * forbid it during the regrid process, this way it's only
172 * executed ones --> saves computation
173 */
174 bool isForbidTransform = amrpbase_p->isForbidTransform();
175
176 if ( !isForbidTransform ) {
177 amrpbase_p->domainMapping();
178 amrpbase_p->setForbidTransform(true);
179 }
180
181 if ( max_level > 0 ) {
182 this->regrid(bunch_mp->getT() * Units::s2ns);
183 }
184
185 if ( !isForbidTransform ) {
186 amrpbase_p->setForbidTransform(false);
187 // map particles back
188 amrpbase_p->domainMapping(true);
189 }
190
191 *gmsg << "* Initialization done." << endl;
192
193 refined_m = true;
194 }
195}
196
197
199 Vector_t maxE = Vector_t(0, 0, 0), minE = Vector_t(0, 0, 0);
200 /* check for maximum / minimum values over all levels
201 * and all components, i.e. Ex, Ey, Ez
202 */
203 for (unsigned int lev = 0; lev < efield_m.size(); ++lev) {
204 for (std::size_t i = 0; i < bunch_mp->Dimension; ++i) {
205 /* calls:
206 * - max(comp, nghost (to search), local)
207 * - min(cmop, nghost, local)
208 */
209 double max = efield_m[lev][i]->max(0, false);
210 maxE[i] = (maxE[i] < max) ? max : maxE[i];
211
212 double min = efield_m[lev][i]->min(0, false);
213 minE[i] = (minE[i] > min) ? min : minE[i];
214 }
215 }
216
217 return VectorPair_t(maxE, minE);
218}
219
220
221double AmrBoxLib::getRho(int /*x*/, int /*y*/, int /*z*/) {
222 //TODO
223 throw OpalException("AmrBoxLib::getRho(x, y, z)", "Not yet Implemented.");
224 return 0.0;
225}
226
227
229 //TODO
230 throw OpalException("AmrBoxLib::computeSelfFields", "Not yet Implemented.");
231}
232
233
235 //TODO
236 throw OpalException("AmrBoxLib::computeSelfFields(int bin)", "Not yet Implemented.");
237}
238
239
241 /*
242 * The potential is not scaled according to domain modification.
243 *
244 */
245 if ( !bunch_mp->hasFieldSolver() )
246 return;
247
248 /*
249 * scatter charges onto grid
250 */
251 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
252
255 bunch_mp->updateLorentzFactor(gamma);
256
257 // map on Amr domain + Lorentz transform
258 double scalefactor = amrpbase_p->domainMapping();
259
260 amrpbase_p->setForbidTransform(true);
261 amrpbase_p->update();
262 amrpbase_p->setForbidTransform(false);
263
264 double invGamma = 1.0 / gamma;
265 int nLevel = finest_level + 1;
266
267 double l0norm = this->solvePoisson_m();
268
269 /* apply scale of electric-field in order to undo the transformation
270 * + undo normalization
271 */
272 for (int i = 0; i <= finestLevel(); ++i) {
273 this->phi_m[i]->mult(scalefactor * l0norm, 0, 1);
274 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
275 this->efield_m[i][j]->mult(scalefactor * scalefactor * l0norm, 0, 1);
276 }
277 }
278
279 for (int i = 0; i <= finest_level; ++i) {
280 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
281 if ( this->efield_m[i][j]->contains_nan(false) )
282 throw OpalException("AmrBoxLib::computeSelfFields_cycl(double gamma) ",
283 "Ef: NANs at level " + std::to_string(i) + ".");
284 }
285 }
286
287 amrpbase_p->gather(bunch_mp->Ef, this->efield_m, bunch_mp->R, 0, finest_level);
288
289 // undo domain change + undo Lorentz transform
290 amrpbase_p->domainMapping(true);
291
293 bunch_mp->Ef *= Vector_t(gamma,
294 1.0,
295 gamma);
296
298 // Relativistic E&M says gamma*v/c^2 = gamma*beta/c = sqrt(gamma*gamma-1)/c
299 // but because we already transformed E_trans into the moving frame we have to
300 // add 1/gamma so we are using the E_trans from the rest frame -DW
301 double betaC = std::sqrt(gamma * gamma - 1.0) * invGamma / Physics::c;
302
304 bunch_mp->Bf(0) = betaC * bunch_mp->Ef(2);
305 bunch_mp->Bf(2) = -betaC * bunch_mp->Ef(0);
306
307 /*
308 * dumping only
309 */
310
311 if ( !(bunch_mp->getLocalTrackStep() % Options::amrYtDumpFreq) ) {
312 AmrYtWriter ytWriter(bunch_mp->getLocalTrackStep());
313
314 AmrIntArray_t rr(finest_level);
315 for (int i = 0; i < finest_level; ++i)
316 rr[i] = this->MaxRefRatio(i);
317
318 double time = bunch_mp->getT(); // in seconds
319
320 // we need to undo coefficient when writing charge density
321 for (int i = 0; i <= finest_level; ++i)
322 this->rho_m[i]->mult(- Physics::epsilon_0 * l0norm, 0, 1);
323
324
325 ytWriter.writeFields(rho_m, phi_m, efield_m, rr, this->geom,
326 nLevel, time, scalefactor);
327
328 ytWriter.writeBunch(bunch_mp, time, gamma);
329 }
330}
331
332
334
335 if ( !bunch_mp->hasFieldSolver() )
336 return;
337
339 double gamma = bunch_mp->getBinGamma(bin);
340
341 /* relativistic factor is always gamma >= 1
342 * --> if no particle --> gamma = 0 --> leave computation
343 */
344 if ( gamma < 1.0 )
345 return;
346
347 /*
348 * scatter charges onto grid
349 */
350 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
351
352 // Lorentz transformation: apply Lorentz factor of that particle bin
353 bunch_mp->updateLorentzFactor(bin);
354
355 // map on Amr domain + Lorentz transform
356 double scalefactor = amrpbase_p->domainMapping();
357
358 amrpbase_p->setForbidTransform(true);
359
360 amrpbase_p->update();
361
362 if ( !(bunch_mp->getLocalTrackStep() % Options::amrRegridFreq) ) {
363 this->regrid(bunch_mp->getT() * Units::s2ns);
364 }
365
366 amrpbase_p->setForbidTransform(false);
367
370 amrpbase_p->scatter(bunch_mp->Q, this->rho_m, bunch_mp->R, 0, finest_level, bunch_mp->Bin, bin);
371
373 // In particle rest frame, the longitudinal length (y for cyclotron) enlarged
374 // calculate Possion equation (with coefficient: -1/(eps))
375 for (int i = 0; i <= finest_level; ++i) {
376 this->rho_m[i]->mult(-1.0 / (Physics::epsilon_0), 0 /*comp*/, 1 /*ncomp*/);
377
378 if ( this->rho_m[i]->contains_nan(false) )
379 throw OpalException("AmrBoxLib::computeSelfFields_cycl(int bin) ",
380 "NANs at level " + std::to_string(i) + ".");
381 }
382
383 // find maximum and normalize each level (faster convergence)
384 double l0norm = 1.0;
385 for (int i = 0; i <= finest_level; ++i)
386 l0norm = std::max(l0norm, this->rho_m[i]->norm0(0));
387
388 for (int i = 0; i <= finest_level; ++i) {
389 this->rho_m[i]->mult(1.0 / l0norm, 0, 1);
390 }
391
392 // mesh scaling for solver
393 meshScaling_m = Vector_t(1.0, 1.0, 1.0);
394
395 PoissonSolver *solver = bunch_mp->getFieldSolver();
396
398
399 for (int i = 0; i <= finest_level; ++i) {
400 phi_m[i]->setVal(0.0, 0, phi_m[i]->nComp(), phi_m[i]->nGrow());
401 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
402 efield_m[i][j]->setVal(0.0, 0, efield_m[i][j]->nComp(), efield_m[i][j]->nGrow());
403 }
404 }
405
406 // in case of binning we reset phi every time
407 solver->solve(rho_m, phi_m, efield_m, 0, finest_level, false);
408
410
411 // make sure ghost cells are filled
412 for (int i = 0; i <= finest_level; ++i) {
413 phi_m[i]->FillBoundary(geom[i].periodicity());
414 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
415 efield_m[i][j]->FillBoundary(geom[i].periodicity());
416 }
417 }
418
419 this->fillPhysbc_m(*(this->phi_m[0]), 0);
420
421
422 /* apply scale of electric-field in order to undo the transformation
423 * + undo normalization
424 */
425 for (int i = 0; i <= finestLevel(); ++i) {
426 this->phi_m[i]->mult(scalefactor * l0norm, 0, 1);
427 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
428 this->efield_m[i][j]->mult(scalefactor * scalefactor * l0norm, 0, 1);
429 }
430 }
431
432 for (int i = 0; i <= finest_level; ++i) {
433 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
434 if ( this->efield_m[i][j]->contains_nan(false) )
435 throw OpalException("AmrBoxLib::computeSelfFields_cycl(double gamma) ",
436 "Ef: NANs at level " + std::to_string(i) + ".");
437 }
438 }
439
440 amrpbase_p->gather(bunch_mp->Eftmp, this->efield_m, bunch_mp->R, 0, finest_level);
441
442 // undo domain change + undo Lorentz transform
443 amrpbase_p->domainMapping(true);
444
446 bunch_mp->Eftmp *= Vector_t(gamma, 1.0, gamma);
447
449 double betaC = std::sqrt(gamma * gamma - 1.0) / gamma / Physics::c;
450
452 bunch_mp->Bf(0) += betaC * bunch_mp->Eftmp(2);
453 bunch_mp->Bf(2) -= betaC * bunch_mp->Eftmp(0);
454
455 bunch_mp->Ef += bunch_mp->Eftmp;
456
457 /*
458 * dumping only
459 */
460 if ( !(bunch_mp->getLocalTrackStep() % Options::amrYtDumpFreq) ) {
461 AmrYtWriter ytWriter(bunch_mp->getLocalTrackStep(), bin);
462
463 int nLevel = finest_level + 1;
464
465 AmrIntArray_t rr(finest_level);
466 for (int i = 0; i < finest_level; ++i)
467 rr[i] = this->MaxRefRatio(i);
468
469 double time = bunch_mp->getT(); // in seconds
470
471 // we need to undo coefficient when writing charge density
472 for (int i = 0; i <= finest_level; ++i)
473 this->rho_m[i]->mult(- Physics::epsilon_0 * l0norm, 0, 1);
474
475
476 ytWriter.writeFields(rho_m, phi_m, efield_m, rr, this->geom,
477 nLevel, time, scalefactor);
478
479 ytWriter.writeBunch(bunch_mp, time, gamma);
480 }
481
482 isPoissonSolved_m = true;
483}
484
485
487 //FIXME What about resizing mesh, i.e. geometry?
488 const AmrReal_t* tmp = this->geom[0].CellSize();
489
490 Vector_t hr;
491 for (int i = 0; i < 3; ++i)
492 hr[i] = tmp[i];
493
494 bunch_mp->setBaseLevelMeshSpacing(hr);
495}
496
497
499 return meshScaling_m;
500}
501
502
504 const Box_t& bx = this->geom[0].Domain();
505
506 const AmrIntVect_t& low = bx.smallEnd();
507 const AmrIntVect_t& high = bx.bigEnd();
508
509 return Vektor<int, 3>(high[0] - low[0] + 1,
510 high[1] - low[1] + 1,
511 high[2] - low[2] + 1);
512}
513
514
515inline const int& AmrBoxLib::maxLevel() const {
516 return this->max_level;
517}
518
519
520inline const int& AmrBoxLib::finestLevel() const {
521 return this->finest_level;
522}
523
524
525inline double AmrBoxLib::getT() const {
526 return bunch_mp->getT();
527}
528
529
531//
532// // copied + modified version of AMReX_Amr.cpp
533// AmrProcMap_t::InitProximityMap();
535//
536// AmrGridContainer_t allBoxes(finest_level + 1);
537// for(unsigned int ilev = 0; ilev < allBoxes.size(); ++ilev) {
538// allBoxes[ilev] = boxArray(ilev);
539// }
540// amrex::Vector<AmrIntArray_t> mLDM;
541//
542// switch ( how ) {
543// case RANK_ZERO:
544// mLDM = AmrProcMap_t::MultiLevelMapRandom(this->ref_ratio,
545// allBoxes,
546// maxGridSize(0),
547// 0/*maxRank*/,
548// 0/*minRank*/);
549// break;
550// case PFC:
551// mLDM = AmrProcMap_t::MultiLevelMapPFC(this->ref_ratio,
552// allBoxes,
553// maxGridSize(0));
554// break;
555// case RANDOM:
556// mLDM = AmrProcMap_t::MultiLevelMapRandom(this->ref_ratio,
557// allBoxes,
558// maxGridSize(0));
559// break;
560// case KNAPSACK:
561// mLDM = AmrProcMap_t::MultiLevelMapKnapSack(this->ref_ratio,
562// allBoxes,
563// maxGridSize(0));
564// break;
565// default:
566// *gmsg << "We didn't redistribute the grids." << endl;
567// return;
568// }
569//
570// for(unsigned int iMap = 0; iMap < mLDM.size(); ++iMap) {
571// AmrField_t::MoveAllFabs(mLDM[iMap]);
572// }
573//
574// /*
575// * particles need to know the BoxArray
576// * and DistributionMapping
577// */
578// for(unsigned int ilev = 0; ilev < allBoxes.size(); ++ilev) {
579// layout_mp->SetParticleBoxArray(ilev, this->grids[ilev]);
580// layout_mp->SetParticleDistributionMap(ilev, this->dmap[ilev]);
581// }
582//
583// for(unsigned int iMap = 0; iMap < mLDM.size(); ++iMap) {
584// rho_m[iMap]->MoveFabs(mLDM[iMap]);
585// phi_m[iMap]->MoveFabs(mLDM[iMap]);
586// efield_m[iMap]->MoveFabs(mLDM[iMap]);
587// }
588}
589
590
591void AmrBoxLib::RemakeLevel (int lev, AmrReal_t /*time*/,
592 const AmrGrid_t& new_grids,
593 const AmrProcMap_t& new_dmap) {
594
595 SetBoxArray(lev, new_grids);
596 SetDistributionMap(lev, new_dmap);
597
598 // #comp #ghosts cells
599 rho_m[lev].reset(new AmrField_t(new_grids, new_dmap, 1, 0));
600 phi_m[lev].reset(new AmrField_t(new_grids, new_dmap, 1, 1));
601 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
602 efield_m[lev][j].reset(new AmrField_t(new_grids, new_dmap, 1, 1));
603 }
604
605 // including nghost = 1
606 rho_m[lev]->setVal(0.0, 0);
607 phi_m[lev]->setVal(0.0, 1);
608
609 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
610 efield_m[lev][j]->setVal(0.0, 1);
611 }
612
613 /*
614 * particles need to know the BoxArray
615 * and DistributionMapping
616 */
617 layout_mp->SetParticleBoxArray(lev, new_grids);
618 layout_mp->SetParticleDistributionMap(lev, new_dmap);
619
620 layout_mp->buildLevelMask(lev, this->nErrorBuf(lev));
621}
622
623
624void AmrBoxLib::MakeNewLevel (int lev, AmrReal_t /*time*/,
625 const AmrGrid_t& new_grids,
626 const AmrProcMap_t& new_dmap) {
627 SetBoxArray(lev, new_grids);
628 SetDistributionMap(lev, new_dmap);
629
630 // #comp #ghosts cells
631 rho_m[lev].reset(new AmrField_t(new_grids, new_dmap, 1, 0));
632 phi_m[lev].reset(new AmrField_t(new_grids, new_dmap, 1, 1));
633
634 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
635 efield_m[lev][j].reset(new AmrField_t(new_grids, new_dmap, 1, 1));
636 }
637
638 // including nghost = 1
639 rho_m[lev]->setVal(0.0, 0);
640 phi_m[lev]->setVal(0.0, 1);
641 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
642 efield_m[lev][j]->setVal(0.0, 1);
643 }
644
645 /*
646 * particles need to know the BoxArray
647 * and DistributionMapping
648 */
649 layout_mp->SetParticleBoxArray(lev, new_grids);
650 layout_mp->SetParticleDistributionMap(lev, new_dmap);
651
652 layout_mp->buildLevelMask(lev, this->nErrorBuf(lev));
653}
654
655
657 rho_m[lev].reset(nullptr);
658 phi_m[lev].reset(nullptr);
659 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
660 efield_m[lev][j].reset(nullptr);
661 }
662 ClearBoxArray(lev);
663 ClearDistributionMap(lev);
664
665 this->isFirstTagging_m[lev] = true;
666}
667
668
670 AmrReal_t time, int ngrow)
671{
672 *gmsg << level2 << "* Start tagging of level " << lev << endl;
673
674 switch ( tagging_m ) {
676 tagForChargeDensity_m(lev, tags, time, ngrow);
677 break;
679 tagForPotentialStrength_m(lev, tags, time, ngrow);
680 break;
682 tagForEfield_m(lev, tags, time, ngrow);
683 break;
685 tagForMomenta_m(lev, tags, time, ngrow);
686 break;
688 tagForMaxNumParticles_m(lev, tags, time, ngrow);
689 break;
691 tagForMinNumParticles_m(lev, tags, time, ngrow);
692 break;
693 default:
694 tagForChargeDensity_m(lev, tags, time, ngrow);
695 break;
696 }
697
698 *gmsg << level2 << "* Finished tagging of level " << lev << endl;
699}
700
701
703 const AmrGrid_t& /*ba*/,
704 const AmrProcMap_t& /*dm*/)
705{
706 throw OpalException("AmrBoxLib::MakeNewLevelFromScratch()", "Shouldn't be called.");
707}
708
709
711 const AmrGrid_t& /*ba*/,
712 const AmrProcMap_t& /*dm*/)
713{
714 throw OpalException("AmrBoxLib::MakeNewLevelFromCoarse()", "Shouldn't be called.");
715}
716
717
718void AmrBoxLib::doRegrid_m(int lbase, double time) {
719 int new_finest = 0;
720 AmrGridContainer_t new_grids(finest_level+2);
721
722 MakeNewGrids(lbase, time, new_finest, new_grids);
723
724 PAssert(new_finest <= finest_level+1);
725
726 for (int lev = lbase+1; lev <= new_finest; ++lev)
727 {
728 if (lev <= finest_level) // an old level
729 {
730 if (new_grids[lev] != grids[lev]) // otherwise nothing
731 {
732 AmrProcMap_t new_dmap(new_grids[lev]);
733 RemakeLevel(lev, time, new_grids[lev], new_dmap);
734 }
735 }
736 else // a new level
737 {
738 AmrProcMap_t new_dmap(new_grids[lev]);
739 MakeNewLevel(lev, time, new_grids[lev], new_dmap);
740 }
741
742 layout_mp->setFinestLevel(new_finest);
743 }
744
745 // now we are safe to delete deprecated levels
746 for (int lev = new_finest+1; lev <= finest_level; ++lev) {
747 ClearLevel(lev);
748 }
749
750 finest_level = new_finest;
751}
752
753
755 /* In case of E-field or potential tagging
756 * in combination of binning we need to make
757 * sure we do not lose accuracy when tagging since
758 * the grid data of the potential or e-field are only
759 * non-zero for the last bin causing the tagging to refine
760 * only there.
761 * So, we need to solve the Poisson problem first assuming
762 * a single bin only
763 */
765 this->solvePoisson_m();
766 }
767}
768
769
770void AmrBoxLib::postRegrid_m(int old_finest) {
771 /* ATTENTION: In this call the bunch has to be updated
772 * since each particle needs to be assigned to its latest
773 * grid and sent to the corresponding MPI-process.
774 *
775 * We need to update the bunch before we clear
776 * levels, otherwise the particle level counter
777 * is invalidated.
778 */
779
780 // redistribute particles and assign correct grid and level
781 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
782 amrpbase_p->update(0, finest_level, true);
783
784 // check if no particles left on higher levels
785 const auto& LocalNumPerLevel = amrpbase_p->getLocalNumPerLevel();
786
787 for (int lev = finest_level+1; lev <= old_finest; ++lev) {
788 if ( LocalNumPerLevel.getLocalNumAtLevel(lev) != 0 ) {
789 throw OpalException("AmrBoxLib::postRegrid_m()",
790 "Still particles on level " + std::to_string(lev) + "!");
791 }
792 layout_mp->ClearParticleBoxArray(lev);
793 layout_mp->ClearParticleDistributionMap(lev);
794 layout_mp->clearLevelMask(lev);
795 }
796
797 if ( bunch_mp->getLocalNum() != LocalNumPerLevel.getLocalNumUpToLevel(finest_level) ) {
798 std::string localnum = std::to_string(bunch_mp->getLocalNum());
799 std::string levelnum = std::to_string(LocalNumPerLevel.getLocalNumUpToLevel(finest_level));
800 throw OpalException("AmrBoxLib::postRegrid_m()",
801 "Number of particles do not agree: " + localnum + " != " + levelnum);
802 }
803
804 PoissonSolver *solver = bunch_mp->getFieldSolver();
805 solver->hasToRegrid();
806}
807
808
810 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
811
813 amrpbase_p->scatter(bunch_mp->Q, this->rho_m, bunch_mp->R, 0, finest_level, bunch_mp->Bin);
814
815 int baseLevel = 0;
816
817 // mesh scaling for solver
818 meshScaling_m = Vector_t(1.0, 1.0, 1.0);
819
820 // charge density is in rho_m
821 // calculate Possion equation (with coefficient: -1/(eps))
822 for (int i = 0; i <= finest_level; ++i) {
823 if ( this->rho_m[i]->contains_nan(false) )
824 throw OpalException("AmrBoxLib::solvePoisson_m() ",
825 "NANs at level " + std::to_string(i) + ".");
826 this->rho_m[i]->mult(-1.0 / Physics::epsilon_0, 0, 1);
827 }
828
829 // find maximum and normalize each level (faster convergence)
830 double l0norm = 1.0;
831 for (int i = 0; i <= finest_level; ++i)
832 l0norm = std::max(l0norm, this->rho_m[i]->norm0(0));
833
834 for (int i = 0; i <= finest_level; ++i) {
835 this->rho_m[i]->mult(1.0 / l0norm, 0, 1);
836 }
837
838 PoissonSolver *solver = bunch_mp->getFieldSolver();
839
841 solver->solve(rho_m, phi_m, efield_m, baseLevel, finest_level);
843
844 // make sure ghost cells are filled
845 for (int i = 0; i <= finest_level; ++i) {
846 phi_m[i]->FillBoundary(geom[i].periodicity());
847 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
848 efield_m[i][j]->FillBoundary(geom[i].periodicity());
849 }
850 }
851
852 this->fillPhysbc_m(*(this->phi_m[0]), 0);
853
854 return l0norm;
855}
856
857
859 AmrReal_t /*time*/, int /*ngrow*/)
860{
861 // we need to update first
862 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
863 amrpbase_p->update(0, lev, true);
864
865 for (int i = lev; i <= finest_level; ++i)
866 rho_m[i]->setVal(0.0, rho_m[i]->nGrow());
867
868 // the new scatter function averages the value also down to the coarsest level
869 amrpbase_p->scatter(bunch_mp->Q, rho_m,
870 bunch_mp->R, 0, lev, bunch_mp->Bin);
871
872 const double& scalefactor = amrpbase_p->getScalingFactor();
873
874 for (int i = lev; i <= finest_level; ++i)
875 rho_m[i]->mult(scalefactor, 0, 1);
876
877 const int clearval = TagBox_t::CLEAR;
878 const int tagval = TagBox_t::SET;
879
880#ifdef _OPENMP
881#pragma omp parallel
882#endif
883 {
884 for (MFIter_t mfi(*rho_m[lev], true); mfi.isValid(); ++mfi) {
885 const Box_t& tilebx = mfi.tilebox();
886
887 TagBox_t& tagfab = tags[mfi];
888 FArrayBox_t& fab = (*rho_m[lev])[mfi];
889
890 const int* tlo = tilebx.loVect();
891 const int* thi = tilebx.hiVect();
892
893 for (int i = tlo[0]; i <= thi[0]; ++i) {
894 for (int j = tlo[1]; j <= thi[1]; ++j) {
895 for (int k = tlo[2]; k <= thi[2]; ++k) {
896
897 amrex::IntVect iv(D_DECL(i,j,k));
898
899 if ( std::abs( fab(iv) ) >= chargedensity_m )
900 tagfab(iv) = tagval;
901 else
902 tagfab(iv) = clearval;
903 }
904 }
905 }
906 }
907 }
908}
909
910
912 AmrReal_t time, int ngrow)
913{
914 /* Tag all cells for refinement
915 * where the value of the potential is higher than 75 percent of the maximum potential
916 * value of this level.
917 */
918 if ( !isPoissonSolved_m || !phi_m[lev]->ok() || this->isFirstTagging_m[lev] ) {
919 *gmsg << level2 << "* Level " << lev << ": We need to perform "
920 << "charge tagging if a new level is created." << endl;
921 this->tagForChargeDensity_m(lev, tags, time, ngrow);
922 this->isFirstTagging_m[lev] = false;
923 return;
924 }
925
926 // tag cells for refinement
927 const int clearval = TagBox_t::CLEAR;
928 const int tagval = TagBox_t::SET;
929
930 AmrReal_t threshold = phi_m[lev]->norm0(0, 0 /*nghost*/, false /*local*/);
931
932 threshold *= scaling_m;
933
934#ifdef _OPENMP
935#pragma omp parallel
936#endif
937 {
938 for (MFIter_t mfi(*phi_m[lev], true); mfi.isValid(); ++mfi) {
939
940 const Box_t& tilebx = mfi.tilebox();
941 TagBox_t& tagfab = tags[mfi];
942 FArrayBox_t& fab = (*phi_m[lev])[mfi];
943
944 const int* tlo = tilebx.loVect();
945 const int* thi = tilebx.hiVect();
946
947 for (int i = tlo[0]; i <= thi[0]; ++i) {
948 for (int j = tlo[1]; j <= thi[1]; ++j) {
949 for (int k = tlo[2]; k <= thi[2]; ++k) {
950
951 amrex::IntVect iv(D_DECL(i,j,k));
952
953 if ( std::abs( fab(iv) ) >= threshold )
954 tagfab(iv) = tagval;
955 else
956 tagfab(iv) = clearval;
957 }
958 }
959 }
960 }
961 }
962}
963
964
966 AmrReal_t time, int ngrow)
967{
968 /* Tag all cells for refinement
969 * where the value of the efield is higher than 75 percent of the maximum efield
970 * value of this level.
971 */
972 if ( !isPoissonSolved_m || !efield_m[lev][0]->ok() || this->isFirstTagging_m[lev] ) {
973 *gmsg << level2 << "* Level " << lev << ": We need to perform "
974 << "charge tagging if a new level is created." << endl;
975 this->tagForChargeDensity_m(lev, tags, time, ngrow);
976 this->isFirstTagging_m[lev] = false;
977 return;
978 }
979
980 // tag cells for refinement
981 const int clearval = TagBox_t::CLEAR;
982 const int tagval = TagBox_t::SET;
983
984 // obtain maximum absolute value
985 amrex::Vector<AmrReal_t> threshold(AMREX_SPACEDIM);
986
987 for (int i = 0; i < 3; ++i) {
988 threshold[i] = efield_m[lev][i]->norm0(0,
989 0 /*nghosts*/,
990 false /*local*/);
991
992 threshold[i] *= scaling_m;
993 }
994
995#ifdef _OPENMP
996#pragma omp parallel
997#endif
998 {
999 for (MFIter_t mfi(this->grids[lev], this->dmap[lev], true); mfi.isValid(); ++mfi) {
1000
1001 const Box_t& tilebx = mfi.tilebox();
1002 TagBox_t& tagfab = tags[mfi];
1003 FArrayBox_t& xfab = (*efield_m[lev][0])[mfi];
1004 FArrayBox_t& yfab = (*efield_m[lev][1])[mfi];
1005 FArrayBox_t& zfab = (*efield_m[lev][2])[mfi];
1006
1007 const int* tlo = tilebx.loVect();
1008 const int* thi = tilebx.hiVect();
1009
1010 for (int i = tlo[0]; i <= thi[0]; ++i) {
1011 for (int j = tlo[1]; j <= thi[1]; ++j) {
1012 for (int k = tlo[2]; k <= thi[2]; ++k) {
1013 AmrIntVect_t iv(i,j,k);
1014 if (std::abs(xfab(iv)) >= threshold[0])
1015 tagfab(iv) = tagval;
1016 else if (std::abs(yfab(iv)) >= threshold[1])
1017 tagfab(iv) = tagval;
1018 else if (std::abs(zfab(iv)) >= threshold[2])
1019 tagfab(iv) = tagval;
1020 else
1021 tagfab(iv) = clearval;
1022 }
1023 }
1024 }
1025 }
1026 }
1027}
1028
1029
1031 AmrReal_t /*time*/, int /*ngrow*/)
1032{
1033 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
1034 // we need to update first
1035 amrpbase_p->update(0, lev, true);
1036 const auto& LocalNumPerLevel = amrpbase_p->getLocalNumPerLevel();
1037
1038 size_t lBegin = LocalNumPerLevel.begin(lev);
1039 size_t lEnd = LocalNumPerLevel.end(lev);
1040
1041 Vector_t pmax = Vector_t(0.0, 0.0, 0.0);
1042 for (size_t i = lBegin; i < lEnd; ++i) {
1043 const Vector_t& tmp = bunch_mp->P[i];
1044 pmax = ( dot(tmp, tmp) > dot(pmax, pmax) ) ? tmp : pmax;
1045 }
1046
1047 double momentum2 = scaling_m * dot(pmax, pmax);
1048
1049 std::map<AmrIntVect_t, bool> cells;
1050 for (size_t i = lBegin; i < lEnd; ++i) {
1051 if ( dot(bunch_mp->P[i], bunch_mp->P[i]) >= momentum2 ) {
1052 AmrIntVect_t iv = layout_mp->Index(bunch_mp->R[i], lev);
1053 cells[iv] = true;
1054 }
1055 }
1056
1057 const int clearval = TagBox_t::CLEAR;
1058 const int tagval = TagBox_t::SET;
1059
1060
1061#ifdef _OPENMP
1062#pragma omp parallel
1063#endif
1064 {
1065 for (MFIter_t mfi(*rho_m[lev], true); mfi.isValid(); ++mfi) {
1066
1067 const Box_t& tilebx = mfi.tilebox();
1068 TagBox_t& tagfab = tags[mfi];
1069
1070 const int* tlo = tilebx.loVect();
1071 const int* thi = tilebx.hiVect();
1072
1073 for (int i = tlo[0]; i <= thi[0]; ++i) {
1074 for (int j = tlo[1]; j <= thi[1]; ++j) {
1075 for (int k = tlo[2]; k <= thi[2]; ++k) {
1076 AmrIntVect_t iv(i, j, k);
1077 if ( cells.find(iv) != cells.end() )
1078 tagfab(iv) = tagval;
1079 else
1080 tagfab(iv) = clearval;
1081 }
1082 }
1083 }
1084 }
1085 }
1086}
1087
1088
1090 AmrReal_t /*time*/, int /*ngrow*/)
1091{
1092 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
1093 // we need to update first
1094 amrpbase_p->update(0, lev, true);
1095 const auto& LocalNumPerLevel = amrpbase_p->getLocalNumPerLevel();
1096
1097 size_t lBegin = LocalNumPerLevel.begin(lev);
1098 size_t lEnd = LocalNumPerLevel.end(lev);
1099
1100 // count particles per cell
1101 std::map<AmrIntVect_t, size_t> cells;
1102 for (size_t i = lBegin; i < lEnd; ++i) {
1103 AmrIntVect_t iv = layout_mp->Index(bunch_mp->R[i], lev);
1104 ++cells[iv];
1105 }
1106
1107 const int clearval = TagBox_t::CLEAR;
1108 const int tagval = TagBox_t::SET;
1109
1110
1111#ifdef _OPENMP
1112#pragma omp parallel
1113#endif
1114 {
1115 for (MFIter_t mfi(*rho_m[lev], true); mfi.isValid(); ++mfi) {
1116
1117 const Box_t& tilebx = mfi.tilebox();
1118 TagBox_t& tagfab = tags[mfi];
1119
1120 const int* tlo = tilebx.loVect();
1121 const int* thi = tilebx.hiVect();
1122
1123 for (int i = tlo[0]; i <= thi[0]; ++i) {
1124 for (int j = tlo[1]; j <= thi[1]; ++j) {
1125 for (int k = tlo[2]; k <= thi[2]; ++k) {
1126 AmrIntVect_t iv(i, j, k);
1127 if ( cells.find(iv) != cells.end() && cells[iv] <= maxNumPart_m )
1128 tagfab(iv) = tagval;
1129 else
1130 tagfab(iv) = clearval;
1131 }
1132 }
1133 }
1134 }
1135 }
1136}
1137
1138
1140 AmrReal_t /*time*/, int /*ngrow*/)
1141{
1142 AmrPartBunch::pbase_t* amrpbase_p = bunch_mp->getAmrParticleBase();
1143 // we need to update first
1144 amrpbase_p->update(0, lev, true);
1145 const auto& LocalNumPerLevel = amrpbase_p->getLocalNumPerLevel();
1146
1147 size_t lBegin = LocalNumPerLevel.begin(lev);
1148 size_t lEnd = LocalNumPerLevel.end(lev);
1149
1150 // count particles per cell
1151 std::map<AmrIntVect_t, size_t> cells;
1152 for (size_t i = lBegin; i < lEnd; ++i) {
1153 AmrIntVect_t iv = layout_mp->Index(bunch_mp->R[i], lev);
1154 ++cells[iv];
1155 }
1156
1157 const int clearval = TagBox_t::CLEAR;
1158 const int tagval = TagBox_t::SET;
1159
1160
1161#ifdef _OPENMP
1162#pragma omp parallel
1163#endif
1164 {
1165 for (MFIter_t mfi(*rho_m[lev], true); mfi.isValid(); ++mfi) {
1166
1167 const Box_t& tilebx = mfi.tilebox();
1168 TagBox_t& tagfab = tags[mfi];
1169
1170 const int* tlo = tilebx.loVect();
1171 const int* thi = tilebx.hiVect();
1172
1173 for (int i = tlo[0]; i <= thi[0]; ++i) {
1174 for (int j = tlo[1]; j <= thi[1]; ++j) {
1175 for (int k = tlo[2]; k <= thi[2]; ++k) {
1176 AmrIntVect_t iv(i, j, k);
1177 if ( cells.find(iv) != cells.end() &&
1178 cells[iv] >= minNumPart_m ) {
1179 tagfab(iv) = tagval;
1180 } else {
1181 tagfab(iv) = clearval;
1182 }
1183 }
1184 }
1185 }
1186 }
1187 }
1188}
1189
1190
1192 // we need to set the AmrCore variable
1193 finest_level = 0;
1194
1195 AmrIntVect_t low(0, 0, 0);
1196 AmrIntVect_t high(nGridPts[0] - 1,
1197 nGridPts[1] - 1,
1198 nGridPts[2] - 1);
1199
1200 const Box_t bx(low, high);
1201 AmrGrid_t ba(bx);
1202 ba.maxSize( this->maxGridSize(0) );
1203
1204 // chop grids to guarantee #grids == #procs on base level
1205 ba = this->MakeBaseGrids();
1206
1207 AmrProcMap_t dmap;
1208 dmap.define(ba, amrex::ParallelDescriptor::NProcs());
1209
1210 this->RemakeLevel(0, 0.0, ba, dmap);
1211
1212 layout_mp->define(this->geom);
1213 layout_mp->define(this->ref_ratio);
1214}
1215
1216
1217void AmrBoxLib::initParmParse_m(const AmrInfo& info, AmrLayout_t* layout_p) {
1218
1219 /*
1220 * All parameters that we set with the OPAL input file
1221 */
1222 amrex::ParmParse pAmr("amr");
1223 pAmr.add("max_grid_size_x", info.maxgrid[0]);
1224 pAmr.add("max_grid_size_y", info.maxgrid[1]);
1225 pAmr.add("max_grid_size_z", info.maxgrid[2]);
1226
1227 pAmr.add("blocking_factor_x", info.bf[0]);
1228 pAmr.add("blocking_factor_y", info.bf[1]);
1229 pAmr.add("blocking_factor_z", info.bf[2]);
1230
1231 const int nratios_vect = info.maxlevel * AMREX_SPACEDIM;
1232
1233 AmrIntArray_t refratio(nratios_vect);
1234
1235 for (int i = 0; i < info.maxlevel; ++i) {
1236 refratio[i * AMREX_SPACEDIM] = info.refratio[0];
1237 refratio[i * AMREX_SPACEDIM + 1] = info.refratio[1];
1238 refratio[i * AMREX_SPACEDIM + 2] = info.refratio[2];
1239 }
1240
1241 pAmr.addarr("ref_ratio_vect", refratio);
1242
1243 amrex::ParmParse pGeom("geometry");
1244 AmrIntArray_t isPeriodic = {
1245 layout_p->Geom(0).isPeriodic(0),
1246 layout_p->Geom(0).isPeriodic(1),
1247 layout_p->Geom(0).isPeriodic(2)
1248 };
1249 pGeom.addarr("is_periodic", isPeriodic);
1250
1251
1252 /*
1253 * "All" other parameters, we moslty take the
1254 * defaults of the code
1255 *
1256 * ATTENTION Be careful with default values since they might
1257 * be changed by the AMReX developers!
1258 *
1259 * Parmparse parameters not to be set because
1260 * alreday given due to constructor
1261 * ----------------------------------
1262 *
1263 * AMReX_AmrMesh:
1264 * - amr.max_level
1265 *
1266 * AMReX_Geometry:
1267 * - geometry.coord_sys
1268 * - geometry.prob_lo
1269 * - geometry.prob_hi
1270 * - geometry.is_periodic
1271 * - geometry.spherical_origin_fix [default: 0]
1272 * (not set because we're using Cartesian coordinates)
1273 *
1274 *
1275 * ParmParse left out
1276 * ------------------
1277 * AMReX_FabArrayBase:
1278 * - fabarray.mfiter_tile_size [default: IntVect(1024000,8,8)]
1279 * - fabarray.mfghostiter_tile_size [default: IntVect(1024000, 8, 8)]
1280 * - fabarray.comm_tile_size [default: IntVect(1024000, 8, 8)]
1281 * - fabarray.maxcomp [default: 25]
1282 * - fabarray.do_async_sends [default: true]
1283 *
1284 * AMReX_MemPool:
1285 * - fab.init_snan [default: 0, if not BL_TESTING or DEBUG enabled]
1286 *
1287 * AMReX_MemProfiler:
1288 * - amrex.memory.log [default: "memlog"]
1289 *
1290 * AMReX_VisMF:
1291 * - vismf.v [verbose, default: 0]
1292 * - vismf.headerversion [default: VisMF::Header::Version_v1 = 1]
1293 * - vismf.groupsets [default: false]
1294 * - vismf.setbuf [default: true]
1295 * - vismf.usesingleread [default: false]
1296 * - vismf.usesinglewrite [default: false]
1297 * - vismf.checkfilepositions [default: false]
1298 * - vismf.usepersistentifstreams [default: true]
1299 * - vismf.usesynchronousreads [default: false]
1300 * - vismf.usedynamicsetselection [default: true]
1301 * - vismf.iobuffersize [default: VisMF::IO_Buffer_Size = 262144 * 8]
1302 *
1303 * AMReX_ParallelDescriptor:
1304 * Needs change of ParallelDescriptor::SetNProcsSidecars in the function
1305 * ParallelDescriptor::StartParallel()
1306 * - amrex.ncolors [default: m_nCommColors = 1]
1307 * - team.size [default: disabled, enable by BL_USE_UPCXX or BL_USE_MPI3]
1308 * - team.reduce [default: disabled, enable by BL_USE_UPCXX or BL_USE_MPI3]
1309 * - mpi.onesided [default: disabled, enable by BL_USE_UPCXX or BL_USE_MPI3]
1310 *
1311 * AMReX:
1312 * Debugging purposes
1313 * - amrex.v [verbose, default: 0]
1314 * - amrex.verbose [default: 0]
1315 * - amrex.fpe_trap_invalid [invalid, default: 0]
1316 * - amrex.fpe_trap_zero [divbyzero, default: 0]
1317 * - amrex.fpe_trap_overflow [overflow, 0]
1318 *
1319 * AMReX_FArrayBox:
1320 * For output only
1321 * - fab.format [default: FABio::FAB_NATIVE]
1322 * For reading in an old FAB
1323 * - fab.ordering
1324 * - fab.initval [default: quiet_NaN() or max()]
1325 * - fab.do_initval [default: false, if not DEBUG or BL_TESTING enabled]
1326 * - fab.init_snan [default: false, if not DEBUG or BL_TESTING enabled]
1327 *
1328 * AMReX_MCMultiGrid:
1329 * We do NOT use this solver
1330 * - mg.maxiter [default: 40]
1331 * - mg.numiter [default: -1]
1332 * - mg.nu_0 [default: 1]
1333 * - mg.nu_1 [default: 2]
1334 * - mg.nu_2 [default: 2]
1335 * - mg.nu_f [default: 8]
1336 * - mg.v [verbose, default: 0]
1337 * - mg.usecg [default: 1]
1338 * - mg.rtol_b [default: 0.01]
1339 * - mg.bot_atol [default: -1.0]
1340 * - mg.nu_b [default: 0]
1341 * - mg.numLevelsMAX [default: 1024]
1342 *
1343 * AMReX_MCLinOp:
1344 * We do NOT use this solver
1345 * - MCLp.harmavg [default: 0]
1346 * - MCLp.v [verbose, default: 0]
1347 * - MCLp.maxorder [default: 2]
1348 *
1349 * AMReX_MCCGSolver:
1350 * We do NOT use this solver
1351 * - cg.maxiter [default: 40]
1352 * - cg.v [verbose, default: 0]
1353 * - cg.isExpert [default: 0]
1354 * - MCCGSolver::def_unstable_criterion [default: 10]
1355 *
1356 * AMReX_LinOp:
1357 * We do NOT use this class
1358 * - Lp.harmavg [default: 0]
1359 * - Lp.v [verbose, default: 0]
1360 * - Lp.maxorder [default: 2]
1361 * - LinOp::LinOp_grow [default: 1]
1362 *
1363 * AMReX_MultiGrid (single level solver):
1364 * We do NOT use this solver, see BoxLibSolvers/FMGPoissonSolver for
1365 * explanation of parameters
1366 * - mg.v [verbose, default: 0]
1367 * - mg.nu_0 [default: 1]
1368 * - mg.nu_1 [default: 2]
1369 * - mg.nu_2 [default: 2]
1370 * - mg.nu_f [default: 8]
1371 * - mg.nu_b [default: 0]
1372 * - mg.usecg [default: 1]
1373 * - mg.rtol_b [default: 0.0001 or 0.01 (if CG_USE_OLD_CONVERGENCE_CRITERIA enabled)]
1374 * - mg.verbose [default: 0]
1375 * - mg.maxiter [default: 40]
1376 * - mg.bot_atol [default: -1.0]
1377 * - mg.maxiter_b [default: 120]
1378 * - mg.numLevelsMAX [default: 1024]
1379 * - mg.smooth_on_cg_unstable [default: 1]
1380 * - mg.use_Anorm_for_convergence [default: 1]
1381 *
1382 * AMReX_CGSolver (single level solver):
1383 * We do NOT use this solver
1384 * - cg.v [verbose, default: 0]
1385 * - cg.SSS [default: SSS_MAX = 4]
1386 * - cg.maxiter [default: 80]
1387 * - cg.verbose [default: 0]
1388 * - cg.variable_SSS [default: true]
1389 * - cg.use_jbb_precond [default: 0]
1390 * - cg.use_jacobi_precond [default: 0]
1391 * - cg.unstable_criterion [default: 10]
1392 * - CGSolver::def_cg_solver [default: BiCGStab]
1393 * - cg.cg_solver [default: CGSolver::def_cg_solver]
1394 *
1395 * AMReX_Amr:
1396 * We do NOT use this class
1397 * - amr.regrid_on_restart [default: 0]
1398 * - amr.use_efficient_regrid [default: 0]
1399 * - amr.plotfile_on_restart [default: 0]
1400 * - amr.checkpoint_on_restart [default: 0]
1401 * - amr.compute_new_dt_on_regrid [default: 0]
1402 * - amr.mffile_nstreams [default: 1]
1403 * - amr.probinit_natonce [default: 32]
1404 * - amr.file_name_digits [default: 5]
1405 * - amr.initial_grid_file
1406 * - amr.regrid_file
1407 * - amr.message_int [default: 10]
1408 * - amr.run_log [default: not parsed]
1409 * - amr.run_log_terse [default: not parsed]
1410 * - amr.grid_log [default: not parsed]
1411 * - amr.data_log [default: not parsed]
1412 * - amr.probin_file [default: not parsed]
1413 * - amr.restart
1414 * - amr.restart_from_plotfile
1415 * - amr.regrid_int [default: 1 at all levels]
1416 * - amr.rebalance_grids [default: 0]
1417 * - amr.nosub [default: not parsed]
1418 * - amr.subcycling_mode [default: "Auto"]
1419 * - amr.subcycling_iterations [default: 0]
1420 * - amr.checkpoint_files_output
1421 * - amr.plot_files_output
1422 * - amr.plot_nfiles
1423 * - amr.checkpoint_nfiles
1424 * - amr.check_file [default: "chk"]
1425 * - amr.check_int [default: -1]
1426 * - amr.check_per [default: -1.0]
1427 * - amr.plot_file [default: "plt"]
1428 * - amr.plot_int [default: -1]
1429 * - amr.plot_per [default: -1.0]
1430 * - amr.small_plot_file [default: "smallplt"]
1431 * - amr.small_plot_int [default: -1]
1432 * - amr.small_plot_per [default: -1.0]
1433 * - amr.write_plotfile_with_checkpoint [default: 1]
1434 * - amr.stream_max_tries [default: 4]
1435 * - amr.abort_on_stream_retry_failure [default: false]
1436 * - amr.precreateDirectories
1437 * - amr.prereadFAHeaders
1438 * - amr.plot_headerversion
1439 * - amr.checkpoint_headerversion
1440 *
1441 * AMReX_SlabStat:
1442 * We do NOT use this class
1443 * - slabstat.boxes
1444 *
1445 * AMReX_AmrLevel:
1446 * We do NOT use this class
1447 * - amr.plot_vars [default: not parsed]
1448 * - amr.derive_plot_vars [default: not parsed]
1449 * - amr.small_plot_vars [default: not parsed]
1450 *
1451 * AMReX_StationData:
1452 * We do NOT use this class
1453 * - StationData.vars [default: not parsed]
1454 * - StationData.coord [default: not parsed]
1455 */
1456
1457 // verbosity
1458 pAmr.add("v", 0);
1459
1460 // # cells required for proper nesting.
1461 pAmr.add("n_proper", 1);
1462
1463 // Grid efficiency. (default: 0.7)
1464 pAmr.add("grid_eff", 0.95);
1465
1466 int nlev = info.maxlevel + 1;
1467
1468 // Buffer cells around each tagged cell.
1469 AmrIntArray_t error_buf(nlev, 4 /*buffer*/);
1470 error_buf[0] = 0;
1471 pAmr.addarr("n_error_buf", error_buf);
1472
1473 // chop up grids to have more grids than the number of procs
1474 pAmr.add("refine_grid_layout", true);
1475
1476
1477 /*
1478 * ParmParse for DistributionMapping
1479 *
1480 * round-robin: FAB i is owned by CPU i%N where N is total number of CPUs
1481 * knapsack: FABs are partitioned across CPUs such that the total volume
1482 * of the Boxes in the underlying BoxArray are as equal across
1483 * CPUs as is possible.
1484 * SFC: Is based on a space filling curve (default)
1485 */
1486 amrex::ParmParse pDmap("DistributionMapping");
1487
1488 // verbosity
1489 pDmap.add("v", 0);
1490
1491 // max_efficiency = 0.9
1492 pDmap.add("efficiency", 0.9);
1493
1494 // SFC = space filling curve
1495 pDmap.add("sfc_threshold", 0);
1496
1497 /* if > 0:
1498 * nworkers = node_size
1499 * nteams = nprocs / node_size
1500 *
1501 * if nwokers * nteams != nprocs:
1502 * nteams = nprocs
1503 * nworkers = 1
1504 */
1505 pDmap.add("node_size", 0);
1506
1507 /* Possibilities:
1508 * - ROUNDROBIN
1509 * - KNAPSACK
1510 * - SFC
1511 * - PFC
1512 * - RRSFC
1513 */
1514 pDmap.add("strategy", "SFC");
1515}
1516
1517
1519 /* Copied from one of the miniapps:
1520 *
1521 * amrex/Src/AmrTask/tutorials/MiniApps/HeatEquation/physbc.cpp
1522 *
1523 * AMReX - Revision:
1524 *
1525 * commit 8174212e898677d6413cf5e4db44148a52b4b732
1526 * Author: Weiqun Zhang <weiqunzhang@lbl.gov>
1527 * Date: Mon Jul 2 10:40:21 2018 -0700
1528 */
1529 if (AmrGeometry_t::isAllPeriodic()) {
1530 return;
1531 }
1532
1533 // Set up BC; see Src/Base/AMReX_BC_TYPES.H for supported types
1534 amrex::Vector<amrex::BCRec> bc(mf.nComp());
1535 for (int n = 0; n < mf.nComp(); ++n) {
1536 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1537 if (AmrGeometry_t::isPeriodic(idim)) {
1538 bc[n].setLo(idim, amrex::BCType::int_dir); // interior
1539 bc[n].setHi(idim, amrex::BCType::int_dir);
1540 } else {
1541 bc[n].setLo(idim, amrex::BCType::foextrap); // first-order extrapolation.
1542 bc[n].setHi(idim, amrex::BCType::foextrap);
1543 }
1544 }
1545 }
1546
1547 amrex::FillDomainBoundary(mf, this->geom[lev], bc);
1548}
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Definition Vector3D.cpp:118
BoxLibLayout< double, 3 > AmrLayout_t
Definition PBunchDefs.h:37
Inform * gmsg
Definition Main.cpp:70
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
T::PETE_Expr_t::PETE_Return_t min(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Inform & level2(Inform &inf)
Definition Inform.cpp:46
Inform & endl(Inform &inf)
Definition Inform.cpp:42
#define PAssert(c)
Definition PAssert.h:102
constexpr double epsilon_0
The permittivity of vacuum in As/Vm.
Definition Physics.h:51
constexpr double c
The velocity of light in m/s.
Definition Physics.h:45
constexpr double s2ns
Definition Units.h:44
int amrRegridFreq
After how many steps the AMR grid hierarchy is updated.
Definition Options.cpp:103
int amrYtDumpFreq
The frequency to dump AMR grid data and particles into file.
Definition Options.cpp:101
ParticleLayout< T, Dim > & getLayout()
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)
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
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
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)
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)
@ MAX_NUM_PARTICLES
max. particles per cell
Definition AmrObject.h:45
@ MIN_NUM_PARTICLES
min. particles per cell
Definition AmrObject.h:44
double scaling_m
Scaling factor for tagging [0, 1].
Definition AmrObject.h:181
IpplTimings::TimerRef amrRegridTimer_m
Definition AmrObject.h:193
std::pair< Vector_t, Vector_t > VectorPair_t
Definition AmrObject.h:35
double chargedensity_m
Tagging value for CHARGE_DENSITY.
Definition AmrObject.h:183
bool refined_m
Only set to true in AmrObject::initFineLevels().
Definition AmrObject.h:189
TaggingCriteria tagging_m
Tagging strategy.
Definition AmrObject.h:179
IpplTimings::TimerRef amrSolveTimer_m
timer for selfField calculation (used in concrete AmrObject classes)
Definition AmrObject.h:192
size_t minNumPart_m
Tagging value for MIN_NUM_PARTICLES.
Definition AmrObject.h:187
size_t maxNumPart_m
Tagging value for MAX_NUM_PARTICLES.
Definition AmrObject.h:185
void writeFields(const amr::AmrScalarFieldContainer_t &rho, const amr::AmrScalarFieldContainer_t &phi, const amr::AmrVectorFieldContainer_t &efield, const amr::AmrIntArray_t &refRatio, const amr::AmrGeomContainer_t &geom, const int &nLevel, const double &time, const double &scale)
void writeBunch(const AmrPartBunch *bunch_p, const double &time, const double &gamma)
void gather(ParticleAttrib< FT > &attrib, AmrVectorFieldContainer_t &f, ParticleAttrib< Vektor< PT, Dim > > &pp, int lbase, int lfine)
void scatter(ParticleAttrib< FT > &attrib, AmrScalarFieldContainer_t &f, ParticleAttrib< Vektor< PT, Dim > > &pp, int lbase, int lfine, const ParticleAttrib< int > &pbin, int bin=-1)
AmrParticle_t pbase_t
virtual void solve(AmrScalarFieldContainer_t &, AmrScalarFieldContainer_t &, AmrVectorFieldContainer_t &, unsigned short, unsigned short, bool=true)
virtual void hasToRegrid()
The base class for all OPAL exceptions.
bool isForbidTransform() const
const double & getScalingFactor() const
const ParticleLevelCounter_t & getLocalNumPerLevel() const
void setForbidTransform(bool forbidTransform)
const double & domainMapping(bool inverse=false)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)
Vektor< double, 3 > Vector_t