The Net Radiation Method System

The GrayDiffuseRadiation syntax invokes the Radiation Transfer Action. It simplifies setting up gray, diffuse radiative exchange problems. It requires the user to provide boundaries, emissivities on these boundaries, and the type of boundary: ADIABATIC, FIXED_TEMPERATURE, VARIABLE_TEMPERATURE; VARIABLE_TEMPERATURE are all boundaries that are not listed in the union of "adiabatic_boundary" and "fixed_temperature_boundary".

The GrayDiffuseRadiation syntax allows to split the boundaries into radiation patches using the "n_patches" parameter. This allows the computation of more accurate view factors between participating boundaries.

Assumed geometry

The GrayDiffuseRadiation syntax assumes that radiation is exchanged in a cavity that is filled by a vacuum or a medium that is transparent to thermal radiation. The cavity must be enclosed by the boundaries specified in parameter "boundary". If the cavity is not enclosed, the result of the computation may be incorrect without MOOSE issuing a warning or error.

By default, MOOSE's raytracing module is used for calculating view factors. The method is discussed in detail in RayTracingViewFactor. The parameter "view_factor_calculator" can be set to "analytical" to force the use of UnobstructedPlanarViewFactor; at this point this is maintained as backward compatibility option and it may be removed in the future.

The "ray_tracing_face_order" parameter is important because it controls the accuracy of the view factor computation. The raytracing computation uses the QGRID quadrature available in libMesh which places quadrature points at uniform distances on the from and to boundary faces. Increasing the "ray_tracing_face_order" increases the accuracy of the view factors. Note, that mesh refinement also enhances the accuracy of view factors, but at the expense of more elements being crossed during the raytracing procedure. Uniformly increasing the number of quadrature points using "ray_tracing_face_order" is essentially equivalent to selectively refining the mesh on the faces in radiative transfer.

Available Actions

  • Heat Transfer App
  • RadiationTransferActionThis action sets up the net radiation calculation between specified sidesets.

Example Input syntax

[Problem]
  kernel_coverage_check = false
[]

[Mesh]
  type = MeshGeneratorMesh

  [./cmg]
    type = CartesianMeshGenerator
    dim = 2
    dx = '1 1.3 1.9'
    ix = '3 3 3'
    dy = '2 1.2 0.9'
    iy = '3 3 3'
    subdomain_id = '0 1 0
                    4 5 2
                    0 3 0'
  [../]

  [./inner_bottom]
    type = SideSetsBetweenSubdomainsGenerator
    input = cmg
    primary_block = 1
    paired_block = 5
    new_boundary = 'inner_bottom'
  [../]

  [./inner_left]
    type = SideSetsBetweenSubdomainsGenerator
    input = inner_bottom
    primary_block = 4
    paired_block = 5
    new_boundary = 'inner_left'
  [../]

  [./inner_right]
    type = SideSetsBetweenSubdomainsGenerator
    input = inner_left
    primary_block = 2
    paired_block = 5
    new_boundary = 'inner_right'
  [../]

  [./inner_top]
    type = SideSetsBetweenSubdomainsGenerator
    input = inner_right
    primary_block = 3
    paired_block = 5
    new_boundary = 'inner_top'
  [../]

  [./rename]
    type = RenameBlockGenerator
    old_block = '1 2 3 4'
    new_block = '0 0 0 0'
    input = inner_top
  [../]
[]

[Variables]
  [./temperature]
    block = 0
  [../]
[]

[Kernels]
  [./heat_conduction]
    type = HeatConduction
    variable = temperature
    block = 0
    diffusion_coefficient = 5
  [../]
[]

[GrayDiffuseRadiation]
  [./cavity]
    boundary = '4 5 6 7'
    emissivity = '0.9 0.8 0.4 1'
    n_patches = '2 2 2 3'
    partitioners = 'centroid centroid centroid centroid'
    centroid_partitioner_directions = 'x y y x'
    temperature = temperature
    adiabatic_boundary = '7'
    fixed_temperature_boundary = '4'
    fixed_boundary_temperatures = '1200'
    view_factor_calculator = analytical
  [../]
[]

[BCs]
  [./left]
    type = DirichletBC
    variable = temperature
    boundary = left
    value = 600
  [../]

  [./right]
    type = DirichletBC
    variable = temperature
    boundary = right
    value = 300
  [../]
[]

[Postprocessors]
  [./average_T_inner_right]
    type = SideAverageValue
    variable = temperature
    boundary = inner_right
  [../]
[]

[Executioner]
  type = Steady
[]

[Outputs]
  exodus = true
[]
(moose/modules/heat_transfer/test/tests/radiation_transfer_action/radiative_transfer_action.i)

Input Parameters

  • boundaryThe boundaries that participate in the radiative exchange.

    C++ Type:std::vector<BoundaryName>

    Unit:(no unit assumed)

    Controllable:No

    Description:The boundaries that participate in the radiative exchange.

  • emissivityEmissivities for each boundary.

    C++ Type:std::vector<FunctionName>

    Unit:(no unit assumed)

    Controllable:No

    Description:Emissivities for each boundary.

  • n_patchesNumber of radiation patches per sideset.

    C++ Type:std::vector<unsigned int>

    Unit:(no unit assumed)

    Controllable:No

    Description:Number of radiation patches per sideset.

  • temperatureThe coupled temperature variable.

    C++ Type:VariableName

    Unit:(no unit assumed)

    Controllable:No

    Description:The coupled temperature variable.

Required Parameters

  • active__all__ If specified only the blocks named will be visited and made active

    Default:__all__

    C++ Type:std::vector<std::string>

    Unit:(no unit assumed)

    Controllable:No

    Description:If specified only the blocks named will be visited and made active

  • adiabatic_boundaryThe adiabatic boundaries that participate in the radiative exchange.

    C++ Type:std::vector<BoundaryName>

    Unit:(no unit assumed)

    Controllable:No

    Description:The adiabatic boundaries that participate in the radiative exchange.

  • azimuthal_quad_order8Order of the azimuthal quadrature per quadrant [azimuthal angle is measured in a plane perpendicular to the normal]. Only used if view_factor_calculator = ray_tracing.

    Default:8

    C++ Type:unsigned int

    Unit:(no unit assumed)

    Controllable:No

    Description:Order of the azimuthal quadrature per quadrant [azimuthal angle is measured in a plane perpendicular to the normal]. Only used if view_factor_calculator = ray_tracing.

  • centroid_partitioner_directionsSpecifies the sort direction if using the centroid partitioner. Available options: x, y, z, radial

    C++ Type:MultiMooseEnum

    Unit:(no unit assumed)

    Options:x, y, z, radial

    Controllable:No

    Description:Specifies the sort direction if using the centroid partitioner. Available options: x, y, z, radial

  • fixed_boundary_temperaturesThe temperatures of the fixed boundary.

    C++ Type:std::vector<FunctionName>

    Unit:(no unit assumed)

    Controllable:No

    Description:The temperatures of the fixed boundary.

  • fixed_temperature_boundaryThe fixed temperature boundaries that participate in the radiative exchange.

    C++ Type:std::vector<BoundaryName>

    Unit:(no unit assumed)

    Controllable:No

    Description:The fixed temperature boundaries that participate in the radiative exchange.

  • inactiveIf specified blocks matching these identifiers will be skipped.

    C++ Type:std::vector<std::string>

    Unit:(no unit assumed)

    Controllable:No

    Description:If specified blocks matching these identifiers will be skipped.

  • normalize_view_factorTrueDetermines if view factors are normalized to sum to one (consistent with their definition).

    Default:True

    C++ Type:bool

    Unit:(no unit assumed)

    Controllable:No

    Description:Determines if view factors are normalized to sum to one (consistent with their definition).

  • partitionersdefaultSpecifies a mesh partitioner to use when preparing the radiation patches.

    Default:default

    C++ Type:MultiMooseEnum

    Unit:(no unit assumed)

    Options:default, metis, parmetis, linear, centroid, hilbert_sfc, morton_sfc, grid

    Controllable:No

    Description:Specifies a mesh partitioner to use when preparing the radiation patches.

  • polar_quad_order16Order of the polar quadrature [polar angle is between ray and normal]. Must be even. Only used if view_factor_calculator = ray_tracing.

    Default:16

    C++ Type:unsigned int

    Unit:(no unit assumed)

    Controllable:No

    Description:Order of the polar quadrature [polar angle is between ray and normal]. Must be even. Only used if view_factor_calculator = ray_tracing.

  • print_view_factor_infoFalseFlag to print information about computed view factors.

    Default:False

    C++ Type:bool

    Unit:(no unit assumed)

    Controllable:No

    Description:Flag to print information about computed view factors.

  • ray_tracing_face_orderCONSTANTThe face quadrature rule order used for ray tracing.

    Default:CONSTANT

    C++ Type:MooseEnum

    Unit:(no unit assumed)

    Options:CONSTANT, FIRST, SECOND, THIRD, FOURTH, FIFTH, SIXTH, SEVENTH, EIGHTH, NINTH, TENTH, ELEVENTH, TWELFTH, THIRTEENTH, FOURTEENTH, FIFTEENTH, SIXTEENTH, SEVENTEENTH, EIGHTTEENTH, NINTEENTH, TWENTIETH

    Controllable:No

    Description:The face quadrature rule order used for ray tracing.

  • ray_tracing_face_typeGRIDThe face quadrature rule type used for ray tracing.

    Default:GRID

    C++ Type:MooseEnum

    Unit:(no unit assumed)

    Options:GAUSS, GRID

    Controllable:No

    Description:The face quadrature rule type used for ray tracing.

  • symmetry_boundaryThe sidesets that represent symmetry lines/planes for the problem. These sidesets do not participate in the radiative exchangeso they should not be listed in the sidesets parameter.

    C++ Type:std::vector<BoundaryName>

    Unit:(no unit assumed)

    Controllable:No

    Description:The sidesets that represent symmetry lines/planes for the problem. These sidesets do not participate in the radiative exchangeso they should not be listed in the sidesets parameter.

  • view_factor_calculatorray_tracingThe view factor calculator being used.

    Default:ray_tracing

    C++ Type:MooseEnum

    Unit:(no unit assumed)

    Options:analytical, ray_tracing

    Controllable:No

    Description:The view factor calculator being used.

Optional Parameters

  • control_tagsAdds user-defined labels for accessing object parameters via control logic.

    C++ Type:std::vector<std::string>

    Unit:(no unit assumed)

    Controllable:No

    Description:Adds user-defined labels for accessing object parameters via control logic.

Advanced Parameters