chombo-discharge
CD_ItoKMCGodunovStepper.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_ItoKMCGodunovStepper_H
13 #define CD_ItoKMCGodunovStepper_H
14 
15 // Our includes
16 #include <CD_ItoKMCStepper.H>
17 #include <CD_Timer.H>
18 #include <CD_PointParticle.H>
19 #include <CD_NamespaceHeader.H>
20 
21 namespace Physics {
22  namespace ItoKMC {
23 
27  template <typename I = ItoSolver, typename C = CdrCTU, typename R = McPhoto, typename F = FieldSolverMultigrid>
28  class ItoKMCGodunovStepper : public ItoKMCStepper<I, C, R, F>
29  {
30  public:
35 
40  ItoKMCGodunovStepper(RefCountedPtr<ItoKMCPhysics>& a_physics);
41 
45  virtual ~ItoKMCGodunovStepper();
46 
51  virtual Real
52  advance(const Real a_dt) override;
53 
57  virtual void
58  allocate() noexcept override;
59 
63  virtual void
64  parseOptions() noexcept override;
65 
69  virtual void
70  parseRuntimeOptions() noexcept override;
71 
79  virtual void
80  preRegrid(const int a_lmin, const int a_oldFinestLevel) noexcept override;
81 
88  virtual void
89  regrid(const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) noexcept override;
90 
91 #ifdef CH_USE_HDF5
97  virtual void
98  writeCheckpointData(HDF5Handle& a_handle, const int a_lvl) const noexcept override;
99 #endif
100 
101 #ifdef CH_USE_HDF5
106  virtual void
107  readCheckpointHeader(HDF5HeaderData& a_header) noexcept override;
108 #endif
109 
110 #ifdef CH_USE_HDF5
115  virtual void
116  writeCheckpointHeader(HDF5HeaderData& a_header) const noexcept override;
117 #endif
118 
119 #ifdef CH_USE_HDF5
125  virtual void
126  readCheckpointData(HDF5Handle& a_handle, const int a_lvl) noexcept override;
127 #endif
128 
132  virtual void
133  postPlot() noexcept override;
134 
135  protected:
139  enum class WhichAlgorithm
140  {
141  EulerMaruyama,
142  };
143 
148 
153 
159 
164 
169 
174 
179 
184  Vector<RefCountedPtr<ParticleContainer<PointParticle>>> m_conductivityParticles;
185 
189  Vector<RefCountedPtr<ParticleContainer<PointParticle>>> m_irregularParticles;
190 
195  Vector<RefCountedPtr<ParticleContainer<PointParticle>>> m_rhoDaggerParticles;
196 
203 
210 
216 
222 
226  Vector<EBAMRCellData> m_cdrDivD;
227 
231  virtual void
232  allocateInternals() noexcept override;
233 
237  virtual void
238  parseAlgorithm() noexcept;
239 
243  virtual void
244  parseCheckpointParticles() noexcept;
245 
249  virtual void
250  setOldPositions() noexcept;
251 
256  virtual void
257  barrier() const noexcept;
258 
264  virtual void
265  remapPointParticles(Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_particles,
266  const SpeciesSubset a_subset) noexcept;
267 
274  virtual void
275  depositPointParticles(const Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_particles,
276  const SpeciesSubset a_subset) noexcept;
277 
283  virtual void
284  clearPointParticles(const Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_particles,
285  const SpeciesSubset a_subset) noexcept;
286 
293  virtual void
294  computeConductivities(const Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_particles) noexcept;
295 
302  virtual void
303  computeCellConductivity(EBAMRCellData& a_conductivityCell,
304  const Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_particles) noexcept;
305 
310  virtual void
311  computeFaceConductivity() noexcept;
312 
316  virtual void
317  computeSemiImplicitRho() noexcept;
318 
322  virtual void
323  setupSemiImplicitPoisson(const Real a_dt) noexcept;
324 
331  virtual void
332  removeCoveredPointParticles(Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_particles,
333  const EBRepresentation a_representation,
334  const Real a_tolerance) const noexcept;
335 
340  virtual void
342  Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_conductivityParticles) noexcept;
343 
348  virtual void
349  advanceEulerMaruyama(const Real a_dt) noexcept;
350 
356  virtual void
357  diffuseParticlesEulerMaruyama(Vector<RefCountedPtr<ParticleContainer<PointParticle>>>& a_rhoDaggerParticles,
358  const Real a_dt) noexcept;
359 
366  virtual void
367  computeDiffusionTermCDR(EBAMRCellData& m_semiImplicitRhoCDR, const Real a_dt) noexcept;
368 
375  virtual void
376  stepEulerMaruyamaParticles(const Real a_dt) noexcept;
377 
383  virtual void
384  stepEulerMaruyamaCDR(const Real a_dt) noexcept;
385 
389  virtual void
390  plotParticles() const noexcept;
391  };
392  } // namespace ItoKMC
393 } // namespace Physics
394 
395 #include <CD_NamespaceFooter.H>
396 
398 
399 #endif
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
Implementation of CD_ItoKMCGodunovStepper.H.
Declaration of an abstract class for integrating the Ito-KMC-Poisson equations.
SpeciesSubset
Enum for differentiating between types of particles.
Definition: CD_ItoKMCStepper.H:40
Declaration of a computational point particle.
Implementation of CD_Timer.H.
Templated class for holding particles on an AMR hierarchy with particle remapping.
Definition: CD_ParticleContainer.H:50
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:1113
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:470
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:1388
virtual void allocate() noexcept override
Allocate storage required for advancing the equations.
Definition: CD_ItoKMCGodunovStepperImplem.H:59
bool m_readCheckpointParticles
If true, then the HDF5 checkpoint file contained particles that we can read.
Definition: CD_ItoKMCGodunovStepper.H:152
virtual Real advance(const Real a_dt) override
Advance the Ito-Poisson-KMC system over a_dt.
Definition: CD_ItoKMCGodunovStepperImplem.H:198
virtual void allocateInternals() noexcept override
Allocate "internal" storage.
Definition: CD_ItoKMCGodunovStepperImplem.H:93
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:1188
virtual void stepEulerMaruyamaCDR(const Real a_dt) noexcept
Step the CDR equations according to the regular Euler-Maruyama scheme.
Definition: CD_ItoKMCGodunovStepperImplem.H:1474
WhichAlgorithm m_algorithm
Which advancement algorithm to use.
Definition: CD_ItoKMCGodunovStepper.H:173
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:752
EBAMRCellData m_semiImplicitRhoCDR
Storage for CDR densities used during the semi-implicit solve.
Definition: CD_ItoKMCGodunovStepper.H:202
virtual void parseAlgorithm() noexcept
Parse advancement algorithm.
Definition: CD_ItoKMCGodunovStepperImplem.H:159
WhichAlgorithm
Simple enum for distinguishing between algorithms.
Definition: CD_ItoKMCGodunovStepper.H:140
virtual void postPlot() noexcept override
Perform post-plot operations.
Definition: CD_ItoKMCGodunovStepperImplem.H:1627
Vector< RefCountedPtr< ParticleContainer< PointParticle > > > m_irregularParticles
Storage for particles that fell inside the EB but should still contribute to the conductivity.
Definition: CD_ItoKMCGodunovStepper.H:189
Vector< EBAMRCellData > m_cdrDivD
Storage for the finite-volume approximation of div(D*grad(phi)) for the CDR equations.
Definition: CD_ItoKMCGodunovStepper.H:226
EBAMRCellData m_scratchSemiImplicitRhoCDR
Scratch storage for CDR contribution to space charge density.
Definition: CD_ItoKMCGodunovStepper.H:215
virtual void parseRuntimeOptions() noexcept override
Parse run-time options.
Definition: CD_ItoKMCGodunovStepperImplem.H:144
virtual void remapPointParticles(Vector< RefCountedPtr< ParticleContainer< PointParticle >>> &a_particles, const SpeciesSubset a_subset) noexcept
Remap the input point particles.
Definition: CD_ItoKMCGodunovStepperImplem.H:646
bool m_canRegridOnRestart
If true, then the class supports regrid-on-restart.
Definition: CD_ItoKMCGodunovStepper.H:158
virtual void computeFaceConductivity() noexcept
Compute the cell-centered conductivity.
Definition: CD_ItoKMCGodunovStepperImplem.H:1046
EBAMRCellData m_semiImplicitConductivityCDR
Storage for conductivity term due to mobile CDR species.
Definition: CD_ItoKMCGodunovStepper.H:209
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:854
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:960
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:163
virtual void setOldPositions() noexcept
Set the starting positions for the ItoSolver particles.
Definition: CD_ItoKMCGodunovStepperImplem.H:610
bool m_smoothConductivity
Smooth the conductivity or not.
Definition: CD_ItoKMCGodunovStepper.H:168
Vector< RefCountedPtr< ParticleContainer< PointParticle > > > m_rhoDaggerParticles
Storage for particles that gave rho^dagger.
Definition: CD_ItoKMCGodunovStepper.H:195
bool m_writeCheckpointParticles
If true, then the particles are checkpointed so we can regrid on checkpoint-restart.
Definition: CD_ItoKMCGodunovStepper.H:147
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:1148
virtual void computeSemiImplicitRho() noexcept
Set up the space charge density for the regrid operation.
Definition: CD_ItoKMCGodunovStepperImplem.H:1076
virtual void computeCellConductivity(EBAMRCellData &a_conductivityCell, const Vector< RefCountedPtr< ParticleContainer< PointParticle >>> &a_particles) noexcept
Compute the cell-centered conductivity.
Definition: CD_ItoKMCGodunovStepperImplem.H:974
virtual void advanceEulerMaruyama(const Real a_dt) noexcept
Advance the particles using the Euler-Maruyama scheme.
Definition: CD_ItoKMCGodunovStepperImplem.H:1242
virtual void preRegrid(const int a_lmin, const int a_oldFinestLevel) noexcept override
Perform pre-regrid operations.
Definition: CD_ItoKMCGodunovStepperImplem.H:431
EBAMRCellData m_scratchSemiImplicitConductivityCDR
Scratch storage for CDR contribution to conductivity.
Definition: CD_ItoKMCGodunovStepper.H:221
virtual void plotParticles() const noexcept
Utility function for plotting the ItoSolver particles. These are written in a particles folder.
Definition: CD_ItoKMCGodunovStepperImplem.H:1641
ItoKMCGodunovStepper()=delete
Disallowed default constructor. Use the full constructor.
Vector< RefCountedPtr< ParticleContainer< PointParticle > > > m_conductivityParticles
Storage for simplified particles that gave us sigma^k.
Definition: CD_ItoKMCGodunovStepper.H:184
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:1328
Timer m_timer
Timer used for run-time logging of routines.
Definition: CD_ItoKMCGodunovStepper.H:178
virtual void parseCheckpointParticles() noexcept
Parse checkpoint-restart functionality.
Definition: CD_ItoKMCGodunovStepperImplem.H:184
virtual void barrier() const noexcept
Set an MPI barrier if using debug mode.
Definition: CD_ItoKMCGodunovStepperImplem.H:115
virtual void parseOptions() noexcept override
Parse options.
Definition: CD_ItoKMCGodunovStepperImplem.H:129
virtual void stepEulerMaruyamaParticles(const Real a_dt) noexcept
Step the particles according to the regular Euler-Maruyama scheme.
Definition: CD_ItoKMCGodunovStepperImplem.H:1428
Base time stepper class that advances the Ito-KMC-Poisson system of equations. If you want a differen...
Definition: CD_ItoKMCStepper.H:60
A particle class that only has a position and a weight.
Definition: CD_PointParticle.H:29
Class which is used for run-time monitoring of events.
Definition: CD_Timer.H:31
Name containing various physics models for running chombo-discharge code.
Definition: CD_AdvectionDiffusion.H:15