Structure containing the radiation module parameters. More...
#include <cs_rad_transfer.h>
Data Fields | |
cs_rad_transfer_model_t | type |
int | nrphas |
int | iimpar |
int | verbosity |
int | imgrey |
int | imoadf |
int | iwrp1t |
int | imfsck |
int | imrcfsk |
double | xnp1mx |
int | idiver |
int | i_quadrature |
int | ndirec |
int | ndirs |
cs_real_3_t * | vect_s |
cs_real_t * | angsol |
int | restart |
int | nwsgg |
cs_real_t * | wq |
int | itpimp |
int | ipgrno |
int | iprefl |
int | ifgrno |
int | ifrefl |
int | itpt1d |
int | ifinfe |
int | atmo_model |
int | atmo_dr_id |
int | atmo_dr_o3_id |
int | atmo_df_id |
int | atmo_df_o3_id |
int | atmo_ir_id |
bool | dispersion |
cs_real_t | dispersion_coeff |
cs_time_control_t | time_control |
Structure containing the radiation module parameters.
cs_real_t* angsol |
Weight of the solid angle.
int atmo_df_id |
Atmospheric radiation model: id of the Diffuse Solar band or -1 if not activated (automatically computed)
int atmo_df_o3_id |
Atmospheric radiation model: id of the Diffuse Solar O3 band (SUV) or -1 if not activated (automatically computed)
int atmo_dr_id |
Atmospheric radiation model: id of the Direct Solar band or -1 if not activated (automatically computed)
int atmo_dr_o3_id |
Atmospheric radiation model: id of the Direct Solar O3 band or -1 if not activated (automatically computed)
int atmo_ir_id |
Atmospheric radiation model: id of the InfraRed band or -1 if not activated (automatically computed)
int atmo_model |
Atmospheric radiation model:
bool dispersion |
add dispersion (through diffusion)
cs_real_t dispersion_coeff |
dispersion coefficient. The dispersion coefficient leading to the best precision may depend on the chosen quadrature, and has been observed to be 3 for 128 directions (T4) and 5 for 32 directions (T2) on a (cube with point source) test case; the default value of 1 already improves precision in both cases.
int i_quadrature |
Index of the quadrature and number of directions for a single octant.
Sn quadrature (n(n+2) directions)
int idiver |
Explicit radiative source term computation mode -1: no renormalization 0: Semi-analytic (mandatory if transparent) 1: Conservative 2: Corrected semi-analytic (to be conservative) REMARK: if transparent, idiver = -1 automatically in DOM
Indicates the method used to calculate the radiative source term:
int ifgrno |
int ifinfe |
Modeling of an infinite extrusion for open boundaries.
int ifrefl |
int iimpar |
Verbosity level in the log concerning the calculation of the wall temperatures:
int imfsck |
FSCK model (0: off, 1: on)
FSCK model:
int imgrey |
Absorption coefficient computation 0: without the Grey Body model 1: with the Grey Body model
When gas or coal combustion is activated, imgrey indicates whether the absorption coefficient shall be calculated "automatically" (=1) using Modak's model or read from the data file (=0).
int imoadf |
ADF model: 0: not used 1: with wavelength interval of a 8 2: with wavelength interval of a 50
ADF model:
int imrcfsk |
RCFSK model (0: off, 1: on)
RCFSK model:
int ipgrno |
int iprefl |
int itpimp |
int itpt1d |
See CS_BOUNDARY_RAD_WALL_GRAY_1D_T
int iwrp1t |
P1 model transparency warnings counter.
int ndirec |
Number of directions for the angular discretisation of the radiation propagation with the DOM model.
No other possible value, because of the way the directions are calculated.
The calculation with 32 directions may break the symmetry of physically axi-symmetric cases (but the cost in CPU time is much lower than with 128 directions).
Useful if and only if the radiation module is activated with the DOM method.
int ndirs |
For the Tn quadrature, ndirec squared
int nrphas |
Phase which radiates (bulk by default, but may be coal class or fuel droplets phase).
int nwsgg |
Spectral radiation models (ADF and FSCK).
Number of ETRs to solve.
int restart |
Indicates whether the radiation variables should be initialized or read from a restart file.
cs_time_control_t time_control |
Determines at which time steps the variables are updated Also, in order to have proper initialization of the variables, the radiation module should always be called at the first time step of a calculation (restart or not).
model activation and type
cs_real_3_t* vect_s |
Direction vectors of angular values of the quadrature sx, sy, sz.
int verbosity |
Radiance resolution verbosity
Verbosity level in the log concerning the calculation of the radiative transfer equation:
cs_real_t* wq |
Weights of the Gaussian quadrature
double xnp1mx |
For the P-1 model, percentage of cells for which we allow the optical thickness to exceed unity, although this should be avoided. (more precisely, where