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