Home > TOSSH > TOSSH_code > calculation_functions > calc_BasicSet.m

calc_BasicSet

PURPOSE ^

calc_BasicSet calculates basic set of signatures.

SYNOPSIS ^

function [results] = calc_BasicSet(Q_mat, t_mat, varargin)

DESCRIPTION ^

calc_BasicSet calculates basic set of signatures.
   The basic set of signatures are designed to cover the five components 
   of a natural streamflow regime as defined by Poff et al. (1997) and 
   Richter et al. (1996): magnitude, frequency, duration, timing and rate
   of change. As Poff et al. state, these components "can be used to 
   characterize the entire range of flows and specific hydrologic 
   phenomena, such as floods or low flows, that are critical to the 
   integrity of river ecosystems".

   INPUT
   Q_mat: streamflow [mm/timestep] matrix (cell array)
   t_mat: time [Matlab datenum] matrix (cell array)
   OPTIONAL
   start_water_year: first month of water year, default = 10 (October)
   plot_results: whether to plot results, default = false

   OUTPUT
   results: struc array with all results (each signature for each time
       series and associated error strings)

   EXAMPLE
   % load example data
   data = load('example/example_data/33029_daily.mat');
   % create consistent cell arrays
   Q_mat = {data.Q};
   t_mat = {data.t};
   results = calc_BasicSet(Q_mat,t_mat);

   References
   Poff, N.L., Allan, J.D., Bain, M.B., Karr, J.R., Prestegaard, K.L.,
   Richter, B.D., Sparks, R.E. and Stromberg, J.C., 1997. The natural flow
   regime. BioScience, 47(11), pp.769-784.
   Richter, B.D., Baumgartner, J.V., Powell, J. and Braun, D.P., 1996. A
   method for assessing hydrologic alteration within ecosystems.

   Copyright (C) 2020
   This software is distributed under the GNU Public License Version 3.
   See <https://www.gnu.org/licenses/gpl-3.0.en.html> for details.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [results] = calc_BasicSet(Q_mat, t_mat, varargin)
0002 %calc_BasicSet calculates basic set of signatures.
0003 %   The basic set of signatures are designed to cover the five components
0004 %   of a natural streamflow regime as defined by Poff et al. (1997) and
0005 %   Richter et al. (1996): magnitude, frequency, duration, timing and rate
0006 %   of change. As Poff et al. state, these components "can be used to
0007 %   characterize the entire range of flows and specific hydrologic
0008 %   phenomena, such as floods or low flows, that are critical to the
0009 %   integrity of river ecosystems".
0010 %
0011 %   INPUT
0012 %   Q_mat: streamflow [mm/timestep] matrix (cell array)
0013 %   t_mat: time [Matlab datenum] matrix (cell array)
0014 %   OPTIONAL
0015 %   start_water_year: first month of water year, default = 10 (October)
0016 %   plot_results: whether to plot results, default = false
0017 %
0018 %   OUTPUT
0019 %   results: struc array with all results (each signature for each time
0020 %       series and associated error strings)
0021 %
0022 %   EXAMPLE
0023 %   % load example data
0024 %   data = load('example/example_data/33029_daily.mat');
0025 %   % create consistent cell arrays
0026 %   Q_mat = {data.Q};
0027 %   t_mat = {data.t};
0028 %   results = calc_BasicSet(Q_mat,t_mat);
0029 %
0030 %   References
0031 %   Poff, N.L., Allan, J.D., Bain, M.B., Karr, J.R., Prestegaard, K.L.,
0032 %   Richter, B.D., Sparks, R.E. and Stromberg, J.C., 1997. The natural flow
0033 %   regime. BioScience, 47(11), pp.769-784.
0034 %   Richter, B.D., Baumgartner, J.V., Powell, J. and Braun, D.P., 1996. A
0035 %   method for assessing hydrologic alteration within ecosystems.
0036 %
0037 %   Copyright (C) 2020
0038 %   This software is distributed under the GNU Public License Version 3.
0039 %   See <https://www.gnu.org/licenses/gpl-3.0.en.html> for details.
0040 
0041 % check input parameters
0042 if nargin < 2
0043     error('Not enough input arguments.')
0044 end
0045 
0046 ip = inputParser;
0047 ip.CaseSensitive = true;
0048 
0049 % required input arguments
0050 % Please input time series as a cell array of the following format:
0051 % {x_1; x_2; ...; x_n}, where each entry (1, 2, ..., n) corresponds to one
0052 % time series, e.g. from one catchment. For one catchment only, please
0053 % input {x}. Example: {Q_1; Q_2; ...; Q_n} for streamflow.
0054 addRequired(ip, 'Q_mat', @(Q_mat) iscell(Q_mat))
0055 addRequired(ip, 't_mat', @(t_mat) iscell(t_mat))
0056 
0057 % optional input arguments
0058 addParameter(ip, 'start_water_year', 10, @isnumeric) % when does the water year start? Default: 10
0059 addParameter(ip, 'plot_results', false, @islogical) % whether to plot results
0060 
0061 parse(ip, Q_mat, t_mat, varargin{:})
0062 start_water_year = ip.Results.start_water_year;
0063 plot_results = ip.Results.plot_results;
0064 
0065 % initialise arrays
0066 Q_mean = NaN(size(Q_mat,1),1);
0067 Q_mean_error_str = strings(size(Q_mat,1),1);
0068 Q5 = NaN(size(Q_mat,1),1);
0069 Q5_error_str = strings(size(Q_mat,1),1);
0070 Q95 = NaN(size(Q_mat,1),1);
0071 Q95_error_str = strings(size(Q_mat,1),1);
0072 Q_mean_monthly = NaN(size(Q_mat,1),12);
0073 Q_mean_monthly_error_str = strings(size(Q_mat,1),1);
0074 Q_7_day_min = NaN(size(Q_mat,1),1);
0075 Q_7_day_min_error_str = strings(size(Q_mat,1),1);
0076 BFI = NaN(size(Q_mat,1),1);
0077 BFI_error_str = strings(size(Q_mat,1),1);
0078 CoV = NaN(size(Q_mat,1),1);
0079 CoV_error_str = strings(size(Q_mat,1),1);
0080 x_Q_frequency = NaN(size(Q_mat,1),1);
0081 x_Q_frequency_error_str = strings(size(Q_mat,1),1);
0082 x_Q_duration = NaN(size(Q_mat,1),1);
0083 x_Q_duration_error_str = strings(size(Q_mat,1),1);
0084 HFD_mean = NaN(size(Q_mat,1),1);
0085 HFD_mean_error_str = strings(size(Q_mat,1),1);
0086 HFI_mean = NaN(size(Q_mat,1),1);
0087 HFI_mean_error_str = strings(size(Q_mat,1),1);
0088 AC1 = NaN(size(Q_mat,1),1);
0089 AC1_error_str = strings(size(Q_mat,1),1);
0090 FDC_slope = NaN(size(Q_mat,1),1);
0091 FDC_slope_error_str = strings(size(Q_mat,1),1);
0092 BaseflowRecessionK = NaN(size(Q_mat,1),1);
0093 BaseflowRecessionK_error_str = strings(size(Q_mat,1),1);
0094 
0095 % loop over all catchments
0096 for i = 1:size(Q_mat,1)
0097     
0098     [Q_mean(i),~,Q_mean_error_str(i)] = sig_Q_mean(Q_mat{i},t_mat{i});
0099     [Q5(i),~,Q5_error_str(i)] = sig_x_percentile(Q_mat{i},t_mat{i},[5]);
0100     [Q95(i),~,Q95_error_str(i)] = sig_x_percentile(Q_mat{i},t_mat{i},[95]);
0101     [Q_mean_monthly(i,:),~,Q_mean_monthly_error_str(i)] = sig_Q_mean_monthly(Q_mat{i},t_mat{i},[1:12]);
0102     [Q_7_day_min(i),~,Q_7_day_min_error_str(i)] = sig_Q_n_day_min(Q_mat{i},t_mat{i},7);
0103     [BFI(i),~,BFI_error_str(i)] = sig_BFI(Q_mat{i},t_mat{i});
0104     [CoV(i),~,CoV_error_str(i)] = sig_Q_CoV(Q_mat{i},t_mat{i});
0105     [x_Q_frequency(i),~,x_Q_frequency_error_str(i)] = sig_x_Q_frequency(Q_mat{i},t_mat{i},'low');
0106     [x_Q_duration(i),~,x_Q_duration_error_str(i)] = sig_x_Q_duration(Q_mat{i},t_mat{i},'low');
0107     [HFD_mean(i),~,HFD_mean_error_str(i)] = sig_HFD_mean(Q_mat{i},t_mat{i},'start_water_year',start_water_year);
0108     [HFI_mean(i),~,HFI_mean_error_str(i)] = sig_HFI_mean(Q_mat{i},t_mat{i},'start_water_year',start_water_year);
0109     [AC1(i),~,AC1_error_str(i)] = sig_Autocorrelation(Q_mat{i},t_mat{i},'lag',1);
0110     [FDC_slope(i),~,FDC_slope_error_str(i)] = sig_FDC_slope(Q_mat{i},t_mat{i});
0111     [BaseflowRecessionK(i),~,BaseflowRecessionK_error_str(i)] = sig_BaseflowRecessionK(...
0112         Q_mat{i},t_mat{i},'eps',0.001*median(Q_mat{i},'omitnan'));
0113     
0114 end
0115 
0116 % add results to struct array
0117 results.Q_mean = Q_mean;
0118 results.Q_mean_error_str = Q_mean_error_str;
0119 results.Q5 = Q5;
0120 results.Q5_error_str = Q5_error_str;
0121 results.Q95 = Q95;
0122 results.Q95_error_str = Q95_error_str;
0123 results.Q_mean_monthly = Q_mean_monthly;
0124 results.Q_mean_monthly_error_str = Q_mean_monthly_error_str;
0125 results.Q_7_day_min = Q_7_day_min;
0126 results.Q_7_day_min_error_str = Q_7_day_min_error_str;
0127 results.BFI = BFI;
0128 results.BFI_error_str = BFI_error_str;
0129 results.CoV = CoV;
0130 results.CoV_error_str = CoV_error_str;
0131 results.x_Q_frequency = x_Q_frequency;
0132 results.x_Q_frequency_error_str = x_Q_frequency_error_str;
0133 results.x_Q_duration = x_Q_duration;
0134 results.x_Q_duration_error_str = x_Q_duration_error_str;
0135 results.HFD_mean = HFD_mean;
0136 results.HFD_mean_error_str = HFD_mean_error_str;
0137 results.HFI_mean = HFI_mean;
0138 results.HFI_mean_error_str = HFI_mean_error_str;
0139 results.AC1 = AC1;
0140 results.AC1_error_str = AC1_error_str;
0141 results.FDC_slope = FDC_slope;
0142 results.FDC_slope_error_str = FDC_slope_error_str;
0143 results.BaseflowRecessionK = BaseflowRecessionK;
0144 results.BaseflowRecessionK_error_str = BaseflowRecessionK_error_str;
0145 
0146 end

Generated on Tue 02-Feb-2021 09:27:04 by m2html © 2005