pyts.specModels package¶
Submodules¶
pyts.specModels.api module¶
PyTurbSim turbulence ‘spectra models’ API.
Available spectral models¶
nwtc.smooth()(alias: smooth)- The nwtc ‘smooth’ spectral model.
nwtc.nwtcup()(alias: nwtcup)- The nwtc ‘upwind’ spectral model.
hydro.tidal(alias: tidal)- The hydro ‘tidal’ spectral model.
hydro.river(alias: river)- The hydro ‘river’ spectral model.
iec.ieckai(alias: ieckai)- The IEC ‘Kaimal’ spectral model.
iec.iecvkm(alias: iecvkm)- The IEC ‘Von-Karman’ spectral model.
specModelBase- This is the base class for spectral models. To create a new spectral model, subclass this class or subclass and modify an existing spectal model.
specObj- This is the ‘spectral object’ class. All spectral model __call__
methods must take a
tsrunas input and return this class.
Example usage¶
>>> import pyts.specModels.api as sm
Create a smooth spectral model with friction velocity Ustar=1m/s and Richardson number Ri=0.6:
>>> my_spec_model=sm.smooth(1.0,0.6)
Now set a pyts.main.tsrun instance to use this spectral model:
>>> tsrun.spec=my_spec_model(tsrun)
Notes¶
For a description of the difference between ‘spectral models’ (e.g. ‘my_spec_model’ in example above) and the spectral array they output (tsrun.spec), see the Code Framework section of the PyTurbSim documentation.
pyts.specModels.base module¶
This is the turbulence spectrum package’s base module.
-
class
pyts.specModels.base.specModelBase[source]¶ Bases:
pyts.base.modelBaseA base class for TurbSim spectral models.
-
pow2_3= 0.6666667¶
-
pow5_3= 1.6666666¶
-
-
class
pyts.specModels.base.specObj(tsrun)[source]¶ Bases:
pyts.base.gridProps,pyts.base.calcObjSpectral objects contain the array (self.array) of turbulence spectra values for a specific PyTurbSim run. This class defines various shortcuts to the data.
Parameters: tsrun : tsrun type
The PyTurbSim run object in which the spectra will be used.
-
Suu¶ This is the u-component of the TKE spectrum.
-
Svv¶ This is the v-component of the TKE spectrum.
-
Sww¶ This is the w-component of the TKE spectrum.
-
flat¶ Return the partially-flattened spectral array. This flattens the spatial dimensions of the spectral array (component and frequency dimensions are retained).
For example, for a 5 x 4 spatial grid (nz=5,ny=4) with 1000 frequency values, the size of the array in this spectral object will be 3 x 5 x 4 x 1000. Performing this flatten operation will return a 3 x 20 x 1000 shaped spectral array. The ordering (row-major or column-major) is defined in the spatial grid class.
This is used for input into the functions that calculate the coherence.
-
tke¶ This is the component-wise turbulent kinetic energy.
-
pyts.specModels.hydro module¶
This module contains the turbulence models for the aquatic environment.
-
class
pyts.specModels.hydro.river(Ustar, Zref)[source]¶ Bases:
pyts.specModels.hydro.tidalRiver turbulence spectral model.
This model is based on measurements from the East River, in New York City. It is identical to the
tidalmodel, but uses different values forcoef.-
coef= array([[1.057, 3.432], [0.351, 0.546], [0.265, 0.341]], dtype=float32)¶
-
-
class
pyts.specModels.hydro.tidal(Ustar, Zref)[source]¶ Bases:
pyts.specModels.base.specModelBaseTidal Channel spectral model.
The tidal spectral model is based on measurements from Admiralty Inlet, in Puget Sound, WA.
Parameters: Ustar : float
bottom boundary friction velocity [m/s].
Zref : float
Reference height at which u’w’ stress nears zero [m].
See also
coef- The ‘fit coefficients’
Notes
This model is similar to
NWTC_stable, but uses different values forcoef, and does not support new fit-coefficients. The form of this model is:\[S_k(f) = \frac{\sigma_k^2 \mathrm{coef}[k,0]}{1+\mathrm{coef}[k,1](f/\hat{f})^{5/3}} \qquad k=0,1,2\ (u,v,w)\]Where,
\(\hat{f}=\frac{\partial \bar{u}}{\partial z}\)
\(\bar{u}\) is the mean velocity from the
profObj.\(\sigma_k^2 = U_{*}^2 \alpha_k exp(-2 z/Zref)\)
\(\alpha_k = (4.4,2.25,0.9)\)
-
__call__(tsrun)[source]¶ Create the spectral object for
tsrun.Parameters: tsrun :
tsrunA TurbSim run object.
Returns: out :
specObjAn IEC spectral object for the grid in
tsrun.
-
coef= array([[1.21, 4.3 ], [0.33, 0.5 ], [0.23, 0.26]], dtype=float32)¶
pyts.specModels.iec module¶
This module contains the IEC turbulence models.
See the Original-TurbSim user manual for more info on IEC spectral models.
-
class
pyts.specModels.iec.IECKai(IECwindtype, IECstandard, IECedition, IECturbc, ETMc=None)[source]¶ Bases:
pyts.specModels.iec.iecbaseIEC Kaimal spectral model.
Notes
The form of this model is,
\[S_k(f) = \frac{4 \sigma_k^2/\hat{f}_k}{(1+6 f/\hat{f}_k )^{5/3}} \qquad \mathrm{for}\ k=u,v,w\]Where,
-
class
pyts.specModels.iec.IECVKm(IECwindtype, IECstandard, IECedition, IECturbc, ETMc=None)[source]¶ Bases:
pyts.specModels.iec.iecbaseIEC Von-Karman spectral model
Notes
The form of this model is,
\[ \begin{align}\begin{aligned}S_u(f) = \frac{4 \sigma^2/\hat{f}}{(1+71(f/\hat{f})^2)^{5/6}}\\S_v(f) = S_w(f) = (1+189(f/\hat{f})^2\frac{2\sigma^2/\hat{f}} {(1+71 (f/\hat{f})^2)^{11/6}}\end{aligned}\end{align} \]Where,
-
class
pyts.specModels.iec.iecbase(IECwindtype, IECstandard, IECedition, IECturbc, ETMc=None)[source]¶ Bases:
pyts.specModels.base.specModelBaseThis is a base class for the IEC spectral models (IECKAI and IECVKM).
Parameters: IECwindtype : str {‘NTM’=normal,
‘xETM’=extreme turbulence, ‘xEWM1’=extreme 1-year wind, ‘xEWM50’=extreme 50-year wind, where x=wind turbine class 1, 2, or 3}
IECstandard : int {1}
Currently this only supports IECstandard==1. IECstandard == 2 and 3 correspond to small and offshore wind, respectively, and have not yet been implemented here.
IECedition : int {2,3}
This is the ‘edition’ number of the -1 standard.
IECturbc : str | float {‘A’,’B’,’C’}
The string form correspodes to the IEC turbulence ‘categories’. If a float is provided, it specifies a specific value of Turbulence Intensity.
ETMc : float, optional (2.0)
ETMc specifies the value of the ‘c’ parameter in the equation for \(\sigma_u\). It is only used when IECwindtype`=`xETM.
Notes
For further details on the IEC spectral models see the Original-TurbSim user manual.
-
IEC_Sigma(uhub)[source]¶ Calculate the default value of the standard deviation of u-component wind speed, \(\sigma\) or \(\sigma_u\).
Notes
For input
IECturbca numeric value, it simply specifies the turbulence intensity. That is,\[\sigma = \mathrm{IECturbc}\times u_{hub}\]Otherwise only IECversion == 1 is supported. In that case:
For IECedition == 2 the windtype must be ‘NTM’, in that case:
\[\sigma = Ti (\frac{15+\beta u_{hub}}{\beta + 1})\]Where,
For IECedition == 3, Ti=(0.16,0.14,0.12) for IECturbc=(a,b,c), respectively. In this case, different formulations are used for different
IECwindtype.For IECwindtype == ‘NTM’:
\[\sigma = Ti(0.75u_{hub} + 5.6)\]For IECwindtype == ‘xETM’, the value of ‘x’ in that string provides another input variable, Vref=(50,42.5,37.5) for x=1,2,3.
\[\sigma = \mathrm{ETMc} \times Ti (0.072(0.2 V_{ref}/\mathrm{ETMc}+3)(u_{hub}/\mathrm{ETMc}-4)+10)\]For IECwindtype == ‘xEWM1’ | ‘xEWM50’, the same value of Vref apply and
\[\sigma = 0.11 V_{ref}\]
-
pyts.specModels.nwtc module¶
This module contains the nwtc spectral models.
-
class
pyts.specModels.nwtc.NWTC_stable(Ustar, zL, coef=None)[source]¶ Bases:
pyts.specModels.nwtc.genNWTCThe NWTC ‘stable’ spectral model.
Parameters: Ustar : float
friction velocity [m/s].
zL : float
The z/L stability parameter [non-dimensional]
coef : array_like((3,2),dtype=float)
spectral coefficients for this model.
See also
s_coef- These are the hard-wired \(\mathrm{scoef}\) coefficients.
Notes
The specific form of this model is,
\[S_k(f) = \frac{U_{*}^2 A_k \hat{f}^{-1}\gamma}{1+B_k(f/\hat{f})^{5/3}} \qquad k=0,1,2\ (u,v,w)\]Where,
\(\gamma=(\phi_E/\phi_M)^{2/3}\)
\(\hat{f}=\frac{\bar{u}\phi_M}{z}\)
\(\phi_E=(1+2.5 (z/L)^{0.6})^{1.5}\)
\(\phi_M=1+4.7*z/L\)
\(\bar{u}\), is the mean velocity at each grid-point, (taken from the profile model)
\(A_k=\mathrm{scoef}[k,0] \mathrm{coef}[k,0]\)
\(B_k=\mathrm{scoef}[k,1] \mathrm{coef}[k,1]^{5/3}\)
-
model(z, u, comp)[source]¶ Calculate the spectral model for height z, velocity u, and velocity component comp.
-
s_coef= array([[ 79. , 263. ], [ 13. , 32. ], [ 3.5, 8.6]])¶
-
class
pyts.specModels.nwtc.NWTC_unstable(Ustar, zL, ZI, p_coefs=None, f_coefs=None)[source]¶ Bases:
pyts.specModels.nwtc.genNWTCThe NWTC ‘unstable’ spectral model.
\[S_k(f) = U_\mathrm{star}^2 G_k(f,\bar{u},z,ZI) \qquad k = u, v, w\]Where \(G_k\) is a function that depends on the frequency, \(f\), the mean velocity \(\bar{u}\), the height \(z\), and the mixing layer depth \(ZI\). The exact form of \(G_k\) can be found in this class’s ‘model’ method.
Parameters: Ustar : float
The friction velocity (at the bottom boundary).
zL : float
The ratio of the HubHt to the Monin-Obhukov length.
ZI : float
The friction boundary layer height.
p_coefs : array_like(3,2), optional
Fit coefficients for this model.
f_coefs : array_like(3,2), optional
Fit coefficients for this model.
-
L¶
-
model(z, u, comp)[source]¶ Computes the spectrum for this ‘unstable’ spectral model.
Parameters: z : array_like (nz)
Height above the surface [m].
u : array_like (nz,ny)
Mean velocity [m/s].
comp : int {0,1,2}
Index (u,v,w) of the spectrum to compute.
Returns: spec :
specObjThe spectral object which contains the ‘array’ (property) of spectra at each point in the grid.
Notes
The form of the u-compenent spectrum is:
\[S_u(f) = U_\mathrm{star}^2 \left( p_{u,1}\frac{\alpha}{1+( F_{u,1} \hat{f})^{5/3}} + p_{u,2}\frac{\beta}{( \delta_u + F_{u,2} f' )^{5/3} } \right )\]The form of the v-compenent spectrum is
\[S_v(f) = U_\mathrm{star}^2 \left( p_{v,1}\frac{\alpha}{(1 + F_{v,1} \hat{f})^{5/3}} + p_{v,2}\frac{\beta}{( \delta_v + F_{v,2} f' )^{5/3} } \right )\]The form of the w-compenent spectrum is:
\[S_w(f) = U_\mathrm{star}^2 \left( p_{w,1}\frac{\alpha}{(1 + F_{w,1} \hat{f})^{5/3}} \gamma + p_{w,2}\frac{\beta}{ 1 + F_{w,2} {f'} ^{5/3} } \right )\]Where,
\(\hat{f} = f ZI/\bar{u}\)
\(f' = f z/ \bar{u}\)
\(\alpha = ZI^{5/3}/(-\bar{u}L^{2/3})\)
\(\beta = z(1-z/ZI)^2/\bar{u}\)
\(\delta_u = 1+15 z/ZI\)
\(\delta_v = 1+2.8 z/ZI\)
\(\gamma = \frac{f'^2+0.09(z/ZI)^2}{f'^2+0.0225}\)
-
-
class
pyts.specModels.nwtc.genNWTC[source]¶ Bases:
pyts.specModels.base.specModelBaseAn abstract base class for NWTC spectral models.
-
pyts.specModels.nwtc.nwtcup(Ustar, Ri, ZI=None)[source]¶ Compute the ‘nwtcup’ spectral model.
Parameters: Ustar : float
The bottom-boundary friction velocity [m/s].
Ri : float
The Richardson number stability parameter.
ZI : float, optional
mixing layer depth [m]. Only needed for Ri<0.
Returns: specModel :
NWTC_stable(Ri>0), orNWTC_unstable(Ri<0)
-
pyts.specModels.nwtc.smooth(Ustar, Ri, ZI=None)[source]¶ Compute the ‘smooth’ spectral model.
Parameters: Ustar : float
The bottom-boundary friction velocity [m/s].
Ri : float
The Richardson number stability parameter.
ZI : float, optional
mixing layer depth [m]. Only needed for Ri<0.
Returns: specModel :
NWTC_stable(Ri>0), orNWTC_unstable(Ri<0)
Module contents¶
This is the PyTurbSim ‘spectral models’ package. This package contains pre-defined spectral models that can be implemented in a PyTurbSim run.
The wind spectral models are more fully described in the Original-TurbSim user manual.
For more information and to use this module import the api package, e.g.
>>> import pyts.specModels.api as sm