ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
|
#include <ABLMost.H>
Public Types | |
enum class | FluxCalcType { MOENG = 0 , DONELAN , CUSTOM } |
enum class | ThetaCalcType { ADIABATIC = 0 , HEAT_FLUX , SURFACE_TEMPERATURE } |
enum class | RoughCalcType { CONSTANT = 0 , CHARNOCK , MODIFIED_CHARNOCK } |
Public Member Functions | |
ABLMost (const amrex::Vector< amrex::Geometry > &geom, bool &use_exp_most, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &z_phys_nd, amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab >>> &sst_lev, amrex::Vector< amrex::Vector< std::unique_ptr< amrex::iMultiFab >>> &lmask_lev, amrex::Vector< amrex::Vector< amrex::MultiFab * >> lsm_data, amrex::Vector< amrex::Vector< amrex::MultiFab * >> lsm_flux, amrex::Real start_bdy_time=0.0, amrex::Real bdy_time_interval=0.0) | |
void | update_fluxes (const int &lev, const amrex::Real &time, int max_iters=25) |
template<typename FluxIter > | |
void | compute_fluxes (const int &lev, const int &max_iters, const FluxIter &most_flux) |
void | impose_most_bcs (const int &lev, const amrex::Vector< amrex::MultiFab * > &mfs, amrex::MultiFab *xzmom_flux, amrex::MultiFab *zxmom_flux, amrex::MultiFab *yzmom_flux, amrex::MultiFab *zymom_flux, amrex::MultiFab *heat_flux, amrex::MultiFab *eddyDiffs, amrex::MultiFab *z_phys) |
template<typename FluxCalc > | |
void | compute_most_bcs (const int &lev, const amrex::Vector< amrex::MultiFab * > &mfs, amrex::MultiFab *xzmom_flux, amrex::MultiFab *zxmom_flux, amrex::MultiFab *yzmom_flux, amrex::MultiFab *zymom_flux, amrex::MultiFab *heat_flux, amrex::MultiFab *eddyDiffs, amrex::MultiFab *z_phys, const amrex::Real &dz_no_terrain, const FluxCalc &flux_comp) |
void | time_interp_sst (const int &lev, const amrex::Real &time) |
void | get_lsm_tsurf (const int &lev) |
void | update_surf_temp (const amrex::Real &time) |
void | update_mac_ptrs (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim) |
const amrex::MultiFab * | get_u_star (const int &lev) |
const amrex::MultiFab * | get_t_star (const int &lev) |
const amrex::MultiFab * | get_olen (const int &lev) |
const amrex::MultiFab * | get_mac_avg (const int &lev, int comp) |
const amrex::MultiFab * | get_t_surf (const int &lev) |
amrex::Real | get_zref () |
const amrex::FArrayBox * | get_z0 (const int &lev) |
template<typename FluxCalc > | |
void | compute_most_bcs (const int &lev, const Vector< MultiFab * > &mfs, MultiFab *xzmom_flux, MultiFab *zxmom_flux, MultiFab *yzmom_flux, MultiFab *zymom_flux, MultiFab *heat_flux, MultiFab *eddyDiffs, MultiFab *z_phys, const Real &dz_no_terrain, const FluxCalc &flux_comp) |
Public Attributes | |
FluxCalcType | flux_type {FluxCalcType::MOENG} |
ThetaCalcType | theta_type {ThetaCalcType::ADIABATIC} |
RoughCalcType | rough_type {RoughCalcType::CONSTANT} |
Private Attributes | |
bool | use_moisture |
bool | m_exp_most = false |
amrex::Real | z0_const |
amrex::Real | surf_temp |
amrex::Real | surf_heating_rate {0} |
amrex::Real | surf_temp_flux {0} |
amrex::Real | custom_ustar {0} |
amrex::Real | custom_tstar {0} |
amrex::Real | custom_qstar {0} |
amrex::Real | cnk_a {0.0185} |
amrex::Real | depth {30.0} |
amrex::Real | m_start_bdy_time |
amrex::Real | m_bdy_time_interval |
amrex::Vector< amrex::Geometry > | m_geom |
amrex::Vector< amrex::FArrayBox > | z_0 |
MOSTAverage | m_ma |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | u_star |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | t_star |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | q_star |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | olen |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | t_surf |
amrex::Vector< amrex::Vector< amrex::MultiFab * > > | m_sst_lev |
amrex::Vector< amrex::Vector< amrex::iMultiFab * > > | m_lmask_lev |
amrex::Vector< amrex::Vector< amrex::MultiFab * > > | m_lsm_data_lev |
amrex::Vector< amrex::Vector< amrex::MultiFab * > > | m_lsm_flux_lev |
Monin-Obukhov surface layer profile
van der Laan, P., Kelly, M. C., & Sørensen, N. N. (2017). A new k-epsilon model consistent with Monin-Obukhov similarity theory. Wind Energy, 20(3), 479–489. https://doi.org/10.1002/we.2017
Consistent with Dyer (1974) formulation from page 57, Chapter 2, Modeling the vertical ABL structure in Modelling of Atmospheric Flow Fields, Demetri P Lalas and Corrado F Ratto, January 1996, https://doi.org/10.1142/2975.
|
strong |
|
strong |
|
strong |
|
inlineexplicit |
void ABLMost::compute_fluxes | ( | const int & | lev, |
const int & | max_iters, | ||
const FluxIter & | most_flux | ||
) |
Function to compute the fluxes (u^star and t^star) for Monin Obukhov similarity theory.
[in] | lev | Current level |
[in] | max_iters | maximum iterations to use |
[in] | most_flux | structure to iteratively compute ustar and tstar |
void ABLMost::compute_most_bcs | ( | const int & | lev, |
const amrex::Vector< amrex::MultiFab * > & | mfs, | ||
amrex::MultiFab * | xzmom_flux, | ||
amrex::MultiFab * | zxmom_flux, | ||
amrex::MultiFab * | yzmom_flux, | ||
amrex::MultiFab * | zymom_flux, | ||
amrex::MultiFab * | heat_flux, | ||
amrex::MultiFab * | eddyDiffs, | ||
amrex::MultiFab * | z_phys, | ||
const amrex::Real & | dz_no_terrain, | ||
const FluxCalc & | flux_comp | ||
) |
void ABLMost::compute_most_bcs | ( | const int & | lev, |
const Vector< MultiFab * > & | mfs, | ||
MultiFab * | xzmom_flux, | ||
MultiFab * | zxmom_flux, | ||
MultiFab * | yzmom_flux, | ||
MultiFab * | zymom_flux, | ||
MultiFab * | heat_flux, | ||
MultiFab * | eddyDiffs, | ||
MultiFab * | z_phys, | ||
const Real & | dz_no_terrain, | ||
const FluxCalc & | flux_comp | ||
) |
Function to calculate MOST fluxes for populating ghost cells.
[in] | lev | Current level |
[in,out] | mfs | Multifabs to populate |
[in] | eddyDiffs | Diffusion coefficients from turbulence model |
[in] | flux_comp | structure to compute fluxes |
void ABLMost::get_lsm_tsurf | ( | const int & | lev | ) |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
void ABLMost::impose_most_bcs | ( | const int & | lev, |
const amrex::Vector< amrex::MultiFab * > & | mfs, | ||
amrex::MultiFab * | xzmom_flux, | ||
amrex::MultiFab * | zxmom_flux, | ||
amrex::MultiFab * | yzmom_flux, | ||
amrex::MultiFab * | zymom_flux, | ||
amrex::MultiFab * | heat_flux, | ||
amrex::MultiFab * | eddyDiffs, | ||
amrex::MultiFab * | z_phys | ||
) |
Wrapper to impose Monin Obukhov similarity theory fluxes by populating ghost cells.
[in] | lev | Current level |
[in,out] | mfs | Multifabs to populate |
[in] | eddyDiffs | Diffusion coefficients from turbulence model |
void ABLMost::time_interp_sst | ( | const int & | lev, |
const amrex::Real & | time | ||
) |
void ABLMost::update_fluxes | ( | const int & | lev, |
const amrex::Real & | time, | ||
int | max_iters = 25 |
||
) |
Wrapper to update ustar and tstar for Monin Obukhov similarity theory.
[in] | lev | Current level |
[in] | max_iters | maximum iterations to use |
|
inline |
|
inline |
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost().
FluxCalcType ABLMost::flux_type {FluxCalcType::MOENG} |
Referenced by ABLMost().
|
private |
|
private |
|
private |
Referenced by ABLMost(), and update_surf_temp().
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost().
|
private |
Referenced by get_mac_avg(), get_zref(), and update_mac_ptrs().
|
private |
Referenced by ABLMost().
|
private |
|
private |
Referenced by ABLMost(), and get_olen().
|
private |
Referenced by ABLMost().
RoughCalcType ABLMost::rough_type {RoughCalcType::CONSTANT} |
Referenced by ABLMost().
|
private |
Referenced by ABLMost(), and update_surf_temp().
|
private |
Referenced by ABLMost(), and update_surf_temp().
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost(), and get_t_star().
|
private |
Referenced by ABLMost(), get_t_surf(), and update_surf_temp().
ThetaCalcType ABLMost::theta_type {ThetaCalcType::ADIABATIC} |
Referenced by ABLMost().
|
private |
Referenced by ABLMost(), and get_u_star().
|
private |
Referenced by ABLMost().
|
private |
Referenced by ABLMost().