chombo-discharge
Loading...
Searching...
No Matches
CD_RadiativeTransferStepperImplem.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// Std includes
13#include <math.h>
14
15// Chombo includes
16#include <CH_Timer.H>
17#include <ParmParse.H>
18#include <PolyGeom.H>
19
20// Our includes
23#include <CD_NamespaceHeader.H>
24
25using namespace Physics::RadiativeTransfer;
26
27template <typename T>
29{
30 CH_TIME("RadiativeTransferStepper::RadiativeTransferStepper");
31
32 // Get shit for setting the source term
33 ParmParse pp("RadiativeTransferStepper");
35
36 pp.get("verbosity", m_verbosity);
37 pp.get("realm", m_realm);
38 pp.get("dt", m_dt);
39
40 m_forceDt = -1.0;
41}
42
43template <typename T>
45{
46 CH_TIME("RadiativeTransferStepper::~RadiativeTransferStepper");
47 if (m_verbosity > 5) {
48 pout() << "RadiativeTransferStepper::~RadiativeTransferStepper" << endl;
49 }
50}
51
52template <typename T>
53void
55{
56 CH_TIME("RadiativeTransferStepper::setupSolvers");
57 if (m_verbosity > 5) {
58 pout() << "RadiativeTransferStepper::setupSolvers" << endl;
59 }
60
62
63 // Solver setup
64 m_solver = RefCountedPtr<RtSolver>(new T());
65 m_solver->setVerbosity(m_verbosity);
66 m_solver->setRtSpecies(m_species);
67 m_solver->parseOptions();
68 m_solver->setPhase(phase::gas);
69 m_solver->setAmr(m_amr);
70 m_solver->setComputationalGeometry(m_computationalGeometry);
71 m_solver->setRealm(m_realm);
72}
73
74template <typename T>
75void
77{
78 CH_TIME("RadiativeTransferStepper::parseRuntimeOptions");
79 if (m_verbosity > 5) {
80 pout() << "RadiativeTransferStepper::parseRuntimeOptions" << endl;
81 }
82
83 ParmParse pp("RadiativeTransferStepper");
84
85 pp.get("verbosity", m_verbosity);
86 pp.get("dt", m_dt);
87
88 m_solver->parseRuntimeOptions();
89}
90
91template <typename T>
92void
94{
95 CH_TIME("RadiativeTransferStepper::allocate");
96 if (m_verbosity > 5) {
97 pout() << "RadiativeTransferStepper::allocate" << endl;
98 }
99
100 m_solver->allocate();
101}
102
103template <typename T>
104void
106{
107 CH_TIME("RadiativeTransferStepper::initialData");
108 if (m_verbosity > 5) {
109 pout() << "RadiativeTransferStepper::initialData" << endl;
110 }
111
112 m_solver->initialData();
113 this->setGaussianSource();
114
115 if (m_solver->isStationary()) {
116 m_solver->advance(0.0, false);
117 }
118}
119
120template <typename T>
121void
123{
124 CH_TIME("RadiativeTransferStepper::postInitialize");
125 if (m_verbosity > 5) {
126 pout() << "RadiativeTransferStepper::postInitialize" << endl;
127 }
128}
129
130template <typename T>
131void
133{
134 CH_TIME("RadiativeTransferStepper::forceDt");
135 if (m_verbosity > 5) {
136 pout() << "RadiativeTransferStepper::forceDt" << endl;
137 }
138
139 m_forceDt = a_dt;
140}
141
142#ifdef CH_USE_HDF5
143template <typename T>
144void
146{
147 CH_TIME("RadiativeTransferStepper::writeCheckpointData");
148 if (m_verbosity > 5) {
149 pout() << "RadiativeTransferStepper::writeCheckpointData" << endl;
150 }
151
152 m_solver->writeCheckpointLevel(a_handle, a_lvl);
153}
154#endif
155
156#ifdef CH_USE_HDF5
157template <typename T>
158void
160{
161 CH_TIME("RadiativeTransferStepper::readCheckpointData");
162 if (m_verbosity > 5) {
163 pout() << "RadiativeTransferStepper::readCheckpointData" << endl;
164 }
165
166 m_solver->readCheckpointLevel(a_handle, a_lvl);
167}
168#endif
169
170template <typename T>
171void
173{
174 CH_TIME("RadiativeTransferStepper::postCheckpointSetup");
175 if (m_verbosity > 5) {
176 pout() << "RadiativeTransferStepper::postCheckpointSetup" << endl;
177 }
178
179 this->setGaussianSource();
180}
181
182template <typename T>
183int
185{
186 CH_TIME("RadiativeTransferStepper::getNumberOfPlotVariables");
187 if (m_verbosity > 5) {
188 pout() << "RadiativeTransferStepper::getNumberOfPlotVariables" << endl;
189 }
190
191 return m_solver->getNumberOfPlotVariables();
192}
193
194template <typename T>
197{
198 CH_TIME("RadiativeTransferStepper::getPlotVariableNames");
199 if (m_verbosity > 5) {
200 pout() << "RadiativeTransferStepper::getPlotVariableNames" << endl;
201 }
202
203 return m_solver->getPlotVariableNames();
204}
205
206template <typename T>
207void
209 int& a_icomp,
211 const int a_level) const
212{
213 CH_TIME("RadiativeTransferStepper::writePlotData");
214 if (m_verbosity > 5) {
215 pout() << "RadiativeTransferStepper::writePlotData" << endl;
216 }
217
218 m_solver->writePlotData(a_output, a_icomp, a_outputRealm, a_level);
219}
220
221template <typename T>
222Real
224{
225 CH_TIME("RadiativeTransferStepper::computeDt");
226 if (m_verbosity > 5) {
227 pout() << "RadiativeTransferStepper::computeDt" << endl;
228 }
229
230 Real dt = m_dt;
231
232 if (m_forceDt > 0.0) {
233 dt = m_forceDt;
234 }
235
236 return dt;
237}
238
239template <typename T>
240Real
242{
243 CH_TIME("RadiativeTransferStepper::advance");
244 if (m_verbosity > 5) {
245 pout() << "RadiativeTransferStepper::advance" << endl;
246 }
247
248 m_solver->advance(a_dt);
249
250 return a_dt;
251}
252
253template <typename T>
254void
256{
257 CH_TIME("RadiativeTransferStepper::synchronizeSolverTimes");
258 if (m_verbosity > 5) {
259 pout() << "RadiativeTransferStepper::synchronizeSolverTimes" << endl;
260 }
261
263 m_time = a_time;
264 m_dt = a_dt;
265
266 m_solver->setTime(m_timeStep, m_time, m_dt);
267}
268
269template <typename T>
270void
272{
273 CH_TIME("RadiativeTransferStepper::printStepReport");
274 if (m_verbosity > 5) {
275 pout() << "RadiativeTransferStepper::printStepReport" << endl;
276 }
277}
278
279template <typename T>
280void
282{
283 CH_TIME("RadiativeTransferStepper::registerRealms");
284 if (m_verbosity > 5) {
285 pout() << "RadiativeTransferStepper::registerRealms" << endl;
286 }
287
288 m_amr->registerRealm(m_realm);
289}
290
291template <typename T>
292void
294{
295 CH_TIME("RadiativeTransferStepper::registerOperators");
296 if (m_verbosity > 5) {
297 pout() << "RadiativeTransferStepper::registerOperators" << endl;
298 }
299
300 m_solver->registerOperators();
301}
302
303template <typename T>
304void
306{
307 CH_TIME("RadiativeTransferStepper::preRegrid");
308 if (m_verbosity > 5) {
309 pout() << "RadiativeTransferStepper::preRegrid" << endl;
310 }
311
312 m_solver->preRegrid(a_base, a_oldFinestLevel);
313}
314
315template <typename T>
316void
318{
319 CH_TIME("RadiativeTransferStepper::regrid");
320 if (m_verbosity > 5) {
321 pout() << "RadiativeTransferStepper::regrid" << endl;
322 }
323
324 m_solver->regrid(a_lmin, a_oldFinestLevel, a_newFinestLevel);
325}
326
327template <typename T>
328void
330{
331 CH_TIME("RadiativeTransferStepper::postRegrid");
332 if (m_verbosity > 5) {
333 pout() << "RadiativeTransferStepper::postRegrid" << endl;
334 }
335
336 this->setGaussianSource();
337}
338
339template <typename T>
340const EBAMRCellData&
342{
343 CH_TIME("RadiativeTransferStepper::getPhi");
344 if (m_verbosity > 5) {
345 pout() << "RadiativeTransferStepper::getPhi" << endl;
346 }
347
348 return m_solver->getPhi();
349}
350
351template <typename T>
352void
354{
355 CH_TIME("RadiativeTransferStepper::setGaussianSource");
356 if (m_verbosity > 5) {
357 pout() << "RadiativeTransferStepper::setGaussianSource" << endl;
358 }
359
360 ParmParse pp("RadiativeTransferStepper");
361
366
367 pp.get("blob_amplitude", blobAmp);
368 pp.get("blob_radius", blobRad);
369 pp.getarr("blob_center", v, 0, SpaceDim);
370 blobCenter = RealVect(D_DECL(v[0], v[1], v[2]));
371
374 const Real dist2 = dist.dotProduct(dist);
375 return blobAmp * exp(-dist2 / (2 * blobRad * blobRad));
376 };
377
378 m_solver->setSource(sourceFunction);
379}
380
381#include <CD_NamespaceFooter.H>
Simple species for Physics/RadiativeTransfer.
TimeStepper class for only evolving radiative transfer modules.
Implementation of RtSpecies for usage in RadiativeTransfer module.
Definition CD_RadiativeTransferSpecies.H:28
Implementation of TimeStepper for solving for a single radiative transfer species....
Definition CD_RadiativeTransferStepper.H:33
void registerOperators() override
Register operators – calls the solver registration routine.
Definition CD_RadiativeTransferStepperImplem.H:293
void postCheckpointSetup() override
Perform post-checkpoint setup routines (sets the source in the solver)
Definition CD_RadiativeTransferStepperImplem.H:172
void registerRealms() override
Register realms.
Definition CD_RadiativeTransferStepperImplem.H:281
void regrid(const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) override
Regrid function. Calls the solver function.
Definition CD_RadiativeTransferStepperImplem.H:317
void allocate() override
Allocate necessary memory for solvers.
Definition CD_RadiativeTransferStepperImplem.H:93
RadiativeTransferStepper()
Constructor. Reads a couple of input options.
Definition CD_RadiativeTransferStepperImplem.H:28
void writePlotData(LevelData< EBCellFAB > &a_output, int &a_icomp, const std::string a_outputRealm, const int a_level) const override
Write plot data to output holder.
Definition CD_RadiativeTransferStepperImplem.H:208
void postInitialize() override
Post-initialization functionality – sets the source term.
Definition CD_RadiativeTransferStepperImplem.H:122
Real advance(const Real a_dt) override
Advancement method. Calls the solver function.
Definition CD_RadiativeTransferStepperImplem.H:241
Vector< std::string > getPlotVariableNames() const override
Get plot variable names.
Definition CD_RadiativeTransferStepperImplem.H:196
int getNumberOfPlotVariables() const override
Get number of plot variables for this physics module.
Definition CD_RadiativeTransferStepperImplem.H:184
void preRegrid(const int a_lmin, const int a_oldFinestLevel) override
Perform pre-regrid operations – calls the solver function. operation takes a lot of memory.
Definition CD_RadiativeTransferStepperImplem.H:305
void parseRuntimeOptions() override
Parse runtime options.
Definition CD_RadiativeTransferStepperImplem.H:76
Real computeDt() override
Compute a time step to be used by Driver.
Definition CD_RadiativeTransferStepperImplem.H:223
void setGaussianSource()
For setting a Gaussian source in the radiative transfer equation solver.
Definition CD_RadiativeTransferStepperImplem.H:353
void forceDt(const Real a_dt)
Force usage of a time step.
Definition CD_RadiativeTransferStepperImplem.H:132
void printStepReport() override
Print a step report (does nothing)
Definition CD_RadiativeTransferStepperImplem.H:271
void setupSolvers() override
Set up the solver.
Definition CD_RadiativeTransferStepperImplem.H:54
void postRegrid() override
Post-regrid function. Sets a Gaussian source.
Definition CD_RadiativeTransferStepperImplem.H:329
void initialData() override
Fill simulation with initial data.
Definition CD_RadiativeTransferStepperImplem.H:105
void synchronizeSolverTimes(const int a_step, const Real a_time, const Real a_dt) override
Synchronzie solver times and time steps.
Definition CD_RadiativeTransferStepperImplem.H:255
virtual ~RadiativeTransferStepper()
Destructor (does nothing)
Definition CD_RadiativeTransferStepperImplem.H:44
const EBAMRCellData & getPhi() const
Get the solver solution.
Definition CD_RadiativeTransferStepperImplem.H:341
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
TracerParticleSolver()
Default constructor.
Definition CD_TracerParticleSolverImplem.H:25
virtual void setTime(const int a_step, const Real a_time, const Real a_dt)
Set the time for this solver.
Definition CD_TracerParticleSolverImplem.H:256
std::string m_realm
Realm where this solver lives.
Definition CD_TracerParticleSolver.H:345
int m_verbosity
Verbosity level.
Definition CD_TracerParticleSolver.H:380
RefCountedPtr< AmrMesh > m_amr
Handle to AMR mesh.
Definition CD_TracerParticleSolver.H:320
RefCountedPtr< ComputationalGeometry > m_computationalGeometry
Handle to computational geometry.
Definition CD_TracerParticleSolver.H:325
virtual void parseRuntimeOptions()
Parse solver run-time options.
Definition CD_TracerParticleSolverImplem.H:76
Real m_dt
Time step.
Definition CD_TracerParticleSolver.H:365
int m_timeStep
Time step.
Definition CD_TracerParticleSolver.H:375
Real m_time
Time.
Definition CD_TracerParticleSolver.H:370
Namespace for encapsulating the radiative transfer physics module.
Definition CD_RadiativeTransferSpecies.H:21