chombo-discharge
Loading...
Searching...
No Matches
CD_ItoKMCGodunovStepperImplem.H
Go to the documentation of this file.
1/* chombo-discharge
2 * Copyright © 2021 SINTEF Energy Research.
3 * Please refer to Copyright.txt and LICENSE in the chombo-discharge root directory.
4 */
5
12#ifndef CD_ItoKMCGodunovStepperImplem_H
13#define CD_ItoKMCGodunovStepperImplem_H
14
15// Chombo includes
16#include <ParmParse.H>
17
18// Our includes
20#include <CD_Timer.H>
21#include <CD_ParallelOps.H>
22#include <CD_DataOps.H>
23#include <CD_Units.H>
24#include <CD_Photon.H>
25#include <CD_DischargeIO.H>
26#include <CD_NamespaceHeader.H>
27
28using namespace Physics::ItoKMC;
29
30template <typename I, typename C, typename R, typename F>
32 : ItoKMCStepper<I, C, R, F>(a_physics)
33{
34 CH_TIME("ItoKMCGodunovStepper::ItoKMCGodunovStepper");
35
36 this->m_name = "ItoKMCGodunovStepper";
37 this->m_prevDt = 0.0;
38 this->m_writeCheckpointParticles = false;
39 this->m_readCheckpointParticles = false;
40 this->m_extendConductivityEB = false;
41 this->m_canRegridOnRestart = true;
42 this->m_prevDt = 0.0;
43 this->m_maxFieldAbort = std::numeric_limits<Real>::max();
44
45 if (a_parseOptions) {
46 this->parseOptions();
47 }
48}
49
50template <typename I, typename C, typename R, typename F>
52{
53 CH_TIME("ItoKMCGodunovStepper::~ItoKMCGodunovStepper");
54 if (this->m_verbosity > 5) {
55 pout() << "ItoKMCGodunovStepper::~ItoKMCGodunovStepper" << endl;
56 }
57}
58
59template <typename I, typename C, typename R, typename F>
60void
62{
63 CH_TIME("ItoKMCGodunovStepper::registerOperators");
64 if (this->m_verbosity > 5) {
65 pout() << "ItoKMCGodunovStepper::registerOperators" << endl;
66 }
67
69
70 // Register this because we must be able to deposit particles inside dielectrics (due to particle diffusion across the EB).
71 (this->m_amr)->registerOperator(s_particle_mesh, this->m_particleRealm, phase::solid);
73
74template <typename I, typename C, typename R, typename F>
75void
77{
78 CH_TIME("ItoKMCGodunovStepper::allocate");
79 if (this->m_verbosity > 5) {
80 pout() << "ItoKMCGodunovStepper::allocate" << endl;
81 }
82
84
85 // Now allocate for the conductivity particles and rho^dagger particles. This is only done in the 'allocate' routine
86 // and not in 'allocateInternals' because that would discard the particles during regrids. That has definitely never
87 // happen, and there's no way I've spent countless hours tracking down such a bug.
88 const int numItoSpecies = this->m_physics->getNumItoSpecies();
89
90 m_conductivityParticles.resize(numItoSpecies);
91 m_irregularParticles.resize(numItoSpecies);
92 m_rhoDaggerParticles.resize(numItoSpecies);
93
94 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
95 const int idx = solverIt.index();
96
97 m_conductivityParticles[idx] = RefCountedPtr<ParticleContainer<PointParticle>>(
101
102 (this->m_amr)->allocate(*m_conductivityParticles[idx], this->m_particleRealm);
103 (this->m_amr)->allocate(*m_irregularParticles[idx], this->m_particleRealm);
104 (this->m_amr)->allocate(*m_rhoDaggerParticles[idx], this->m_particleRealm);
105 }
106}
107
108template <typename I, typename C, typename R, typename F>
109void
111{
112 CH_TIME("ItoKMCGodunovStepper::allocateInternals");
113 if (this->m_verbosity > 5) {
114 pout() << this->m_name + "::allocateInternals" << endl;
115 }
116
118
119 const int numCdrSpecies = this->m_physics->getNumCdrSpecies();
120
121 m_cdrDivD.resize(numCdrSpecies);
122 for (int i = 0; i < numCdrSpecies; i++) {
123 this->m_amr->allocate(m_cdrDivD[i], this->m_fluidRealm, this->m_plasmaPhase, 1);
124 }
125
126 this->m_amr->allocate(m_semiImplicitRhoCDR, this->m_fluidRealm, this->m_plasmaPhase, 1);
127 this->m_amr->allocate(m_semiImplicitConductivityCDR, this->m_fluidRealm, this->m_plasmaPhase, 1);
128}
129
130template <typename I, typename C, typename R, typename F>
131void
133{
134 CH_TIME("ItoKMCGodunovStepper::barrier");
135 if (this->m_verbosity > 5) {
136 pout() << this->m_name + "::barrier" << endl;
137 }
138
139 if ((this->m_profile)) {
141 }
142}
143
144template <typename I, typename C, typename R, typename F>
145void
148 CH_TIME("ItoKMCGodunovStepper::parseOptions");
149 if (this->m_verbosity > 5) {
150 pout() << this->m_name + "::parseOptions" << endl;
151 }
152
154
155 this->parseAlgorithm();
156 this->parseFiltering();
157 this->parseCheckpointParticles();
158 this->parseSecondaryEmissionSpecification();
159}
160
161template <typename I, typename C, typename R, typename F>
162void
164{
165 CH_TIME("ItoKMCGodunovStepper::parseRuntimeOptions");
166 if (this->m_verbosity > 5) {
167 pout() << this->m_name + "::parseRuntimeOptions" << endl;
168 }
169
171
172 this->parseAlgorithm();
173 this->parseFiltering();
174 this->parseCheckpointParticles();
175 this->parseSecondaryEmissionSpecification();
176}
177
178template <typename I, typename C, typename R, typename F>
179void
181{
182 CH_TIME("ItoKMCGodunovStepper::parseAlgorithm");
183 if (this->m_verbosity > 5) {
184 pout() << this->m_name + "::parseAlgorithm" << endl;
185 }
186
187 ParmParse pp(this->m_name.c_str());
189
190 pp.get("extend_conductivity", m_extendConductivityEB);
191 pp.get("algorithm", str);
192 pp.get("abort_max_field", m_maxFieldAbort);
193
194 // Get algorithm
195 if (str == "euler_maruyama") {
196 m_algorithm = WhichAlgorithm::EulerMaruyama;
197 }
198 else {
199 MayDay::Abort("ItoKMCGodunovStepper::parseAlgorithm - unknown algorithm requested");
200 }
201}
202
203template <typename I, typename C, typename R, typename F>
204void
206{
207 CH_TIME("ItoKMCGodunovStepper::parseFiltering");
208 if (this->m_verbosity > 5) {
209 pout() << this->m_name + "::parseFiltering" << endl;
210 }
211
212 ParmParse pp(this->m_name.c_str());
214
215 m_rhoFilterNum = -1;
216 m_rhoFilterMaxStride = 1;
217 m_rhoFilterAlpha = 0.5;
218
219 m_condFilterNum = -1;
220 m_condFilterMaxStride = 1;
221 m_condFilterAlpha = 0.5;
222
223 pp.get("rho_filter_num", m_rhoFilterNum);
224 pp.get("rho_filter_max_stride", m_rhoFilterMaxStride);
225 pp.get("rho_filter_alpha", m_rhoFilterAlpha);
226
227 pp.get("cond_filter_num", m_condFilterNum);
228 pp.get("cond_filter_max_stride", m_condFilterMaxStride);
229 pp.get("cond_filter_alpha", m_condFilterAlpha);
230
232 MayDay::Abort("ItoKMCGodunovStepper::parseFiltering -- cannot have alpha <= 0 or alpha >= 1 for rho_filter");
233 }
235 MayDay::Abort("ItoKMCGodunovStepper::parseFiltering -- cannot have alpha <= 0 or alpha >= 1 for cond_filter");
236 }
237}
238
239template <typename I, typename C, typename R, typename F>
240void
242{
243 CH_TIME("ItoKMCGodunovStepper::parseCheckpointParticles");
244 if (this->m_verbosity > 5) {
245 pout() << this->m_name + "::parseCheckpointParticles" << endl;
246 }
247
248 ParmParse pp(this->m_name.c_str());
249
250 pp.query("checkpoint_particles", m_writeCheckpointParticles);
251}
252
253template <typename I, typename C, typename R, typename F>
254void
256{
257 CH_TIME("ItoKMCGodunovStepper::parseSecondaryEmissionSpecifiation");
258 if (this->m_verbosity > 5) {
259 pout() << this->m_name + "::parseSecondaryEmissionSpecification" << endl;
260 }
261
262 ParmParse pp(this->m_name.c_str());
263
265
266 pp.query("secondary_emission", str);
267
268 if (str == "before_reactions") {
269 m_emitSecondaryParticlesBeforeReactions = true;
270 }
271 else if (str == "after_reactions") {
272 m_emitSecondaryParticlesBeforeReactions = false;
273 }
274 else {
276
277 err = "ItoKMCGodunovStepper::parseSecondaryEmissionSpecification - expected 'before_reactions' or 'after_reactions'";
278 err += "but got" + str;
279
280 MayDay::Abort(err.c_str());
282}
283
284template <typename I, typename C, typename R, typename F>
285Real
288 CH_TIME("ItoKMCGodunovStepper::computeDt");
289 if (this->m_verbosity > 5) {
290 pout() << this->m_name + "::computeDt" << endl;
291 }
292
294
295 if ((this->m_maxReducedField > m_maxFieldAbort) && (m_maxFieldAbort > 0.0)) {
296 pout() << this->m_name + " stopping because maximum field is too high (" << this->m_maxReducedField << ")" << endl;
297
298 this->m_keepGoing = false;
300
301 return dt;
302}
303
304template <typename I, typename C, typename R, typename F>
307{
308 CH_TIME("ItoKMCGodunovStepper::advance");
309 if (this->m_verbosity > 5) {
310 pout() << this->m_name + "::advance" << endl;
312
313 auto debugCharge = [this](const std::string a_message) -> void {
314 const Real Qtot = this->computeTotalCharge();
315 if (procID() == 0) {
316 std::cout << a_message << " = " << Qtot << "\n";
317 }
318 };
319
320 // Special flag for telling the class that we have the necessary things in place for doing a regrid. This is
321 // an if-but-maybe situation where the user chose not to checkpoint the particles we need for regrids, yet tries
322 // to restart a simulation and regrid without all the prerequisites being in place. This flag is set to true
323 // because these requirements are checked during the regrid routine.
324 m_canRegridOnRestart = true;
325
326 m_timer = Timer("ItoKMCGodunovStepper::advance");
327
328 // debugCharge("advance");
329
330 // Previous time step is needed when regridding.
331 this->m_prevDt = a_dt;
332
333 // Done only so we can plot the absorbed photons (advanceReactionNetwork absorbs them)
334 m_timer.startEvent("Deposit photons");
335 for (auto solverIt = (this->m_rte)->iterator(); solverIt.ok(); ++solverIt) {
337
338 EBAMRCellData& phi = solver->getPhi();
339 ParticleContainer<Photon>& photons = solver->getBulkPhotons();
340
341 solver->depositPhotons<Photon, const Real&, &Photon::weight>(phi, photons, DepositionType::NGP);
342 }
343 m_timer.stopEvent("Deposit photons");
344
345 // ====== BEGIN TRANSPORT STEP ======
346 // Semi-implicitly advance the particles and the field.
347 switch (m_algorithm) {
348 case WhichAlgorithm::EulerMaruyama: {
349 this->advanceEulerMaruyama(a_dt);
350
351 break;
352 }
353 default: {
354 MayDay::Abort("ItoKMCGodunovStepper::advance - logic bust");
356 break;
357 }
358 }
359
360 // Do intersection test and remove particles that struck the EB or domain. Transfer them to appropriate containers. Then recompute the number of particles per cell.
361 this->barrier();
362 m_timer.startEvent("EB/Particle intersection");
363 if (m_extendConductivityEB) {
365 // TLDR: This hook does a special intersection routine where instead of transferring the particles that were intersected, they
366 // are put in a separate data container. We want to do this in order to reduce artificial gradients in the particle
367 // densities near the EB. The way we do this is that we flag the particles during the intersection; particles that are flagged
368 // are transferred to a container which is added to the conductivity-related particles. The particles are later removed.
369 auto setFlag = [](ItoParticle& p) -> void {
370 p.tmpReal() = 1.0;
371 };
372 auto nonDeletionModifier = [](ItoParticle& p) -> void {
373 p.tmpReal() = -1.0;
374 };
375 auto removeCrit = [](const ItoParticle& p) -> bool {
376 return p.tmpReal() < 0.0;
377 };
379 // Set flag for identifying which particles were intersected by not removed.
380 for (auto it = this->m_ito->iterator(); it.ok(); ++it) {
381 ParticleOps::setData<ItoParticle>(it()->getParticles(ItoSolver::WhichContainer::Bulk), setFlag);
382 }
383
384 // Intersect particles, but don't remove them.
385 const bool deleteParticles = false;
386 this->intersectParticles(SpeciesSubset::AllMobileOrDiffusive, deleteParticles, nonDeletionModifier);
387
388 // Particles that were not removed are copied to a separate data container.
389 for (auto it = this->m_ito->iterator(); it.ok(); ++it) {
391 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
392
393 const int idx = it.index();
394 const int Z = species->getChargeNumber();
395
396 ParticleContainer<PointParticle>& irregParticles = *m_irregularParticles[idx];
397 ParticleContainer<ItoParticle>& ebParticles = solver->getParticles(ItoSolver::WhichContainer::EB);
398 ParticleContainer<ItoParticle>& bulkParticles = it()->getParticles(ItoSolver::WhichContainer::Bulk);
399
400 irregParticles.clearParticles();
401
402 if (Z != 0 && solver->isMobile()) {
403 for (int lvl = 0; lvl <= (this->m_amr)->getFinestLevel(); lvl++) {
404 const DisjointBoxLayout& dbl = (this->m_amr)->getGrids((this->m_particleRealm))[lvl];
405 const DataIterator& dit = dbl.dataIterator();
406
407 const int nbox = dit.size();
408
409#pragma omp parallel for schedule(runtime)
410 for (int mybox = 0; mybox < nbox; mybox++) {
412
415
417 const ItoParticle& p = lit();
418
419 if (p.tmpReal() < 0.0) {
420 const RealVect& pos = p.position();
421 const Real& weight = p.weight();
422 const Real& mobility = p.mobility();
423
424 pointParticles.add(PointParticle(pos, weight * mobility));
425 }
426 }
428 }
429 }
430 }
431
432 // Finally, delete the original particles
433 for (auto it = this->m_ito->iterator(); it.ok(); ++it) {
434 ParticleContainer<ItoParticle>& bulkParticles = it()->getParticles(ItoSolver::WhichContainer::Bulk);
435 ParticleOps::removeParticles<ItoParticle>(bulkParticles, removeCrit);
436 }
438 else {
439 const bool deleteParticles = true;
440
441 this->intersectParticles(SpeciesSubset::AllMobileOrDiffusive, deleteParticles);
442 }
443
444 // The intersection tests above may not have caught all particles -- do a cleanup sweep where particles on the wrong side of
445 // the EB are put on the EB.
446 for (auto it = this->m_ito->iterator(); it.ok(); ++it) {
448
449 ParticleContainer<ItoParticle>& ebParticles = solver->getParticles(ItoSolver::WhichContainer::EB);
450 ParticleContainer<ItoParticle>& bulkParticles = solver->getParticles(ItoSolver::WhichContainer::Bulk);
451
452 this->m_amr->transferIrregularParticles(ebParticles, bulkParticles, this->m_plasmaPhase);
453 }
454 m_timer.stopEvent("EB/Particle intersection");
455 // ====== END TRANSPORT STEP ======
456
457 // Photon transport
458 this->barrier();
459 m_timer.startEvent("Photon transport");
460 this->advancePhotons(a_dt);
461 m_timer.stopEvent("Photon transport");
462
463 // Compute the gradients of the various species densities - this is used in the KMC kernels.
464 if ((this->m_physics)->needGradients()) {
465 m_timer.startEvent("Gradient calculation");
466 (this->m_ito)->depositParticles();
467 this->computeDensityGradients();
468 m_timer.stopEvent("Gradient calculation");
469 }
470
471 // Resolve secondary emission. We have filled the relevant particles in the transport step -- this can be done
472 // either before or after the reactions.
473 if (m_emitSecondaryParticlesBeforeReactions) {
474 this->barrier();
475 m_timer.startEvent("EB particle injection");
476 this->fillSecondaryEmissionEB(a_dt);
477 this->resolveSecondaryEmissionEB(a_dt);
478 m_timer.stopEvent("EB particle injection");
479 }
480
481 // Sort the particles and photons per cell so we can call reaction algorithms
482 this->barrier();
483 m_timer.startEvent("Sort by cell");
484 (this->m_ito)->organizeParticlesByCell(ItoSolver::WhichContainer::Bulk);
485 this->sortPhotonsByCell(McPhoto::WhichContainer::Bulk);
486 this->sortPhotonsByCell(McPhoto::WhichContainer::Source);
487 m_timer.stopEvent("Sort by cell");
488
489 // Run the Kinetic Monte Carlo reaction kernels.
490 this->barrier();
491 m_timer.startEvent("Reaction network");
492 this->advanceReactionNetwork(a_dt);
493 m_timer.stopEvent("Reaction network");
494
495 // Sort particles per patch.
496 this->barrier();
497 m_timer.startEvent("Sort by patch");
498 (this->m_ito)->organizeParticlesByPatch(ItoSolver::WhichContainer::Bulk);
499 this->sortPhotonsByPatch(McPhoto::WhichContainer::Bulk);
500 this->sortPhotonsByPatch(McPhoto::WhichContainer::Source);
501 m_timer.stopEvent("Sort by patch");
502
503 // Resolve secondary emission. We have filled the relevant particles in the transport step -- this can be done
504 // either before or after the reactions.
505 if (!m_emitSecondaryParticlesBeforeReactions) {
506 this->barrier();
507 m_timer.startEvent("EB particle injection");
508 this->fillSecondaryEmissionEB(a_dt);
509 this->resolveSecondaryEmissionEB(a_dt);
510 m_timer.stopEvent("EB particle injection");
511 }
512
513 // Remove particles that are inside the EB -- this is not a part of the algorithm, just a safety measure to make sure
514 // we don't do things with particles that lie inside the EB.
515 this->barrier();
516 m_timer.startEvent("Remove covered");
517 this->removeCoveredParticles(SpeciesSubset::AllMobileOrDiffusive, EBRepresentation::Discrete, this->m_toleranceEB);
518 m_timer.stopEvent("Remove covered");
519
520 // Clear BC data holders.
521 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
522 solverIt()->clear(ItoSolver::WhichContainer::EB);
523 solverIt()->clear(ItoSolver::WhichContainer::Domain);
524 }
525
526 // Prepare for the next time step
527 this->barrier();
528 m_timer.startEvent("Post-compute v");
529 this->computeDriftVelocities();
530 m_timer.stopEvent("Post-compute v");
531
532 this->barrier();
533 m_timer.startEvent("Post-compute D");
534 this->computeDiffusionCoefficients();
535 m_timer.stopEvent("Post-compute D");
536
537 this->computePhysicsDt();
538
539 if ((this->m_profile)) {
540 m_timer.eventReport(pout(), false);
541 }
542
543 m_timer.clear();
544
545 // Compute the maximum field (in Townsend).
546 this->m_maxReducedField = this->computeMaxReducedElectricField(this->m_plasmaPhase);
547
548 return a_dt;
549}
550
551template <typename I, typename C, typename R, typename F>
552void
554{
555 CH_TIME("ItoKMCGodunovStepper::preRegrid");
556 if (this->m_verbosity > 5) {
557 pout() << "ItoKMCGodunovStepper::preRegrid" << endl;
558 }
559
560 const int numItoSpecies = (this->m_physics)->getNumItoSpecies();
561 const int numCdrSpecies = (this->m_physics)->getNumCdrSpecies();
562 const int numPlasmaSpecies = (this->m_physics)->getNumPlasmaSpecies();
563 const int numPhotonSpecies = (this->m_physics)->getNumPhotonSpecies();
564
566
567 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
568 const int idx = solverIt.index();
569
570 m_conductivityParticles[idx]->preRegrid(a_lmin);
571 m_irregularParticles[idx]->preRegrid(a_lmin);
572 m_rhoDaggerParticles[idx]->preRegrid(a_lmin);
573 }
574
575 this->m_amr->allocate(m_scratchSemiImplicitRhoCDR, this->m_fluidRealm, this->m_plasmaPhase, 1);
576 this->m_amr->allocate(m_scratchSemiImplicitConductivityCDR, this->m_fluidRealm, this->m_plasmaPhase, 1);
577
578 DataOps::copy(m_scratchSemiImplicitRhoCDR, m_semiImplicitRhoCDR);
579 DataOps::copy(m_scratchSemiImplicitConductivityCDR, m_semiImplicitConductivityCDR);
580
581 // Release unecessary storage.
582 for (int i = 0; i < numCdrSpecies; i++) {
583 m_cdrDivD[i].clear();
584 }
585
586 m_semiImplicitRhoCDR.clear();
587 m_semiImplicitConductivityCDR.clear();
588}
589
590template <typename I, typename C, typename R, typename F>
591void
593 const int a_oldFinestLevel,
594 const int a_newFinestLevel) noexcept
595{
596 CH_TIME("ItoKMCGodunovStepper::regrid");
597 if (this->m_verbosity > 5) {
598 pout() << "ItoKMCGodunovStepper::regrid" << endl;
599 }
600
601 m_timer = Timer("ItoKMCGodunovStepper::regrid");
602
603 // A special flag for aborting the simulation if the user did NOT put checkpoint particles in the checkpoint
604 // file but still try to regrid on restart.
605 if (!m_canRegridOnRestart) {
606 const std::string baseErr = "ItoKMCGodunovStepper::regrid -- can't regrid because";
607 const std::string err1 = "checkpoint file does not contain particles. Set Driver.initial_regrids=0";
608
609 pout() << baseErr + err1 << endl;
610
611 MayDay::Error((baseErr + err1).c_str());
612 }
613
614 // Regrid solvers
615 m_timer.startEvent("Regrid ItoSolver");
616 (this->m_ito)->regrid(a_lmin, a_oldFinestLevel, a_newFinestLevel);
617 if (this->m_timeStep == 0) {
618 // Necessary because first time step is not semi-implicit, so it will use the densities from the
619 // Ito and CDr solvers. Howver, ItoSolver does not re-deposit the particles during regrid, so we
620 // enforce it here.
621 (this->m_ito)->depositParticles();
622 }
623 m_timer.stopEvent("Regrid ItoSolver");
624
625 m_timer.startEvent("Regrid CdrSolver");
626 (this->m_cdr)->regrid(a_lmin, a_oldFinestLevel, a_newFinestLevel);
627 m_timer.stopEvent("Regrid CdrSolver");
628
629 m_timer.startEvent("Regrid FieldSolver");
630 (this->m_fieldSolver)->regrid(a_lmin, a_oldFinestLevel, a_newFinestLevel);
631 m_timer.stopEvent("Regrid FieldSolver");
632
633 m_timer.startEvent("Regrid RTE");
634 (this->m_rte)->regrid(a_lmin, a_oldFinestLevel, a_newFinestLevel);
635 m_timer.stopEvent("Regrid RTE");
636
637 m_timer.startEvent("Regrid SurfaceODESolver");
638 this->m_sigmaSolver->regrid(a_lmin, a_oldFinestLevel, a_newFinestLevel);
639 m_timer.stopEvent("Regrid SurfaceODESolver");
640
641 // Allocate internal memory for ItoKMCGodunovStepper now....
642 m_timer.startEvent("Allocate internals");
643 this->allocateInternals();
644 m_timer.stopEvent("Allocate internals");
645
646 // We need to remap/regrid the stored data required for the semi-implicit update as well.
647 m_timer.startEvent("Remap algorithm-particles");
648 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
649 const int idx = solverIt.index();
650 (this->m_amr)->remapToNewGrids(*m_rhoDaggerParticles[idx], a_lmin, a_newFinestLevel);
651 (this->m_amr)->remapToNewGrids(*m_conductivityParticles[idx], a_lmin, a_newFinestLevel);
652 (this->m_amr)->remapToNewGrids(*m_irregularParticles[idx], a_lmin, a_newFinestLevel);
653 }
654 m_timer.stopEvent("Remap algorithm-particles");
655
656 // Also regrid the space charge density contribution from the CDR equations
657 this->m_amr->interpToNewGrids(m_semiImplicitRhoCDR,
658 m_scratchSemiImplicitRhoCDR,
659 this->m_plasmaPhase,
660 a_lmin,
663 EBCoarseToFineInterp::Type::ConservativePWC);
664
665 this->m_amr->interpToNewGrids(m_semiImplicitConductivityCDR,
666 m_scratchSemiImplicitConductivityCDR,
667 this->m_plasmaPhase,
668 a_lmin,
671 EBCoarseToFineInterp::Type::ConservativePWC);
672
673 // Set up the field solver with standard coefficients or with
674 // modified coefficients if we are reusing data from the last time step.
675 m_timer.startEvent("Setup field solver");
676 (this->m_fieldSolver)->setupSolver();
677 this->computeConductivities(m_conductivityParticles);
678 this->setupSemiImplicitPoisson(this->m_prevDt);
679 m_timer.stopEvent("Setup field solver");
680
681 // Solve the Poisson equation.
682 m_timer.startEvent("Solve Poisson");
683 if (this->m_timeStep == 0) {
684 this->computeSpaceChargeDensity();
685 }
686 else {
687 this->depositPointParticles(m_rhoDaggerParticles, SpeciesSubset::All);
688 this->computeSemiImplicitRho();
689 }
690
691 const bool converged = this->solvePoisson();
692
693 if (!converged) {
694 const std::string errMsg = "ItoKMCGodunovStepper::regrid - Poisson solve did not converge after regrid";
695
696 pout() << errMsg << endl;
697
698 if (this->m_abortOnFailure) {
699 MayDay::Error(errMsg.c_str());
700 }
701 }
702 m_timer.stopEvent("Solve Poisson");
703
704 // Regrid superparticles.
705 if (this->m_regridSuperparticles) {
706 m_timer.startEvent("Make superparticles");
707 (this->m_ito)->organizeParticlesByCell(ItoSolver::WhichContainer::Bulk);
708 (this->m_ito)->makeSuperparticles(ItoSolver::WhichContainer::Bulk, (this->m_particlesPerCell));
709 (this->m_ito)->organizeParticlesByPatch(ItoSolver::WhichContainer::Bulk);
710 m_timer.stopEvent("Make superparticles");
711 }
712
713 // Now let Ihe ito solver deposit its actual particles... In the above it deposit m_rhoDaggerParticles.
714 m_timer.startEvent("Deposit particles");
715 (this->m_ito)->depositParticles();
716 m_timer.stopEvent("Deposit particles");
717
718 // Recompute new velocities and diffusion coefficients
719 m_timer.startEvent("Prepare next step");
720 this->computeDiffusionCoefficients();
721 this->computeDriftVelocities();
722 m_timer.stopEvent("Prepare next step");
723
724 m_timer.eventReport(pout(), false);
725
726 // No reason to have these lying around.
727 m_scratchSemiImplicitRhoCDR.clear();
728 m_scratchSemiImplicitConductivityCDR.clear();
729
730 // Fill the neutral density on the mesh
731 this->fillNeutralDensity();
732}
733
734template <typename I, typename C, typename R, typename F>
735void
737{
738 CH_TIME("ItoKMCGodunovStepper::setOldPositions");
739 if (this->m_verbosity > 5) {
740 pout() << this->m_name + "::setOldPositions" << endl;
741 }
742
743 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
745
746 for (int lvl = 0; lvl <= (this->m_amr)->getFinestLevel(); lvl++) {
747 const DisjointBoxLayout& dbl = (this->m_amr)->getGrids((this->m_particleRealm))[lvl];
748 const DataIterator& dit = dbl.dataIterator();
749
750 ParticleData<ItoParticle>& particles = solver->getParticles(ItoSolver::WhichContainer::Bulk)[lvl];
751
752 const int nbox = dit.size();
753
754#pragma omp parallel for schedule(runtime)
755 for (int mybox = 0; mybox < nbox; mybox++) {
756 const DataIndex& din = dit[mybox];
757
759
761 ItoParticle& p = lit();
762
763 p.oldPosition() = p.position();
764 }
765 }
766 }
767 }
768}
769
770template <typename I, typename C, typename R, typename F>
771void
774 const SpeciesSubset a_subset) noexcept
775{
776 CH_TIME("ItoKMCGodunovStepper::remapPointParticles");
777 if (this->m_verbosity > 5) {
778 pout() << this->m_name + "::remapPointParticles" << endl;
779 }
780
781 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
783 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
784
785 const int idx = solverIt.index();
786
787 const bool mobile = solver->isMobile();
788 const bool diffusive = solver->isDiffusive();
789 const bool charged = species->getChargeNumber() != 0;
790
791 switch (a_subset) {
792 case SpeciesSubset::All: {
794
795 break;
796 }
797 case SpeciesSubset::AllMobile: {
798 if (mobile) {
800 }
801
802 break;
803 }
804 case SpeciesSubset::AllDiffusive: {
805 if (diffusive) {
807 }
808
809 break;
810 }
811 case SpeciesSubset::AllMobileOrDiffusive: {
812 if (mobile || diffusive) {
814 }
815
816 break;
817 }
818 case SpeciesSubset::AllMobileAndDiffusive: {
819 if (mobile && diffusive) {
821 }
822
823 break;
824 }
825 case SpeciesSubset::Charged: {
826 if (charged) {
828 }
829
830 break;
831 }
832 case SpeciesSubset::ChargedMobile: {
833 if (charged && mobile) {
835 }
836
837 break;
838 }
839 case SpeciesSubset::ChargedDiffusive: {
840 if (charged && diffusive) {
842 }
843
844 break;
845 }
846 case SpeciesSubset::ChargedMobileOrDiffusive: {
847 if (charged && (mobile || diffusive)) {
849 }
850
851 break;
852 }
853 case SpeciesSubset::ChargedMobileAndDiffusive: {
854 if (charged && (mobile && diffusive)) {
856 }
857
858 break;
859 }
860 case SpeciesSubset::Stationary: {
861 if (!mobile && !diffusive) {
863 }
864
865 break;
866 }
867 default: {
868 MayDay::Abort("ItoKMCGodunovStepper::remapPointParticles - logic bust");
869
870 break;
871 }
872 }
873 }
874}
875
876template <typename I, typename C, typename R, typename F>
877void
880 const SpeciesSubset a_subset) noexcept
881{
882 CH_TIME("ItoKMCGodunovStepper::depositPointParticles");
883 if (this->m_verbosity > 5) {
884 pout() << this->m_name + "::depositPointParticles" << endl;
885 }
886
887 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
889 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
890
891 const int idx = solverIt.index();
892
893 const bool mobile = solver->isMobile();
894 const bool diffusive = solver->isDiffusive();
895 const bool charged = species->getChargeNumber() != 0;
896
897 switch (a_subset) {
898 case SpeciesSubset::All: {
900
901 break;
902 }
903 case SpeciesSubset::AllMobile: {
904 if (mobile) {
906 *a_particles[idx]);
907 }
908
909 break;
910 }
911 case SpeciesSubset::AllDiffusive: {
912 if (diffusive) {
914 *a_particles[idx]);
915 }
916
917 break;
918 }
919 case SpeciesSubset::AllMobileOrDiffusive: {
920 if (mobile || diffusive) {
922 *a_particles[idx]);
923 }
924 break;
925 }
926 case SpeciesSubset::AllMobileAndDiffusive: {
927 if (mobile && diffusive) {
929 *a_particles[idx]);
930 }
931 break;
932 }
933 case SpeciesSubset::Charged: {
934 if (charged) {
936 *a_particles[idx]);
937 }
938 break;
939 }
940 case SpeciesSubset::ChargedMobile: {
941 if (charged && mobile) {
943 *a_particles[idx]);
944 }
945 break;
946 }
947 case SpeciesSubset::ChargedDiffusive: {
948 if (charged && diffusive) {
950 *a_particles[idx]);
951 }
952
953 break;
954 }
955 case SpeciesSubset::ChargedMobileOrDiffusive: {
956 if (charged && (mobile || diffusive)) {
958 *a_particles[idx]);
959 }
960
961 break;
962 }
963 case SpeciesSubset::ChargedMobileAndDiffusive: {
964 if (charged && (mobile && diffusive)) {
966 *a_particles[idx]);
967 }
968
969 break;
970 }
971 case SpeciesSubset::Stationary: {
972 if (!mobile && !diffusive) {
974 *a_particles[idx]);
975 }
976
977 break;
978 }
979 default: {
980 MayDay::Abort("ItoKMCGodunovStepper::depositPointParticles - logic bust");
981
982 break;
983 }
984 }
985 }
986}
987
988template <typename I, typename C, typename R, typename F>
989void
992 const SpeciesSubset a_subset) noexcept
993{
994 CH_TIME("ItoKMCGodunovStepper::clearPointParticles");
995 if (this->m_verbosity > 5) {
996 pout() << this->m_name + "::clearPointParticles" << endl;
997 }
998
999 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
1001 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
1002
1003 const int idx = solverIt.index();
1004
1005 const bool mobile = solver->isMobile();
1006 const bool diffusive = solver->isDiffusive();
1007 const bool charged = species->getChargeNumber() != 0;
1008
1009 switch (a_subset) {
1010 case SpeciesSubset::All: {
1011 a_particles[idx]->clearParticles();
1012
1013 break;
1014 }
1015 case SpeciesSubset::AllMobile: {
1016 if (mobile) {
1017 a_particles[idx]->clearParticles();
1018 }
1019
1020 break;
1021 }
1022 case SpeciesSubset::AllDiffusive: {
1023 if (diffusive) {
1024 a_particles[idx]->clearParticles();
1025 }
1026
1027 break;
1028 }
1029 case SpeciesSubset::AllMobileOrDiffusive: {
1030 if (mobile || diffusive) {
1031 a_particles[idx]->clearParticles();
1032 }
1033
1034 break;
1035 }
1036 case SpeciesSubset::AllMobileAndDiffusive: {
1037 if (mobile && diffusive) {
1038 a_particles[idx]->clearParticles();
1039 }
1040
1041 break;
1042 }
1043 case SpeciesSubset::Charged: {
1044 if (charged) {
1045 a_particles[idx]->clearParticles();
1046 }
1047
1048 break;
1049 }
1050 case SpeciesSubset::ChargedMobile: {
1051 if (charged && mobile) {
1052 a_particles[idx]->clearParticles();
1053 }
1054
1055 break;
1056 }
1057 case SpeciesSubset::ChargedDiffusive: {
1058 if (charged && diffusive) {
1059 a_particles[idx]->clearParticles();
1060 }
1061
1062 break;
1063 }
1064 case SpeciesSubset::ChargedMobileOrDiffusive: {
1065 if (charged && (mobile || diffusive)) {
1066 a_particles[idx]->clearParticles();
1067 }
1068
1069 break;
1070 }
1071 case SpeciesSubset::ChargedMobileAndDiffusive: {
1072 if (charged && (mobile && diffusive)) {
1073 a_particles[idx]->clearParticles();
1074 }
1075
1076 break;
1077 }
1078 case SpeciesSubset::Stationary: {
1079 if (!mobile && !diffusive) {
1080 a_particles[idx]->clearParticles();
1081 }
1082
1083 break;
1084 }
1085 default: {
1086 MayDay::Abort("ItoKMCGodunovStepper::clearPointParticles - logic bust");
1087
1088 break;
1089 }
1090 }
1091 }
1092}
1093
1094template <typename I, typename C, typename R, typename F>
1095void
1098{
1099 CH_TIME("ItoKMCGodunovStepper::computeConductivities");
1100 if (this->m_verbosity > 5) {
1101 pout() << this->m_name + "::computeConductivities" << endl;
1102 }
1103
1104 this->computeCellConductivity((this->m_conductivityCell), a_particles);
1105 this->computeFaceConductivity();
1106}
1107
1108template <typename I, typename C, typename R, typename F>
1109void
1113{
1114 CH_TIME("ItoKMCGodunovStepper::computeCellConductivity(EBAMRCellData, PointParticle");
1115 if (this->m_verbosity > 5) {
1116 pout() << this->m_name + "::computeCellConductivity(EBAMRCellData, PointParticle)" << endl;
1117 }
1118
1120
1121 // Contribution from Ito solvers.
1122 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
1124 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
1125
1126 const int idx = solverIt.index();
1127 const int Z = species->getChargeNumber();
1128
1129 if (Z != 0 && solver->isMobile()) {
1130 // Deposit on the particle realm.
1131 DataOps::setValue(this->m_particleScratch1, 0.0);
1132 solver->depositParticles<PointParticle, const Real&, &PointParticle::weight>(this->m_particleScratch1,
1133 *a_particles[idx]);
1134
1135 // Copy to fluid realm and add to total conductivity.
1136 (this->m_amr)->copyData(this->m_fluidScratch1, this->m_particleScratch1);
1137 DataOps::incr(a_conductivityCell, this->m_fluidScratch1, 1.0 * std::abs(Z));
1138 }
1139 }
1140
1141 // Contribution from CDR solvers.
1142 DataOps::setValue(m_semiImplicitConductivityCDR, 0.0);
1143 for (auto solverIt = (this->m_cdr)->iterator(); solverIt.ok(); ++solverIt) {
1145 const RefCountedPtr<CdrSpecies>& species = solver->getSpecies();
1146
1147 const int index = solverIt.index();
1148 const int Z = species->getChargeNumber();
1149
1150 if (Z != 0 && solver->isMobile()) {
1151 const EBAMRCellData& phi = solver->getPhi();
1152 const EBAMRCellData& mu = this->m_cdrMobilities[index];
1153
1154 DataOps::copy(this->m_fluidScratch1, phi);
1155 DataOps::multiply(this->m_fluidScratch1, mu);
1156
1157 DataOps::incr(a_conductivityCell, this->m_fluidScratch1, 1.0 * std::abs(Z));
1158 }
1159 }
1160
1161 // Conductivity is mobility * weight * Q
1163
1164 // Coarsen, update ghost cells and interpolate to centroids. User can ask for bi/tri-linear filtering.
1165 (this->m_amr)->arithmeticAverage(a_conductivityCell, this->m_fluidRealm, (this->m_plasmaPhase));
1166 (this->m_amr)->interpGhostPwl(a_conductivityCell, this->m_fluidRealm, (this->m_plasmaPhase));
1167
1168 // User can choose to filter the conductivity. Need to sync after each smoothing.
1169 if (m_condFilterNum > 0 && m_condFilterMaxStride > 0) {
1170 for (int i = 0; i < m_condFilterNum; i++) {
1171 for (int curStride = 1; curStride <= m_condFilterMaxStride; curStride++) {
1172 DataOps::filterSmooth(a_conductivityCell, m_condFilterAlpha, curStride, true);
1173
1174 (this->m_amr)->arithmeticAverage(a_conductivityCell, this->m_fluidRealm, (this->m_plasmaPhase));
1175 (this->m_amr)->interpGhostPwl(a_conductivityCell, this->m_fluidRealm, (this->m_plasmaPhase));
1176 }
1177 }
1178 }
1179
1180 (this->m_amr)->interpToCentroids(a_conductivityCell, this->m_fluidRealm, (this->m_plasmaPhase));
1181}
1182
1183template <typename I, typename C, typename R, typename F>
1184void
1186{
1187 CH_TIME("ItoKMCGodunovStepper::computeFaceConductivity");
1188 if (this->m_verbosity > 5) {
1189 pout() << this->m_name + "::computeFaceConductivity" << endl;
1190 }
1191
1192 DataOps::setValue((this->m_conductivityFace), 0.0);
1193 DataOps::setValue((this->m_conductivityEB), 0.0);
1194
1195 // Average the cell-centered conductivity to faces. Note that this includes one "ghost face", which we need
1196 // because the multigrid solver will interpolate face-centered conductivities to face centroids.
1197 const Average average = Average::Arithmetic;
1198 const int tanGhost = 1;
1199 const Interval interv(0, 0);
1200
1201 DataOps::averageCellToFace((this->m_conductivityFace),
1202 (this->m_conductivityCell),
1203 (this->m_amr)->getDomains(),
1204 tanGhost,
1205 interv,
1206 interv,
1207 average);
1208
1209 // Set the EB conductivity.
1210 DataOps::incr((this->m_conductivityEB), (this->m_conductivityCell), 1.0);
1211}
1212
1213template <typename I, typename C, typename R, typename F>
1214void
1216{
1217 CH_TIMERS("ItoKMCGodunovStepper::computeSemiImplicitRho");
1218 CH_TIMER("ItoKMCGodunovStepper::computeSemiImplicitRho::plasma_phase", t1);
1219 CH_TIMER("ItoKMCGodunovStepper::computeSemiImplicitRho::solid_phase", t2);
1220 CH_TIMER("ItoKMCGodunovStepper::computeSemiImplicitRho::filter", t3);
1221 if (this->m_verbosity > 5) {
1222 pout() << this->m_name + "::computeSemiImplicitRho" << endl;
1223 }
1224
1225 // Soft requirement?
1226 CH_assert(this->m_plasmaPhase == phase::gas);
1227
1228 const RefCountedPtr<MultiFluidIndexSpace>& mfis = (this->m_computationalGeometry)->getMfIndexSpace();
1229 const Vector<Dielectric>& dielectrics = (this->m_computationalGeometry)->getDielectrics();
1230
1231 const bool hasDielectrics = (mfis->numPhases() > 1) && (dielectrics.size() > 0);
1232
1233 MFAMRCellData& rho = this->m_fieldSolver->getRho();
1234 EBAMRCellData rhoGas = (this->m_amr)->alias(phase::gas, rho);
1236
1237 if (hasDielectrics) {
1238 rhoSolid = (this->m_amr)->alias(phase::solid, rho);
1239 }
1240
1241 DataOps::setValue(rho, 0.0);
1242
1243 // Contribution from Ito solvers to the gas-side space charge density.
1244 CH_START(t1);
1245 for (auto solverIt = this->m_ito->iterator(); solverIt.ok(); ++solverIt) {
1247 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
1248 const int Z = species->getChargeNumber();
1249
1250 if (Z != 0) {
1251 (this->m_amr)->copyData(this->m_fluidScratch1, solver->getPhi());
1252
1253 DataOps::incr(rhoGas, this->m_fluidScratch1, 1.0 * Z * Units::Qe);
1254 }
1255 }
1256
1257 // Contribution from CDR equations to the gas-side space charge density.
1258 DataOps::incr(rhoGas, m_semiImplicitRhoCDR, 1.0);
1259 CH_STOP(t1);
1260
1261 // Contribution from gas-side particles that diffused into EBs. This might be necessary near dielectric EBs as a comparatively small
1262 // correction in the space charge density. There should be no contribution from the CDR solvers because the diffusion-only divergence
1263 // term was computed using homogeneous Neumann boundary conditions.
1264 if (hasDielectrics) {
1265 CH_START(t2);
1266
1269
1270 (this->m_amr)->allocate(particleScratch, this->m_particleRealm, phase::solid, 1);
1271 (this->m_amr)->allocate(fluidScratch, this->m_fluidRealm, phase::solid, 1);
1272
1273 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
1275 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
1276 const int Z = species->getChargeNumber();
1277
1278 if (Z != 0) {
1281
1282 (this->m_amr)
1285 this->m_particleRealm,
1286 phase::solid,
1287 *m_rhoDaggerParticles[solverIt.index()],
1288 DepositionType::CIC,
1289 CoarseFineDeposition::Halo,
1290 true);
1291
1292 (this->m_amr)->copyData(fluidScratch, particleScratch);
1293
1295 }
1296 }
1297
1298 CH_STOP(t2);
1299 }
1300
1301 // Sync across levels
1302 this->m_amr->arithmeticAverage(rho, this->m_fluidRealm);
1303 this->m_amr->interpGhostPwl(rho, this->m_fluidRealm);
1304
1305 // User can choose to filter the space charge density. Need to sync after each smoothing.
1306 if (m_rhoFilterNum > 0 && m_rhoFilterMaxStride > 0) {
1307 CH_START(t3);
1308 for (int i = 0; i < m_rhoFilterNum; i++) {
1309 for (int curStride = 1; curStride <= m_rhoFilterMaxStride; curStride++) {
1310
1311 DataOps::filterSmooth(rhoGas, m_rhoFilterAlpha, curStride, true);
1312
1313 this->m_amr->arithmeticAverage(rhoGas, this->m_fluidRealm, this->m_plasmaPhase);
1314 this->m_amr->interpGhost(rhoGas, this->m_fluidRealm, this->m_plasmaPhase);
1315 }
1316 }
1317 CH_STOP(t3);
1318 }
1319
1320 // Put data on centroids.
1321 this->m_amr->interpToCentroids(rhoGas, this->m_fluidRealm, phase::gas);
1322 if (hasDielectrics) {
1323 this->m_amr->interpToCentroids(rhoSolid, this->m_fluidRealm, phase::solid);
1324 }
1325}
1326
1327template <typename I, typename C, typename R, typename F>
1328void
1330{
1331 CH_TIME("ItoKMCGodunovStepper::setupSemiImplicitPoisson");
1332 if (this->m_verbosity > 5) {
1333 pout() << this->m_name + "::setupSemiImplicitPoisson" << endl;
1334 }
1335
1336 // Set coefficients as usual
1337 (this->m_fieldSolver)->setPermittivities();
1338
1339 // Get the permittivities
1340 MFAMRCellData& permCell = (this->m_fieldSolver)->getPermittivityCell();
1341 MFAMRFluxData& permFace = (this->m_fieldSolver)->getPermittivityFace();
1342 MFAMRIVData& permEB = (this->m_fieldSolver)->getPermittivityEB();
1343
1344 // Get handles to the gas-phase permittivities.
1345 EBAMRFluxData permFaceGas = (this->m_amr)->alias((this->m_plasmaPhase), permFace);
1346 EBAMRIVData permEBGas = (this->m_amr)->alias((this->m_plasmaPhase), permEB);
1347
1348 // Increment the field solver permittivities by a_factor*sigma. After this, the "permittivities" are
1349 // given by epsr + a_factor*sigma
1350 DataOps::incr(permFaceGas, (this->m_conductivityFace), a_dt / Units::eps0);
1351 DataOps::incr(permEBGas, (this->m_conductivityEB), a_dt / Units::eps0);
1352
1353 // Coarsen coefficients.
1354 (this->m_amr)->arithmeticAverage(permFaceGas, this->m_fluidRealm, (this->m_plasmaPhase));
1355 (this->m_amr)->arithmeticAverage(permEBGas, this->m_fluidRealm, (this->m_plasmaPhase));
1356
1357 // Set up the solver with the "permittivities"
1358 (this->m_fieldSolver)->setSolverPermittivities(permCell, permFace, permEB);
1359}
1360
1361template <typename I, typename C, typename R, typename F>
1362void
1366 const Real a_tolerance) const noexcept
1367{
1368 CH_TIME("ItoKMCGodunovStepper::removeCoveredPointParticles");
1369 if (this->m_verbosity > 5) {
1370 pout() << this->m_name + "::removeCoveredPointParticles" << endl;
1371 }
1372
1373 for (int i = 0; i < a_particles.size(); i++) {
1374 if (a_particles[i] != nullptr) {
1376
1377 switch (a_representation) {
1378 case EBRepresentation::Discrete: {
1379 (this->m_amr)->removeCoveredParticlesDiscrete(particles, (this->m_plasmaPhase), a_tolerance);
1380
1381 break;
1382 }
1383 case EBRepresentation::ImplicitFunction: {
1384 (this->m_amr)->removeCoveredParticlesIF(particles, (this->m_plasmaPhase), a_tolerance);
1385
1386 break;
1387 }
1388 case EBRepresentation::Voxel: {
1389 (this->m_amr)->removeCoveredParticlesVoxels(particles, (this->m_plasmaPhase));
1390
1391 break;
1392 }
1393 default: {
1394 MayDay::Error("ItoKMCGodunovStepper::removeCoveredParticles - logic bust");
1395 }
1396 }
1397 }
1398 }
1399}
1400
1401template <typename I, typename C, typename R, typename F>
1402void
1405{
1406 CH_TIME("ItoKMCGodunovStepper::copyConductivityParticles");
1407 if (this->m_verbosity > 5) {
1408 pout() << this->m_name + "::copyConductivityParticles" << endl;
1409 }
1410
1411 // Clear particles first.
1412 this->clearPointParticles(a_conductivityParticles, SpeciesSubset::All);
1413
1414 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
1416 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
1417
1418 const int idx = solverIt.index();
1419 const int Z = species->getChargeNumber();
1420
1421 if (Z != 0 && solver->isMobile()) {
1422 const ParticleContainer<ItoParticle>& solverParticles = solver->getParticles(ItoSolver::WhichContainer::Bulk);
1423
1424 for (int lvl = 0; lvl <= (this->m_amr)->getFinestLevel(); lvl++) {
1425 const DisjointBoxLayout& dbl = (this->m_amr)->getGrids((this->m_particleRealm))[lvl];
1426 const DataIterator& dit = dbl.dataIterator();
1427
1428 const int nbox = dit.size();
1429
1430#pragma omp parallel for schedule(runtime)
1431 for (int mybox = 0; mybox < nbox; mybox++) {
1432 const DataIndex& din = dit[mybox];
1433
1435
1437 List<PointParticle>& irregParticles = (*m_irregularParticles[idx])[lvl][din].listItems();
1438
1440 const ItoParticle& p = lit();
1441 const RealVect& pos = p.position();
1442 const Real& weight = p.weight();
1443 const Real& mobility = p.mobility();
1444
1445 pointParticles.add(PointParticle(pos, weight * mobility));
1446 }
1447
1449 }
1450 }
1451 }
1452 }
1453}
1454
1455template <typename I, typename C, typename R, typename F>
1456bool
1458{
1459 CH_TIME("ItoKMCGodunovStepper::solvePoisson()");
1460 if (this->m_verbosity > 5) {
1461 pout() << this->m_name + "::solvePoisson()" << endl;
1462 }
1463
1464 // Solve the Poisson equation and compute the cell-centered electric field.
1465 MFAMRCellData& phi = this->m_fieldSolver->getPotential();
1466 MFAMRCellData& rho = this->m_fieldSolver->getRho();
1467 EBAMRIVData& sigma = this->m_sigmaSolver->getPhi();
1468
1469 const bool converged = (this->m_fieldSolver)->solve(phi, rho, sigma, false);
1470
1471 (this->m_fieldSolver)->computeElectricField();
1472
1473 // Copy the electric field to appropriate data holders and perform center-to-centroid
1474 // interpolation.
1476 (this->m_amr)->allocatePointer(E, this->m_fluidRealm);
1477 (this->m_amr)->alias(E, this->m_plasmaPhase, (this->m_fieldSolver)->getElectricField());
1478
1479 // Fluid realm
1480 (this->m_amr)->copyData(this->m_electricFieldFluid, E);
1481 (this->m_amr)->conservativeAverage(this->m_electricFieldFluid, this->m_fluidRealm, this->m_plasmaPhase);
1482 (this->m_amr)->interpGhostPwl(this->m_electricFieldFluid, this->m_fluidRealm, this->m_plasmaPhase);
1483 (this->m_amr)->interpToCentroids(this->m_electricFieldFluid, this->m_fluidRealm, this->m_plasmaPhase);
1484
1485 // Particle realm
1486 (this->m_amr)->copyData(this->m_electricFieldParticle, E);
1487 (this->m_amr)->conservativeAverage(this->m_electricFieldParticle, this->m_particleRealm, this->m_plasmaPhase);
1488 (this->m_amr)->interpGhostPwl(this->m_electricFieldParticle, this->m_particleRealm, this->m_plasmaPhase);
1489 (this->m_amr)->interpToCentroids(this->m_electricFieldParticle, this->m_particleRealm, this->m_plasmaPhase);
1490
1491 return converged;
1492}
1493
1494template <typename I, typename C, typename R, typename F>
1495void
1497{
1498 CH_TIME("ItoKMCGodunovStepper::advanceEulerMaruyama");
1499 if (this->m_verbosity > 5) {
1500 pout() << this->m_name + "::advanceEulerMaruyama" << endl;
1501 }
1502
1503 // Store X^k positions.
1504 this->setOldPositions();
1505
1506 // Diffuse the particles. This copies onto m_rhoDaggerParticles and stores the hop on the full particles. We need
1507 // to remap the particles species that made a diffusion hop.
1508 this->barrier();
1509 m_timer.startEvent("Diffuse particles");
1510 this->diffuseParticlesEulerMaruyama(m_rhoDaggerParticles, a_dt);
1511 this->remapPointParticles(m_rhoDaggerParticles, SpeciesSubset::ChargedDiffusive);
1512 m_timer.stopEvent("Diffuse particles");
1513
1514 // Perform the diffusive CDR advance.
1515 this->barrier();
1516 m_timer.startEvent("Diffuse CDR");
1517 this->computeDiffusionTermCDR(m_semiImplicitRhoCDR, a_dt);
1518 m_timer.stopEvent("Diffuse CDR");
1519
1520 // Compute the conductivity on the mesh. This deposits q_e * Z * w * mu on the mesh.
1521 this->barrier();
1522 m_timer.startEvent("Compute conductivities");
1523 this->copyConductivityParticles(m_conductivityParticles);
1524 this->computeConductivities(m_conductivityParticles);
1525 m_timer.stopEvent("Compute conductivities");
1526
1527 // Set up the semi-implicit Poisson solver with the computed conductivities.
1528 this->barrier();
1529 m_timer.startEvent("Setup Poisson");
1530 this->setupSemiImplicitPoisson(a_dt);
1531 m_timer.stopEvent("Setup Poisson");
1532
1533 // Compute space charge density arising from the new particle positions X^k + sqrt(2*D*dt)*W. Only need to
1534 // do the diffusive and charged species.
1535 this->barrier();
1536 m_timer.startEvent("Deposit point particles");
1537 this->depositPointParticles(m_rhoDaggerParticles, SpeciesSubset::Charged);
1538 this->computeSemiImplicitRho();
1539 m_timer.stopEvent("Deposit point particles");
1540
1541 // Solve the semi-implicit Poisson equation.
1542 this->barrier();
1543 m_timer.startEvent("Solve Poisson");
1544 const bool converged = this->solvePoisson();
1545 if (!converged) {
1546 const std::string errMsg = "ItoKMCGodunovStepper::advanceEulerMaruyama - Poisson solve did not converge";
1547
1548 pout() << errMsg << endl;
1549
1550 if (this->m_abortOnFailure) {
1551 MayDay::Error(errMsg.c_str());
1552 }
1553 }
1554 m_timer.stopEvent("Solve Poisson");
1555
1556 // Recompute velocities with the new electric field. This interpolates the velocities to the current particle positions, i.e.
1557 // we compute V^(k+1)(X^k) = mu^k * E^(k+1)(X^k)
1558 this->barrier();
1559 m_timer.startEvent("Step-compute v");
1560#if 1 // This is what the algorithm says.
1561 this->setCdrVelocityFunctions();
1562 this->setItoVelocityFunctions();
1563 (this->m_ito)->interpolateVelocities();
1564 this->multiplyCdrVelocitiesByMobilities();
1565#else // Have to use this for LEA - need to debug.
1566 this->computeDriftVelocities();
1567#endif
1568 m_timer.stopEvent("Step-compute v");
1569
1570 // Finalize the Euler-Maruyama update.
1571 this->barrier();
1572 m_timer.startEvent("Euler-Maruyama step");
1573 this->stepEulerMaruyamaParticles(a_dt);
1574 this->remapParticles(SpeciesSubset::AllMobileOrDiffusive);
1575 this->stepEulerMaruyamaCDR(a_dt);
1576 m_timer.stopEvent("Euler-Maruyama step");
1577}
1578
1579template <typename I, typename C, typename R, typename F>
1580void
1583 const Real a_dt) noexcept
1584{
1585 CH_TIME("ItoKMCGodunovStepper::diffuseParticlesEulerMaruyama");
1586 if (this->m_verbosity > 5) {
1587 pout() << this->m_name + "::diffuseParticlesEulerMaruyama" << endl;
1588 }
1589
1590 this->clearPointParticles(a_rhoDaggerParticles, SpeciesSubset::All);
1591
1592 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
1594 const RefCountedPtr<ItoSpecies>& species = solver->getSpecies();
1595
1596 const int idx = solverIt.index();
1597
1598 const bool mobile = solver->isMobile();
1599 const bool diffusive = solver->isDiffusive();
1600 const int Z = species->getChargeNumber();
1601
1602 const auto& diffusionFunction = (this->m_physics)->getItoDiffusionFunctions()[idx];
1603
1604 for (int lvl = 0; lvl <= (this->m_amr)->getFinestLevel(); lvl++) {
1605 const DisjointBoxLayout& dbl = (this->m_amr)->getGrids((this->m_particleRealm))[lvl];
1606 const DataIterator& dit = dbl.dataIterator();
1607
1608 ParticleData<ItoParticle>& particles = solver->getParticles(ItoSolver::WhichContainer::Bulk)[lvl];
1609
1610 const int nbox = dit.size();
1611
1612#pragma omp parallel for schedule(runtime)
1613 for (int mybox = 0; mybox < nbox; mybox++) {
1614 const DataIndex& din = dit[mybox];
1615
1618
1620 ItoParticle& p = lit();
1621 const Real& weight = p.weight();
1622 const RealVect& pos = p.position();
1623
1624 // Compute a particle hop and store it on the run-time storage.
1625 RealVect& hop = p.tmpVect();
1626
1628
1629 if (Z != 0) {
1630 pointParticles.add(PointParticle(pos + hop, weight));
1631 }
1632 }
1633 }
1634 }
1635 }
1636}
1637
1638template <typename I, typename C, typename R, typename F>
1639void
1641{
1642 CH_TIME("ItoKMCGodunovStepper::diffuseCDREulerMaruyama");
1643 if (this->m_verbosity > 5) {
1644 pout() << this->m_name + "::diffuseCDREulerMaruyama" << endl;
1645 }
1646
1648
1649 for (auto solverIt = this->m_cdr->iterator(); solverIt.ok(); ++solverIt) {
1651 const RefCountedPtr<CdrSpecies>& species = solver->getSpecies();
1652
1653 const int index = solverIt.index();
1654 const int Z = species->getChargeNumber();
1655
1656 const EBAMRCellData& phi = solver->getPhi();
1657
1658 // Compute finite-volume approximation to div(D*grad(phi))
1659 if (solver->isDiffusive()) {
1660 solver->computeDivD(m_cdrDivD[index], solver->getPhi(), false, false, false);
1661 }
1662
1663 // Update the space charge density arising from the update phi^dagger = phi^k + dt * div(D*grad(phi^k))
1664 if (Z != 0) {
1666 if (solver->isDiffusive()) {
1667 DataOps::incr(a_semiImplicitRhoCDR, m_cdrDivD[index], 1.0 * Z * a_dt);
1668 }
1669 }
1670 }
1671
1673
1674 this->m_amr->arithmeticAverage(a_semiImplicitRhoCDR, this->m_fluidRealm, this->m_plasmaPhase);
1675 this->m_amr->interpGhostPwl(a_semiImplicitRhoCDR, this->m_fluidRealm, this->m_plasmaPhase);
1676}
1677
1678template <typename I, typename C, typename R, typename F>
1679void
1681{
1682 CH_TIME("ItoKMCGodunovStepper::stepEulerMaruyamaParticles");
1683 if (this->m_verbosity > 5) {
1684 pout() << this->m_name + "::stepEulerMaruyamaParticles" << endl;
1685 }
1686
1687 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
1689
1690 const bool mobile = solver->isMobile();
1691 const bool diffusive = solver->isDiffusive();
1692
1693 const Real f = mobile ? a_dt : 0.0;
1694 const Real g = diffusive ? 1.0 : 0.0;
1695
1696 if (mobile || diffusive) {
1697 for (int lvl = 0; lvl <= (this->m_amr)->getFinestLevel(); lvl++) {
1698 const DisjointBoxLayout& dbl = (this->m_amr)->getGrids((this->m_particleRealm))[lvl];
1699 const DataIterator& dit = dbl.dataIterator();
1700
1701 ParticleData<ItoParticle>& particles = solver->getParticles(ItoSolver::WhichContainer::Bulk)[lvl];
1702
1703 const int nbox = dit.size();
1704
1705#pragma omp parallel for schedule(runtime)
1706 for (int mybox = 0; mybox < nbox; mybox++) {
1707 const DataIndex& din = dit[mybox];
1708
1710
1712 ItoParticle& p = lit();
1713
1714 // Add in the diffusion hop and advective contribution.
1715 const RealVect& hop = p.tmpVect();
1716 p.position() = p.oldPosition() + f * p.velocity() + g * hop;
1717 }
1718 }
1719 }
1720 }
1721 }
1722}
1723
1724template <typename I, typename C, typename R, typename F>
1725void
1727{
1728 CH_TIME("ItoKMCGodunovStepper::stepEulerMaruyamaCDR");
1729 if (this->m_verbosity > 5) {
1730 pout() << this->m_name + "::stepEulerMaruyamaCDR" << endl;
1731 }
1732
1733 for (auto solverIt = (this->m_cdr)->iterator(); solverIt.ok(); ++solverIt) {
1735
1736 const int index = solverIt.index();
1737
1738 EBAMRCellData& phi = solver->getPhi();
1739
1740 this->m_amr->conservativeAverage(phi, this->m_fluidRealm, this->m_plasmaPhase);
1741 this->m_amr->interpGhostPwl(phi, this->m_fluidRealm, this->m_plasmaPhase);
1742
1743 // Add in the advective term -- BC comes later.
1744 if (solver->isMobile()) {
1745 DataOps::setValue(solver->getEbFlux(), 0.0);
1746
1747 // Compute the FV approximation to the advective term and do the Euler advance. If the underlying
1748 // CDR solver is a CTU solver, we also add in the transverse terms.
1749 solver->computeDivF(this->m_fluidScratch1, phi, 0.5 * a_dt, false, true, true);
1750
1751 DataOps::incr(phi, this->m_fluidScratch1, -a_dt);
1752 }
1753
1754 // Add in the diffusion term.
1755 if (solver->isDiffusive()) {
1756 DataOps::incr(phi, this->m_cdrDivD[index], a_dt);
1757 }
1758
1759 DataOps::floor(phi, 0.0);
1760 }
1761
1762 this->coarsenCDRSolvers();
1763}
1764
1765#ifdef CH_USE_HDF5
1766template <typename I, typename C, typename R, typename F>
1767void
1769{
1770 CH_TIME("ItoKMCGodunovStepper::writeCheckpointHeader");
1771 if (this->m_verbosity > 5) {
1772 pout() << this->m_name + "::writeCheckpointHeader" << endl;
1773 }
1774
1775 a_header.m_real["prev_dt"] = this->m_prevDt;
1776 a_header.m_real["physics_dt"] = this->m_physicsDt;
1777 a_header.m_int["checkpoint_particles"] = m_writeCheckpointParticles ? 1 : 0;
1778}
1779#endif
1780
1781#ifdef CH_USE_HDF5
1782template <typename I, typename C, typename R, typename F>
1783void
1785{
1786 CH_TIME("ItoKMCGodunovStepper::readCheckpointHeader");
1787 if (this->m_verbosity > 5) {
1788 pout() << this->m_name + "::readCheckpointHeader" << endl;
1789 }
1790
1791 this->m_prevDt = a_header.m_real["prev_dt"];
1792 this->m_physicsDt = a_header.m_real["physics_dt"];
1793
1794 m_readCheckpointParticles = (a_header.m_int["checkpoint_particles"] != 0) ? true : false;
1795 m_canRegridOnRestart = m_readCheckpointParticles;
1796}
1797#endif
1798
1799#ifdef CH_USE_HDF5
1800template <typename I, typename C, typename R, typename F>
1801void
1803{
1804 CH_TIME("ItoKMCGodunovStepper::writeCheckpointData");
1805 if (this->m_verbosity > 5) {
1806 pout() << this->m_name + "::writeCheckpointData" << endl;
1807 }
1808
1810
1811 // Write the point-particles.
1812 if (m_writeCheckpointParticles) {
1813 for (int i = 0; i < (this->m_physics)->getNumItoSpecies(); i++) {
1814 const std::string identifierSigma = "ItoKMCGodunovStepper::conductivityParticles_" + std::to_string(i);
1815 const std::string identifierRho = "ItoKMCGodunovStepper::spaceChargeParticles_" + std::to_string(i);
1816
1817 const ParticleContainer<PointParticle>& conductivityParticles = *m_conductivityParticles[i];
1818 const ParticleContainer<PointParticle>& rhoDaggerParticles = *m_rhoDaggerParticles[i];
1819
1820 DischargeIO::writeCheckParticlesToHDF(a_handle, conductivityParticles[a_lvl], identifierSigma);
1821 DischargeIO::writeCheckParticlesToHDF(a_handle, rhoDaggerParticles[a_lvl], identifierRho);
1822 }
1823 }
1824
1825 // Write data required for the semi-implicit regrid.
1826 if (this->m_physics->getNumCdrSpecies() > 0) {
1827 write(a_handle, *m_semiImplicitRhoCDR[a_lvl], "ItoKMCGodunovStepper::semiImplicitRhoCDR");
1828 write(a_handle, *m_semiImplicitConductivityCDR[a_lvl], "ItoKMCGodunovStepper::semiImplicitConductivityCDR");
1829 }
1830}
1831#endif
1832
1833#ifdef CH_USE_HDF5
1834template <typename I, typename C, typename R, typename F>
1835void
1837{
1838 CH_TIME("ItoKMCGodunovStepper::readCheckpointData");
1839 if (this->m_verbosity > 5) {
1840 pout() << this->m_name + "::readCheckpointData" << endl;
1841 }
1842
1844
1845 // Write the point-particles.
1846 if (m_readCheckpointParticles) {
1847 for (int i = 0; i < (this->m_physics)->getNumItoSpecies(); i++) {
1848 const std::string identifierSigma = "ItoKMCGodunovStepper::conductivityParticles_" + std::to_string(i);
1849 const std::string identifierRho = "ItoKMCGodunovStepper::spaceChargeParticles_" + std::to_string(i);
1850
1851 ParticleContainer<PointParticle>& conductivityParticles = *m_conductivityParticles[i];
1852 ParticleContainer<PointParticle>& rhoDaggerParticles = *m_rhoDaggerParticles[i];
1853
1854 DischargeIO::readCheckParticlesFromHDF(a_handle, conductivityParticles[a_lvl], identifierSigma);
1855 DischargeIO::readCheckParticlesFromHDF(a_handle, rhoDaggerParticles[a_lvl], identifierRho);
1856 }
1857 }
1858
1859 // Read in data that is required for the semi-implicit regrid.
1860 if (this->m_physics->getNumCdrSpecies() > 0) {
1861 const Interval interv(0, 0);
1862
1864 *m_semiImplicitRhoCDR[a_lvl],
1865 "ItoKMCGodunovStepper::semiImplicitRhoCDR",
1866 this->m_amr->getGrids(this->m_fluidRealm)[a_lvl],
1867 interv,
1868 false);
1869
1871 *m_semiImplicitConductivityCDR[a_lvl],
1872 "ItoKMCGodunovStepper::semiImplicitConductivityCDR",
1873 this->m_amr->getGrids(this->m_fluidRealm)[a_lvl],
1874 interv,
1875 false);
1876 }
1877}
1878#endif
1879
1880template <typename I, typename C, typename R, typename F>
1881void
1883{
1884 CH_TIME("ItoKMCGodunovStepper::postPlot");
1885 if (this->m_verbosity > 5) {
1886 pout() << this->m_name + "::postPlot" << endl;
1887 }
1888
1889 this->m_physicsPlotVariables.clear();
1890
1891 this->plotParticles();
1892}
1893
1894template <typename I, typename C, typename R, typename F>
1895void
1897{
1898 CH_TIME("ItoKMCGodunovStepper::plotParticles");
1899 if (this->m_verbosity > 2) {
1900 pout() << this->m_name + "::plotParticles" << endl;
1901 }
1902
1903 bool plotParticles = false;
1904
1905 ParmParse pp(this->m_name.c_str());
1906
1907 pp.query("plot_particles", plotParticles);
1908
1909 if (plotParticles) {
1910
1911 for (auto solverIt = (this->m_ito)->iterator(); solverIt.ok(); ++solverIt) {
1913 const ParticleContainer<ItoParticle>& particles = solver->getParticles(ItoSolver::WhichContainer::Bulk);
1914
1915 // Create the output folder
1916 std::string cmd = "mkdir -p particles/" + solver->getName();
1917 int success = 0;
1918 if (procID() == 0) {
1919 success = system(cmd.c_str());
1920 }
1921
1922 if (success != 0) {
1923 MayDay::Error("ItoKMCGodunovStepper::plotParticles - could not create 'particles' directory");
1924 }
1925
1926 // Set plot file name
1927 const std::string prefix = "./particles/" + solver->getName() + "/" + solver->getName();
1928 char fileChar[1000];
1929 sprintf(fileChar, "%s.step%07d.%dd.h5part", prefix.c_str(), this->m_timeStep, SpaceDim);
1930
1931 // Plot the particles
1936 this->m_amr->getProbLo(),
1937 this->m_time);
1938 }
1939 }
1940}
1941
1942#include <CD_NamespaceFooter.H>
1943
1944#endif
Average
Various averaging methods.
Definition CD_Average.H:24
Agglomeration of useful data operations.
Silly, but useful functions that override standard Chombo HDF5 IO.
EBRepresentation
Enum for putting some logic into how we think about EBs. This is just a simply supporting class for v...
Definition CD_EBRepresentation.H:22
Declaration of a class which uses a semi-implicit Godunov method for Ito plasma equations.
SpeciesSubset
Enum for differentiating between types of particles.
Definition CD_ItoKMCStepper.H:40
Agglomeration of basic MPI reductions.
Declaration of a photon class for particle methods.
Implementation of CD_Timer.H.
Declaration of various useful units.
static void scale(MFAMRCellData &a_lhs, const Real &a_scale) noexcept
Scale data by factor.
Definition CD_DataOps.cpp:2398
static void floor(EBAMRCellData &a_lhs, const Real a_value)
Floor values in data holder. This sets all values below a_value to a_value.
Definition CD_DataOps.cpp:1394
static void filterSmooth(EBAMRCellData &a_data, const Real a_alpha, const int a_stride, const bool a_zeroEB) noexcept
Apply a convolved filter phi = alpha * phi_i + 0.5*(1-alpha) * [phi_(i+s) + phi_(i-s)] in each direct...
Definition CD_DataOps.cpp:660
static void multiply(EBAMRCellData &a_lhs, const EBAMRCellData &a_rhs)
Multiply data holder by another data holder.
Definition CD_DataOps.cpp:2156
static void incr(MFAMRCellData &a_lhs, const MFAMRCellData &a_rhs, const Real a_scale) noexcept
Function which increments data in the form a_lhs = a_lhs + a_rhs*a_scale for all components.
Definition CD_DataOps.cpp:801
static void setValue(LevelData< MFInterfaceFAB< T > > &a_lhs, const T &a_value)
Set value in an MFInterfaceFAB data holder.
Definition CD_DataOpsImplem.H:23
static void averageCellToFace(EBAMRFluxData &a_faceData, const EBAMRCellData &a_cellData, const Vector< ProblemDomain > &a_domains)
Average all components of the cell-centered data to faces.
Definition CD_DataOps.cpp:153
static void copy(MFAMRCellData &a_dst, const MFAMRCellData &a_src)
Copy data from one data holder to another.
Definition CD_DataOps.cpp:1132
RealVect & position()
Get the particle position.
Definition CD_GenericParticleImplem.H:49
A particle class for use with ItoSolvers, i.e. drifting Brownian walkers.
Definition CD_ItoParticle.H:40
Real & mobility()
Get mobility coefficient.
Definition CD_ItoParticleImplem.H:83
RealVect & oldPosition()
Get the old particle position.
Definition CD_ItoParticleImplem.H:119
RealVect & tmpVect()
Return scratch RealVect storage.
Definition CD_ItoParticleImplem.H:173
static std::vector< std::string > s_vectVariables
Naming convention for vector fields.
Definition CD_ItoParticle.H:17
Real & weight()
Get particle weight.
Definition CD_ItoParticleImplem.H:71
static std::vector< std::string > s_realVariables
Naming convention for scalar fields.
Definition CD_ItoParticle.H:16
RealVect & velocity()
Get the particle velocity.
Definition CD_ItoParticleImplem.H:131
Real & tmpReal()
Return scratch scalar storage.
Definition CD_ItoParticleImplem.H:161
Templated class for holding particles on an AMR hierarchy with particle remapping.
Definition CD_ParticleContainer.H:51
Particle class for usage with Monte Carlo radiative transfer.
Definition CD_Photon.H:29
Real & weight()
Get photon weight.
Definition CD_PhotonImplem.H:40
Implementation of ItoKMCStepper that uses a semi-implicit split-step formalism for advancing the Ito-...
Definition CD_ItoKMCGodunovStepper.H:29
virtual void setupSemiImplicitPoisson(const Real a_dt) noexcept
Set up the semi-implicit Poisson solver.
Definition CD_ItoKMCGodunovStepperImplem.H:1329
virtual void regrid(const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) noexcept override
Regrid methods – puts all data on the new mesh.
Definition CD_ItoKMCGodunovStepperImplem.H:592
virtual void removeCoveredPointParticles(Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_particles, const EBRepresentation a_representation, const Real a_tolerance) const noexcept
Remove covered particles.
Definition CD_ItoKMCGodunovStepperImplem.H:1363
virtual void computeDiffusionTermCDR(EBAMRCellData &m_semiImplicitRhoCDR, const Real a_dt) noexcept
Compute the diffusion term for the CDR equations as well as the resulting CDR-contributions to the sp...
Definition CD_ItoKMCGodunovStepperImplem.H:1640
virtual void allocate() noexcept override
Allocate storage required for advancing the equations.
Definition CD_ItoKMCGodunovStepperImplem.H:76
bool m_readCheckpointParticles
If true, then the HDF5 checkpoint file contained particles that we can read.
Definition CD_ItoKMCGodunovStepper.H:166
virtual Real advance(const Real a_dt) override
Advance the Ito-Poisson-KMC system over a_dt.
Definition CD_ItoKMCGodunovStepperImplem.H:306
virtual void depositPointParticles(const Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_particles, const SpeciesSubset a_subset) noexcept
Deposit the input point particles on the mesh.
Definition CD_ItoKMCGodunovStepperImplem.H:878
virtual void allocateInternals() noexcept override
Allocate "internal" storage.
Definition CD_ItoKMCGodunovStepperImplem.H:110
Real m_maxFieldAbort
Limit for maximum field abort.
Definition CD_ItoKMCGodunovStepper.H:227
virtual void computeConductivities(const Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_particles) noexcept
Compute all conductivities (cell, face, and EB) from the input point particles.
Definition CD_ItoKMCGodunovStepperImplem.H:1096
virtual void stepEulerMaruyamaCDR(const Real a_dt) noexcept
Step the CDR equations according to the regular Euler-Maruyama scheme.
Definition CD_ItoKMCGodunovStepperImplem.H:1726
virtual void diffuseParticlesEulerMaruyama(Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_rhoDaggerParticles, const Real a_dt) noexcept
Perform the diffusive Ito advance in the Euler-Maruyama step.
Definition CD_ItoKMCGodunovStepperImplem.H:1581
virtual void parseAlgorithm() noexcept
Parse advancement algorithm.
Definition CD_ItoKMCGodunovStepperImplem.H:180
virtual void postPlot() noexcept override
Perform post-plot operations.
Definition CD_ItoKMCGodunovStepperImplem.H:1882
virtual Real computeDt() override
Compute a time step used for the advance method.
Definition CD_ItoKMCGodunovStepperImplem.H:286
virtual void parseSecondaryEmissionSpecification() noexcept
Parse when secondary particles are emitted.
Definition CD_ItoKMCGodunovStepperImplem.H:255
virtual void parseFiltering() noexcept
Parse filter settings.
Definition CD_ItoKMCGodunovStepperImplem.H:205
virtual void parseRuntimeOptions() noexcept override
Parse run-time options.
Definition CD_ItoKMCGodunovStepperImplem.H:163
virtual bool solvePoisson() noexcept override
Solve the electrostatic problem.
Definition CD_ItoKMCGodunovStepperImplem.H:1457
bool m_canRegridOnRestart
If true, then the class supports regrid-on-restart.
Definition CD_ItoKMCGodunovStepper.H:172
virtual void computeFaceConductivity() noexcept
Compute the cell-centered conductivity.
Definition CD_ItoKMCGodunovStepperImplem.H:1185
virtual ~ItoKMCGodunovStepper()
Destructor. Does nothing.
Definition CD_ItoKMCGodunovStepperImplem.H:51
bool m_extendConductivityEB
For achieving a slightly smoother gradient in the conductivity near the EB.
Definition CD_ItoKMCGodunovStepper.H:177
virtual void setOldPositions() noexcept
Set the starting positions for the ItoSolver particles.
Definition CD_ItoKMCGodunovStepperImplem.H:736
bool m_writeCheckpointParticles
If true, then the particles are checkpointed so we can regrid on checkpoint-restart.
Definition CD_ItoKMCGodunovStepper.H:161
virtual void remapPointParticles(Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_particles, const SpeciesSubset a_subset) noexcept
Remap the input point particles.
Definition CD_ItoKMCGodunovStepperImplem.H:772
virtual void computeSemiImplicitRho() noexcept
Set up the space charge density for the regrid operation.
Definition CD_ItoKMCGodunovStepperImplem.H:1215
virtual void advanceEulerMaruyama(const Real a_dt) noexcept
Advance the particles using the Euler-Maruyama scheme.
Definition CD_ItoKMCGodunovStepperImplem.H:1496
virtual void preRegrid(const int a_lmin, const int a_oldFinestLevel) noexcept override
Perform pre-regrid operations.
Definition CD_ItoKMCGodunovStepperImplem.H:553
virtual void registerOperators() noexcept override
Register operators used for the simulation.
Definition CD_ItoKMCGodunovStepperImplem.H:61
virtual void plotParticles() const noexcept
Utility function for plotting the ItoSolver particles. These are written in a particles folder.
Definition CD_ItoKMCGodunovStepperImplem.H:1896
ItoKMCGodunovStepper()=delete
Disallowed default constructor. Use the full constructor.
virtual void computeCellConductivity(EBAMRCellData &a_conductivityCell, const Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_particles) noexcept
Compute the cell-centered conductivity.
Definition CD_ItoKMCGodunovStepperImplem.H:1110
virtual void copyConductivityParticles(Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_conductivityParticles) noexcept
Copy particles from the ItoSolver into PointParticles whose weight are ItoParticle::m_weight * ItoPar...
Definition CD_ItoKMCGodunovStepperImplem.H:1403
virtual void clearPointParticles(const Vector< RefCountedPtr< ParticleContainer< PointParticle > > > &a_particles, const SpeciesSubset a_subset) noexcept
Clear the input particle data holders.
Definition CD_ItoKMCGodunovStepperImplem.H:990
virtual void parseCheckpointParticles() noexcept
Parse checkpoint-restart functionality.
Definition CD_ItoKMCGodunovStepperImplem.H:241
virtual void barrier() const noexcept
Set an MPI barrier if using debug mode.
Definition CD_ItoKMCGodunovStepperImplem.H:132
virtual void parseOptions() noexcept override
Parse options.
Definition CD_ItoKMCGodunovStepperImplem.H:146
virtual void stepEulerMaruyamaParticles(const Real a_dt) noexcept
Step the particles according to the regular Euler-Maruyama scheme.
Definition CD_ItoKMCGodunovStepperImplem.H:1680
Base time stepper class that advances the Ito-KMC-Poisson system of equations. If you want a differen...
Definition CD_ItoKMCStepper.H:60
virtual void parseRuntimeOptions() noexcept override
Parse runtime configurable options.
Definition CD_ItoKMCStepperImplem.H:109
std::string m_name
Time stepper name.
Definition CD_ItoKMCStepper.H:367
virtual void registerOperators() noexcept override
Register operators used for the simulation.
Definition CD_ItoKMCStepperImplem.H:1550
virtual void parseOptions() noexcept
Parse options.
Definition CD_ItoKMCStepperImplem.H:89
virtual void preRegrid(const int a_lmin, const int a_oldFinestLevel) noexcept override
Perform pre-regrid operations - storing relevant data from the old grids.
Definition CD_ItoKMCStepperImplem.H:1599
Real m_prevDt
Previous time step.
Definition CD_ItoKMCStepper.H:506
virtual Real computeDt() override
Compute a time step used for the advance method.
Definition CD_ItoKMCStepperImplem.H:1421
virtual void allocateInternals() noexcept
Allocate "internal" storage.
Definition CD_ItoKMCStepperImplem.H:535
virtual void allocate() noexcept override
Allocate storage for solvers.
Definition CD_ItoKMCStepperImplem.H:517
A particle class that only has a position and a weight.
Definition CD_PointParticle.H:29
Real & weight()
Get weight.
Definition CD_PointParticleImplem.H:38
Class which is used for run-time monitoring of events.
Definition CD_Timer.H:31
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
virtual void remap()
Remap particles.
Definition CD_TracerParticleSolverImplem.H:338
void depositParticles(EBAMRCellData &a_phi, const ParticleContainer< T > &a_particles, const DepositionType a_baseDeposition, const CoarseFineDeposition a_coarseFineDeposition) const noexcept
Generic particle deposition method for putting a scalar field onto the mesh.
Definition CD_TracerParticleSolverImplem.H:735
TracerParticleSolver()
Default constructor.
Definition CD_TracerParticleSolverImplem.H:25
virtual void preRegrid(const int a_lbase, const int a_oldFinestLevel)
Perform pre-regrid operations.
Definition CD_TracerParticleSolverImplem.H:306
virtual ParticleContainer< P > & getParticles()
Get all particles.
Definition CD_TracerParticleSolverImplem.H:662
virtual void allocate()
Allocate storage for this solver.
Definition CD_TracerParticleSolverImplem.H:194
int m_timeStep
Time step.
Definition CD_TracerParticleSolver.H:375
Real m_time
Time.
Definition CD_TracerParticleSolver.H:370
void writeH5Part(const std::string a_filename, const ParticleContainer< GenericParticle< M, N > > &a_particles, const std::vector< std::string > a_realVars=std::vector< std::string >(), const std::vector< std::string > a_vectVars=std::vector< std::string >(), const RealVect a_shift=RealVect::Zero, const Real a_time=0.0) noexcept
A shameless copy of Chombo's writeEBHDF5 but including the lower-left corner of the physical domain a...
Definition CD_DischargeIOImplem.H:29
void barrier() noexcept
MPI barrier.
Definition CD_ParallelOpsImplem.H:25
constexpr Real eps0
Permittivity of free space.
Definition CD_Units.H:29
constexpr Real Qe
Elementary charge.
Definition CD_Units.H:34