Pre-processing Utilities

PreProcessDriver

class PreProcessDriver

A driver that runs all preprocessor tasks.

This class is responsible for reading the input file, parsing the user-requested list of tasks, initializing the task instances, executing them, and finally writing out the updated Exodus database with changed inputs.

Public Functions

PreProcessDriver(stk::ParallelMachine &comm, const std::string filename)

Parameters
  • comm: MPI Communicator reference

  • filename: Name of the YAML input file

void run()

Run all tasks and output the updated Exodus database.

PreProcessingTask

class PreProcessingTask

An abstract implementation of a PreProcessingTask.

This class defines the interface for a pre-processing task and contains the infrastructure to allow concrete implementations of pre-processing tasks to register themselves for automatic runtime discovery. Derived classes must implement two methods:

  • initialize - Perform actions on STK MetaData before processing BulkData

  • run - All actions on BulkData and other operations on mesh after it has been loaded from the disk.

For automatic class registration, the derived classes must implement a constructor that takes two arguments: a CFDMesh reference, and a const reference to YAML::Node that contains the inputs necessary for the concrete task implementation. It is the derived class’ responsibility to process the input dictionary and perform error checking. No STK mesh manipulations must occur in the constructor.

Subclassed by sierra::nalu::ABLFields, sierra::nalu::BdyIOPlanes, sierra::nalu::ChannelFields, sierra::nalu::HITFields, sierra::nalu::InflowHistory, sierra::nalu::NDTW2D, sierra::nalu::NestedRefinement, sierra::nalu::RotateMesh, sierra::nalu::SamplingPlanes, sierra::nalu::TranslateMesh

Public Functions

PreProcessingTask(CFDMesh &mesh)

Parameters

virtual void initialize() = 0

Initialize the STK MetaData instance.

This method handles the registration and creation of new parts and fields. All subclasses must implement this method.

virtual void run() = 0

Process the STK BulkData instance.

This method handles the creating of new entities, manipulating coordinates, and populating fields.

Public Static Functions

PreProcessingTask *create(CFDMesh &mesh, const YAML::Node &node, std::string lookup)

Runtime creation of concrete task instance.

Protected Attributes

CFDMesh &mesh_

Reference to the CFDMesh instance.

ABLFields

class ABLFields : public sierra::nalu::PreProcessingTask

Initialize velocity and temperature fields for ABL simulations.

This task is activated by using the init_abl_fields task in the preprocessing input file. It requires a section init_abl_fields in the nalu_preprocess section with the following parameters:

init_abl_fields:
  fluid_parts: [Unspecified-2-HEX]

  temperature:
    heights: [    0, 650.0, 750.0, 10750.0]
    values:  [280.0, 280.0, 288.0,   318.0]

  velocity:
    heights: [0.0, 10.0, 30.0, 70.0, 100.0, 650.0, 10000.0]
    values:
      - [ 0.0, 0.0, 0.0]
      - [4.81947, -4.81947, 0.0]
      - [5.63845, -5.63845, 0.0]
      - [6.36396, -6.36396, 0.0]
      - [6.69663, -6.69663, 0.0]
      - [8.74957, -8.74957, 0.0]
      - [8.74957, -8.74957, 0.0]

The sections temperature and velocity are optional, allowing the user to initialize only the temperature or the velocity as desired. The heights are in meters, the temperature is the potential temperature in Kelvin, and the velocity is the actual vector in m/s. Currently, the code does not include the ability to automatically convert (mangitude, direction) to velocity vectors.

Public Functions

ABLFields(CFDMesh &mesh, const YAML::Node &node)

Parameters

void initialize()

Declare velocity and temperature fields and register them for output.

void run()

Initialize the velocity and/or temperature fields by linear interpolation.

Private Functions

void load(const YAML::Node &abl)

Parse the YAML file and initialize parameters.

void load_velocity_info(const YAML::Node &abl)

Helper function to parse and initialize velocity inputs.

void load_temperature_info(const YAML::Node &abl)

Helper function to parse and initialize temperature inputs.

void init_velocity_field()

Initialize the velocity field through linear interpolation.

void init_temperature_field()

Intialize the temperature field through linear interpolation.

void perturb_velocity_field()

Add perturbations to velocity field.

void perturb_temperature_field()

Add perturbations to temperature field.

Private Members

stk::mesh::MetaData &meta_

STK Metadata object.

stk::mesh::BulkData &bulk_

STK Bulkdata object.

stk::mesh::PartVector fluid_parts_

Parts of the fluid mesh where velocity/temperature is initialized.

std::vector<double> vHeights_

List of heights where velocity is defined.

Array2D<double> velocity_

List of velocity (3-d components) at the user-defined heights.

std::vector<double> THeights_

List of heights where temperature is defined.

std::vector<double> TValues_

List of temperatures (K) at user-defined heights (THeights_)

std::vector<std::string> periodicParts_

List of periodic parts.

double deltaU_ = {1.0}

Velocity perturbation amplitude for Ux.

double deltaV_ = {1.0}

Velocity perturbation amplitude for Uy.

double Uperiods_ = {4.0}

Number of periods for Ux.

double Vperiods_ = {4.0}

Number of periods for Uy.

double zRefHeight_ = {50.0}

Reference height for velocity perturbations.

double thetaAmplitude_

Amplitude of temperature perturbations.

double thetaGaussMean_ = {0.0}

Mean for the Gaussian random number generator.

double thetaGaussVar_ = {1.0}

Variance of the Gaussian random number generator.

double thetaCutoffHt_

Cutoff height for temperature fluctuations.

int ndim_

Dimensionality of the mesh.

bool doVelocity_

Flag indicating whether velocity is initialized.

bool doTemperature_

Flag indicating whether temperature is initialized.

bool perturbU_ = {false}

Flag indicating whether velocity perturbations are added during initialization.

bool perturbT_ = {false}

Flag indicating whether temperature perturbations are added.

BdyIOPlanes

class BdyIOPlanes : public sierra::nalu::PreProcessingTask

Extract boundary planes for I/O mesh.

Given an ABL precursor mesh, this utility extracts the specified boundaries and creates a new IO Transfer mesh for use with ABL precursor simulations.

Public Functions

BdyIOPlanes(CFDMesh &mesh, const YAML::Node &node)

Parameters

void initialize()

Register boundary parts and attach coordinates to the parts.

The parts are created as SHELL elements to as needed by the Nalu Transfer class.

void run()

Copy user specified boundaries and save the IO Transfer mesh.

Private Functions

void load(const YAML::Node &node)

Parse user inputs from the YAML file.

void create_boundary(const std::string bdyName)

Copy the boundary from Fluid mesh to the IO Xfer mesh.

Private Members

CFDMesh &mesh_

Original mesh DB information.

CFDMesh iomesh_

IO Mesh db STK meta and bulk data.

std::vector<std::string> bdyNames_

User specified list of boundaries to be extracted.

std::string output_db_ = {""}

Name of the I/O db where the boundaries are written out.

SamplingPlanes

class SamplingPlanes : public sierra::nalu::PreProcessingTask

Generate 2-D grids/planes for data sampling.

Currently only generates horizontal planes at user-defined heights.

Requires a section generate_planes in the input file within the nalu_preprocess section:

generate_planes:
  fluid_part: Unspecified-2-hex

  heights: [ 70.0 ]
  part_name_format: "zplane_%06.1f"

  dx: 12.0
  dy: 12.0

With the above input definition, it will use the bounding box of the fluid_part to determine the bounding box of the plane to be generated. This will provide coordinate axis aligned sapling planes in x and y directions. Alternately, the user can specify boundary_type to be quad_vertices and provide the vertices of the quadrilateral that is used to generate the sampling plane as shown below:

generate_planes:
  boundary_type: quad_vertices
  fluid_part: Unspecified-2-hex

  heights: [ 50.0, 70.0, 90.0 ]
  part_name_format: "zplane_%06.1f"

  nx: 25  # Number of divisions along (1-2) and (4-3) vertices
  ny: 25  # Number of divisions along (1-4) and (2-3) vertices
  vertices:
    - [250.0, 0.0]
    - [500.0, -250.0]
    - [750.0, 0.0]
    - [500.0, 250.0]

part_name_format is a printf-like format specification that takes one argument - the height as a floating point number. The user can use this to tailor how the nodesets or the shell parts are named in the output Exodus file.

Public Types

enum PlaneBoundaryType

Sampling Plane boundary type.

Values:

BOUND_BOX = 0

Use bounding box of the fluid mesh defined.

QUAD_VERTICES

Use user-defined vertex list for plane boundary.

Public Functions

void initialize()

Initialize the STK MetaData instance.

This method handles the registration and creation of new parts and fields. All subclasses must implement this method.

void run()

Process the STK BulkData instance.

This method handles the creating of new entities, manipulating coordinates, and populating fields.

Private Functions

void calc_bounding_box()

Use fluid Realm mesh to estimate the x-y bounding box for the sampling planes.

void generate_zplane(const double zh)

Generate entities and update coordinates for a given sampling plane.

Private Members

stk::mesh::MetaData &meta_

STK Metadata object.

stk::mesh::BulkData &bulk_

STK Bulkdata object.

std::vector<double> heights_

Heights where the averaging planes are generated.

std::array<std::array<double, 3>, 2> bBox_

Bounding box of the original mesh.

std::string name_format_

Format specification for the part name.

std::vector<std::string> fluidPartNames_

Fluid realm parts (to determine mesh bounding box)

stk::mesh::PartVector fluidParts_

Parts of the fluid mesh (to determine mesh bounding box)

double dx_

Spatial resolution in x and y directions.

double dy_

Spatial resolution in x and y directions.

size_t nx_

Number of nodes in x and y directions.

size_t mx_

Number of elements in x and y directions.

int ndim_

Dimensionality of the mesh.

PlaneBoundaryType bdyType_ = {BOUND_BOX}

User defined selection of plane boundary type.

NestedRefinement

class NestedRefinement : public sierra::nalu::PreProcessingTask

Tag regions in mesh for refinement with Percept mesh_adapt utility.

This utility creates a field turbine_refinement_field that is populated with an indicator value between [0, 1] that can be used with the Percept mesh_adapt utility to locally refine regions of interest.

A typical use of this utility is to refine an ABL mesh around turbines, especially for use with actuator line wind farm simulations.

Public Functions

void initialize()

Initialize the refinement field and register to parts.

void run()

Perform search and tag elements with appropriate values for subsequent refinement.

Private Functions

void load(const YAML::Node &node)

Parse the YAML file and initialize the necessary parameters.

void process_inputs()

Process input data and populate necessary data structures for subsequent use.

double compute_refine_fraction(Vec3D &point)

Estimate the refinement fraction [0,1] for a given element, indicated by the element mid point.

void write_percept_inputs()

Write out the input files that can be used with Percept.

Private Members

std::vector<std::string> fluidPartNames_

Partnames for the ABL mesh.

stk::mesh::PartVector fluidParts_

Parts of the ABL mesh where refinement is performed.

std::vector<double> turbineDia_

List of turbine diameters for the turbines in the wind farm [numTurbines].

std::vector<double> turbineHt_

List of turbine tower heights for the turbines in wind farm [numTurbines].

std::vector<Vec3D> turbineLocs_

List of turbine pad locations [numTurbines, 3].

std::vector<std::vector<double>> refineLevels_

List of refinement levels [numLevels, 3].

std::vector<TrMat> boxAxes_

Transformation matrices for each turbine [numTurbines].

std::vector<Vec3D> corners_

The minimum corners for each refinement box [numTurbines * numLevels].

std::vector<Vec3D> boxLengths_

The dimensions of each box [numTurbines * numLevels].

std::string refineFieldName_ = {"turbine_refinement_field"}

Field name used in the Exodus mesh for the error indicator field.

std::string perceptFilePrefix_ = {"adapt"}

Prefix for the input file name.

double searchTol_ = {10.0}

Search tolerance used when searching for box inclusion.

double windAngle_ = {270.0}

Compass direction of the wind (in degrees)

size_t numTurbines_

The number of turbines in the wind farm.

size_t numLevels_

The number of refinement levels.

bool writePercept_ = {true}

Write input files for use with subsequent percept run.

ChannelFields

class ChannelFields : public sierra::nalu::PreProcessingTask

Initialize velocity fields for channel flow simulations.

This task is activated by using the init_channel_fields task in the preprocessing input file. It requires a section init_channel_fields in the nalu_preprocess section with the following parameters:

init_channel_fields:
  fluid_parts: [Unspecified-2-HEX]

  velocity:
    Re_tau : 550
    viscosity : 0.0000157

The user specified the friction Reynolds number, Re_tau, and the kinematic viscosity (in m^2/s). The velocity field is initialized to a Reichardt function, with an imposed sinusoidal perturbation and random perturbation in the wall parallel directions.

Public Functions

void initialize()

Declare velocity fields and register them for output.

void run()

Initialize the velocity fields by linear interpolation.

RotateMesh

class RotateMesh : public sierra::nalu::PreProcessingTask

Rotate a mesh.

rotate_mesh:
  mesh_parts:
    - unspecified-2-hex

  angle: 45.0
  origin: [500.0, 0.0, 0.0]
  axis: [0.0, 0.0, 1.0]

Public Functions

void initialize()

Initialize the STK MetaData instance.

This method handles the registration and creation of new parts and fields. All subclasses must implement this method.

void run()

Process the STK BulkData instance.

This method handles the creating of new entities, manipulating coordinates, and populating fields.

Private Members

stk::mesh::MetaData &meta_

STK Metadata object.

stk::mesh::BulkData &bulk_

STK Bulkdata object.

std::vector<std::string> meshPartNames_

Part names of the mesh that needs to be rotated.

stk::mesh::PartVector meshParts_

Parts of the mesh that need to be rotated.

double angle_

Angle of rotation.

std::vector<double> origin_

Point about which rotation is performed.

std::vector<double> axis_

Axis around which the rotation is performed.

int ndim_

Dimensionality of the mesh.

NDTW2D

class NDTW2D : public sierra::nalu::PreProcessingTask

2-D Nearest distance to wall calculator

Calculates a new field NDTW containing the wall distance for 2-D airfoil-like applications used in RANS wall models.

Public Functions

void initialize()

Initialize the NDTW field and register for output.

void run()

Calculate wall distance and update NDTW field.