AddFieldTransferAction

Overview

The AddFieldTransferAction is responsible for performing data transfers between NekRS and MOOSE of field variables (quantities in NekRS defined at the nodal points, distributed in space throughout a mesh). This is done with the [FieldTransfers] block in a Cardinal input file. Examples of field transfers which can be added include:

Creating Variables

When writing data to_nek, Cardinal will read from an AuxVariable; when reading data from_nek, Cardinal will write that data into an AuxVariable. Therefore, this object also automatically creates the needed AuxVariable using the name of the object. The order of the auxvariable is matched to the order of the NekRSMesh.

Passing Data

The order of the internal NekRS mesh may not necessarily match that of NekRSMesh. Therefore, this object will interpolate between the NekRS GLL points and the NekRSMesh using Vandermonde matrices.

Usrwrk Array

When direction = to_nek, data is "sent" into NekRS by writing into the nrs->usrwrk scratch space array and also copying from this array into its corresponding array on the host (nrs->o_usrwrk) so that the data is accessible both on host and device within NekRS. For field transfers, each field occupies a single "slot" in the nrs->usrwrk array. For transfers into NekRS, you are required to list which slot in the array to write. Note that Cardinal automatically allocates this array, but you can control the size of it by setting n_usrwrk_slots in NekRSProblem.

warningwarning

Allocation of nrs->usrwrk and nrs->o_usrwrk is done automatically. If you attempt to run a NekRS input file that accesses bc->usrwrk in the .oudf file without a Cardinal executable (e.g. like nrsmpi case 4), then that scratch space will have to be manually allocated in the .udf file, or else your input will seg fault. This will not be typically encountered by most users, but if you really do want to run the NekRS input files intended for a Cardinal case with the NekRS executable (perhaps for debugging), we recommend simply replacing bc->usrwrk by a dummy value, such as bc->flux = 0.0 for the boundary heat flux use case. This just replaces a value that normally comes from MOOSE by a fixed value. All other aspects of the NekRS case files should not require modification.

Example Input File Syntax

As an example, a NekBoundaryFlux and a NekVolumetricSource are added in the [Problem/FieldTransfers] block. These objects automatically create auxvariables with names avg_flux and heat_source, respectively; these variables will be read from in order to write data into NekRS. A NekFieldVariable is also added, which will create a variable named temp; this auxiliary variable will receive the NekRS internal field variable, temperature.

[Problem<<<{"href": "../../syntax/Problem/index.html"}>>>]
  type = NekRSProblem
  casename<<<{"description": "Case name for the NekRS input files; this is <case> in <case>.par, <case>.udf, <case>.oudf, and <case>.re2."}>>> = 'pyramid'

  [FieldTransfers<<<{"href": "../../syntax/Problem/FieldTransfers/index.html"}>>>]
    [avg_flux]
      type = NekBoundaryFlux<<<{"description": "Reads/writes boundary flux data between NekRS and MOOSE."}>>>
      usrwrk_slot<<<{"description": "When 'direction = to_nek', the slot(s) in the usrwrk array to write the incoming data; provide one entry for each quantity being passed"}>>> = 0
      direction<<<{"description": "Direction in which to send data"}>>> = to_nek
      postprocessor_to_conserve<<<{"description": "Name of the postprocessor/vectorpostprocessor containing the integral(s) used to ensure conservation; defaults to the name of the object plus '_integral'"}>>> = flux_integral
    []
    [heat_source]
      type = NekVolumetricSource<<<{"description": "Reads/writes volumetric source data between NekRS and MOOSE."}>>>
      usrwrk_slot<<<{"description": "When 'direction = to_nek', the slot(s) in the usrwrk array to write the incoming data; provide one entry for each quantity being passed"}>>> = 1
      direction<<<{"description": "Direction in which to send data"}>>> = to_nek
      postprocessor_to_conserve<<<{"description": "Name of the postprocessor/vectorpostprocessor containing the integral(s) used to ensure conservation; defaults to the name of the object plus '_integral'"}>>> = source_integral
    []
    [temp]
      type = NekFieldVariable<<<{"description": "Reads/writes volumetric field data between NekRS and MOOSE."}>>>
      field<<<{"description": "NekRS field variable to read/write; defaults to the name of the object"}>>> = temperature
      direction<<<{"description": "Direction in which to send data"}>>> = from_nek
    []
  []
[]
(test/tests/conduction/boundary_and_volume/prism/nek.i)

When running this case, Cardinal will print out a table showing a summary of what data is occupying the various slots in the usrwrk arrays and how they can be accessed. This table for the above example is listed below, which shows that heat flux values will exist in the first slice and volumetric heat source values will exist in the second slice.


-----------------------------------------------------------------------------------------------------------
| Slot | MOOSE quantity |          How to Access (.oudf)          |         How to Access (.udf)          |
-----------------------------------------------------------------------------------------------------------
|    0 | avg_flux       | bc->usrwrk[0 * bc->fieldOffset+bc->idM] | nrs->usrwrk[0 * nrs->fieldOffset + n] |
|    1 | heat_source    | bc->usrwrk[1 * bc->fieldOffset+bc->idM] | nrs->usrwrk[1 * nrs->fieldOffset + n] |
|    2 | unused         | bc->usrwrk[2 * bc->fieldOffset+bc->idM] | nrs->usrwrk[2 * nrs->fieldOffset + n] |
|    3 | unused         | bc->usrwrk[3 * bc->fieldOffset+bc->idM] | nrs->usrwrk[3 * nrs->fieldOffset + n] |
|    4 | unused         | bc->usrwrk[4 * bc->fieldOffset+bc->idM] | nrs->usrwrk[4 * nrs->fieldOffset + n] |
|    5 | unused         | bc->usrwrk[5 * bc->fieldOffset+bc->idM] | nrs->usrwrk[5 * nrs->fieldOffset + n] |
|    6 | unused         | bc->usrwrk[6 * bc->fieldOffset+bc->idM] | nrs->usrwrk[6 * nrs->fieldOffset + n] |
-----------------------------------------------------------------------------------------------------------

For instance, the heat flux is then used in a boundary condition on the device in the .oudf file.

void scalarNeumannConditions(bcData *bc)
(test/tests/conduction/boundary_and_volume/prism/pyramid.oudf)

Available FieldTransfer Objects

  • Cardinal App
  • NekBoundaryFluxReads/writes boundary flux data between NekRS and MOOSE.
  • NekFieldVariableReads/writes volumetric field data between NekRS and MOOSE.
  • NekMeshDeformationReads/writes mesh deformation between NekRS and MOOSE.
  • NekVolumetricSourceReads/writes volumetric source data between NekRS and MOOSE.

Input Parameters

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

    Default:__all__

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

    Controllable:No

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

  • inactiveIf specified blocks matching these identifiers will be skipped.

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

    Controllable:No

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

Optional Parameters

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

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

    Controllable:No

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

Advanced Parameters