pyTOPKAPI
logo Sponsor logos

Table Of Contents

Previous topic

Tutorial

Next topic

Bug Reports

This Page

PyTOPKAPI API Documentation

Main module

Core logic of PyTOPKAPI.

The run function in this module contains the logic to run a TOPKAPI simulation based on the parameters specified in an INI file.

pytopkapi.model.run(ini_file='TOPKAPI.ini')
Run the model with the set-up defined by ini_file.

ODE module

This module contains the routines used to solve the Ordinary Differential Equations in TOPKAPI.

pytopkapi.ode.solve_storage_eq(I, b, alpha, Vt0, Dt, solve_method=1)

Compute the final volume in a TOPKAPI water store for one time-step

This function calculates the volume in a generic store at the end of the defined time-step, given the initial volume, inflow and ODE parameters. This is done by solving the storage ODE . Todini’s quasi-analytical solution can be used as a speed up.

Parameters:

I : scalar

Constant inflow rate to store during the time-step ().

b : scalar

The b parameter in the ODE.

alpha : scalar

The exponent parameter in the ODE.

Vt0 : scalar

The water volume at the beginning of the current time-step ().

Dt : scalar

The time-step over which to solve the ODE ().

solve_method : int

Flag specifying whether to try using the quasi-analytical solution to speed up the solution of the ODE. A value of 1 (default) means the qas will be used, a value of 0 specifies use of the full RKF formulation.

Returns:

Vt1 : scalar

The estimated volume at the end of the current time-step.

pytopkapi.ode.qas(a, b, alpha, V0, delta_t, derivative=0)

Quasi-analytical solution for volume in a TOPKAPI water store

This function calculates the volume in a generic store at the end of the defined time-step, given the initial volume, inflow and ODE parameters. This is done using Todini’s quasi-analytical solution of the storage ODE .

Parameters:

a : scalar

Constant inflow rate to store during the time-step ().

b : scalar

The b parameter in the ODE.

alpha : scalar

The exponent parameter in the ODE.

V0 : scalar

The water volume at the beginning of the current time-step ().

Dt : scalar

The time-step over which to solve the ODE ().

Returns:

V1 : scalar

The estimated volume at the end of the current time-step.

class pytopkapi.ode.RKF(min_step=1.0000000000000001e-09, max_step=3600, min_tol=9.9999999999999995e-08, max_tol=0.001, init_time_step=3600)

Adaptive Runge-Kutte Fehlberg solver

Methods

compute_error(f, x, t, delta_t)
getnewdelta_t(f, x, t, max_delta_t=0)
step(f, x, t, delta_t=0)

Fluxes module

Functions required by the TOPKAPI model for the management of input and output of cells.

Note

The subroutines solving the differential equation are not in this module (see ode for more information)

pytopkapi.fluxes.Qout_computing(V_t0, V_t1_prim, a, Dt)

Compute the outflow Qout from a generic water store.

Calculate the mean outflow rate during the current time-step from a generic water store. The outflow is calculated as the difference between the inflow and rate of change in storage during the time-step.

Parameters:

V_t0 : scalar

Volume of water at the start of the current time-step ()

V_t1_prim : scalar

Volume of water at the end of the current time-step ()

a : scalar

The inflow rate to the water store during the time-step ()

Dt : scalar

The length of the current time-step ()

Returns:

Qout : scalar

The outflow rate from the water store during the time-step ()

pytopkapi.fluxes.Qout_computing2(V_t0, V_t1, b, alpha)

Compute the output flows Qout from the computed water volumes:

Returns:

Qout : scalar

The outflow rate from the water store during the time-step ()

pytopkapi.fluxes.flow_partitioning(Lambda, Qs_out, Qo_out, W, X, Xc)

Partition the outflows from soil and overland stores.

Calculates the correct partitioning of outflows from the soil and overland stores into contributions going to the downstream cell and/or the channel store of the current cell (if this exists).

Parameters:

Lambda : int

Switch indicating whether the current cell contains a channel. A value of 1 indicates a channel cell, 0 indicates no channel

Qs_out : scalar

Outflow from the soil store during the current time-step ()

Qo_out : scalar

Outflow from the overland store during the current time-step ()

W : scalar

Width of the channel ()

X : scalar

The lateral dimension of the grid-cell ()

Xc : scalar

The length of the channel cell, this can be different from X if the channel runs along the cell diagonal ()

Returns:

Q_to_next_cell : scalar

Combined outflow from soil and overland to the downstream cell ()

Q_to_channel : scalar

Combined outflow from soil and overland to the channel store ()

Q_to_channel_sub : scalar

Outflow from soil store to the channel store ()

pytopkapi.fluxes.initial_volume_channel(ar_Q, ar_W, X, ar_n_c)

Compute initial channel volume.

The initial volume ar_Vc_t0 is computed for all the channel cells, from a known initial flow.

pytopkapi.fluxes.initial_volume_soil(ar_pVs_t0, ar_Vsm)
initialize_volume_soil Compute the intial content of water (volume of water) from the initial saturated values given in the parameter file.
pytopkapi.fluxes.input_channel(ar_Qc_out, Q_to_channel, ar_cell_up)

Compute the total inflow to the channel of a channel cell.

Calculate the total inflow to the channel as the sum of inflows from upstream cells and inflows from the soil and overland stores in the current cell.

Parameters:

ar_Qc_out : (N,) Numpy array

Array of channel outflows from each cell in the catchment for the current time-step ()

Q_to_channel : scalar

Combined outflow from soil and overland to the channel store ()

ar_cell_up : list of int

List of integer indices into ar_Qc_out. The indices point to the cells upstream of the current cell.

Returns:

a_c : scalar

The total inflow to the channel store during the current time-step ()

Qc_cell_up : scalar

The contribution to the total inflow to the channel store from upstream cells during the current time-step ()

pytopkapi.fluxes.input_soil(P, Dt, X, ar_Q_to_next_cell, ar_cell_up)

Compute the total input to a soil store.

Calculate the total rate of input to a single soil store. This comprises the sum of rainfall input, subsurface contribution from upstream cells and the overland contribution from upstream cells.

Parameters:

P : scalar

Precipitation input to the cell during the current time-step ()

Dt : scalar

The length of the current time-step in seconds

X : scalar

The lateral dimension of the grid-cell ()

ar_Q_to_next_cell : (N,) Numpy array

The total contribution from each cell to it’s downstream neighbour as a result of subsurface and overland fluxes calculated during the previous timestep ()

ar_cell_up : list of int

List of integer indices into ar_Q_to_next_cell. The indices point to the cells upstream of the current cell

Returns:

a_s : scalar

The input flux to the soil store during the current time-step in

pytopkapi.fluxes.manning_depth(Q, W, X, n)

Compute Manning depth for flow Q.

Compute with the manning equation (high width hypothesis) the water depth h corresponding to the flow Q for the channel characteristics: width W, manning n. The volume is then returned V=hWX

pytopkapi.fluxes.output_soil(Vs_t0, Vs_t1_prim, Vsm, a_s, b_s, alpha_s, Dt)

Compute the outflow from and volume in a soil store.

Calculate the outflow rate and final volume in a soil store after transferring volume greater than the saturated soil moisture content to the overland store of the same model cell.

Parameters:

Vs_t0 : scalar

Volume of water in the soil store at the beginning of the current time- step ()

Vs_t1_prim : scalar

Volume of water in the soil store at the end of the current time-step as a result of combined input to and drainage from the non-linear reservoir ()

Vsm : scalar

Volume of soil moisture for the store under saturated conditions ()

a_s : scalar

The input flux to the soil store during the current time-step ()

b_s : scalar

The constant term of the non-linear differential equation

alpha_s : scalar

The dimensionless pore-size distribution parameter for the soil store

Dt : scalar

The length of the current time-step ()

Returns:

Qs_out : scalar

Rate of flow out of the soil store during the time-step ()

Vs_out : scalar

Volume remaining in the soil store at the end of the time-step ()

Evaporation module

Evaporation and evapotranspiration routines.

A selection of routines that can be used to compute the evapotranspiration losses in PyTOPKAPI.

pytopkapi.evap.evapor_channel(Vc0, ETo, W, X)

Compute the evaporation from a channel store.

Calculate the evaporation loss from a channel store. The function attempts to remove the demand from the open water potential evaporation from the channel store, if there isn’t sufficient water in the store then the total volume in the channel is removed.

Parameters:

Vc0 : scalar

Volume in the channel store before evaporation removal ()

ETo : scalar

The potential evaporation from an open water surface ()

W : scalar

Width of the channel ()

X : scalar

The length of the channel cell, this can be different from the cell dimension if the channel runs along the cell diagonal ()

Returns:

ET_channel : scalar

The actual depth of water removed from the channel store ()

Vc1 : scalar

Volume in the channel store after evaporation removal ()

pytopkapi.evap.evapot_soil_Liu_and_Todini(Vs0, Vsm, kc, ETr, X)
The evapotranspiration taken up from the soil:
  • at the rate ks*ETc ks being the soil saturation rate
  • without overland runoff infiltration
pytopkapi.evap.evapot_soil_Liu_and_Todini_ETc(Vs0, Vsm, kc, ETr, X)
The evapotranspiration taken up from the soil:
  • at the reference crop rate ETc (potential evapotranspiration rate)
  • without overland runoff infiltration
pytopkapi.evap.evapot_soil_overland(Vo0, Vs0, Vsm, kc, ETr, X)

Compute the evapotranspiration from a model cell.

The evapotranspiration loss is first taken from the overland store, if the storage in the overland store cannot satisfy the ET demand then the water is extracted from the soil store. In cases where the ET demand is greater than the available water, both the soil and overland store are totally drained.

Parameters:

Vo0 : scalar

Volume in the overland store before ET removal ()

Vs0 : scalar

Volume in the soil store before ET removal ()

Vsm : scalar

Volume of soil moisture for the store under saturated conditions ()

kc : scalar

Dimensionless crop co-efficient for the model cell

ETr : scalar

Reference crop ET for the cell during the current time-step ()

X : scalar

The lateral dimension of the grid-cell ()

Returns:

ETa : scalar

The actual ET removed from the soil store, calculated as with the relative saturation of the soil store before ET is removed ()

Vs1 : scalar

Volume in the soil store after ET removal ()

Vo1 : scalar

Volume in the overland store after ET removal ()

pytopkapi.evap.intercept_rain_ET(P, ETr, kc)
The evapotranspiration is taken from the precipitation:
  • at the reference crop rate ETc

Pretreatment module

Functions required to compute the intrinsic TOPKAPI parameters from the physical parameters

pytopkapi.pretreatment.compute_cell_param(X, ar_Xc, Dt, alpha_s, alpha_o, alpha_c, nb_cell, A_thres, W_max, W_min, ar_lambda, ar_tan_beta, ar_tan_beta_channel, ar_L, ar_Ks, ar_theta_r, ar_theta_s, ar_n_o, ar_n_c, ar_A_drained)

Compute model parameters from physical parameters.

This function uses the physically based parameters and constraints for each model cell to compute, the saturated soil moisture volume, channel width and constant terms for the differential equations of the soil, overland and channel stores.

Parameters:

X : scalar

The lateral dimension of the grid-cell ()

ar_Xc : (N,) float array

The length of the channel in a cell, this can be different from the lateral dimension of the grid cell if the channel runs along the cell diagonal ()

Dt : scalar

The length of the current time-step in seconds

alpha_s : scalar

The dimensionless pore-size distribution parameter for the soil store

alpha_o : scalar

Power co-efficient for Mannings Equation applied to the overland store the value is typically 5/3

alpha_c : scalar

Power co-efficient for Mannings Equation applied to the channel store the value is typically 5/3

nb_cell : scalar

The number of cells in the catchment

A_thres : scalar

The minimum area of upstream contributing cells required before a cell is considered to initiate a river channel ()

W_min : scalar

The minimum width of a channel ()

W_max : scalar

The maximum width of a channel ()

ar_lambda : (N,) int array

Switch indicating whether the current cell contains a channel. A value of 1 indicates a channel cell, 0 indicates no channel

ar_tan_beta : (N,) float array

The tangent of the surface slope for each cell. The surface slope affects processes in the overland and soil stores

ar_tan_beta_channel : (N,) float array

The tangent of the channel slope for each cell. This is conceptually different from the surface slope and affects the channel store

ar_L : (N,) float array

The depth of the soil store in each cell ()

ar_Ks : (N,) float array

The saturated hydraulic conductivity of each cell ()

ar_theta_r : (N,) float array

The residual soil moisture content for each cell

ar_theta_s : (N,) float array

The saturated soil moisture content for each cell

ar_n_o : (N,) float array

Manning’s roughness coefficient for overland flows in each cell

ar_n_c : (N,) float array

Manning’s roughness coefficient for channel flows in each cell

ar_A_drained : (N,) float array

The total drained area associated with each cell ()

Returns:

ar_Vsm : (N,) float array

The saturated moisture volume of the soil store for each cell ()

ar_b_s : (N,) float array

The constant term of the non differential equation for each soil store –>

ar_b_o : (N,) float array

The constant term of the non differential equation for each overland store –>

ar_W : (N,) float array

The channel width for each cell ()

ar_b_c : (N,) float array

The constant term of the non differential equation for each channel store –>

pytopkapi.pretreatment.direct_up_cell(ar_cell_label, ar_cell_down, ar_label_sort)

Calculate the upstream cells for each cell.

This function calculates the immediate upstream cells contributing flow to each cell in the catchment.

Parameters:

ar_cell_label : (N,) int array

Numbers labelling each cell

ar_cell_down : (N,) int array

The label (from ar_cell_label) associated with the cell downstream of the current cell

ar_label_sort : (N,) int array

Sorted array of cell labels

Returns:

li_cell_up : List of (M,) int arrays

A list of arrays containing the upstream cell labels for each cell in a catchment

pytopkapi.pretreatment.drained_area(ar_label_sort, li_cell_up, X)

Compute the drained area for each cell.

This function calculates the total area drained for each cell in the catchment, as the sum of it’s area and the upstream area.

Parameters:

ar_label_sort : (N,) int array

Sorted array of cell labels

li_cell_up : List of (M,) int arrays

A list of arrays containing the upstream cell labels for each cell in a catchment

X : scalar

The lateral dimension of the grid-cell ()

Returns:

ar_A_drained : (N,) float array

The total drained area associated with each cell ()

pytopkapi.pretreatment.read_cell_parameters(file_name)

Read the spatially variable cell parameters from file.

Read the file containing the physical parameters of each cell. This information governs the distributed behaviour of the model.

Parameters:

file_name : string

Full name of the parameter file (incl. path)

Returns:

ar_cell_label : (N,) int array

Numbers labelling each cell

ar_coorx : (N,) float array

The x co-ordinate of the centre of each cell (). This is the Longitude expressed in metres using an appropriate map projection

ar_coory : (N,) float array

The y co-ordinate of the centre of each cell (). This is the Latitude expressed in metres using an appropriate map projection

ar_lambda : (N,) int array

Switch indicating whether the current cell contains a channel. A value of 1 indicates a channel cell, 0 indicates no channel

ar_Xc : (N,) float array

The length of the channel in a cell, this can be different from the lateral dimension of the grid cell if the channel runs along the cell diagonal ()

ar_dam : (N,) int array

Switch indicating whether the current cell contains a dam. A value of 1 indicates a dam cell, 0 indicates no dam. The switch currently has no influence in the model

ar_tan_beta : (N,) float array

The tangent of the surface slope for each cell. The surface slope affects processes in the overland and soil stores

ar_tan_beta_channel : (N,) float array

The tangent of the channel slope for each cell. This is conceptually different from the surface slope and affects the channel store

ar_L : (N,) float array

The depth of the soil store in each cell ()

ar_Ks : (N,) float array

The saturated hydraulic conductivity of each cell ()

ar_theta_r : (N,) float array

The residual soil moisture content for each cell

ar_theta_s : (N,) float array

The saturated soil moisture content for each cell

ar_n_o : (N,) float array

Manning’s roughness coefficient for overland flows in each cell

ar_n_c : (N,) float array

Manning’s roughness coefficient for channel flows in each cell

ar_cell_down : (N,) int array

The label (from ar_cell_label) associated with the cell downstream of the current cell

ar_pVs_t0 : (N,) float array

The initial saturation of each soil store (%)

ar_Vo_t0 : (N,) float array

The initial volume of water in each overland store ()

ar_Qc_t0 : (N,) float array

The initial channel discharge for each cell, if applicable ()

ar_kc : (N,) float array

The crop co-efficient for each cell

pytopkapi.pretreatment.read_column_input(file_name, nb_cell)

read_column_input

Read the file containing data in column format: Cell1 Cell2 Cell3 ... 1.3 4.3 5.2 ... 2.3 5.6 4.2 ...

Return a matrix mat_out(nrow, ncol)

pytopkapi.pretreatment.read_global_parameters(file_name)

Read global model parameters from file.

Read the file that specifies the parameters of the model, which are globally applied to all cells in the model.

Parameters:

file_name : string

Full name of the parameter file (incl. path)

Returns:

X : scalar

The lateral dimension of the grid-cell ()

Dt : scalar

The length of the current time-step in seconds

alpha_s : scalar

The dimensionless pore-size distribution parameter for the soil store

alpha_o : scalar

Power co-efficient for Mannings Equation applied to the overland store the value is typically 5/3

alpha_c : scalar

Power co-efficient for Mannings Equation applied to the channel store the value is typically 5/3

A_thres : scalar

The minimum area of upstream contributing cells required before a cell is considered to initiate a river channel ()

W_min : scalar

The minimum width of a channel ()

W_max : scalar

The maximum width of a channel ()

pytopkapi.pretreatment.sort_cell(ar_cell_label, ar_cell_down)

Determine a suitable computation order for the cells.

Sort the cells into a valid computation order based on the distance to the outlet of the catchment. This is necessary as each cell depends on the contribution from upstream cells.

Parameters:

ar_cell_label : (N,) int array

Numbers labelling each cell

ar_cell_down : (N,) int array

The label (from ar_cell_label) associated with the cell downstream of the current cell

Returns:

ar_label_sort : (N,) int array

Sorted array of cell labels

Utils module

pytopkapi.utils.CRange(ar_x)
Returns the range of an array
pytopkapi.utils.distance(x1, y1, x2, y2)
Compute the distance between two points
pytopkapi.utils.f_axe(p, xc)
Returns the value in the array xc associated to a relative value inside [0,1]
pytopkapi.utils.find_cell_coordinates(ar_cell_label, Xoutlet, Youtlet, ar_coorx, ar_coory, ar_lambda, channel=True)

Find the label of the cell closest to (Xoutlet, Youtlet).

Find the label of the model cell containing the specified location. The co-ordinates of the location must be given in the same co-ordinate system as that specifying the model catchment.

Parameters:

ar_cell_label : (N,) int array

Numbers labelling each cell.

Xoutlet : float

The x co-ordinate of a point. This is the Longitude expressed in metres using the same projection as ar_coorx.

Youtlet : float

The y co-ordinate of a point. This is the Longitude expressed in metres using the same projection as ar_coory.

ar_coorx : (N,) float array

The x co-ordinate of the centre of each cell (m). This is the Longitude expressed in metres using a Transverse Mercator projection, but any appropriate projection can be used.

ar_coory : (N,) float array

The y co-ordinate of the centre of each cell (m). This is the Latitude expressed in metres using a Transverse Mercator projection, but any appropriate projection can be used.

ar_lambda : (N,) int array

Switch indicating whether the current cell contains a channel. A value of 1 indicates a channel cell, 0 indicates no channel.

channel : boolean (default=True)

Allows cells with or without channels to be chosen.

Returns:

cell_outlet : int

The label for the cell closest to the defined location.

pytopkapi.utils.find_dist_max(ar_coorx, ar_coory)
Compute the maximum distance between several points defined by their coordinates ar_coorx and ar_coory
pytopkapi.utils.mov_avg(ar_float, period)
period is a multiple of 2
pytopkapi.utils.string(integer, len_str_out)
From a given integer, return an string of length len_str_out completed by zero Example: ut.string(1,3)–>‘001’

Arcfltgrid module

Utilities for reading and plotting ArcGIS binary files.

This module contains some simple functions to make it easier to read and plot the data contained in the binary grid files produced by ArcGIS.

pytopkapi.arcfltgrid.plot(bin_name, fig_name, title='Raster Plot')
Create a plot of the data in an ArcGIS binary file.
pytopkapi.arcfltgrid.read(bingrid_name)

Read the data field and headers from an ArcGIS binary grid

This function reads the header and data from the ArcGIS binary data files produced by the “Raster to Float” tool in ArcGIS 9.1

pytopkapi.arcfltgrid.read_bin(fname)

Read data from a ArcGIS binary file into an array.

Read the data from a binary file created by ArcGIS into a Numpy array. The file is expected to be in binary format with floating point precision. e.g. “.flt” extension.

pytopkapi.arcfltgrid.read_headers(bingrid_name)

Read the ascii headers of the ArcGIS binary grid file

The headers have the following format:

ncols 62 nrows 121 xllcorner -288595.47161281 yllcorner -3158065.5722693 cellsize 1000 NODATA_value -9999 byteorder LSBFIRST

Creative Commons License Copyright 2008-2010, Theo Vischel, Scott Sinclair & Geoff Pegram
Last updated on Nov 01, 2010 | Created using Sphinx