Grey, absorbing-emitting, non-scattering Discrete Ordinates Method (DOM) radiation solver model.
[Physical Assumptions] - Grey (single wavenumber/band) absorbing-emitting medium. - Non-scattering; no scattering source is included in the transport equations.
[Boundary Condition Assumptions & Limitations] - Wall/open-boundary inflow currently assumes black diffuse behavior, epsilon = 1.0: I_in = (sigma/pi)*T_w^4. - Input emissivity values below 1 are detected and warned about, but diffuse-gray reflection is not implemented in this local independent-ordinate solve. - Symmetry boundaries are simplified/skipped. A true DOM symmetry condition should reflect the ordinate direction through the symmetry plane.
[X-1 Contract] - This module computes and populates only the unreduced local partial absorption-like term kappaG_local in source%qrad. - The radiation driver is responsible for performing the MPI reduction and then locally subtracting the isotropic emission term 4kappasigmaT^4 on all ranks.
This module solves the discretized finite-volume grey DOM transport equations: s_m . grad(I_m) + kappaI_m = kappaI_b where I_b = sigma/pi*T^4. Angular decomposition is used to distribute directions over MPI processes.
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| real(kind=rk), | private, | parameter | :: | pi | = | 3.141592653589793_rk | |
| real(kind=rk), | private, | parameter | :: | sigma | = | 5.670374419e-8_rk |
Stefan-Boltzmann constant [W/(m^2*K^4)] |
Returns first-octant normalized weights. Full-space weights are w_ijk*pi/2.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| character(len=*), | intent(in) | :: | quad_type | |||
| integer, | intent(in) | :: | i | |||
| integer, | intent(in) | :: | j | |||
| integer, | intent(in) | :: | k |
Main entry point for the DOM radiation model.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(mesh_t), | intent(in) | :: | mesh | |||
| type(bc_set_t), | intent(in) | :: | bc | |||
| type(radiation_context_t), | intent(in) | :: | context | |||
| type(radiation_state_t), | intent(in) | :: | state | |||
| type(radiation_source_t), | intent(inout) | :: | source | |||
| type(case_params_t), | intent(in) | :: | params |
Sparse matrix-vector product for the assembled DOM transport system.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(mesh_t), | intent(in) | :: | mesh | |||
| real(kind=rk), | intent(in) | :: | diag(:) | |||
| real(kind=rk), | intent(in) | :: | coeff(:,:) | |||
| integer, | intent(in) | :: | nb(:,:) | |||
| real(kind=rk), | intent(in) | :: | x(:) | |||
| real(kind=rk), | intent(out) | :: | ax(:) |
Symmetric level-style quadrature lookup/generator.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| character(len=*), | intent(in) | :: | quad_type | |||
| integer, | intent(out) | :: | n_dirs | |||
| real(kind=rk), | intent(out), | allocatable | :: | s_x(:) | ||
| real(kind=rk), | intent(out), | allocatable | :: | s_y(:) | ||
| real(kind=rk), | intent(out), | allocatable | :: | s_z(:) | ||
| real(kind=rk), | intent(out), | allocatable | :: | w(:) |
Jacobi-preconditioned BiCGStab solver for the unsymmetric DOM transport system.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(mesh_t), | intent(in) | :: | mesh | |||
| real(kind=rk), | intent(in) | :: | diag(:) | |||
| real(kind=rk), | intent(in) | :: | coeff(:,:) | |||
| integer, | intent(in) | :: | nb(:,:) | |||
| real(kind=rk), | intent(in) | :: | rhs(:) | |||
| real(kind=rk), | intent(inout) | :: | x(:) | |||
| integer, | intent(in) | :: | max_iter | |||
| real(kind=rk), | intent(in) | :: | tol |
Runtime quadrature consistency checks.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| integer, | intent(in) | :: | n_dirs | |||
| real(kind=rk), | intent(in) | :: | s_x(:) | |||
| real(kind=rk), | intent(in) | :: | s_y(:) | |||
| real(kind=rk), | intent(in) | :: | s_z(:) | |||
| real(kind=rk), | intent(in) | :: | w(:) |