4.12.1. MPMD Radiation
This section outlines the use of Multiple-Program-Multiple-Data (MPMD) coupling to perform radiation solves with Sierra/PMR. MPMD coupling allows the radiation solve to occur on an entirely different set of processors using an executable other than Fuego. In the past, PMR was solved with Syrinx (compiled into the Fuego executable). That has been deprecated and removed in favor of MPMD coupling - first with Sierra/Nalu and Scefire, and now with Sierra/PMR.
4.12.1.1. Input Deck Conversion from Nalu to Sierra/PMR
The new Sierra/PMR capability maintains backwards compatibility with Nalu input decks, so conversion from Sierra/Nalu to Sierra/PMR for the PMR solve requires simply changing the name of the executable in your launch scripts from nalu to pmr.
4.12.1.2. Input Deck Conversion from Syrinx to Sierra/PMR
The Conversion of a Fuego/Syrinx input deck to a Fuego MPMD input deck is straightforward. The basic steps to convert a Syrinx case to a Fuego-PMR case are:
Delete the Syrinx region, PMR material, PMR solver, and PMR mesh (finite element) blocks from the Fuego input deck.
Delete all transfer blocks to or from the PMR region and all calls to them from solution control (see Transfers).
Add the
USE MPMD RADIATIONcommand to the Fuego input deck in the Fuego region block to activate the MPMD coupling.Add the
SKIP STEPS FOR PMRcommand in the Fuego procedure or time blocks if you want to run the PMR solve less frequently than every time step.Specify emissivity, transmissivity, radiation environment temperature, and radiation boundary temperature in the Fuego boundary condition blocks. If any of these are omitted defaults will be used. Default conditions are emissivity = 1, transmissivity = 0, and boundary temperature equal to the local or wall temperature depending on the boundary type.
Copy the template PMR input file from PMR Input Deck and update the name of the mesh file.
(Optional) Update the quadrature and numerical scheme in the PMR input deck per your problem requirements.
(Optional) If using a particle region, transfer absorption and radiation sources from the particle region to the Fuego region, and transfer scalar flux from the Fuego region to the particle region (see Particle Radiation Terms).
Run your case using an MPMD launch command (see Running MPMD Jobs).
4.12.1.2.1. Transfers
In prior Fuego-Syrinx simulations the user would specifically call out which fields to transfer between Fuego and Syrinx in the input file, or call out a pre-defined transfer that would automatically define fields to transfer. With MPMD PMR the transfers are always determined automatically, so nothing is required in the input deck with regards to transfers to or from the PMR code.
The Fuego to PMR transfers are:
Radiative source (
rad_source,). This is calculated by Fuego either from the temperature field or combustion model (EDC or mixture fraction) and is sent on all volume nodes. For the non-reacting laminar case, this term is
(4.36)
and for reacting cases it is either tabulated or otherwise modified to include sub-grid effects.
Absorption coefficient (
absorption,). Calculated using the current property model and sent on all volume nodes.
Emissivity. Calculated in Fuego and sent on all surface nodes.
Transmissivity. Calculated in Fuego and sent on all surface nodes.
Boundary Source. Calculated in Fuego using the emissivity (
), transmissivity (
), radiation boundary temperature (
), radiation environment temperature (
), and band fractions (
,
) if doing a spectral calculation. This is sent on all surface nodes.
(4.37)
Boundary Beam Source. Calculated in Fuego on all surface nodes only if the beam model is in use.
The PMR to Fuego transfers are:
Scalar flux (
scalar_flux,). Calculated by the PMR code and sent back on all volume nodes, this is used to apply the radiative source term to the enthalpy equation.
(4.38)
This source term is linearized using the laminar form of
to give
(4.39)
Radiative Flux (
radiative_flux). Calculated by the PMR code this flux vector is sent back on all volume nodes and is primarily used for post-processing.Incident Flux (
incident_flux). The incident flux on all surfaces is calculated the PMR code and sent back to Fuego.
4.12.1.2.2. PMR Input Deck
The corresponding PMR input deck for Fuego/PMR MPMD calculations is shown in the following sections. The PMR input file is YAML formatted, so spaces and indentation matter. The complete user manual for Sierra/PMR can be found in the PMR user manual at https://pmr.sandia.gov, so details about the input options will be omitted here for brevity.
simulation_settings:
mesh: mesh.g
method: edge_upwind
quadrature:
type: thurgood
order: 4
4.12.1.2.3. Running MPMD Jobs
The command to run in MPMD mode is different from what was used to run Fuego-Syrinx cases. To run Fuego-Syrinx cases you simply ran a parallel job of Fuego using something like
mpirun -np 100 fuego -i fire.i
In this case one executable would use 100 cores, so both Fuego and Syrinx were taking turns using the same CPU resources. With MPMD runs you are launching two separate MPI jobs with two different codes that can communicate. An example MPMD launch command would look like
mpirun -np 100 fuego -i fire.i : -np 100 pmr -i pmr.i
However, there is no requirement any more that the two codes use the same number of cores, so depending on the mesh and computational costs you may choose a different allocation per code. For example, if your PMR solve is very expensive you may allocate more cores to PMR than Fuego:
mpirun -np 50 fuego -i fire.i : -np 150 pmr -i pmr.i
Special care must be taken when submitting MPMD jobs on the HPCs or any queued environment. By default, the two MPMD codes cannot share cores so to launch the case above on an HPC you would need to request an allocation 200 cores. This is unnecessarily wasteful though since Sierra/PMR would not be using its 150 cores while Fuego runs, and Fuego would not be using its 50 cores while Sierra/PMR runs. To get around this, you must enable oversubscription. Since Fuego and Sierra/PMR run sequentially (never executing at the same time) you can allow them to share resources. To get an allocation of 100 cores and use them all for both codes you must add additional mpi flags:
mpiexec --oversubscribe --bind-to core:overload-allowed -np 100 \
fuego -i fire.i : --bind-to core:overload-allowed -np 100 pmr -i pmr.i
Keep in mind that the specific command to use can be platform dependent. A more complete example submission script on an HPC may look like
#!/bin/bash
#SBATCH --nodes=10
#SBATCH --time=48:00:00
#SBATCH --account=xxxxxxxx
#SBATCH --job-name=fire
#SBATCH --partition=batch
nodes=$SLURM_JOB_NUM_NODES
cores=36
module load sierra
export OMPI_MCA_rmaps_base_oversubscribe=1
mpiexec --oversubscribe \
--bind-to core:overload-allowed \
--npernode $cores --n $(($cores*$nodes)) fuego -i fire.i : \
--bind-to core:overload-allowed \
--npernode $cores --n $(($cores*$nodes)) pmr -i pmr.i
On vortex or Sierra, an explicit resource file is required to use MPMD. An example fuego_pmr.erf explicit resource file is
overlapping-rs : allow
oversubscribe-cpu : allow
app 0 : fuego -i fire.i
app 1 : pmr -i pmr.i
4 : {host: *; cpu: *; gpu: * } : app 0
4 : {host: *; cpu: *; gpu: * } : app 1
where both apps are set to use 4 ranks on each node and allowed to use both the CPU and GPU. Using this file, a submission script on vortex could look like
#!/bin/bash
#BSUB -nnodes 4
module load sierra
export CUDA_LAUNCH_BLOCKING=0
export CUDA_MANAGED_FORCE_DEVICE_ALLOC=1
export TPETRA_ASSUME_GPU_AWARE_MPI=1
rm -f stderr-b stdout-b
jsrun -M "-gpu" -k stderr-b -o stdout-b --erf_input=fuego_pmr.erf
Contact sierra-help@sandia.gov if you need more help or encounter issues running MPMD jobs.
4.12.1.3. Spectral Radiation Transport
Beta Capability
Spectral radiation transport is a beta feature, with limited testing.
Fuego now supports spectral radiation calculations through MPMD coupling with Sierra/PMR. Spectral radiation transport is instantiated in an MPMD Fuego input deck through a spectral model section similar to the following:
BEGIN RAD TRANSPORT SPECTRAL MODEL SPECIFICATION spectralRadModel
Spectral Band Model = LINEAR
Spectral Subband Model = AVERAGE
Minimum Frequency = 1.0e13
Maximum Frequency = 2.0e13
Number Spectral Bands = 2
END RAD TRANSPORT SPECTRAL MODEL SPECIFICATION
The spectral model utilized in Fuego is based on the concept of spectral bands, which means that the user must describe frequency bands over which radiation absorptivities are averaged. In this context, a richer set of data is available to the analyst, including radiation transport terms (banded absorptivity, scalar flux, and radiative flux) than is available through the standard gray (nonbanded) radiation transport description that has been available to Fuego analysts for many years.
Spectral transport makes use of spectral data files that the user must provide. Files are named as: Species-data.dat
Where “Species” is the name of the chemical species whose spectral transport properties are being described. Each species present (all species for EDC simulations) or those output through the use of output variables in mixture fraction simulations must have all necessary spectral data files available in the working directory when a spectral transport calculation is performed. Spectral data files have the following structure:
For each temperature of the M temperatures specified in the spectral data files, there are
frequencies (
through
) specified with corresponding absorptivity values (
through
).
Temperatures within these files are required to be in non-decreasing order and frequencies must be ordered in increasing value within a specific temperature.
Absorptivities are expected to be in units of
, which is a number density weighted absorptivity.
Unit conversion within the Fuego input deck is required if the desired units of the simulation are not cgs.
The ideal gas law and knowledge of the mass fraction of each species is used to calculate the number density.
For the spectral radiation model specified here, the name of the particular spectral model is given as spectralRadModel, but could be any valid string.
For this case, a LINEAR “Spectral Band Model” is used, where the set of spectral bands, in this case 2 as specified by “Number of Spectral Bands” are equally divided along the linear spectral frequency range with “Minimum Frequency = 1.0e13” Hz and “Maximum Frequency = 2.0e13” Hz.
One could alternatively choose a LOGARITHMIC Spectral Band Model where bands are divided up equally in logarithmic frequency space.
In this case, we have chosen an AVERAGE “Spectral Subband Model”, which indicates that averaging of absorptivities for frequencies within a spectral band will follow standard (unweighted) averaging over the particular spectral band.
A PLANCK_AVERAGE Spectral Subband Model is also available where integrations are weighted by their Planck black-body spectrum contribution.
PLANCK_AVERAGE_WITH_REFERENCE_TEMPERATURE is a final choice that a user can select where the band contributions are weighted by the Planck black body-spectrum at a specified user reference temperature.
If this final option is used, a user must also specify a reference temperature, which can be done through:
Planck Subband Model Reference Temperature = 2000.0
Here the Reference temperature is set to 2000.0 K, but could be any nonzero temperature.
4.12.1.3.1. Spectral properties on surfaces
Beta Capability
Spectral surface properties is a beta feature, with limited testing.
Fuego can now perform spectral radiation calculations with banded emissivities defined on surfaces. This is a beta feature, and still requires additional testing. A command similar to the one shown below is used to define banded emissivities on surfaces:
BEGIN WALL BOUNDARY CONDITION ON SURFACE SURFACE_1
INTERFACE BOUNDARY
EMISSIVITY SPECTRAL FILE NAME = SURFACE_1_SPECTRAL.DATA
TRANSPARENT BAND EMISSIVITY = 1.0
END WALL BOUNDARY CONDITION ON SURFACE SURFACE_1
Similar to the spectral radiation model described in Section
Spectral Radiation Transport, a spectral data file containing the emissivity bands must
be provided for each surface (provided by the EMISSIVITY SPECTRAL FILE NAME = command in the above). These files have a structure similar to the species
spectral files:
Similar to the species spectral files, each temperature of the
temperatures specified in the spectral file have
frequencies defined,
each with a corresponding emissivity value (
) through
. Temperatures within these files must be defined in a
non-decreasing order, with frequencies within each temperature group ordered in
increasing value.
In addition to defining a spectral emissivity file, an additional emissivity for the transparent band can be defined as a general expression; this value is defaulted to 1. In the above example, the transparent band emissivity is defined as a constant (of 1).
4.12.1.4. Particle Radiation Terms
Particle Radiation contributions are also still available in Fuego simulations with Lagrangian particle types that have radiation contributions (heated, CPD, evaporating, general chemistry, wildfire). Previously, the user would specify transfers of necessary thermal/radiation fields between Fuego and Syrinx, which would resemble the following:
BEGIN TRANSFER pmr_to_particle
COPY VOLUME NODES FROM pmr_region TO particle_region
SEND BLOCK block_1 TO block_1
SEND FIELD scalar_flux STATE none TO scalar_flux STATE new
END TRANSFER pmr_to_particle
BEGIN TRANSFER particle_to_pmr
COPY VOLUME NODES FROM particle_region to pmr_region
send field particle_absorption_coeff state new to particle_absorption_coeff state none
SEND FIELD particle_rte_source STATE new TO particle_rte_source STATE none
END TRANSFER particle_to_pmr
For MPMD Radiation transport simulations involving Lagrangian particles, the user now removes the PMR (Syrinx) to particle transfers and instead adds these to the already-present set of fluid to particle transfers as is seen in the following:
BEGIN TRANSFER fluid_to_particle
COPY VOLUME NODES FROM fluid_region TO particle_region
SEND BLOCK block_1 TO block_1
#...
SEND FIELD scalar_flux STATE none TO scalar_flux STATE new
END TRANSFER fluid_to_particle
BEGIN TRANSFER particle_to_fluid
COPY VOLUME NODES FROM particle_region TO fluid_region
SEND BLOCK block_1 TO block_1
#...
send field particle_absorption_coeff state new to particle_absorption_coeff state none
send field particle_rte_source state new to particle_rte_source state none
END TRANSFER particle_to_fluid
4.12.1.5. Beam (Directed Flux)
Fuego users are able to now specify beam or directed flux boundary conditions for use in MPMD simulations with PMR. A beam boundary condition is specified in the Fuego input deck in the following manner:
#CONVERGING/DIVERGING BEAM
BEGIN BEAM RADIATION BOUNDARY SPECIFICATION
BEAM RADIATION BOUNDARY MODEL = FOCUSED
FOCAL POINT = FX FY FZ
FOCAL POWER = FP
FOCUSED BEAM TYPE = DIVERGING (OR CONVERGING)
CONSTRAIN TO SURFACES = sideA sideB ...
BEAM TEMPERATURE = BT
END BEAM RADIATION BOUNDARY SPECIFICATION
In the above case, a focused beam boundary BEAM RADIATION BOUNDARY MODEL = FOCUSED is defined where the FOCAL POINT of the beam, either sourcing from or arriving at this point are (FX, FY, FZ).
The user specifies whether the beam diverges from this point or converges to this focal point.
For the DIVERGING beam, the focal point needs to lie outside the simulation domain, whereas for the CONVERGING case, the focal point could lie in or outside the domain.
The focal power is specified as FP which is the integrated power of the beam assuming the intensity is spherically isotropic.
For a converging case, the power should be properly scaled to respect the actual solid angle over which the beam converges to the focal point.
For instance if the converging beam actually occupies a solid angle of , the
FOCAL POWER should be scaled up by a factor of since if the beam were incident from all solid angles, the total beam power would increase by this factor.
Beams can be constrained to arrive through the domain through only a fixed set of external surfaces if so desired as indicated through
CONSTRAIN TO SURFACES = sideA sideB....
The user also specifies the BEAM TEMPERATURE for cases using spectral radiation transport, since the beam power must then be parceled out among the defined spectral bands.
Another “directed” or beam-type boundary condition is available, where the user specifies not a focused beam but rather a plane wave as is seen below:
#PLANE WAVE
BEGIN BEAM RADIATION BOUNDARY SPECIFICATION
BEAM RADIATION BOUNDARY MODEL = PLANE_WAVE
PLANE WAVE BEAM DIRECTION = WX WY WZ
PLANE WAVE INTENSITY = PWI
CONSTRAIN TO SURFACES = sideA sideB ...
BEAM TEMPERATURE = 3000.0
END BEAM RADIATION BOUNDARY SPECIFICATION
In this case, the BEAM RADIATION BOUNDARY MODEL = PLANE_WAVE, indicates a plane-wave type source is to be used.
The user also specifies the direction for the plane wave as PLANE WAVE BEAM DIRECTION = WX WY WZ.
The direction must be non-zero.
The intensity of the plane-wave in units of Power per area is defined through PLANE WAVE INTENSITY = PWI.
The ability to constrain the beam to only certain external surfaces is identical to that for the focused beam above as is the beam temperature.
One should keep in mind that for both beam types (focused and plane wave), beam vectors at external surfaces are only calculated for nodes on surfaces where the inner (dot) product of outwardly directed area vector and the beam flux vector is less than 0, indicating that a user specified beam must enter rather than exit the simulation domain. Of course, given the physics internal to the domain of simulation, a radiative flux vector can exit the simulation, but we disallow the possibility of a user specifying an outwardly directed beam that does not result from internal physics.