smrf.utils.wind package

Submodules

smrf.utils.wind.model module

class smrf.utils.wind.model.wind_model(x, y, dem, nthreads=1)[source]

Bases: object

Estimating wind speed and direction is complex terrain can be difficult due to the interaction of the local topography with the wind. The methods described here follow the work developed by Winstral and Marks (2002) and Winstral et al. (2009) [19] [20] which parameterizes the terrain based on the upwind direction. The underlying method calulates the maximum upwind slope (maxus) within a search distance to determine if a cell is sheltered or exposed.

The azimuth A is the direction of the prevailing wind for which the maxus value will be calculated within a maximum search distance dmax. The maxus (Sx) parameter can then be estimated as the maximum value of the slope from the cell of interest to all of the grid cells along the search vector. The efficiency in selection of the maximum value can be increased by using the techniques from the horizon functio which calculates the horizon for each pixel. Therefore, less calculations can be performed. Negative Sx values indicate an exposed pixel location (shelter pixel was lower) and positive Sx values indicate a sheltered pixel (shelter pixel was higher).

After all the upwind direction are calculated, the average Sx over a window is calculated. The average Sx accounts for larger lanscape obsticles that may be adjacent to the upwind direction and affect the flow. A window size in degrees takes the average of all Sx.

Parameters
  • x – array of x locations

  • y – array of y locations

  • dem – matrix of the dem elevation values

  • nthread – number of threads to use for maxus calculation

bresenham(start, end)[source]

Python implementation of the Bresenham algorthim to find all the pixels that a line between start and end interscet

Parameters
  • start – list of start point

  • end – list of end point

Returns

Array path of all points between start and end

find_maxus(index)[source]

Calculate the maxus given the start and end point

Parameters

index – index to a point in the array

Returns

maxus value for the point

hord(x, y, z)[source]

Calculate the horizon pixel for all z This mimics the simple algorthim from Dozier 1981 but was adapated for use in finding the maximum upwind slope

Works backwards from the end but looks forwards for the horizon

Parameters
  • x – x locations for the points

  • y – y locations for the points

  • z – elevations for the points

Returns

array of the horizon index for each point

ismember(a, b)[source]
maxus(dmax, inc=5, inst=2, out_file='smrf_maxus.nc')[source]

Calculate the maxus values

Parameters
  • dmax – length of outlying upwind search vector (meters)

  • inc – increment between direction calculations (degrees)

  • inst – Anemometer height (meters)

  • out_file – NetCDF file for output results

Returns

None, outputs maxus array straight to file

maxus_angle(angle, dmax)[source]

Calculate the maxus for a single direction for a search distance dmax

Note

This will produce different results than the original maxus program. The differences are due to:

  1. Using dtype=double for the elevations

  2. Using different type of search method to find the endpoints.

However, if the elevations are rounded to integers, the cardinal directions will reproduce the original results.

Parameters
  • angle – middle upwind direction around which to run model (degrees)

  • dmax – length of outlying upwind search vector (meters)

Returns

array of maximum upwind slope values within dmax

Return type

maxus

output(ptype, index)[source]

Output the data into the out file that has previously been initialized.

Parameters
  • ptype – type of calculation that will be saved, either ‘maxus’ or ‘tbreak’

  • index – index into the file for where to place the output

output_init(ptype, filename, ex_att=None)[source]

Initialize a NetCDF file for outputing the maxus values or tbreak

Parameters
  • ptype – type of calculation that will be saved, either ‘maxus’ or ‘tbreak’

  • filename – filename to save the output into

  • ex_att – extra attributes to add

tbreak(dmax, sepdist, inc=5, inst=2, out_file='smrf_tbreak.nc')[source]

Calculate the topobreak values

Parameters
  • dmax – length of outlying upwind search vector (meters)

  • sepdist – length of local max upwind slope search vector (meters)

  • angle – middle upwind direction around which to run model (degrees)

  • inc – increment between direction calculations (degrees)

  • inst – Anemometer height (meters)

  • out_file – NetCDF file for output results

Returns

None, outputs maxus array straight to file

windower(maxus_file, window_width, wtype)[source]

Take the maxus output and average over the window width

Parameters
  • maxus_file – location of the previously calculated maxus values

  • window_width – window width about the wind direction

  • wtype – type of wind calculation ‘maxus’ or ‘tbreak’

Returns

New file containing the windowed values

smrf.utils.wind.wind_c module

Cython wrapper to the underlying C code

20160816

smrf.utils.wind.wind_c.call_maxus()

Call the function maxus_grid in calc_wind.c which will iterate over the grid within the C code

Parameters
  • - [nsta x nsta] matrix of distances between stations (ad) –

  • - [ngrid x nsta] matrix of distances between grid points and stations (dgrid) –

  • - [nsta] array of station elevations (elevations) –

  • weights (return) –

  • - number of threads to use in parallel processing (nthreads) –

Out:

weights changed in place

20160222 Scott Havens

smrf.utils.wind.wind_c module

Cython wrapper to the underlying C code

20160816

smrf.utils.wind.wind_c.call_maxus()

Call the function maxus_grid in calc_wind.c which will iterate over the grid within the C code

Parameters
  • - [nsta x nsta] matrix of distances between stations (ad) –

  • - [ngrid x nsta] matrix of distances between grid points and stations (dgrid) –

  • - [nsta] array of station elevations (elevations) –

  • weights (return) –

  • - number of threads to use in parallel processing (nthreads) –

Out:

weights changed in place

20160222 Scott Havens

Module contents