chombo-discharge
CD_ItoKMCStepper.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_ItoKMCStepper_H
13 #define CD_ItoKMCStepper_H
14 
15 // Std includes
16 #include <functional>
17 
18 // Our includes
19 #include <CD_TimeStepper.H>
20 #include <CD_ItoKMCPhysics.H>
21 #include <CD_ItoLayout.H>
22 #include <CD_CdrLayout.H>
23 #include <CD_PointParticle.H>
24 #include <CD_RtLayout.H>
25 #include <CD_McPhoto.H>
26 #include "CD_FieldSolver.H"
27 #include "CD_ItoSolver.H"
28 #include "CD_CdrCTU.H"
30 #include "CD_SurfaceODESolver.H"
31 #include <CD_NamespaceHeader.H>
32 
33 namespace Physics {
34  namespace ItoKMC {
35 
39  enum class SpeciesSubset
40  {
41  All,
42  AllMobile,
43  AllDiffusive,
44  AllMobileOrDiffusive,
45  AllMobileAndDiffusive,
46  Charged,
47  ChargedMobile,
48  ChargedDiffusive,
49  ChargedMobileOrDiffusive,
50  ChargedMobileAndDiffusive,
51  Stationary,
52  };
53 
58  template <typename I = ItoSolver, typename C = CdrCTU, typename R = McPhoto, typename F = FieldSolverMultigrid>
59  class ItoKMCStepper : public TimeStepper
60  {
61  public:
62  static_assert(std::is_base_of<ItoSolver, I>::value, "I must derive from ItoSolver");
63  static_assert(std::is_base_of<CdrSolver, C>::value, "C must derive from CdrSolver");
64  static_assert(std::is_base_of<McPhoto, R>::value, "R must derive from McPhoto");
65  static_assert(std::is_base_of<FieldSolver, FieldSolverMultigrid>::value, "F must derive from FieldSolver");
66 
70  ItoKMCStepper() noexcept;
71 
76  ItoKMCStepper(RefCountedPtr<ItoKMCPhysics>& a_physics) noexcept;
77 
81  virtual ~ItoKMCStepper() noexcept;
82 
86  virtual void
87  parseOptions() noexcept;
88 
92  virtual void
93  parseRuntimeOptions() noexcept override;
94 
98  virtual void
99  allocate() noexcept override;
100 
104  virtual void
105  setupSolvers() noexcept override;
106 
110  virtual void
111  initialData() noexcept override;
112 
117  virtual void
118  postCheckpointSetup() noexcept override;
119 
123  virtual void
124  postInitialize() noexcept override;
125 
126 #ifdef CH_USE_HDF5
131  virtual void
132  writeCheckpointHeader(HDF5HeaderData& a_header) const noexcept override;
133 #endif
134 
135 #ifdef CH_USE_HDF5
140  virtual void
141  readCheckpointHeader(HDF5HeaderData& a_header) noexcept override;
142 #endif
143 
144 #ifdef CH_USE_HDF5
150  virtual void
151  writeCheckpointData(HDF5Handle& a_handle, const int a_lvl) const noexcept override;
152 #endif
153 
154 #ifdef CH_USE_HDF5
160  virtual void
161  readCheckpointData(HDF5Handle& a_handle, const int a_lvl) noexcept override;
162 #endif
163 
167  virtual int
168  getNumberOfPlotVariables() const noexcept override;
169 
173  virtual Vector<std::string>
174  getPlotVariableNames() const noexcept override;
175 
183  virtual void
184  writePlotData(LevelData<EBCellFAB>& a_output,
185  int& a_icomp,
186  const std::string a_outputRealm,
187  const int a_level) const noexcept override;
188 
198  virtual Vector<long int>
199  getCheckpointLoads(const std::string a_realm, const int a_level) const override;
200 
206  virtual Real
207  advance(const Real a_dt) override = 0;
208 
212  virtual Real
213  computeDt() override;
214 
221  virtual void
222  synchronizeSolverTimes(const int a_step, const Real a_time, const Real a_dt) noexcept override;
223 
227  virtual void
228  printStepReport() noexcept override;
229 
233  virtual void
234  registerRealms() noexcept override;
235 
239  virtual void
240  registerOperators() noexcept override;
241 
245  virtual void
246  prePlot() noexcept override;
247 
251  virtual void
252  postPlot() noexcept override;
253 
259  virtual void
260  preRegrid(const int a_lmin, const int a_oldFinestLevel) noexcept override;
261 
268  virtual void
269  regrid(const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) noexcept override;
270 
274  virtual void
275  postRegrid() noexcept override;
276 
281  virtual bool
282  loadBalanceThisRealm(const std::string a_realm) const override;
283 
297  virtual void
298  loadBalanceBoxes(Vector<Vector<int>>& a_procs,
299  Vector<Vector<Box>>& a_boxes,
300  const std::string a_realm,
301  const Vector<DisjointBoxLayout>& a_grids,
302  const int a_lmin,
303  const int a_finestLevel) override;
304 
310  virtual void
311  setVoltage(const std::function<Real(const Real a_time)>& a_voltage) noexcept;
312 
317  virtual Real
318  getTime() const noexcept;
319 
325  virtual void
326  computeElectricField(EBAMRCellData& a_electricField, const phase::which_phase a_phase) const noexcept;
327 
328  protected:
332  enum class TimeCode
333  {
334  AdvectionIto,
335  DiffusionIto,
336  AdvectionDiffusionIto,
337  AdvectionDiffusionCDR,
338  RelaxationTime,
339  Hardcap,
340  Physics
341  };
342 
347 
352 
357  std::string m_fluidRealm;
358 
362  std::string m_particleRealm;
363 
367  std::string m_name;
368 
373  phase::which_phase m_plasmaPhase;
374 
379  RefCountedPtr<ItoKMCPhysics> m_physics;
380 
384  RefCountedPtr<ItoLayout<ItoSolver>> m_ito;
385 
389  RefCountedPtr<CdrLayout<CdrSolver>> m_cdr;
390 
394  RefCountedPtr<RtLayout<McPhoto>> m_rte;
395 
399  RefCountedPtr<FieldSolver> m_fieldSolver;
400 
404  RefCountedPtr<SurfaceODESolver<1>> m_sigmaSolver;
405 
411  Vector<ParticleContainer<ItoParticle>> m_secondaryParticles;
412 
418  Vector<ParticleContainer<Photon>> m_secondaryPhotons;
419 
425  Vector<EBAMRIVData> m_cdrFluxes;
426 
430  Vector<EBAMRIVData> m_cdrFluxesExtrap;
431 
435  std::function<Real(const Real a_time)> m_voltage;
436 
441 
446 
451 
455  bool m_profile;
456 
461 
466 
471 
476 
481 
486 
490  Vector<int> m_particlesPerCell;
491 
497 
501  Real m_prevDt;
502 
507 
512 
519 
524 
529 
534 
539 
544 
549 
554 
559 
565 
569  Real m_minDt;
570 
574  Real m_maxDt;
575 
580 
585 
590 
595 
600 
605 
609  Vector<EBAMRCellData> m_loadBalancePPC;
610 
615  Vector<EBAMRCellData> m_cdrMobilities;
616 
621  Vector<EBAMRCellData> m_fluidGradPhiCDR;
622 
627  Vector<EBAMRCellData> m_fluidGradPhiIto;
628 
632  Vector<EBAMRCellData> m_fluidPhiIto;
633 
639  mutable Vector<ParticleContainer<PointParticle>> m_cdrPhotoiProducts;
640 
644  Vector<int> m_loadBalanceIndices;
645 
651 
656 
662 
668 
674 
680 
686 
692 
698 
704 
710 
716 
722 
727 
732 
737  Vector<EBAMRCellData> m_energySources;
738 
743 
748 
753 
758 
763 
768 
773 
777  virtual void
778  allocateInternals() noexcept;
779 
783  virtual void
784  setupIto() noexcept;
785 
789  virtual void
790  setupCdr() noexcept;
791 
795  virtual void
796  setupRadiativeTransfer() noexcept;
797 
801  virtual void
802  setupPoisson() noexcept;
803 
807  virtual void
808  setupSigma() noexcept;
809 
819  virtual void
821  const SpeciesSubset a_speciesSubset,
822  const bool a_delete,
823  const std::function<void(ItoParticle&)> a_nonDeletionModifier = [](ItoParticle&) -> void {
824  return;
825  }) noexcept;
826 
837  virtual void
839  const SpeciesSubset a_speciesSubset,
840  const ItoSolver::WhichContainer a_containerBulk,
841  const ItoSolver::WhichContainer a_containerEB,
842  const ItoSolver::WhichContainer a_containerDomain,
843  const bool a_delete,
844  const std::function<void(ItoParticle&)> a_nonDeletionModifier = [](ItoParticle&) -> void {
845  return;
846  }) noexcept;
847 
855  virtual void
857  const EBRepresentation a_representation,
858  const Real a_tolerance) noexcept;
859 
868  virtual void
870  const ItoSolver::WhichContainer a_container,
871  const EBRepresentation a_representation,
872  const Real a_tolerance) noexcept;
873 
881  virtual void
882  transferCoveredParticles(const SpeciesSubset a_speciesSubset,
883  const EBRepresentation a_representation,
884  const Real a_tolerance) noexcept;
885 
895  virtual void
896  transferCoveredParticles(const SpeciesSubset a_speciesSubset,
897  const ItoSolver::WhichContainer a_containerFrom,
898  const ItoSolver::WhichContainer a_containerTo,
899  const EBRepresentation a_representation,
900  const Real a_tolerance) noexcept;
901 
912  virtual void
913  writeData(LevelData<EBCellFAB>& a_output,
914  int& a_comp,
915  const EBAMRCellData& a_data,
916  const std::string a_outputRealm,
917  const int a_level,
918  const bool a_interpToCentroids,
919  const bool a_interpGhost) const noexcept;
920 
928  virtual void
929  writeNumberOfParticlesPerPatch(LevelData<EBCellFAB>& a_output,
930  int& a_icomp,
931  const std::string a_outputRealm,
932  const int a_level) const noexcept;
933 
941  virtual void
942  getMaxMinItoDensity(Real& a_maxDensity,
943  Real& a_minDensity,
944  std::string& a_maxSolver,
945  std::string& a_minSolver) const noexcept;
946 
954  virtual void
955  getMaxMinCDRDensity(Real& a_maxDensity,
956  Real& a_minDensity,
957  std::string& a_maxSolver,
958  std::string& a_minSolver) const noexcept;
959 
969  virtual void
970  getParticleStatistics(Real& a_avgParticles,
971  Real& a_sigma,
972  Real& a_minParticles,
973  Real& a_maxParticles,
974  int& a_minRank,
975  int& a_maxRank);
976 
986  virtual void
987  loadBalanceParticleRealm(Vector<Vector<int>>& a_procs,
988  Vector<Vector<Box>>& a_boxes,
989  const std::string a_realm,
990  const Vector<DisjointBoxLayout>& a_grids,
991  const int a_lmin,
992  const int a_finestLevel) noexcept;
993 
1003  virtual void
1004  loadBalanceFluidRealm(Vector<Vector<int>>& a_procs,
1005  Vector<Vector<Box>>& a_boxes,
1006  const std::string a_realm,
1007  const Vector<DisjointBoxLayout>& a_grids,
1008  const int a_lmin,
1009  const int a_finestLevel) noexcept;
1010 
1015  virtual Vector<RefCountedPtr<ItoSolver>>
1016  getLoadBalanceSolvers() const noexcept;
1017 
1022  virtual Real
1023  computeMaxElectricField(const phase::which_phase a_phase) const noexcept;
1024 
1029  virtual void
1030  computeSpaceChargeDensity() noexcept;
1031 
1038  virtual void
1040  const Vector<EBAMRCellData*>& a_itoDensities,
1041  const Vector<EBAMRCellData*>& a_cdrDensities) noexcept;
1042 
1048  virtual void
1049  computeConductivityCell(EBAMRCellData& a_conductivity) noexcept;
1050 
1057  virtual void
1058  computeConductivityCell(EBAMRCellData& a_conductivity,
1059  const Vector<ParticleContainer<ItoParticle>*>& a_particles) noexcept;
1060 
1066  virtual void
1067  computeDensityGradients() noexcept;
1068 
1073  virtual void
1074  computeCurrentDensity(EBAMRCellData& a_J) noexcept;
1075 
1079  virtual Real
1080  computeRelaxationTime() noexcept;
1081 
1088  virtual bool
1089  solvePoisson() noexcept;
1090 
1095  virtual void
1096  depositParticles(const SpeciesSubset a_speciesSubset) noexcept;
1097 
1103  virtual void
1104  depositParticles(const SpeciesSubset a_speciesSubset, const ItoSolver::WhichContainer a_container) noexcept;
1105 
1111  virtual void
1112  remapParticles(const SpeciesSubset a_speciesSubset) noexcept;
1113 
1119  virtual void
1120  remapParticles(const SpeciesSubset a_speciesSubset, const ItoSolver::WhichContainer a_container) noexcept;
1121 
1126  virtual void
1127  computeDriftVelocities() noexcept;
1128 
1133  virtual void
1134  setItoVelocityFunctions() noexcept;
1135 
1140  virtual void
1141  setCdrVelocityFunctions() noexcept;
1142 
1146  virtual void
1148 
1154  virtual void
1155  computeMobilities() noexcept;
1156 
1164  virtual void
1165  computeMobilities(Vector<EBAMRCellData*>& a_itoMobilities,
1166  Vector<EBAMRCellData>& a_cdrMobilities,
1167  const EBAMRCellData& a_electricField,
1168  const Real a_time) noexcept;
1169 
1178  virtual void
1179  computeMobilities(Vector<LevelData<EBCellFAB>*>& a_itoMobilities,
1180  Vector<LevelData<EBCellFAB>*>& a_cdrMobilities,
1181  const LevelData<EBCellFAB>& a_E,
1182  const int a_level,
1183  const Real a_time) noexcept;
1184 
1195  virtual void
1196  computeMobilities(Vector<EBCellFAB*>& a_itoMobilities,
1197  Vector<EBCellFAB*>& a_cdrMobilities,
1198  const EBCellFAB& a_electricField,
1199  const int a_level,
1200  const DataIndex a_dit,
1201  const Box a_cellBox,
1202  const Real a_time) noexcept;
1203 
1208  virtual void
1209  computeDiffusionCoefficients() noexcept;
1210 
1218  virtual void
1219  computeDiffusionCoefficients(Vector<EBAMRCellData*>& a_itoDiffusionCoefficients,
1220  Vector<EBAMRCellData*>& a_cdrDiffusionCoefficients,
1221  const EBAMRCellData& a_electricField,
1222  const Real a_time) noexcept;
1223 
1232  virtual void
1233  computeDiffusionCoefficients(Vector<LevelData<EBCellFAB>*>& a_itoDiffusionCoefficients,
1234  Vector<LevelData<EBCellFAB>*>& a_cdrDiffusionCoefficients,
1235  const LevelData<EBCellFAB>& a_electricField,
1236  const int a_level,
1237  const Real a_time) noexcept;
1238 
1251  virtual void
1252  computeDiffusionCoefficients(Vector<EBCellFAB*>& a_itoDiffusionCoefficients,
1253  Vector<EBCellFAB*>& a_cdrDiffusionCoefficients,
1254  const EBCellFAB& a_E,
1255  const int a_level,
1256  const DataIndex a_dit,
1257  const Box a_box,
1258  const Real a_time) noexcept;
1259 
1263  virtual void
1265 
1270  virtual void
1271  getPhysicalParticlesPerCell(EBAMRCellData& a_ppc) const noexcept;
1272 
1278  virtual void
1280 
1287  virtual void
1288  computeReactiveItoParticlesPerCell(LevelData<EBCellFAB>& a_ppc, const int a_level) noexcept;
1289 
1299  virtual void
1300  computeReactiveItoParticlesPerCell(EBCellFAB& a_ppc,
1301  const int a_level,
1302  const DataIndex a_dit,
1303  const Box a_box,
1304  const EBISBox& a_ebisbox) noexcept;
1305 
1312  virtual void
1314 
1321  virtual void
1322  computeReactiveCdrParticlesPerCell(LevelData<EBCellFAB>& a_ppc, const int a_level) noexcept;
1323 
1333  virtual void
1334  computeReactiveCdrParticlesPerCell(EBCellFAB& a_ppc,
1335  const int a_level,
1336  const DataIndex a_dit,
1337  const Box a_box,
1338  const EBISBox& a_ebisbox) noexcept;
1339 
1345  virtual void
1346  computeReactiveMeanEnergiesPerCell(EBAMRCellData& a_meanEnergies) noexcept;
1347 
1354  virtual void
1355  computeReactiveMeanEnergiesPerCell(LevelData<EBCellFAB>& a_meanEnergies, const int a_level) noexcept;
1356 
1366  virtual void
1367  computeReactiveMeanEnergiesPerCell(EBCellFAB& a_meanEnergies,
1368  const int a_level,
1369  const DataIndex a_dit,
1370  const Box a_box,
1371  const EBISBox& a_ebisbox) noexcept;
1372 
1378  virtual void
1379  advanceReactionNetwork(const Real a_dt) noexcept;
1380 
1387  virtual void
1388  advanceReactionNetwork(const EBAMRCellData& a_E, const Real a_dt) noexcept;
1389 
1400  inline void
1401  advanceReactionNetwork(LevelData<EBCellFAB>& a_particlesPerCell,
1402  LevelData<EBCellFAB>& a_newPhotonsPerCell,
1403  const LevelData<EBCellFAB>& a_electricField,
1404  const int a_level,
1405  const Real a_dt) const noexcept;
1406 
1420  inline void
1421  advanceReactionNetwork(EBCellFAB& a_particlesPerCell,
1422  EBCellFAB& a_newPhotonsPerCell,
1423  const EBCellFAB& a_electricField,
1424  const int a_level,
1425  const DataIndex a_dit,
1426  const Box a_box,
1427  const Real a_dx,
1428  const Real a_dt) const noexcept;
1429 
1437  inline void
1438  reconcileParticles(const EBAMRCellData& a_newParticlesPerCell,
1439  const EBAMRCellData& a_oldParticlesPerCell,
1440  const EBAMRCellData& a_newPhotonsPerCell) const noexcept;
1441 
1450  inline void
1451  reconcileParticles(const LevelData<EBCellFAB>& a_newParticlesPerCell,
1452  const LevelData<EBCellFAB>& a_oldParticlesPerCell,
1453  const LevelData<EBCellFAB>& a_newPhotonsPerCell,
1454  const int a_level) const noexcept;
1455 
1467  inline void
1468  reconcileParticles(const EBCellFAB& a_newParticlesPerCell,
1469  const EBCellFAB& a_oldParticlesPerCell,
1470  const EBCellFAB& a_newPhotonsPerCell,
1471  const int a_level,
1472  const DataIndex a_dit,
1473  const Box a_box,
1474  const Real a_dx) const noexcept;
1475 
1482  virtual void
1483  reconcilePhotoionization() noexcept;
1484 
1492  virtual void
1493  reconcileCdrDensities(const EBAMRCellData& a_newParticlesPerCell,
1494  const EBAMRCellData& a_oldParticlesPerCell,
1495  const Real a_dt) noexcept;
1496 
1505  virtual void
1506  reconcileCdrDensities(const LevelData<EBCellFAB>& a_newParticlesPerCell,
1507  const LevelData<EBCellFAB>& a_oldParticlesPerCell,
1508  const int a_level,
1509  const Real a_dt) noexcept;
1510 
1522  virtual void
1523  reconcileCdrDensities(const EBCellFAB& a_newParticlesPerCell,
1524  const EBCellFAB& a_oldParticlesPerCell,
1525  const int a_level,
1526  const DataIndex a_dit,
1527  const Box a_box,
1528  const Real a_dx,
1529  const Real a_dt) noexcept;
1530 
1534  virtual void
1535  coarsenCDRSolvers() noexcept;
1536 
1542  virtual void
1543  fillSecondaryEmissionEB(const Real a_dt) noexcept;
1544 
1557  virtual void
1558  fillSecondaryEmissionEB(Vector<ParticleContainer<ItoParticle>>& a_secondaryParticles,
1559  Vector<EBAMRIVData>& a_cdrFluxes,
1560  Vector<ParticleContainer<Photon>>& a_secondaryPhotons,
1561  Vector<ParticleContainer<ItoParticle>*>& a_primaryParticles,
1562  Vector<EBAMRIVData>& a_cdrFluxesExtrap,
1563  Vector<ParticleContainer<Photon>*>& a_primaryPhotons,
1564  const EBAMRCellData& a_electricField,
1565  const Real a_dt) noexcept;
1566 
1572  virtual void
1573  resolveSecondaryEmissionEB(const Real a_dt) noexcept;
1574 
1587  virtual void
1588  resolveSecondaryEmissionEB(Vector<ParticleContainer<ItoParticle>*>& a_secondaryParticles,
1589  Vector<ParticleContainer<ItoParticle>*>& a_primaryParticles,
1590  Vector<EBAMRIVData*>& a_cdrFluxes,
1591  EBAMRIVData& a_surfaceChargeDensity,
1592  const Real a_dt) noexcept;
1593 
1599  virtual Real
1600  computePhysicsDt() noexcept;
1601 
1608  virtual Real
1609  computePhysicsDt(const EBAMRCellData& a_electricField) noexcept;
1610 
1618  virtual Real
1619  computePhysicsDt(const LevelData<EBCellFAB>& a_electricField,
1620  const LevelData<EBCellFAB>& a_particlesPerCell,
1621  const int a_level) noexcept;
1622 
1632  virtual Real
1633  computePhysicsDt(const EBCellFAB& a_electricField,
1634  const EBCellFAB& a_particlesPercCell,
1635  const int a_level,
1636  const DataIndex a_dit,
1637  const Box a_box) noexcept;
1638 
1642  virtual Real
1643  computeTotalCharge() const noexcept;
1644 
1648  virtual Real
1649  computeQplus() const noexcept;
1650 
1654  virtual Real
1655  computeQminu() const noexcept;
1656 
1660  virtual Real
1661  computeQsurf() const noexcept;
1662 
1667  virtual void
1668  advancePhotons(const Real a_dt) noexcept;
1669 
1674  virtual void
1675  sortPhotonsByCell(const McPhoto::WhichContainer a_which) noexcept;
1676 
1681  virtual void
1682  sortPhotonsByPatch(const McPhoto::WhichContainer a_which) noexcept;
1683 
1688  virtual void
1689  postCheckpointPoisson() noexcept;
1690 
1696  virtual void
1697  computeEdotJSource(const Real a_dt) noexcept;
1698 
1702  virtual void
1703  initialSigma() noexcept;
1704 
1709  virtual void
1710  parseVerbosity() noexcept;
1711 
1715  virtual void
1716  parseExitOnFailure() noexcept;
1717 
1721  virtual void
1722  parseRedistributeCDR() noexcept;
1723 
1728  virtual void
1729  parseSuperParticles() noexcept;
1730 
1735  virtual void
1736  parseDualGrid() noexcept;
1737 
1741  virtual void
1742  parseLoadBalance() noexcept;
1743 
1747  virtual void
1748  parseTimeStepRestrictions() noexcept;
1749 
1753  virtual void
1754  parseParametersEB() noexcept;
1755 
1759  virtual void
1760  parsePlotVariables() noexcept;
1761 
1766  virtual void
1767  computePhysicsPlotVariables(EBAMRCellData& a_physicsPlotVars) noexcept;
1768  };
1769  } // namespace ItoKMC
1770 } // namespace Physics
1771 
1772 #include <CD_NamespaceFooter.H>
1773 
1774 #include <CD_ItoKMCStepperImplem.H>
1775 
1776 #endif
BoxSorting
Enum for sorting boxes.
Definition: CD_BoxSorting.H:21
Declaration of a class which implements CdrMultigrid using MUSCL for advection.
Declaration of a class that holds a set of CdrSolvers (to cut down on typing).
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 FieldSolverMultigrid.
Contains declaration of a base electrostatics solver class.
Main file for describing Ito-based plasma physics.
Implementation of CD_ItoKMCStepper.H.
SpeciesSubset
Enum for differentiating between types of particles.
Definition: CD_ItoKMCStepper.H:40
Declaration of a class that holds a set of ItoSolvers.
Declaration of solver class for Ito diffusion.
Declaration of a radiative transfer solver which uses Monte Carlo sampling of computational or real p...
Declaration of a computational point particle.
Declaration of a class that holds a set of RtSolvers.
Declaration of a cut-cell ODE solver.
Declaration of main (abstract) time stepper class.
A particle class for use with ItoSolvers, i.e. drifting Brownian walkers.
Definition: CD_ItoParticle.H:40
Base class for Ito diffusion particle models.
Definition: CD_ItoSolver.H:35
WhichContainer
Enum class for distinguishing various types of particle containers.
Definition: CD_ItoSolver.H:49
Radiative tranfer equation solver using Monte-Carlo simulation.
Definition: CD_McPhoto.H:39
Templated class for holding particles on an AMR hierarchy with particle remapping.
Definition: CD_ParticleContainer.H:50
Particle class for usage with Monte Carlo radiative transfer.
Definition: CD_Photon.H:29
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 loadBalanceBoxes(Vector< Vector< int >> &a_procs, Vector< Vector< Box >> &a_boxes, const std::string a_realm, const Vector< DisjointBoxLayout > &a_grids, const int a_lmin, const int a_finestLevel) override
Load balance grid boxes for a specified realm.
Definition: CD_ItoKMCStepperImplem.H:5509
virtual void transferCoveredParticles(const SpeciesSubset a_speciesSubset, const EBRepresentation a_representation, const Real a_tolerance) noexcept
Transfer covered particles (i.e., particles inside the EB) from the ItoSolver bulk container to EB co...
Definition: CD_ItoKMCStepperImplem.H:2335
Real m_relaxationTime
The relaxation time eps0/sigma.
Definition: CD_ItoKMCStepper.H:599
virtual void setupCdr() noexcept
Set up the CDR solvers.
Definition: CD_ItoKMCStepperImplem.H:433
bool m_abortOnFailure
Flag for abandoning simulation of Poisson solver fails.
Definition: CD_ItoKMCStepper.H:440
virtual void getParticleStatistics(Real &a_avgParticles, Real &a_sigma, Real &a_minParticles, Real &a_maxParticles, int &a_minRank, int &a_maxRank)
Compute some particle statistics.
Definition: CD_ItoKMCStepperImplem.H:1339
RefCountedPtr< SurfaceODESolver< 1 > > m_sigmaSolver
Surface charge solver.
Definition: CD_ItoKMCStepper.H:404
virtual void computePhysicsPlotVariables(EBAMRCellData &a_physicsPlotVars) noexcept
Compute physics plot variables.
Definition: CD_ItoKMCStepperImplem.H:5954
Vector< EBAMRCellData > m_fluidPhiIto
For holding the Ito species densities on the fluid realm.
Definition: CD_ItoKMCStepper.H:632
virtual void computeReactiveMeanEnergiesPerCell(EBAMRCellData &a_meanEnergies) noexcept
Compute the mean particle energy in all grid cells.
Definition: CD_ItoKMCStepperImplem.H:3603
EBAMRCellData m_conductivityCell
Cell-centered conductivity.
Definition: CD_ItoKMCStepper.H:661
EBAMRCellData m_particleOldItoPPC
For holding the previous number of particles per cell for all species.
Definition: CD_ItoKMCStepper.H:697
Real m_particleAdvectionDt
The particle advective time step.
Definition: CD_ItoKMCStepper.H:579
Real m_minParticleAdvectionCFL
Minimum CFL-like time step for particle advection.
Definition: CD_ItoKMCStepper.H:528
Vector< EBAMRIVData > m_cdrFluxes
CDR fluxes for CDR BCs.
Definition: CD_ItoKMCStepper.H:425
virtual void parseRuntimeOptions() noexcept override
Parse runtime configurable options.
Definition: CD_ItoKMCStepperImplem.H:106
std::string m_name
Time stepper name.
Definition: CD_ItoKMCStepper.H:367
phase::which_phase m_plasmaPhase
Phase where we solve for the plasma.
Definition: CD_ItoKMCStepper.H:373
bool m_profile
Profile kernels or not.
Definition: CD_ItoKMCStepper.H:455
virtual void computeElectricField(EBAMRCellData &a_electricField, const phase::which_phase a_phase) const noexcept
Recompute the electric field onto the specified data holder.
Definition: CD_ItoKMCStepperImplem.H:1740
virtual void removeCoveredParticles(const SpeciesSubset a_which, const EBRepresentation a_representation, const Real a_tolerance) noexcept
Remove covered particles (i.e., particles inside the EB)
Definition: CD_ItoKMCStepperImplem.H:2216
Vector< ParticleContainer< PointParticle > > m_cdrPhotoiProducts
Photoionization products to be put in the CDR equations.
Definition: CD_ItoKMCStepper.H:639
BoxSorting m_boxSort
Box sorting method when using dual-grid with particle load balancing.
Definition: CD_ItoKMCStepper.H:346
Real m_fluidAdvectionDiffusionCFL
CFL-like time step for fluid advection-diffusion.
Definition: CD_ItoKMCStepper.H:558
virtual void setupRadiativeTransfer() noexcept
Set up the radiative transfer solver.
Definition: CD_ItoKMCStepperImplem.H:452
Real m_physicsDt
The physics-based time step.
Definition: CD_ItoKMCStepper.H:604
EBAMRCellData m_electricFieldParticle
Storage for holding the plasma phase electric field on the particle realm.
Definition: CD_ItoKMCStepper.H:731
virtual void registerRealms() noexcept override
Register realms used for the simulation.
Definition: CD_ItoKMCStepperImplem.H:1485
Real m_particleAdvectionDiffusionDt
The advection-diffusion time step.
Definition: CD_ItoKMCStepper.H:589
virtual Vector< RefCountedPtr< ItoSolver > > getLoadBalanceSolvers() const noexcept
Get the solvers used for load balancing.
Definition: CD_ItoKMCStepperImplem.H:5453
virtual void computeSpaceChargeDensity() noexcept
Compute the space charge. Calls the other version.
Definition: CD_ItoKMCStepperImplem.H:1767
bool m_plotCurrentDensity
Plot conductivity or not.
Definition: CD_ItoKMCStepper.H:480
Real m_maxShrinkDt
Maximum permissible time step shrinkage.
Definition: CD_ItoKMCStepper.H:511
virtual void setVoltage(const std::function< Real(const Real a_time)> &a_voltage) noexcept
Set voltage used for the simulation.
Definition: CD_ItoKMCStepperImplem.H:1700
virtual void advanceReactionNetwork(const Real a_dt) noexcept
Chemistry advance over time a_dt.
Definition: CD_ItoKMCStepperImplem.H:3738
virtual Real getTime() const noexcept
Get current simulation time.
Definition: CD_ItoKMCStepperImplem.H:1755
virtual int getNumberOfPlotVariables() const noexcept override
Get number of plot variables for the output file.
Definition: CD_ItoKMCStepperImplem.H:864
Vector< ParticleContainer< Photon > > m_secondaryPhotons
List of secondary photons injected through the EB.
Definition: CD_ItoKMCStepper.H:418
Real m_loadPerCell
The "background" load per cell when using particle load balancing.
Definition: CD_ItoKMCStepper.H:518
EBAMRCellData m_physicsPlotVariables
Storage for physics plot variables.
Definition: CD_ItoKMCStepper.H:655
virtual void remapParticles(const SpeciesSubset a_speciesSubset) noexcept
Remap a subset of ItoSolver particles.
Definition: CD_ItoKMCStepperImplem.H:2459
bool m_redistributeCDR
Flag for abandoning simulation of Poisson solver fails.
Definition: CD_ItoKMCStepper.H:445
virtual void parseVerbosity() noexcept
Parse chattiness.
Definition: CD_ItoKMCStepperImplem.H:133
virtual void computeReactiveCdrParticlesPerCell(EBAMRCellData &a_ppc) noexcept
Compute the number of reactive particles per cell for the CDR solvers.
Definition: CD_ItoKMCStepperImplem.H:3498
RefCountedPtr< FieldSolver > m_fieldSolver
Field solver.
Definition: CD_ItoKMCStepper.H:399
virtual void registerOperators() noexcept override
Register operators used for the simulation.
Definition: CD_ItoKMCStepperImplem.H:1499
EBAMRCellData m_particleScratch1
Scratch storage on the particle realm with 1 component.
Definition: CD_ItoKMCStepper.H:757
EBAMRCellData m_EdotJ
Storage for holdnig E*J on the fluid realm. 1 component.
Definition: CD_ItoKMCStepper.H:742
virtual void averageDiffusionCoefficientsCellToFace() noexcept
Average cell-centered diffusion coefficient to faces.
Definition: CD_ItoKMCStepperImplem.H:3314
virtual void parsePlotVariables() noexcept
Parse plot variables.
Definition: CD_ItoKMCStepperImplem.H:176
virtual void parseSuperParticles() noexcept
Parse the desired number of particles per cell.
Definition: CD_ItoKMCStepperImplem.H:212
virtual void writeData(LevelData< EBCellFAB > &a_output, int &a_comp, const EBAMRCellData &a_data, const std::string a_outputRealm, const int a_level, const bool a_interpToCentroids, const bool a_interpGhost) const noexcept
Write data to output. Convenience function.
Definition: CD_ItoKMCStepperImplem.H:1022
Vector< ParticleContainer< ItoParticle > > m_secondaryParticles
List of secondary particles injected through the EB.
Definition: CD_ItoKMCStepper.H:411
virtual void computeDensityGradients() noexcept
Compute grad(phi) and phi for both CDR and Ito species and put the result on the fluid realm.
Definition: CD_ItoKMCStepperImplem.H:1902
virtual void computeEdotJSource(const Real a_dt) noexcept
Compute the energy source term for the various plasma species.
Definition: CD_ItoKMCStepperImplem.H:5781
EBAMRCellData m_particleItoPPC
For holding the number of physical particles per cell for all Ito species.
Definition: CD_ItoKMCStepper.H:679
virtual bool solvePoisson() noexcept
Solve the electrostatic problem.
Definition: CD_ItoKMCStepperImplem.H:1994
virtual Real computeQminu() const noexcept
Compute negative charge.
Definition: CD_ItoKMCStepperImplem.H:5318
virtual void multiplyCdrVelocitiesByMobilities() noexcept
Multiply CDR solver velocities by mobilities.
Definition: CD_ItoKMCStepperImplem.H:2754
virtual void parseRedistributeCDR() noexcept
Parse CDR mass redistribution when assigning reactive products.
Definition: CD_ItoKMCStepperImplem.H:162
virtual void computeCurrentDensity(EBAMRCellData &a_J) noexcept
Compute the current density.
Definition: CD_ItoKMCStepperImplem.H:1942
TimeCode
An enum for encapsulating how time steps were restricted.
Definition: CD_ItoKMCStepper.H:333
Vector< int > m_particlesPerCell
Target number of particles per cell when squishing ItoParticle's into superparticles.
Definition: CD_ItoKMCStepper.H:490
Real m_toleranceEB
Accepted tolerance (relative to dx) for EB intersection.
Definition: CD_ItoKMCStepper.H:523
virtual void writeNumberOfParticlesPerPatch(LevelData< EBCellFAB > &a_output, int &a_icomp, const std::string a_outputRealm, const int a_level) const noexcept
Write number of particles per patch to output holder.
Definition: CD_ItoKMCStepperImplem.H:1089
EBAMRCellData m_fluidScratchD
Scratch storage on the fluid realm with SpaceDim components.
Definition: CD_ItoKMCStepper.H:752
bool m_loadBalanceFluid
Load balance fluid realm or not.
Definition: CD_ItoKMCStepper.H:470
virtual void parseTimeStepRestrictions() noexcept
Parse time step restrictions.
Definition: CD_ItoKMCStepperImplem.H:309
EBAMRIVData m_particleScratchEB
Scratch storage for EB-only data on the particle realm. One component.
Definition: CD_ItoKMCStepper.H:772
virtual void fillSecondaryEmissionEB(const Real a_dt) noexcept
Resolve particle injection at EBs.
Definition: CD_ItoKMCStepperImplem.H:4651
EBAMRCellData m_particleYPC
For holding the number of generated photons per cell.
Definition: CD_ItoKMCStepper.H:685
virtual void initialSigma() noexcept
Fill surface charge solver with initial data taken from the physics interface.
Definition: CD_ItoKMCStepperImplem.H:679
ItoKMCStepper() noexcept
Default constructor. Sets default options.
Definition: CD_ItoKMCStepperImplem.H:35
EBAMRCellData m_currentDensity
Storage for current density.
Definition: CD_ItoKMCStepper.H:650
virtual void printStepReport() noexcept override
Print a step report. Used by Driver for user monitoring of simulation.
Definition: CD_ItoKMCStepperImplem.H:1154
virtual void depositParticles(const SpeciesSubset a_speciesSubset) noexcept
Deposit a subset of the ItoSolver particles on the mesh.
Definition: CD_ItoKMCStepperImplem.H:2574
Vector< EBAMRCellData > m_cdrMobilities
For holding the mobilities for the CDR species.
Definition: CD_ItoKMCStepper.H:615
virtual Real computePhysicsDt() noexcept
Compute a maximum time step from the physics interface.
Definition: CD_ItoKMCStepperImplem.H:5075
bool m_regridSuperparticles
Do or do not superparticle merging/splitting on regrids.
Definition: CD_ItoKMCStepper.H:460
virtual void getMaxMinItoDensity(Real &a_maxDensity, Real &a_minDensity, std::string &a_maxSolver, std::string &a_minSolver) const noexcept
Get maximum density of the Ito species.
Definition: CD_ItoKMCStepperImplem.H:1273
virtual void setupSigma() noexcept
Set up the surface charge solver.
Definition: CD_ItoKMCStepperImplem.H:489
virtual void loadBalanceFluidRealm(Vector< Vector< int >> &a_procs, Vector< Vector< Box >> &a_boxes, const std::string a_realm, const Vector< DisjointBoxLayout > &a_grids, const int a_lmin, const int a_finestLevel) noexcept
Routine called by loadBalanceBoxes and used for particle-based load balancing.
Definition: CD_ItoKMCStepperImplem.H:5684
Vector< EBAMRCellData > m_loadBalancePPC
For holding the number of computational particles per cell when load balancing.
Definition: CD_ItoKMCStepper.H:609
virtual void setupSolvers() noexcept override
Set up solvers.
Definition: CD_ItoKMCStepperImplem.H:398
EBAMRCellData m_fluidPPC
For holding the number of particles per cell for all species.
Definition: CD_ItoKMCStepper.H:715
virtual void parseOptions() noexcept
Parse options.
Definition: CD_ItoKMCStepperImplem.H:86
virtual void advancePhotons(const Real a_dt) noexcept
Photon advancement routine.
Definition: CD_ItoKMCStepperImplem.H:5374
Real m_minDt
Minimum permitted time step.
Definition: CD_ItoKMCStepper.H:569
int m_mergeInterval
How often to merge superparticles.
Definition: CD_ItoKMCStepper.H:496
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:1547
virtual void loadBalanceParticleRealm(Vector< Vector< int >> &a_procs, Vector< Vector< Box >> &a_boxes, const std::string a_realm, const Vector< DisjointBoxLayout > &a_grids, const int a_lmin, const int a_finestLevel) noexcept
Routine called by loadBalanceBoxes and used for particle-based load balancing.
Definition: CD_ItoKMCStepperImplem.H:5531
EBAMRIVData m_fluidScratchEB
Scratch storage for EB-only data on the fluid realm. One component.
Definition: CD_ItoKMCStepper.H:767
virtual void parseDualGrid() noexcept
Parse dual or single realm calculations.
Definition: CD_ItoKMCStepperImplem.H:236
virtual void reconcilePhotoionization() noexcept
Reconcile the results from photoionization reactions.
Definition: CD_ItoKMCStepperImplem.H:4426
EBAMRCellData m_fluidScratch1
Scratch storage on the fluid realm having 1 component.
Definition: CD_ItoKMCStepper.H:747
EBAMRCellData m_fluidOldCdrPPC
For holding the previous number of physical particles per cell for all CDR species.
Definition: CD_ItoKMCStepper.H:709
virtual Vector< std::string > getPlotVariableNames() const noexcept override
Get plot variable names.
Definition: CD_ItoKMCStepperImplem.H:917
std::string m_fluidRealm
Realm used for the fluid part (i.e., electrostatic) part of the simulation.
Definition: CD_ItoKMCStepper.H:357
virtual void sortPhotonsByCell(const McPhoto::WhichContainer a_which) noexcept
Sort photons by cells.
Definition: CD_ItoKMCStepperImplem.H:5425
virtual Real computeQplus() const noexcept
Compute positive charge.
Definition: CD_ItoKMCStepperImplem.H:5274
virtual void parseLoadBalance() noexcept
Parse load balancing.
Definition: CD_ItoKMCStepperImplem.H:259
virtual void computeDriftVelocities() noexcept
Compute ItoSolver velocities.
Definition: CD_ItoKMCStepperImplem.H:2780
virtual void getMaxMinCDRDensity(Real &a_maxDensity, Real &a_minDensity, std::string &a_maxSolver, std::string &a_minSolver) const noexcept
Get maximum density of the CDr species.
Definition: CD_ItoKMCStepperImplem.H:1306
EBAMRCellData m_electricFieldFluid
Storage for holding the plasma phase electric field on the fluid realm.
Definition: CD_ItoKMCStepper.H:726
Real m_prevDt
Previous time step.
Definition: CD_ItoKMCStepper.H:501
std::string m_particleRealm
Realm used for the particle part of the simulation.
Definition: CD_ItoKMCStepper.H:362
EBAMRCellData m_particleScratchD
Scratch storage on the particle realm with SpaceDim components.
Definition: CD_ItoKMCStepper.H:762
Real m_fluidAdvectionDiffusionDt
The advection-diffusion time step for the CDR equations.
Definition: CD_ItoKMCStepper.H:594
Vector< EBAMRCellData > m_fluidGradPhiCDR
For holding the gradient of all CDR species densities.
Definition: CD_ItoKMCStepper.H:621
virtual void setupPoisson() noexcept
Set up the electrostatic field solver.
Definition: CD_ItoKMCStepperImplem.H:472
Real m_relaxTimeFactor
Factor proportional to the dielectric relaxation time dtRelax = eps0/sigma.
Definition: CD_ItoKMCStepper.H:564
RefCountedPtr< CdrLayout< CdrSolver > > m_cdr
CDR solvers.
Definition: CD_ItoKMCStepper.H:389
virtual void setupIto() noexcept
Set up the Ito particle solvers.
Definition: CD_ItoKMCStepperImplem.H:414
virtual Real computeQsurf() const noexcept
Compute surface charge.
Definition: CD_ItoKMCStepperImplem.H:5362
virtual Real computeDt() override
Compute a time step used for the advance method.
Definition: CD_ItoKMCStepperImplem.H:1371
virtual void synchronizeSolverTimes(const int a_step, const Real a_time, const Real a_dt) noexcept override
Synchronize solver times for all the solvers.
Definition: CD_ItoKMCStepperImplem.H:1135
virtual Real computeTotalCharge() const noexcept
Compute total charge.
Definition: CD_ItoKMCStepperImplem.H:5254
virtual void resolveSecondaryEmissionEB(const Real a_dt) noexcept
Resolve secondary emission at the EB.
Definition: CD_ItoKMCStepperImplem.H:4924
Real m_maxDt
Maximum permitted time step.
Definition: CD_ItoKMCStepper.H:574
TimeCode m_timeCode
Time code for understanding how the time step was restricted.
Definition: CD_ItoKMCStepper.H:351
void reconcileParticles(const EBAMRCellData &a_newParticlesPerCell, const EBAMRCellData &a_oldParticlesPerCell, const EBAMRCellData &a_newPhotonsPerCell) const noexcept
Reconcile particles. At the bottom, this will call the physics interface for particle reconciliation.
Definition: CD_ItoKMCStepperImplem.H:4094
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_ItoKMCStepperImplem.H:1645
Real m_minParticleAdvectionDiffusionCFL
Maximum CFL-like time step for particle advection-diffusion.
Definition: CD_ItoKMCStepper.H:548
virtual Vector< long int > getCheckpointLoads(const std::string a_realm, const int a_level) const override
Get computational loads to be checkpointed.
Definition: CD_ItoKMCStepperImplem.H:5732
EBAMRCellData m_fluidCdrPPC
For holding the number of physical particles per cell for all CDR species.
Definition: CD_ItoKMCStepper.H:703
EBAMRFluxData m_conductivityFace
Face-centered conductivity.
Definition: CD_ItoKMCStepper.H:667
virtual void computeDiffusionCoefficients() noexcept
Compute mesh-based diffusion coefficients for LFA coupling.
Definition: CD_ItoKMCStepperImplem.H:3030
virtual void allocateInternals() noexcept
Allocate "internal" storage.
Definition: CD_ItoKMCStepperImplem.H:524
bool m_plotParticlesPerPatch
Plot number of particles per patch or not.
Definition: CD_ItoKMCStepper.H:485
EBAMRIVData m_conductivityEB
EB-centered conductivity.
Definition: CD_ItoKMCStepper.H:673
Vector< EBAMRIVData > m_cdrFluxesExtrap
Extrapolated CDR fluxes.
Definition: CD_ItoKMCStepper.H:430
Real m_maxParticleAdvectionCFL
Maximum CFL-like time step for particle advection.
Definition: CD_ItoKMCStepper.H:533
virtual Real computeMaxElectricField(const phase::which_phase a_phase) const noexcept
Compute the maximum electric field (norm)
Definition: CD_ItoKMCStepperImplem.H:1712
virtual void setCdrVelocityFunctions() noexcept
Set the Cdr velocities to be sgn(charge) * E.
Definition: CD_ItoKMCStepperImplem.H:2720
virtual void postRegrid() noexcept override
Perform post-regrid operations.
Definition: CD_ItoKMCStepperImplem.H:1687
virtual Real computeRelaxationTime() noexcept
Compute the dielectric relaxation time.
Definition: CD_ItoKMCStepperImplem.H:1962
virtual void parseParametersEB() noexcept
Parse parameters related to how we treat particle-EB interaction.
Definition: CD_ItoKMCStepperImplem.H:382
virtual Real advance(const Real a_dt) override=0
Advancement method. Needs to be implemented by subclasses.
virtual void initialData() noexcept override
Fill solvers with initial data.
Definition: CD_ItoKMCStepperImplem.H:644
Real m_maxParticleDiffusionCFL
Maximum CFL-like time step for particle diffusion.
Definition: CD_ItoKMCStepper.H:543
virtual void postPlot() noexcept override
Perform post-plot operations.
Definition: CD_ItoKMCStepperImplem.H:1535
virtual ~ItoKMCStepper() noexcept
Destructor.
Definition: CD_ItoKMCStepperImplem.H:79
RefCountedPtr< ItoKMCPhysics > m_physics
Implementation of ItoKMCPhysics.
Definition: CD_ItoKMCStepper.H:379
virtual void computeReactiveItoParticlesPerCell(EBAMRCellData &a_ppc) noexcept
Compute the number of reactive particles per cell.
Definition: CD_ItoKMCStepperImplem.H:3378
virtual void writePlotData(LevelData< EBCellFAB > &a_output, int &a_icomp, const std::string a_outputRealm, const int a_level) const noexcept override
Write plot data to output holder.
Definition: CD_ItoKMCStepperImplem.H:968
bool m_dualGrid
Using dual grid or not.
Definition: CD_ItoKMCStepper.H:450
EBAMRCellData m_fluidYPC
For holding the number of generated photons per cell.
Definition: CD_ItoKMCStepper.H:721
virtual void allocate() noexcept override
Allocate storage for solvers.
Definition: CD_ItoKMCStepperImplem.H:506
RefCountedPtr< RtLayout< McPhoto > > m_rte
Radiative transfer solvers.
Definition: CD_ItoKMCStepper.H:394
bool m_plotConductivity
Plot conductivity or not.
Definition: CD_ItoKMCStepper.H:475
virtual void parseExitOnFailure() noexcept
Parse exit on failure.
Definition: CD_ItoKMCStepperImplem.H:148
RefCountedPtr< ItoLayout< ItoSolver > > m_ito
Ito solvers.
Definition: CD_ItoKMCStepper.H:384
Vector< int > m_loadBalanceIndices
Solver indices used when load-balancing the particle solvers.
Definition: CD_ItoKMCStepper.H:644
EBAMRCellData m_particleEPS
For holding the mean particle energy.
Definition: CD_ItoKMCStepper.H:691
std::function< Real(const Real a_time)> m_voltage
Voltage curve on the electrodes used in the simulation.
Definition: CD_ItoKMCStepper.H:435
virtual void reconcileCdrDensities(const EBAMRCellData &a_newParticlesPerCell, const EBAMRCellData &a_oldParticlesPerCell, const Real a_dt) noexcept
Reconcile the CDR densities after the reaction network.
Definition: CD_ItoKMCStepperImplem.H:4476
virtual void computeConductivityCell(EBAMRCellData &a_conductivity) noexcept
Compute the cell-centered conductiivty.
Definition: CD_ItoKMCStepperImplem.H:1834
Real m_minParticleDiffusionCFL
Minium CFL-like time step for particle diffusion.
Definition: CD_ItoKMCStepper.H:538
virtual void postCheckpointPoisson() noexcept
Do some post-checkpoint operations for the electrostatic part.
Definition: CD_ItoKMCStepperImplem.H:747
virtual bool loadBalanceThisRealm(const std::string a_realm) const override
Load balancing query for a specified realm. If this returns true for a_realm, load balancing routines...
Definition: CD_ItoKMCStepperImplem.H:5488
virtual void prePlot() noexcept override
Perform pre-plot operations.
Definition: CD_ItoKMCStepperImplem.H:1515
virtual void postInitialize() noexcept override
Post-initialization operations. Default does nothing.
Definition: CD_ItoKMCStepperImplem.H:634
Real m_maxParticleAdvectionDiffusionCFL
Maximum CFL-like time step for particle advection-diffusion.
Definition: CD_ItoKMCStepper.H:553
virtual void coarsenCDRSolvers() noexcept
Coarsen data for CDR solvers.
Definition: CD_ItoKMCStepperImplem.H:4625
virtual void postCheckpointSetup() noexcept override
Perform post-checkpoint operations.
Definition: CD_ItoKMCStepperImplem.H:728
Real m_particleDiffusionDt
The particle diffusive time step.
Definition: CD_ItoKMCStepper.H:584
Vector< EBAMRCellData > m_energySources
Storage for holding the energy sources for each species.
Definition: CD_ItoKMCStepper.H:737
virtual void intersectParticles(const SpeciesSubset a_speciesSubset, const bool a_delete, const std::function< void(ItoParticle &)> a_nonDeletionModifier=[](ItoParticle &) -> void { return;}) noexcept
Intersect a subset of the particles with the domain and embedded boundary.
Definition: CD_ItoKMCStepperImplem.H:2033
Real m_maxGrowthDt
Maximum permissible time step growth.
Definition: CD_ItoKMCStepper.H:506
virtual void computeMobilities() noexcept
Compute mesh-based mobilities for LFA coupling.
Definition: CD_ItoKMCStepperImplem.H:2805
virtual void setItoVelocityFunctions() noexcept
Set the Ito velocity functions. This is sgn(charge) * E.
Definition: CD_ItoKMCStepperImplem.H:2689
Vector< EBAMRCellData > m_fluidGradPhiIto
For holding the gradient of all Ito species densities.
Definition: CD_ItoKMCStepper.H:627
bool m_loadBalanceParticles
Load balance particle realm or not.
Definition: CD_ItoKMCStepper.H:465
virtual void sortPhotonsByPatch(const McPhoto::WhichContainer a_which) noexcept
Sort photons by patch.
Definition: CD_ItoKMCStepperImplem.H:5439
virtual void getPhysicalParticlesPerCell(EBAMRCellData &a_ppc) const noexcept
Get the physical number of particles per cell.
Definition: CD_ItoKMCStepperImplem.H:3356
Base class for advancing equations.
Definition: CD_TimeStepper.H:30
Name containing various physics models for running chombo-discharge code.
Definition: CD_AdvectionDiffusion.H:15
phase names
Definition: CD_MultiFluidIndexSpace.H:27