-
Notifications
You must be signed in to change notification settings - Fork 83
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
wall model refactor #335
Merged
Merged
wall model refactor #335
Changes from all commits
Commits
Show all changes
4 commits
Select commit
Hold shift + click to select a range
File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,164 @@ | ||
#ifndef ShearStress_H | ||
#define ShearStress_H | ||
|
||
#include "amr-wind/wind_energy/MOData.H" | ||
|
||
/** | ||
* \defgroup shear_stress Shear Stress Wall Models | ||
* | ||
* ShearStress contains functions to compute velocity and temperature shear | ||
* stress wall models the default is the Moeng wall model specifying the wall | ||
* model is done through the input file using ABL.wall_shear_stress_type options | ||
* include "constant", "local", "Schumann", and "Moeng" | ||
* | ||
* \ingroup we_abl | ||
*/ | ||
|
||
struct ShearStressConstant | ||
{ | ||
explicit ShearStressConstant(const amr_wind::MOData& mo) | ||
: utau2(mo.utau * mo.utau) | ||
, u_mean(mo.vel_mean[0]) | ||
, v_mean(mo.vel_mean[1]) | ||
, wspd_mean(mo.vmag_mean) | ||
, theta_mean(mo.theta_mean) | ||
, theta_surface(mo.surf_temp) | ||
, term1(mo.utau * mo.kappa / mo.phi_h()) | ||
{} | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_x(amrex::Real /* u */, amrex::Real /* wspd */) const | ||
{ | ||
return u_mean / wspd_mean * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_y(amrex::Real /* u */, amrex::Real /* wspd */) const | ||
{ | ||
return v_mean / wspd_mean * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_theta(amrex::Real /* wspd */, amrex::Real /* theta */) const | ||
{ | ||
return term1 * (theta_mean - theta_surface); | ||
}; | ||
|
||
amrex::Real utau2; | ||
amrex::Real u_mean; | ||
amrex::Real v_mean; | ||
amrex::Real wspd_mean; | ||
amrex::Real theta_mean; | ||
amrex::Real theta_surface; | ||
amrex::Real term1; | ||
}; | ||
|
||
struct ShearStressLocal | ||
{ | ||
explicit ShearStressLocal(const amr_wind::MOData& mo) | ||
: utau2(mo.utau * mo.utau) | ||
, theta_surface(mo.surf_temp) | ||
, term1(mo.utau * mo.kappa / mo.phi_h()) | ||
{} | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_x(amrex::Real u, amrex::Real wspd) const | ||
{ | ||
return u / amrex::max(wspd, small_vel) * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_y(amrex::Real v, amrex::Real wspd) const | ||
{ | ||
return v / amrex::max(wspd, small_vel) * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_theta(amrex::Real /* wspd */, amrex::Real theta) const | ||
{ | ||
return term1 * (theta - theta_surface); | ||
}; | ||
|
||
amrex::Real utau2; | ||
amrex::Real theta_surface; | ||
amrex::Real term1; | ||
amrex::Real small_vel{1.0e-6}; | ||
}; | ||
|
||
struct ShearStressSchumann | ||
{ | ||
explicit ShearStressSchumann(const amr_wind::MOData& mo) | ||
: utau2(mo.utau * mo.utau) | ||
, wspd_mean(mo.vmag_mean) | ||
, theta_surface(mo.surf_temp) | ||
, term1(mo.utau * mo.kappa / mo.phi_h()) | ||
{} | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_x(amrex::Real u, amrex::Real /* wspd */) const | ||
{ | ||
return u / wspd_mean * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_y(amrex::Real v, amrex::Real /* wspd */) const | ||
{ | ||
return v / wspd_mean * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_theta(amrex::Real /* wspd */, amrex::Real theta) const | ||
{ | ||
return term1 * (theta - theta_surface); | ||
}; | ||
|
||
amrex::Real utau2; | ||
amrex::Real wspd_mean; | ||
amrex::Real theta_surface; | ||
amrex::Real term1; | ||
}; | ||
|
||
struct ShearStressMoeng | ||
{ | ||
explicit ShearStressMoeng(const amr_wind::MOData& mo) | ||
: utau2(mo.utau * mo.utau) | ||
, u_mean(mo.vel_mean[0]) | ||
, v_mean(mo.vel_mean[1]) | ||
, wspd_mean(mo.vmag_mean) | ||
, theta_surface(mo.surf_temp) | ||
, theta_mean(mo.theta_mean) | ||
, term1(mo.utau * mo.kappa / (mo.vmag_mean * mo.phi_h())) | ||
{} | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_x(amrex::Real u, amrex::Real wspd) const | ||
{ | ||
return ((u - u_mean) * wspd_mean + wspd * u_mean) / | ||
(wspd_mean * wspd_mean) * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_vel_y(amrex::Real v, amrex::Real wspd) const | ||
{ | ||
return ((v - v_mean) * wspd_mean + wspd * v_mean) / | ||
(wspd_mean * wspd_mean) * utau2; | ||
}; | ||
|
||
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | ||
calc_theta(amrex::Real wspd, amrex::Real theta) const | ||
{ | ||
const amrex::Real num1 = (theta - theta_mean) * wspd_mean; | ||
const amrex::Real num2 = (theta_mean - theta_surface) * wspd; | ||
return term1 * (num1 + num2); | ||
}; | ||
|
||
amrex::Real utau2; | ||
amrex::Real u_mean; | ||
amrex::Real v_mean; | ||
amrex::Real wspd_mean; | ||
amrex::Real theta_surface; | ||
amrex::Real theta_mean; | ||
amrex::Real term1; | ||
}; | ||
|
||
#endif /* ShearStress_H */ |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
May be add doxygen docs and assign them to ABL group?