-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathflameSolver.h
executable file
·285 lines (244 loc) · 9.77 KB
/
flameSolver.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
#pragma once
#include "readConfig.h"
#include "chemistry0d.h"
#include "grid.h"
#include "splitSolver.h"
#include "readConfig.h"
#include "perfTimer.h"
#include "integrator.h"
#include "sourceSystem.h"
#include "diffusionSystem.h"
#include "convectionSystem.h"
#include "quasi2d.h"
#include "callback.h"
#include <iostream>
#include <boost/ptr_container/ptr_vector.hpp>
#include <boost/shared_ptr.hpp>
#include "tbb_tools.h"
using std::string;
class ScalarFunction;
//! Class which manages the main integration loop.
//! Contains the split solvers and is responsible for the large-scale time integration.
class FlameSolver : public GridBased, public SplitSolver
{
public:
FlameSolver();
virtual ~FlameSolver();
double P; // GET_RHO_U
double XMDOT;// GET_RHO_U
double XMDT;
int GFAC;
double PDFA;
double PDFB;
//double dt;
int NTS;
double NFL;
int NSIM;
int NSIM_counter=1;
int NTSPSIM;
int NTS_COUNT=1;
double XNU;
double DOM;
double C_lambda;
double Rate;
int NTS_PE;
int NC;
int NCM1;
int print_counter=0;
double flamevelocity,unitmovement;
double movement=0.0;
int NCP1;
int movecount;
int moveflag=0;
double DX;
double TAU;
double XLint;
double XLk;
double Re;
int MTS;
// random number
double random;
// eddy size
int L;
// starting point of triplet map
int M;
int init_flag=1;
int Check_flag=1;
int Print_flag=0;
double Print_counter = 0;
int error_flag=0;
double check_velocity;
int target1,target2;
int sizedatacounter=0;
struct Config
{
double endtime;
double timestep;
double Re_t;
double dom;// cm
double pressure; // dynes/cm2
double velocity; // cm/s
double Th; // K
double cp; // erg/g-K
double kinematic_viscosity; // cm2/s
double D,lambda,r_datas,trip_map,GFAC,FAL,Intlength,Neta,Clambda,Sl;
};
double RHOM,RHOP,SUMYK,SUMX,TDOT,XMDXM;
void VolumeExpansion();
void ReadParameters(Config& config);
void INIT_AllParameters();
void DiffusionVelocityCalculator();
void Random_Number();
void eddyLength();
void BTriplet(double var[]);
void TM();
void CFUEL();
void XRecord();
void PREMIXADV();
void velocityCalculator();
void Find2NearPoints(int local);
void FlamePositionCorrection();
void setOptions(const ConfigOptions& options); //!< Set options read from the configuration file
void initialize(); //!< call to generate profiles and perform one-time setup
void finalize();
//! Load initial temperature, mass fraction and velocity profiles.
//! Profiles are loaded either from an HDF5 "restart" file or from the
//! the ConfigOptions object, populated from the Python "Config" object.
virtual void loadProfile();
//! Determine whether to terminate integration based on reaching steady-state solution.
bool checkTerminationCondition(void);
void writeStateFile(const std::string& fileName="", bool errorFile=false, bool updateDerivatives=true);
void saveTimeSeriesData(const std::string& filename, bool write); //!< Collect info for out.h5 file
//! Compute integral heat release rate [W/m^2].
//! Assumes that #qDot has already been updated.
double getHeatReleaseRate(void);
//! Compute flame consumption speed [m/s]. Computed by integrating the energy equation:
//! \f[ S_c = \frac{\int_{-\infty}^\infty \dot{q}/c_p}{\rho_u(T_b-T_u)} \f]
//! Assumes that #qDot has already been updated.
double getConsumptionSpeed(void);
//! Compute position of flame (heat release centroid) [m]
double getFlamePosition(void);
// Time-series data
dvector timeVector; //!< Time [s] at the end of every ConfigOptions::outputStepInterval timesteps.
dvector heatReleaseRate; //!< Integral heat release rate [W/m^2] at the times in #timeVector.
long int nTotal; //!< total number of timesteps taken
int nRegrid; //!< number of time steps since regridding/adaptation
int nOutput; //!< number of time steps since storing integral flame parameters
int nProfile; //!< number of time steps since saving flame profiles
int nTerminate; //!< number of steps since last checking termination condition
int nCurrentState; //!< number of time steps since profNow.h5 and out.h5 were written
double terminationCondition; //!< Current value to be compared with terminationTolerance
double tOutput; //!< time of next integral flame parameters output (this step)
double tRegrid; //!< time of next regridding
double tProfile; //!< time of next profile output
boost::ptr_vector<SourceSystem> sourceTerms; // One for each grid point
vector<bool> useCVODE; // At each grid point, a flag for which integrator to use
boost::ptr_vector<DiffusionSystem> diffusionTerms; // One for each species (plus T and U)
boost::ptr_vector<TridiagonalIntegrator> diffusionSolvers; // One for each state variable
//! System representing the convection terms and containing their integrators
ConvectionSystemSplit convectionSystem;
// Boundary values
double rhou, rhob, rhoLeft, rhoRight;
double Tleft, Tright;
dvec Yleft, Yright;
void setupStep();
void prepareIntegrators();
int finishStep();
void resizeAuxiliary(); //!< Handle resizing of data structures as grid size changes
void resizeMappedArrays(); //!< update data that shadows SplitSolver arrays
void updateCrossTerms(); //!< calculates values of cross-component terms: jSoret, sumcpj, and jCorr
void updateChemicalProperties(); //!< Update thermodynamic, transport, and kinetic properties
// void updateChemicalProperties(size_t j1, size_t j2); //!< Update thermodynamic, transport, and kinetic properties
void updateBC(); //!< Set boundary condition for left edge of domain
void calculateQdot(); //!< Compute heat release rate using the current temperature and mass fractions
//! Correct the drift of the total mass fractions and reset any negative mass fractions.
void correctMassFractions();
void setDiffusionSolverState(double tInitial);
void setConvectionSolverState(double tInitial);
void setProductionSolverState(double tInitial);
// Steps in the Strang split integration process
void integrateConvectionTerms();
void integrateProductionTerms();
void integrateDiffusionTerms();
void integrateProductionTerms(size_t j1, size_t j2);
void integrateDiffusionTerms(size_t k1, size_t k2);
size_t nSpec; //!< Number of chemical species
size_t nVars; //!< Number of state variables at each grid point (nSpec + 2)
// State variables:
VecMap U; //!< normalized tangential velocity (u*a/u_inf) [1/s]
VecMap T; //!< temperature [K]
MatrixMap Y; //!< species mass fractions, Y(k,j) [-]
// Auxiliary variables:
dvec rho; //!< density [kg/m^3]
dvec rho_old;
dvec T_old;
dvec T_transient;
dvec DCvolume;
dvec position;
dvec U_velocity;// velocity matrix
dvec uniformGrid;
dvec drhodt; //!< time derivative of density [kg/m^3*s]
dvec jCorr; //!< Correction to ensure sum of mass fractions = 1
dvec sumcpj; //!< part of the enthalpy flux term
dvec qDot; //!< Heat release rate [W/m^3]
dmatrix wDot; //!< species production rates [kmol/m^3*s]
dvec Wmx; //!< mixture molecular weight [kg/kmol]
dvec W; //!< species molecular weights [kg/kmol]
dvec mu; //!< dynamic viscosity [Pa*s]
dvec sizedata;
dvec lambda; //!< thermal conductivity [W/m*K]
dvec cp; //!< mixture heat capacity [J/kg*K]
dmatrix cpSpec; //!< species molar heat capacities [J/kmol*K]
dmatrix rhoD; //!< density * diffusivity [kg/m*s]
dmatrix Dkt; //!< thermal diffusivity
dmatrix Dkm;
dmatrix hk; //!< species molar enthalpies [J/kmol]
dmatrix jFick; //!< Fickian mass flux [kg/m^2*s]
dmatrix jSoret; //!< Soret mass flux [kg/m^2*s]
dmatrix dVel;
dmatrix Y_old;
dmatrix Y_transient;
dmatrix YB;
dvec TB;
dmatrix F;
// jCorr is a correction to force the net diffusion mass flux to be zero
// jCorrSystem / jCorrSolver are used to introduce numerical diffusion into
// jCorr to eliminate spatial instabilities
DiffusionSystem jCorrSystem;
TridiagonalIntegrator jCorrSolver;
//! Function which describes strain rate a(t) and its derivative
ScalarFunction* strainfunc;
//! Function which describes a multiplier for the chemical reaction rates
ScalarFunction* rateMultiplierFunction;
LoggerCallback* stateWriter;
LoggerCallback* timeseriesWriter;
IntegratorCallback* heatLossFunction;
double rVcenter; //!< mass flux at centerline [kg/m^2 or kg/m*rad*s]
double rVzero; //!< mass flux at j=0
double tFlamePrev, tFlameNext;
double xFlameTarget, xFlameActual;
double flamePosIntegralError;
//! Cantera data
CanteraGas gas;
tbb::enumerable_thread_specific<CanteraGas> gases;
tbb::mutex gasInitMutex;
tbb::task_scheduler_init tbbTaskSched;
void rollVectorVector(vector<dvector>& vv, const dmatrix& M) const;
void unrollVectorVector(vector<dvector>& vv, dmatrix& M, size_t i) const;
void update_xStag(const double t, const bool updateIntError);
double targetFlamePosition(double t); //!< [m]
void printPerformanceStats(void);
void printPerfString(std::ostream& stats, const std::string& label, const PerfTimer& T);
//! Data for solving quasi-2D method-of-lines problems
boost::shared_ptr<BilinearInterpolator> vzInterp, vrInterp, TInterp;
// Performance Timers
// Just the total time:
PerfTimer totalTimer;
// These add up to the total run time:
PerfTimer setupTimer, reactionTimer, diffusionTimer,
convectionTimer, regridTimer;
// These account for special parts of the code
PerfTimer reactionRatesTimer, transportTimer, thermoTimer;
PerfTimer jacobianTimer;
PerfTimer conductivityTimer, viscosityTimer, diffusivityTimer;
};