Data Transfers in MOOSE

In this tutorial, you will learn how to:

  • Select a Transfer appropriate for your simulation

To access this tutorial,


cd cardinal/tutorials/transfers

The MOOSE Transfer system is used to send data between applications. In Cardinal, these transfers are used to communicate fields between NekRS, OpenMC, and MOOSE. When using Cardinal, we are most often communicating 3-D fields (temperature, heat flux, power, etc.) between our MultiApps, so we typically only use a subset of the Transfers available in MOOSE. When selecting a transfer, these are some considerations to be aware of; this page provides example use cases for all of the applicable field-based transfers in MOOSE.

This tutorial will use two generic MOOSE applications, because the advice on this page is agnostic of whether the physics solve is performed with OpenMC, NekRS, or another MOOSE application - the main requirement in selecting a tranfser for a field variable is an understanding of the source/target meshes for the transfer. We will consider the most generic use case of a transfer between two entirely different meshes that in some regions of space do not overlap one another. In Cardinal, we often encounter this due to NekRS's very high-order solution (many degrees of freedom per element allows fairly coarse elements). For each case, we will walk through all the applicable MOOSE transfers.

We will transfer data between two cylinders, one with a much finer mesh than the other. First, generate the meshes:


cardinal-opt -i mesh1.i --mesh-only
cardinal-opt -i mesh2.i --mesh-only

The two meshes are shown in Figure 1; the wireframe overlay on the right shows that, while both meshes represent a cylinder with the same outer radius, that due to different azimuthal refinement, there is not perfect spatial overlap of the meshes.

Figure 1: Meshes to be used for data transfers in this tutorial

For all cases, we will send a variable from mesh 1 to mesh 2. We simply set this variable to a generic function:

[AuxVariables]
  [u]
  []
[]

[Functions]
  [u]
    type = ParsedFunction
    expression = 'exp(x)*exp(y)'
  []
[]

[AuxKernels]
  [u]
    type = FunctionAux
    variable = u
    function = u
  []
[]
(tutorials/transfers/main.i)

We will then send this variable using a few different transfers to a sub-application.

[Transfers]
  [interpolation]
    type = MultiAppGeometricInterpolationTransfer
    source_variable = u
    variable = u_interpolation
    to_multi_app = sub
  []
  [mesh_function]
    type = MultiAppGeneralFieldShapeEvaluationTransfer
    source_variable = u
    variable = u_mesh_function
    to_multi_app = sub
  []
  [nearest_node]
    type = MultiAppGeneralFieldNearestLocationTransfer
    source_variable = u
    variable = u_nearest_node
    to_multi_app = sub
  []
  [projection]
    type = MultiAppProjectionTransfer
    source_variable = u
    variable = u_projection
    to_multi_app = sub
  []
  [sample]
    type = MultiAppVariableValueSampleTransfer
    source_variable = u
    variable = u_sample
    to_multi_app = sub
  []
[]
(tutorials/transfers/main.i)

MultiAppGeometricInterpolationTransfer link

The MultiAppGeometricInterpolationTransfer transfers the nearest node's source variable to the nearest node on the target mesh using mesh interpolation. Any non-overlapping domains are set values using extrapolation.

Figure 2: Results of a MultiAppGeometricInterpolationTransfer from the main app to the sub app

MultiAppShapeEvaluationTransfer link

The MultiAppShapeEvaluationTransfer queries the finite element solution from the source variable at each receiving point in the receiving application. Any non-overlapping domains will simply be set to a value of zero, so this transfer is only recommended for meshes that share the same outer boundary.

Figure 3: Results of a MultiAppShapeEvaluationTransfer from the main app to the sub app

MultiAppGeneralFieldNearestLocationTransfer link

The MultiAppGeneralFieldNearestLocationTransfer transfers a variable between the target and source domains using a nearest node lookup. For performance, we highly recommend setting the fixed_meshes = true parameter to cache the nearest node matches after the first transfer.

Figure 4: Results of a MultiAppGeneralFieldNearestLocationTransfer from the main app to the sub app

MultiAppProjectionTransfer link

The MultiAppProjectionTransfer transfers a variable between the target and source domains using a projection mapping. For performance, we highly recommend setting the fixed_meshes = true parameter to cache the nearest node matches after the first transfer.

Figure 5: Results of a MultiAppProjectionTransfer from the main app to the sub app

MultiAppVariableValueSampleTransfer link

The MultiAppVariableValueSampleTransfer transfers the value of a variable within the source application at each receiver application's position in the positions parameter on the MultiApp to a field variable. This transfer is mostly meant for multiscale simulations, where you want to send one value from a main application to sub-applications. The value sent to each sub-application is received as a singly-valued variable by evaluating the source variable at the application's position.

Figure 6: Results of a MultiAppVariableValueSampleTransfer from the main app to the sub app