e2dpolar class

class e2dpolar.main.e2dpolar(**kwargs)[source]

Bases: object

A python class for executing and processing 2D polars.

The class relies on pyellipsys for setting up and running the CFD simulations. All other methods are python-based, gathering a number of existing fortran and bash scripts. The class attributes act as basis for all methods, meaning that they are essentially are the inputs, which should be modified by the user accordingly. As this method is inheriting all its functionality from hypgrid2d and ellipsys2d, it also inherits their input interfaces. Which inputs are inherited is indicated herein and for a full list of input options consult the respective tools’ documentation.

main_dir

main directory of polar computations in which the grid and run directories are created

(default: current working directory)

Type:

string

run_dir_ext

run folder name extension ie Re<>e6_<case> + run_dir_ext, add an ‘_’ if you like some space

(default: ‘’)

Type:

string

pf_coords

profile coordinates, first column x and second y. Note that the direction of the data is assumed to run clockwise (with z pointing out of the plane); from the TE on the pressure side over LE to the suction side TE.

(default: None)

Type:

2D-array(float)

pf_te_thick_lim

set a thickness limit for the trailing edge. If original coordinates should not be modified set to None. The float sets the desired limit value followed by a string defining how the limit should be enforced. Possible are ‘equal’, ‘upper’, ‘lower’ ie upper limit etc …

(default: None)

Type:

list or None

pf_pgl_input

inputs to the PGL AirfoilShape() class used for surface meshing

(default: {‘even’: False,

‘dist’: None, ‘dLE’: False, ‘dTE’: ‘ds_TE’, ‘close_te’: True, ‘linear’: False})

Type:

float

pf_norm

normalize coordinates with chord length, where the chord length is determined as the maximum distance between TE and LE by PGL using an optimization. Should have little influence if already normalized.

(default: True)

Type:

bool

pf_rotate_align

derotate aerofoil such that LE and TE are on the horizon. This is important when it comes to the definiton of the AoA as it is computed with respect to the horizontal.

(default: True)

Type:

bool

grid_type

mesh topology (hypgrid2d input format)

(default: ‘omesh-bcdyn’)

Type:

string

grid_ni

no. of grid points along surface, multiple of 32 and grid block size

(default: 512)

Type:

int

grid_nj

no. of grid points in normal direction, multiple of 32 and grid block size

(default: 256)

Type:

int

grid_bsize

block size in which to partition the grid, must be multiple of ni and nj

(default: 128)

Type:

int

grid_domain_height

distance of outer domain boundaries from profile

(default: 45.)

Type:

float

grid_first_cell_size

first grid cell height on profile surface, unused if grid_normal_distribution is specified, however it is used in the computation of y+

(default: 1.e-7)

Type:

float

grid_normal_distribution

define control points for normal mesh growth, (hypgrid2d input) Integer(1), {Real(1), Real(1), Real(1)} …. {Real(1), Real(1), Real(1))}, (Nr, {S(1), ds(1), I(1)} …… ,{S(Nr) ,ds(Nr) ,I(Nr)}) A command allowing a detailed control of the normal spacing. Where Nr is the number of triplets {S, ds, I}. Here S is the normalized curve length s divided by the target height of the domain, ds is the cell size at that curve length in [m], and I is the normalized cell index defined by (K-1)/(Layers-1). At minimum two triplets must be specified, the first must be (0.0 ds(1) 0.0) and the last must be (1.0 ds(Nr) 1.0). If ds=-1 the spacing is found automatically.

(default: [3, 0.0, 1.e-7, 0.0, 1.875e-2, -1, 0.5, 1.0, 2.4, 1.0]

Type:

list

grid_stretching_function

normal stretching function (hypgrid2d input)

(default: ‘tanh’)

Type:

string

grid_volume_blend

volume blend factor (hypgrid2d input)

(default: 0.2)

Type:

float

grid_dissipation_factor

grid dissipation (hypgrid2d input)

(default: 1.0)

Type:

float

grid_wake_contraction

control for pulling grid points into wake (hypgrid2d input)

(default: 0.3)

Type:

float

h2d_add_inputs

add any additional HypGrid2D inputs; more for superusers! These inputs are appended to the H2D input file and thus overwrite any other inputs as it is read in sequentially.

Type:

dict

Re

Reynolds number

(default : 3e6)

Type:

float

aoas

angles-of-attack to simulate

Type:

1D-array(float)

case

flow case to simulate ‘tr’: transitional boundary layer ‘tu’: turbulent boundary layer

(default: ‘tr’)

Type:

string

ti

turbulence intensity at profile (not percent), only important for transitional bl sims

(default: 0.001)

Type:

float

turbulence_model

turbulence model (ellipsys2d input)

(default: ‘komega’)

Type:

string

ko_version

k-omega version (ellipsys2d input)

(default: ‘sst’)

Type:

string

koSST_a1

SST model constant, active in the SST limiter. (ellipsys2d input)

(default: 0.31)

Type:

float

koSST_hellstenlaine

acitvate Hellsten & Laine SST inhibiting function in the viscous sublayer. Recommended as Bradshaw’s relation is not valid in the viscous sublayer and thus the SST limiter should not be active in this region. (ellipsys2d input)

(default: ‘true’)

Type:

string [‘true’/’false’]

ko_katolaunder

Kato-Launder tke production limiter, for simulations with correlation based transition model by default active. (ellipsys2d input)

(default: None)

Type:

string [‘true’/’false’]

transition_model

transition model (ellipsys2d input) Linear stability theory type models: ‘drela’, ‘drela_bp’, ‘orrsommerfeld’, ‘orrsommerfeld_bp’ Correlation function model: ‘gamma-Retheta’

(default: ‘drela’)

Type:

string

tr_correlation_functions

definition of correlation functions, if correlation based transition model used (ellipsys2d input)

(default: ‘menter’)

Type:

string

tr_nfac

critical amplification factor. If not set (default), then computed from TI by bounded Mack’s relation (ellipsys2d input)

(default: None)

Type:

float

tr_start

iteration at which to start the transition model (ellipsys2d input)

(default: 100)

Type:

int

tr_level

grid level at which to start transition model (ellipsys2d input)

(default: 2)

Type:

int

tr_trip

trip boundary layer (ellipsys2d input)

(default: False)

Type:

bool

tr_trip_sss

suction side trip curve location, dimensional

Type:

float

tr_trip_sps

pressure side trip curve location, dimensional

Type:

float

rough

switch for roughness boundary condition (ellipsys2d input)

(default: False)

Type:

bool

roughness_model

roughness model enforcing rough boundary condition (ellipsys2d input)

(default: ‘knopp’)

Type:

string

rough_patches

list holding dictionaries describing rough regions ie. [{‘start’: [location(float), loc. definition(string)],

‘end’: [location(float), loc. definition(string)], ‘height’: [float, definition(string)]]}, ….].

For each rough region define the start and end position and give the definition of the location. This can either be chordwise ‘x’ or in terms of the curve length from the ‘s’ LE. Defined from LE, +ve increasing towards suction side. Also give the roughness height desired for the region. You can choose ‘d0’, the artificial offset from the wall (Knopp 2009) also called roughness length or ‘ks’ for Nikuradse’s sandgrain roughness. For E2D it is always translated into ‘d0’ ie d0=0.03*ks for consistency. Note that all inputs are dimesional.

(default: [{‘start’: [0.01, ‘s’],

‘end’: [0.06, ‘s’], ‘height’: [1e-4, ‘d0’]}]

Type:

list(dict)

rough_else

roughness outside rough region (needed otherwise the model fails). Set extremely small, the result is identical to the one obtained with the turbulent boundary layer model.

(default: {‘height’: [1e-25, ‘d0’]})

Type:

dict

ambient_turb

switching on sustaining terms for ambient k and omega. The k-omega model exhibits large dissipation of omega and tke from the inflow plane, as far from the aerofoil only the destruction terms are active. This leads to a variation in the viscosity ratio (turbulent/molecular) throughout the domain, that depends on the distance from the inlet, which is undesirable. Spalart and Rumsey (2007) showed that sustaining terms can be be implemented without negative impact and by default these are active. The ambient values of omega and tke are computed as recommended by them by default. Inlet and farfield values are equivalent to _amb.

(default: 5*U*c)

Type:

bool

omega_amb

ambient dissipation. The k-omega model exhibits large dissipation of omega and tke from the inflow plane, as far from the aerofoil only the destruction terms are active. This leads to a variation in the viscosity ratio (turbulent/molecular) throughout the domain, that depends on the distance from the inlet, which is undesirable. Spalart and Rumsey (2007) showed that sustaining terms can be be implemented without negative impact and by default these are active. The ambient values of omega and tke are computed as recommended by them by default. Inlet and farfield values are equivalent to _amb.

(default: 5*U*c)

Type:

float

tke_amb

ambient turbulent kinetic energy, see above for description.

(default: 1e-6*U**2)

Type:

float

start_glevel

starting grid level

(default: 4)

Type:

int

mstep

iterations per grid level [mstep[gl=1], mstep[gl=2], …]

(default: [8000, 2000, 2000, 500])

Type:

list(int)

mstepp

number of subiterations

(default: 6)

Type:

int

reslim

convergence limit per grid level for momentum equations [res[gl=1], res[gl=2], …]

(default: [1.e-6, 1.e-5, 1.e-5, 1.e-5])

Type:

list(float)

reslimp

convergence limit for pressure correction

(default: 1.e-1)

Type:

float

diff_scheme

specifies the convective scheme for each grid level [diff[gl=1], diff[gl=2], …]

(default: [‘quick’, ‘suds’, ‘suds’, ‘suds’])

Type:

list(float)

pres_corr

pressure correction scheme

(default: ‘simple’)

Type:

float

interpolation_order

gives the interpolation order for face values

(default: 2)

Type:

int

relaxu

underrelaxation factor for momentum equations

(default: 0.7)

Type:

float

relaxp

underrelaxation factor for pressure correction equation

(default: 0.1)

Type:

float

relaxtr

underrelaxation factor for transition model

(default: 0.5)

Type:

float

relaxturb

underrelaxation factor for turbulence model

(default: 0.5)

Type:

float

turbcrossterms

a logical parameter controlling whether to include the so called crossterms in the k and omega and transport equations for transition.

(default: ‘false’)

Type:

string [‘true’, ‘false’]

e2d_add_inputs

add any additional EllipSys2D inputs; more for superusers! These inputs are appended to the E2D input file and thus overwrite any other inputs as it is read in sequentially.

example: e2d_add_inputs = [

(‘project’, ‘grid’), (‘mstep’, 1000), (‘refpoint’, [30, 30, 30, 20]), (‘katolaunder’, ‘false’)]

(default: None)

Type:

list(tuple)

cpu_time

time needed for computation of each angle-of-attack format: hh:mm:ss

(default: ‘00:10:00’)

Type:

string

nodes

requested number of computational nodes

(default: 1)

Type:

int

procs

requested number of processors

(default: 1)

Type:

int

queue

cluster queue name

(default: ‘workq’)

Type:

string

glevels

grid levels for which to extract polars

(default: np.array([1]))

Type:

1D-array(int)

pkl_output

write pickle output during post-processing

(default: True)

Type:

bool

iter_avg

no. of iterations to average over to get aero coefficients

(default: 100)

Type:

int

pitch_loc

location about which to compute the pitching moment as a fraction of chord.

(default: 0.25)

Type:

float

tr_limit

threshold above which the flow is assumed turbulent and marked as transitioned. The LST models use the intermittency, whereas the correlation model use the turbulence index by Carnes&Coder (2021, 10.1007/s10494-021-00288-5, Eq.(19)) as gamma is unreliable.

(default: LST:0.025, Corr:0.5)

Type:

float

bl_nout

set at which interval to output boundary layer profiles

(default: 0, ie all surface points)

Type:

int

noise_vinf

velocity used by noise model

(default: np.array([25]))

Type:

1D-array(float)

noise_model

noise model selection

(default: 41)

Type:

int

noise_pressure_select

SPL (1) or surface pressure (2)

(default: 1)

Type:

int

noise_octave_band

octave band

(default: 3)

Type:

int

noise_chord

chord used in noise model

(default: 1.)

Type:

float

noise_span

span used in noise model

(default: 1.)

Type:

float

noise_temperature

temparature [C]

(default: 15.)

Type:

float

noise_atmos_pressure

pressure [Pa]

(default: 1013e3)

Type:

float

noise_observer

observer’s location relative to TE and mid-span [X1(Chordwise =0 at TE), X2(Above TE), X3(Spanwise)]

(default: np.array([0.0, 1.0, 0.0]))

Type:

1D-array(float)

check_finish(wait=True)[source]

Method checking whether the CFD simulations have finished, stalling script execution before another mehtod is executed that relies on the simulation outputs.

compute_coeffs(aoa, file='grid.1force', nr_average=100, path='', pitch_loc=0.25)[source]

Method computing aerodynamic coefficients from a steady-state CFD simulation at a single angle-of-attack.

Parameters:
  • aoa (float) – angle-of-attack in degrees

  • file (string) –

    name of the force output file from e2d

    (default: ‘grid.1force’)

  • nr_average (int) –

    iteration number over which to average the aerodynamic coefficients

    (default: 100)

  • path (string) –

    path to file specified previously, total file name is given by path + file

    (default: ‘’)

  • pitch_loc (float) –

    location about which to compute the pitching moment as a fraction of the chord.

    (default: 0.25)

Returns:

  • out (dict) – dictonary holding aerodynamic coefficients and parameters used for computing the normalisation of the forces.

  • out[‘aoa’] (float) – angle-of-attack

  • out[‘rho’] (float) – air density

  • out[‘vel’] (float) – velocity

  • out[‘chord’] (float) – airfoil chord

  • out[‘pitch_loc_%’] (float) – airfoil pitch centre in chord percentage

  • out[‘cl’] (float) – airfoil sectional lift coeffcient

  • out[‘cd’] (float) – airfoil sectional drag coeffcient

  • out[‘clv’] (float) – airfoil viscous lift coeffcient (non-pressure contribution)

  • out[‘cdv’] (float) – airfoil viscous drag coeffcient (non-pressure contribution)

  • out[‘cm’] (float) – airfoil sectional moment coeffcient

ext_bl(xqin, ff=None, prdat=None, varnames=['vtvn'], jmax=None, pf_spline='cubic', output=False, outputname='bl_profiles', afpoints=None, aflimits=None, pflimits=None, interp='grid', scipy_interp_method='linear')[source]

Extract boundary layer profiles @ locations xq

Specify a numpy array of locations, xq, where to extract boundary layer profiles, negative values indicate extraction on the pressure side, positive on the suction side. The script relies on inputting a flowfield directory with the flow data. The extraction is done at height levels similar to those of the mesh itself, however profiles are extracted normal to the surface. The profiles are automatically plotted.

Parameters:
  • xq (1D-array(float)) – specify where to BL profiles where xq < 0 indicates extraction on the pressure side, and xq > 0 on the suction side. The number of profiles is not limited, nor are their locations.

  • ff (dict) –

    output dictionary from the self.ext_flowfield() method, which holds the mesh data. It tries to load it from the attribute self.flowfield if it is not directly specified here, so make sure to run self.ext_flowfield() beforehand.

    (default: taken from self.flowfield)

  • varname (string) –

    select the variable from which the profiles should be extracted. Names are the same as in the flowfield dict except if you specify ‘vtvn’ this routine automatically computes the wall parallel (vt) and normal (vn) components. For this flowfield needs to hold ‘u’ and ‘v’.

    (default: ‘vtvn’)

  • jmax (float) –

    specify the number of normal grid layer to extract. If not specified the routine will try to estimate it from the no. of normal grid layers self.grid_nj (it will try to extract the first 60%).

    (default: None)

  • pf_spline (float) –

    selects the spline type to use for interpolating profile based qunatities. The ‘cubic’ method is more robust than the NaturalCubicSpline used as default by PGL and is used here instead. The difference is negligble.

    (default: ‘cubic’)

  • output (bool) –

    should output be written. The bl_profiles dict is saved as a pickle and the profiles a plotted in <outputname>.pdf.

    (default: False)

  • outputname (string) –

    name of output data and plots, give the full path if you would like to store it in another folder.

    (default: ‘bl_profiles’)

  • afpoints (1D-array(float)) –

    auxillary airfoil coordinates. By default the first normal mesh layer is taken to be the airfoil surface, however if afpoints is specified the wall normal direction vector is computed from afpoints. This is useful when extracting BLs from eroded leading edges for instance.

    (default: None)

  • aflimits (1D-array(float)) –

    specify plot region for airfoil [xmin, xmax, ymin, ymax]

    (default: None)

  • pflimits (1D-array(float)) –

    specify plot limits for BL profiles [xmin, xmax, ymin, ymax]

    (default: None)

  • interp (string) –

    select how to extract the BL profiles. There are two options: ‘grid’: linearly interpolate between chordwise grid node positions, which yields very accuracte results, but where the BL profiles are not normally extracted to the surface but are curving with the mesh. Mesh curvature is only large around the LE and TE, and limited close to the surface, so it is nearly identical to extracting normally. ‘normal’: truely normal extraction from the wall, however the scipy based interpoaltion can lead to interpoaltion errors, which are small, however depending on the quantity you are interested in be siginificant. The ‘grid’ interpolation as default as it is more accurate.

    (default: ‘grid’)

  • scipy_interp_method (string) –

    interpolation is performed with scipy’s griddata when interp=’normal’, this option selects the interpolator to for this operation. Linear is the safest option.

    (default: ‘linear’)

Returns:

  • bl_profiles (dict) – ouput dictionary holding all boundary layer profiles

  • bl_profiles[‘xy’] (3D-array) – profile points in mesh coordinate system, shape: [nr_profiles, nj, 2], here nj is the no. of normal layers extracted

  • bl_profiles[‘yn’] (2D-array) – wall normal BL profile coordinates, shape: [nr_profiles, nj]

  • bl_profiles[‘varname’] (2D-array) – profile values, shape: [nr_profiles, nj, (2 if ‘vtvn’)]

ext_bls(rerun_bl_ext=True, outputname='bl_dat_prof')[source]

Method for extracting the boundary layer profiles and BL quantities like thickness etc along aerofoil. This is the only model chain that is not integrated in a pythonic way and relies on some precompiled Fortran scripts and bash scripts. It first runs the postprocessor from EllipSys2D to get the boundary layer data, creating two files:

grid.BLdata: Integral BL parameters at each surface grid station grid.BLprof: Actual profiles for each station. blthick_cfd.dat: BL thickness parameters and shape factor

Parameters:
  • rerun_bl_ext (bool) –

    should the boundary layer be extracted again from restart file using the noise post processor. This process is slow!

    (default: True)

  • outputname (string) –

    file name used to save noise and boundary layer data

    (default: ‘bl_dat_prof’)

Returns:

  • self.bl (dict) – dictonary holding boundary layer data

  • self.bl[‘gl<grid_level>’][‘bl_data’] (dict) – holds surface data (length no. of surface mesh points unless self.nout > 1)

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘x’] (1D-array) – surface x coordinates

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘y’] (1D-array) – surface y coordinates

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘cf’] (1D-array) – wall friction coefficient

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘delta’] (1D-array) – BL thickness

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘delta_star’] (1D-array) – BL displacement thickness

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘theta’] (1D-array) – BL momentum thickness

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘Hk’] (1D-array) – BL shape parameter

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘U_edge’] (1D-array) – BL edge velocity

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘x_edge’] (1D-array) – BL edge x coordinate

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘y_edge’] (1D-array) – BL edge y coordinate

  • self.bl[‘gl<grid_level>’][‘bl_data’][‘cp’] (1D-array) – pressure coefficient

  • self.bl[‘gl<grid_level>’][‘bl_profiles’] (dict) – holds BL profile quantities, shape: [nr_profiles, vertical stations]

  • self.bl[‘gl<grid_level>’][‘bl_profiles’][‘yn’] (2D-array) – BL profile wall normal coordinate

  • self.bl[‘gl<grid_level>’][‘bl_profiles’][‘vt’] (2D-array) – BL profile wall tangential velocity

  • self.bl[‘gl<grid_level>’][‘bl_profiles’][‘dvtdy’] (2D-array) – BL profile velocity gradient

  • self.bl[‘gl<grid_level>’][‘bl_profiles’][‘tke’] (2D-array) – BL profile turbulent kinetic energy

  • self.bl[‘gl<grid_level>’][‘bl_profiles’][‘omega’] (2D-array) – BL profile specific dissipation

  • self.bl[‘gl<grid_level>’][‘noise’] (dict) – holds noise outputs

ext_flowfield(aoa, ext_fields=['u', 'v', 'p'], glevel=1, write_plot3d=False, plot=True, plt_fields=['u', 'v', 'p'])[source]

Extract vetex-based flowfield variables from restart file

Parameters:
  • aoa (float) – angle-of-attack for which to extract the flowfield variables

  • ext_fields (list(string)) –

    specify which variables to extract (need to be the variable names of e2d).

    (default: [‘u’, ‘v’, ‘p’])

  • glevel (int) –

    grid level

    (default: 1)

  • write_plot3d (bool) –

    write out the standard plot3d files (grid.f, grid.nam etc)

    (default: False)

  • plot (bool) –

    plot the flowfield variables. The plots are placed inside the angle-of-attack folder and are named ‘ff_<varname>_gl<grid_level>_aoa<AoA>.pdf’

    (default: True)

  • plt_fields (list(string)) –

    specify which variables to plot (need to be the variable names of e2d).

    (default: [‘u’, ‘v’, ‘p’])

Returns:

  • self.flowfield (dict) – ouput dictionary holding the flowfield data

    The keys correspond to the desired fields to be extracted and also include the grid coordinates ‘x’, ‘y’. The shape is identical to the layers speficied + 1 as vertices are extracted, so the shape is [self.ni_grid + 1, self.nj_grid + 1]. The nj layers are the wall normal layers.

  • self.flowfield[‘x’] (2D-array(float)) – mesh vertices x coordinates

  • self.flowfield[‘y’] (2D-array(float)) – mesh vertices y coordinates

  • self.flowfield[‘<varname>’] (2D-array(float)) – values of selected <varname> at vertices, the variable naming follows from the EllipSys2D convention.

ext_noise(rerun_bl_ext=True, outputname='bl_noise')[source]

Method for calculating the trailing edge noise emissions from airfoils. This is the only model chain that is not integrated in a pythonic way and relies on some precompiled Fortran scripts and bash scripts. It first runs the postprocessor from EllipSys2D to get the boundary layer data, creating two files:

grid.BLdata: Integral BL parameters at each surface grid station grid.BLprof: Actual profiles for each station. blthick_cfd.dat: BL thickness parameters and shape factor

It also creates the file ‘pp_noise.out’ storing the outputs from the process itself. If the file exists the BL extraction is not performed again unless you switch the flag aero.noise_postpro=True. The default is True so if you do not want to run the BL extraction again set it to False. The BL files act as input to the noise model. In each AoA directory the script launching the noise computations is created as ‘noise_ext.sh’, so one can run the model also directly.

From those BL files two profiles are extracted close to the TE (not at the TE for numerical reasons) and stored in:

u2_shear_l2_x2pres.dat u2_shear_l2_x2suct.dat

The run output of the noise model can be found in ext_noise_V<>.out, the noise model outputs are:

psd.dat: Power spectal density spl.dat: Sound pressure level spltot.dat: Integrated SPL

Parameters:
  • rerun_bl_ext (bool) –

    should the boundary layer be extracted again from restart file using the noise post processor. This process is slow!

    (default: True)

  • outputname (string) –

    file name used to save noise and boundary layer data

    (default: ‘bl_noise’)

Returns:

self.bl – dictonary holding boundary layer and noise data

Return type:

dict

ext_polar(plot=True, rerun=True)[source]

Method for extracting aerodynamic coefficients and transition locations from specified angle-of-attack range (self.aoas) and assembling of output file (grid-<grid_level>.clcd). Additionally a dictionary is returned to the user (self.gldat), which contains all extracted data for each grid level. The function also automatically plots the extracted data in the run directory.

Parameters:

plot (bool) – indicate whether to write out plots

Returns:

  • self.polar (dict) – ouput dictionary holding all coefficient and surface data for each grid level spcified by self.glevels.

    For each grid level the surface data at each anlge-of-attack (‘surface_data’) and the polars are stored (‘aero_data’). So to access the polars use self.gldat[‘gl<grid_level>’][‘aero_data’] and have a look at the keys. And to get the surface data like Cp and Cf you have to use self.gldat[‘gl<grid_level>’][‘surface_data’][‘a<AoA>].

  • self.polar[‘gl<grid_level>’][‘surface_data’][‘a<AoA>’] (dict) – profile surface data, output from self.ext_pr() for each angle-of-attack

  • self.polar[‘gl<grid_level>’][‘aero_data’] (dict) – polar data, output from self.compute_coeffs()

ext_pr(file='grid.1pr', path='')[source]

This method finds the transition locations on the suction and pressure side. It is located by finding the location where the transition factor exceeds the specified threshold (limit).

Python translation of locatetran.f by NNS.

Parameters:
  • file (string) –

    name of the surface data e2d output file

    (default:’grid.1pr’)

  • path (string) –

    path to file

    (default:’’)

Returns:

  • out (dict) – dictonary holding information on surface data

  • out[‘s’] (1D-array(float)) – airfoil running curve length

  • out[‘x’] (1D-array(float)) – airfoil x coordinates

  • out[‘y’] (1D-array(float)) – airfoil y coordinates

  • out[‘cp’] (1D-array(float)) – profile pressure coefficients

  • out[‘cf’] (1D-array(float)) – profile friction coefficients

  • out[‘gamma’] (1D-array(float)) – profile intermittency (only used for transition)

  • out[‘ru_height’] (1D-array(float)) – roughness height in input specified definition

  • out[‘s_stag’] (1D-array(float)) – stagnation point location [x location, y location, s curve location, grid index]

  • out[‘uf’] (1D-array(float)) – friction velocity

  • out[‘nu] (float) – kinematic viscosity

  • out[‘y+’] (1D-array(float)) – first cell centre wall distance in wall coordinates

get_sec(time_str)[source]

Get seconds as float from time string hh:mm:ss

get_stitched_flowfield(var)[source]

Make single block data from hyogrid2d O-mesh Each block has a ghost layer around it so that the number of vertices in each direction is bsize + 3. In hypgrid the blocks are written with the normal index, j, in the inner loop. By convention the grid grows from the pressure side TE over the LE to the suction side TE so clockwise with z pointing out of the plane. Blocks overlap at the corner vertices.

locatetran(file='grid.1pr', path='', limit=0.025)[source]

This method finds the transition locations on the suction and pressure side. It is located by finding the location where the transition factor exceeds the specified threshold (limit).

Python translation of locatetran.f by NNS.

Parameters:
  • file (string) –

    name of the transition data e2d output file

    (default:’grid.1pr’)

  • path (string) –

    path to file specified previously, total file name is given by path + file

    (default:’’)

  • limit (float) –

    amplification factor threshold above which the flow is assumed turbulent.

    (default:0.025)

Returns:

  • out (dict) – dictonary holding information on the suction side transition point (‘stp’) and the one on the pressure side (‘ptp’).

  • out[‘stp’] (1D-array(float)) – suction side transition position [x location, y location, s curve location, grid index]

  • out[‘ptp’] (1D-array(float)) – pressure side transition position [x location, y location, s curve location, grid index]

make_mesh(rerun=False)[source]

Generate a ellipsys2d ready volume grid from profile coordinates.

make_surface_mesh(rerun=True)[source]

PGL-based generation of surface mesh from profile coordinates. This routine uses the inputted self.pf_coords and creates the surface mesh from it. Note that for the CFD simulations it is always assumed that the profile data is normalized by the chord length and that they run from pressure side TE over LE to suction side TE (clockwise). This routine checks your input coordinates and ensures that they adhere to this definition. Furthermore you can open the TE to have a more realistic airfoil profile (it is sometimes also needed by some PGL routines to work.). You can inspect the raw CFD standard confirming profile coordinates in self.af.pf_in and the final mesh is inside self.af.points.

make_volume_mesh(rerun=True)[source]

Create volume grid from surface mesh. Note that the surface mesh needs to be created beforehand (self.af) by using self.make_surface_mesh() or by using the AirfoilShape class from PGL.

print_runinfo(cool_print=None)[source]
rain_impact(u_r=100, d_r=0.003, rho_r=1000.0, aoas=None, dp=None, xy=None, plot=True, xzoom=[-0.15, 0.15], outname='rain_impact')[source]

Compute rain impact on aerofoil leading edge for loaded profile

Input rain parameters and optionally AoA (otherwise self.aoas are used) and surface direction vectors.

[1] M.H. Keegan, D.H. Nash, and M.M. Stack. “On Erosion Issues Associated With the Leading Edge of Wind Turbine Blades”. In: Journal of Physics D: Applied Physics 46.38 (2013)

[2] G. Fiore, G.E.C. Camarinha Fujiwara, and M.S. Selig. “A Damage Assessment for Wind Turbine Blades From Heavy Atmospheric Particles”. In: 53rd AIAA Aerospace Sciences Meeting. American (2015)

Parameters:
  • u_r (float) –

    rain impact velocity, [m/s]

    (default: 100.)

  • d_r (float) –

    droplet diameter, [m]

    (default: 3e-3)

  • rho_r (float) –

    droplet density, [kg/m^3]

    (default: 1e3)

  • aoas (1D-array(float), Optional) –

    angles-of-attack, [deg], for which to compute rain impact

    (default: self.aoas)

  • dp (2D-array(float), Optional) –

    direction vectors between aerofoil coordinates, first column e_x, second e_y.

    (default: self.af.points)

  • xy (2D-array(float), Optional) –

    aerofoil coordinates

    (default: self.af.points)

  • plot (bool) –

    plot impact

    (default: True)

  • xzoom (list(float)) –

    zoom area in x for plotting

    (default: [-0.15, 0.15])

  • outname (string) –

    name for plot and data output file

    (default: ‘rain_impact’)

Returns:

  • out (dict) – dictonary holding rain impact output

  • out[‘aoas’] (1D-array(float)) – angles-of-attack

  • out[‘xy’] (2D-array(float)) – aerofoil coordinates

  • out[‘s’] (1D-array(float)) – aerofoil curve length

  • out[‘f_r’] (2D-array(float)) – impact force, shape: [aoas, s]

  • out[‘theta_r’] (2D-array(float)) – impact angle, shape: [aoas, s]

run_all(rerun=False, plot=True)[source]

Most complete 2D polar executable. Method creating the surface & volume mesh, launching the CFD simulations and extracts polars.

run_extract_polar(rerun=False, plot=True)[source]

Runs CFD simulations and extracts the polars once the simulations are finished.

run_polar(group_sims=True, dry_run=False, rerun=True)[source]

Method preparing - linking files and creating directories - simulations and then launching them. Depending on the runmachine it either runs the simulations locally or submits each AoA simulation as independant job.

Parameters:
  • group_sims (bool) –

    set to ‘True’ if simulations should be grouped into a single batch submission. Like this aoa simulations are launched in sequence as a single batch job.

    (default: True)

  • dry_run (bool) –

    set to ‘True’ if simulations should not be launched and only all folders and links should be created for launching the simulations. Could be useful to check the setup.

    (default: False)

  • rerun (bool) –

    set to ‘True’ if simulations should be recomputed.

    (default: True)

set_grid_input(more_inputs=None)[source]

Set input to hypgrid2d, the volume mesher

set_run_input(more_inputs=None)[source]

Set ellipsys2d inputs according to the selected flow case

update_rundir()[source]

Update the run folder name. This is required after a change in Reynolds number or flow case, otherwise the change in Re or case won’t take effect.