incflo Class Reference
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
|
#include <incflo.H>
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::CFDSim & | sim () |
const amr_wind::SimTime & | time () const |
amr_wind::FieldRepo & | repo () |
const amr_wind::FieldRepo & | repo () const |
amr_wind::pde::PDEBase & | icns () |
const amr_wind::pde::PDEBase & | icns () const |
amr_wind::pde::PDEMgr::TypeVector & | scalar_eqns () |
const amr_wind::pde::PDEMgr::TypeVector & | scalar_eqns () const |
amr_wind::Field & | velocity () const |
amr_wind::Field & | density () const |
amr_wind::Field & | temperature () const |
amr_wind::Field & | pressure () const |
amr_wind::Field & | grad_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::SimTime & | m_time |
amr_wind::FieldRepo & | m_repo |
amr_wind::OversetOps | m_ovst_ops |
std::unique_ptr< amr_wind::RefineCriteriaManager > | m_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()
|
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.
◆ ApplyCorrector()
void incflo::ApplyCorrector | ( | ) |
Corrector step for MOL scheme
-
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}
- 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.
-
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}
- 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
orn+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
-
density Vector of multifabs containing the density field at all levels time Current time scaling_factor Scaling factor \Delta t incremental Flag indicating if it is an incremental projection.
◆ CheckAndSetUpDryRun()
|
private |
◆ CheckForNans()
|
private |
◆ ClearLevel()
|
override |
◆ ComputeDt()
void incflo::ComputeDt | ( | bool | explicit_diffusion | ) |
Estimate the new timestep for adaptive timestepping algorithm
- Parameters
-
explicit_diffusion Flag 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()
|
inline |
◆ do_advance()
void incflo::do_advance | ( | const int | fixed_point_iteration | ) |
◆ ErrorEst()
|
override |
◆ Evolve()
void incflo::Evolve | ( | ) |
Perform time-integration for user-defined time or timesteps.
◆ Evolve_MultiBlock()
void incflo::Evolve_MultiBlock | ( | int | multiblock_step, |
int | max_block_step ) |
◆ Factory()
|
inlineprivatenoexcept |
◆ grad_p()
|
inline |
◆ icns() [1/2]
|
inline |
◆ icns() [2/2]
|
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.
◆ 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()
|
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()
|
override |
◆ MakeNewLevelFromScratch()
|
override |
◆ need_divtau()
|
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()
|
inline |
◆ PrintMaxGp()
|
private |
◆ PrintMaxValues()
|
private |
◆ PrintMaxVel()
|
private |
◆ PrintMaxVelLocations()
|
private |
◆ ReadCheckpointFile()
void incflo::ReadCheckpointFile | ( | ) |
◆ 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()
|
override |
◆ repo() [1/2]
|
inline |
◆ repo() [2/2]
|
inline |
◆ scalar_eqns() [1/2]
|
inline |
◆ scalar_eqns() [2/2]
|
inline |
◆ set_read_erf()
|
inline |
◆ SetMultiBlockPointer()
|
inline |
◆ sim()
|
inline |
◆ temperature()
|
inline |
◆ time()
|
inline |
◆ velocity()
|
inline |
Member Data Documentation
◆ m_cell_count
|
private |
number of cells on all levels including covered cells
◆ m_diff_type
|
private |
◆ m_do_initial_proj
|
private |
◆ m_dry_run
|
private |
◆ m_fixed_point_iterations
|
private |
◆ m_initial_iterations
|
private |
◆ m_mesh_refiner
|
private |
◆ m_ovst_ops
|
private |
◆ m_prescribe_vel
|
private |
◆ m_reconstruct_true_pressure
|
private |
reconstruct true pressure
◆ m_repo
|
private |
◆ m_sim
|
private |
◆ m_time
|
private |
◆ m_use_godunov
|
private |
◆ m_verbose
|
private |
The documentation for this class was generated from the following files:
- /home/runner/work/amr-wind/amr-wind/amr-wind/incflo.H
- /home/runner/work/amr-wind/amr-wind/amr-wind/incflo.cpp
- /home/runner/work/amr-wind/amr-wind/amr-wind/incflo_advance.cpp
- /home/runner/work/amr-wind/amr-wind/amr-wind/incflo_compute_dt.cpp
- /home/runner/work/amr-wind/amr-wind/amr-wind/incflo_regrid.cpp
- /home/runner/work/amr-wind/amr-wind/amr-wind/projection/incflo_apply_nodal_projection.cpp
- /home/runner/work/amr-wind/amr-wind/amr-wind/setup/init.cpp
- /home/runner/work/amr-wind/amr-wind/amr-wind/utilities/diagnostics.cpp
- /home/runner/work/amr-wind/amr-wind/amr-wind/utilities/io.cpp
Generated by