incflo Class Reference

AMR-Wind API: incflo Class Reference
AMR-Wind API v0.1.0
CFD solver for wind plant simulations

#include <incflo.H>

Inheritance diagram for incflo:
[legend]
Collaboration diagram for incflo:
[legend]

Public Member Functions

 incflo ()
 
 ~incflo () override
 
void InitData ()
 
void Evolve ()
 
void Evolve_MultiBlock (int multiblock_step, int max_block_step)
 
void ErrorEst (int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
 
void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray &new_grids, const amrex::DistributionMapping &new_dmap) override
 
void MakeNewLevelFromCoarse (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
 
void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
 
void ClearLevel (int lev) override
 
void init_mesh ()
 
void init_amr_wind_modules ()
 
void prepare_for_time_integration ()
 
bool regrid_and_update ()
 
void pre_advance_stage1 ()
 
void pre_advance_stage2 ()
 
void prepare_time_step ()
 
void do_advance (const int fixed_point_iteration)
 
void advance (const int fixed_point_iteration)
 
void prescribe_advance ()
 
void post_advance_work ()
 
amr_wind::CFDSimsim ()
 
const amr_wind::SimTimetime () const
 
amr_wind::FieldReporepo ()
 
const amr_wind::FieldReporepo () const
 
amr_wind::pde::PDEBaseicns ()
 
const amr_wind::pde::PDEBaseicns () const
 
amr_wind::pde::PDEMgr::TypeVectorscalar_eqns ()
 
const amr_wind::pde::PDEMgr::TypeVectorscalar_eqns () const
 
amr_wind::Fieldvelocity () const
 
amr_wind::Fielddensity () const
 
amr_wind::Fieldtemperature () const
 
amr_wind::Fieldpressure () const
 
amr_wind::Fieldgrad_p () const
 
void ComputeDt (bool explicit_diffusion)
 
void ComputePrescribeDt ()
 
void ApplyPredictor (const bool incremental_projection=false, const int fixed_point_iteration=0)
 
void ApplyCorrector ()
 
void ApplyPrescribeStep ()
 
void ApplyProjection (amrex::Vector< amrex::MultiFab const * > density, amrex::Real time, amrex::Real scaling_factor, bool incremental)
 
void init_physics_and_pde ()
 Initialize Physics instances as well as PDEs (include turbulence models)
 
void ReadCheckpointFile ()
 
void SetMultiBlockPointer (MultiBlockContainer *mbc)
 
void set_read_erf (ReadERFFunction fcn)
 

Private Member Functions

amrex::FabFactory< amrex::FArrayBox > const & Factory (int lev) const noexcept
 
bool need_divtau () const
 
void CheckAndSetUpDryRun ()
 
void ReadParameters ()
 
void InitialProjection ()
 
void InitialIterations ()
 
void PrintMaxValues (const std::string &header)
 
void PrintMaxVelLocations (const std::string &header)
 
void PrintMaxVel (int lev) const
 
void PrintMaxGp (int lev) const
 
void CheckForNans (int lev) const
 

Private Attributes

amr_wind::CFDSim m_sim
 
amr_wind::SimTimem_time
 
amr_wind::FieldRepom_repo
 
amr_wind::OversetOps m_ovst_ops
 
std::unique_ptr< amr_wind::RefineCriteriaManagerm_mesh_refiner
 
int m_verbose = 0
 
bool m_do_initial_proj = true
 
int m_initial_iterations = 3
 
bool m_use_godunov = true
 
bool m_prescribe_vel = false
 
bool m_dry_run = false
 
int m_fixed_point_iterations {1}
 
amrex::Long m_cell_count {-1}
 number of cells on all levels including covered cells
 
DiffusionType m_diff_type = DiffusionType::Crank_Nicolson
 
bool m_reconstruct_true_pressure {false}
 reconstruct true pressure
 

Detailed Description

AMR-Wind driver class

Constructor & Destructor Documentation

◆ incflo()

incflo::incflo ( )

◆ ~incflo()

incflo::~incflo ( )
overridedefault

Member Function Documentation

◆ advance()

void incflo::advance ( const int fixed_point_iteration)

Advance simulation state by one timestep

Performs the following actions at a given timestep

  • Compute \Delta t
  • Advance all computational fields to new timestate in preparation for time integration
  • Call pre-advance work for all registered physics modules
  • For Godunov scheme, advance to new time state
  • For MOL scheme, call predictor corrector steps
  • Perform any post-advance work

Much of the heavy-lifting is done by incflo::ApplyPredictor and incflo::ApplyCorrector. Please refer to the documentation of those methods for detailed information on the various equations being solved.

Here is the call graph for this function:

◆ ApplyCorrector()

void incflo::ApplyCorrector ( )

Corrector step for MOL scheme

  1. Solve transport equation for momentum and scalars

    \begin{align} \left[1 - \kappa \frac{\Delta t}{\rho} \nabla \cdot \left( \mu \nabla \right)\right] u^{*} &= u^n - \Delta t C_u + (1 - \kappa) \frac{\Delta t}{\rho} \nabla \cdot \left( \mu \nabla\right) u^{n} + \frac{\Delta t}{\rho} \left( S_u - \nabla(p + p_0)\right) \\ \end{align}

    where

    \begin{align} \kappa = \begin{cases} 0 & \text{Explicit} \\ 0.5 & \text{Crank-Nicolson} \\ 1 & \text{Implicit} \end{cases} \end{align}

  2. Apply projection

◆ ApplyPredictor()

void incflo::ApplyPredictor ( const bool incremental_projection = false,
const int fixed_point_iteration = 0 )

Apply predictor step

For Godunov, this completes the timestep. For MOL, this is the first part of the predictor/corrector within a timestep.

  1. Solve transport equation for momentum and scalars

    \begin{align} \left[1 - \kappa \frac{\Delta t}{\rho^{n+1/2}} \nabla \cdot \left( \mu \nabla \right)\right] u^{*} &= u^n - \Delta t (u \cdot \nabla) u + (1 - \kappa) \frac{\Delta t}{\rho^n} \nabla \cdot \left( \mu^{n} \nabla\right) u^{n} + \frac{\Delta t}{\rho^{n+1/2}} \left( S_u - \nabla(p + p_0)\right) \\ \end{align}

    where

    \begin{align} \kappa = \begin{cases} 0 & \text{Explicit} \\ 0.5 & \text{Crank-Nicolson} \\ 1 & \text{Implicit} \end{cases} \end{align}

  2. Apply projection

◆ ApplyPrescribeStep()

void incflo::ApplyPrescribeStep ( )

◆ ApplyProjection()

void incflo::ApplyProjection ( amrex::Vector< amrex::MultiFab const * > density,
amrex::Real time,
amrex::Real scaling_factor,
bool incremental )

Perform nodal projection

Computes the following decomposition:

\begin{align} u + \frac{\Delta t}{\rho} \nabla \phi &= u^{*} & \nabla \cdot u = 0 \end{align}

where u^{*} is a non-div-free velocity field, stored by components in u, v, and w. The resulting div-free velocity field, u, overwrites the value of u* in u, v, and w.

\phi is an auxiliary function related to the pressure p by the relation:

p^{\mathrm{new}} = \phi

except in the initial projection when

p^{\mathrm{new}} = p^{\mathrm{old}} + \phi

Velocity is updated using the following relation

\begin{align} u^{**} &= u^{*} + \frac{\Delta t}{\rho} \nabla p^{\mathrm{old}} \\ \nabla \cdot \left( \frac{\Delta t}{\rho} \nabla \phi \right) &= \nabla \cdot u^{**} \\ u^{\mathrm{new}} &= u^{**} - \frac{\Delta t}{\rho} \nabla p^\mathrm{new} \end{align}

Notes:

  • scaling_factor equals \Delta t except when called during initial projection, when it is 1.0
  • \rho in the above expressions is either at state n+1 or n+1/2 depending on whether this method was called from incflo::ApplyPredictor or incflo::ApplyCorrector.
  • If incremental == true, then the pressure term is not added to u^{**} and the update is in delta-form.

Please consult AMReX Linear Solvers documentation for more information on the nodal projection operator.

Parameters
densityVector of multifabs containing the density field at all levels
timeCurrent time
scaling_factorScaling factor \Delta t
incrementalFlag indicating if it is an incremental projection.

◆ CheckAndSetUpDryRun()

void incflo::CheckAndSetUpDryRun ( )
private

◆ CheckForNans()

void incflo::CheckForNans ( int lev) const
private

◆ ClearLevel()

void incflo::ClearLevel ( int lev)
override

◆ ComputeDt()

void incflo::ComputeDt ( bool explicit_diffusion)

Estimate the new timestep for adaptive timestepping algorithm

Parameters
explicit_diffusionFlag indicating whether user has selected explicit treatment of diffusion term.

Compute new \Delta t by using the formula derived in "A Boundary Condition Capturing Method for Multiphase Incompressible Flow" by Kang et al. (JCP).

\begin{align} \text{CFL} &= \frac{\Delta t}{2} \left[\left(C+V \right) + \sqrt{\left(C+V\right)^2 + 4 \left(\frac{|F_x|}{\Delta x} + \frac{|F_y|}{\Delta y} + \frac{|F_z|}{\Delta z} \right)} \right] && \\ C &= \text{max} \left(\frac{|U_x|}{\Delta x} ,\ \frac{|U_y|}{\Delta y} ,\ \frac{|U_z|}{\Delta z}\right) && \text{Convection} \\ V &= 2 \left(\frac{\mu}{\rho}\right)_\mathrm{max} \left[\frac{1}{\Delta x^2} + \frac{1}{\Delta y^2} + \frac{1}{\Delta z^2} \right] && \text{Diffusion} \end{align}

F_x, F_y, F_z are acceleration due to forcing terms. By default, C is always used for computing CFL. The term V is only used when the user has chosen implicit or Crank-Nicolson scheme for diffusion, and contributions from forcing term when time.use_force_cfl is true (default is true).

◆ ComputePrescribeDt()

void incflo::ComputePrescribeDt ( )

◆ density()

amr_wind::Field & incflo::density ( ) const
inline

◆ do_advance()

void incflo::do_advance ( const int fixed_point_iteration)

◆ ErrorEst()

void incflo::ErrorEst ( int lev,
amrex::TagBoxArray & tags,
amrex::Real time,
int ngrow )
override

◆ Evolve()

void incflo::Evolve ( )

Perform time-integration for user-defined time or timesteps.

Here is the call graph for this function:

◆ Evolve_MultiBlock()

void incflo::Evolve_MultiBlock ( int multiblock_step,
int max_block_step )

◆ Factory()

amrex::FabFactory< amrex::FArrayBox > const & incflo::Factory ( int lev) const
inlineprivatenoexcept

◆ grad_p()

amr_wind::Field & incflo::grad_p ( ) const
inline

◆ icns() [1/2]

amr_wind::pde::PDEBase & incflo::icns ( )
inline

◆ icns() [2/2]

const amr_wind::pde::PDEBase & incflo::icns ( ) const
inline

◆ init_amr_wind_modules()

void incflo::init_amr_wind_modules ( )

Initialize AMR-Wind data structures after mesh has been created.

Modules initialized:

  • Registered Physics models classes
  • Registered PDE systems
  • Registered post-processing classes

◆ init_mesh()

void incflo::init_mesh ( )

Initialize AMR mesh data structure.

Calls the AmrMesh/AmrCore functions to initialize the mesh. For restarted simulations, it loads the checkpoint file.

◆ init_physics_and_pde()

void incflo::init_physics_and_pde ( )

Initialize Physics instances as well as PDEs (include turbulence models)

◆ InitData()

void incflo::InitData ( )

Perform all initialization actions for AMR-Wind.

This is a wrapper method that calls the following methods to perform the actual work.

Here is the call graph for this function:

◆ InitialIterations()

void incflo::InitialIterations ( )
private

Perform initial pressure iterations

Performs a user-defined number of iterations to compute the pressure based on the initial conditions. This method is only invoked for new simulations and skipped for restarted simulations from a checkpoint file.

◆ InitialProjection()

void incflo::InitialProjection ( )
private

Ensure initial velocity field is divergence-free.

Performs a projection step to ensure that the user-provided initial velocity field is divergence free. This method is only invoked for new simulations. For restarted simulations using a checkpoint file, this is not necessary.

◆ MakeNewLevelFromCoarse()

void incflo::MakeNewLevelFromCoarse ( int lev,
amrex::Real time,
const amrex::BoxArray & ba,
const amrex::DistributionMapping & dm )
override

◆ MakeNewLevelFromScratch()

void incflo::MakeNewLevelFromScratch ( int lev,
amrex::Real time,
const amrex::BoxArray & new_grids,
const amrex::DistributionMapping & new_dmap )
override

◆ need_divtau()

bool incflo::need_divtau ( ) const
inlineprivate

◆ post_advance_work()

void incflo::post_advance_work ( )

Perform actions after a timestep

Performs any necessary I/O before advancing to next timestep

◆ pre_advance_stage1()

void incflo::pre_advance_stage1 ( )

◆ pre_advance_stage2()

void incflo::pre_advance_stage2 ( )

◆ prepare_for_time_integration()

void incflo::prepare_for_time_integration ( )

Initialize flow-field before performing time-integration.

This method calls the incflo::InitialProjection step to ensure that the velocity field is divergence free. Then it performs a user-defined number of initial iterations to compute p^{n - 1/2} necessary for advancing into the first timestep. These actions are only performed when starting a new simulation, but not for a restarted simulation.

◆ prepare_time_step()

void incflo::prepare_time_step ( )

◆ prescribe_advance()

void incflo::prescribe_advance ( )

◆ pressure()

amr_wind::Field & incflo::pressure ( ) const
inline

◆ PrintMaxGp()

void incflo::PrintMaxGp ( int lev) const
private

◆ PrintMaxValues()

void incflo::PrintMaxValues ( const std::string & header)
private

◆ PrintMaxVel()

void incflo::PrintMaxVel ( int lev) const
private

◆ PrintMaxVelLocations()

void incflo::PrintMaxVelLocations ( const std::string & header)
private

◆ ReadCheckpointFile()

void incflo::ReadCheckpointFile ( )

◆ ReadParameters()

void incflo::ReadParameters ( )
private

Parse the input file and populate parameters

◆ regrid_and_update()

bool incflo::regrid_and_update ( )

Perform regrid actions at a given timestep.

Returns
Flag indicating if regrid was performed

◆ RemakeLevel()

void incflo::RemakeLevel ( int lev,
amrex::Real time,
const amrex::BoxArray & ba,
const amrex::DistributionMapping & dm )
override

◆ repo() [1/2]

amr_wind::FieldRepo & incflo::repo ( )
inline

◆ repo() [2/2]

const amr_wind::FieldRepo & incflo::repo ( ) const
inline

◆ scalar_eqns() [1/2]

amr_wind::pde::PDEMgr::TypeVector & incflo::scalar_eqns ( )
inline

◆ scalar_eqns() [2/2]

const amr_wind::pde::PDEMgr::TypeVector & incflo::scalar_eqns ( ) const
inline

◆ set_read_erf()

void incflo::set_read_erf ( ReadERFFunction fcn)
inline

◆ SetMultiBlockPointer()

void incflo::SetMultiBlockPointer ( MultiBlockContainer * mbc)
inline

◆ sim()

amr_wind::CFDSim & incflo::sim ( )
inline

◆ temperature()

amr_wind::Field & incflo::temperature ( ) const
inline

◆ time()

const amr_wind::SimTime & incflo::time ( ) const
inline

◆ velocity()

amr_wind::Field & incflo::velocity ( ) const
inline

Member Data Documentation

◆ m_cell_count

amrex::Long incflo::m_cell_count {-1}
private

number of cells on all levels including covered cells

◆ m_diff_type

DiffusionType incflo::m_diff_type = DiffusionType::Crank_Nicolson
private

◆ m_do_initial_proj

bool incflo::m_do_initial_proj = true
private

◆ m_dry_run

bool incflo::m_dry_run = false
private

◆ m_fixed_point_iterations

int incflo::m_fixed_point_iterations {1}
private

◆ m_initial_iterations

int incflo::m_initial_iterations = 3
private

◆ m_mesh_refiner

std::unique_ptr<amr_wind::RefineCriteriaManager> incflo::m_mesh_refiner
private

◆ m_ovst_ops

amr_wind::OversetOps incflo::m_ovst_ops
private

◆ m_prescribe_vel

bool incflo::m_prescribe_vel = false
private

◆ m_reconstruct_true_pressure

bool incflo::m_reconstruct_true_pressure {false}
private

reconstruct true pressure

◆ m_repo

amr_wind::FieldRepo& incflo::m_repo
private

◆ m_sim

amr_wind::CFDSim incflo::m_sim
private

◆ m_time

amr_wind::SimTime& incflo::m_time
private

◆ m_use_godunov

bool incflo::m_use_godunov = true
private

◆ m_verbose

int incflo::m_verbose = 0
private

The documentation for this class was generated from the following files: