Home > TOSSH_development > 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)

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)

   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)
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 %
0015 %   OUTPUT
0016 %   results: struc array with all results (each signature for each time
0017 %       series and associated error strings)
0018 %
0019 %   EXAMPLE
0020 %   % load example data
0021 %   data = load('example/example_data/33029_daily.mat');
0022 %   % create consistent cell arrays
0023 %   Q_mat = {data.Q};
0024 %   t_mat = {data.t};
0025 %   results = calc_BasicSet(Q_mat,t_mat);
0026 %
0027 %   References
0028 %   Poff, N.L., Allan, J.D., Bain, M.B., Karr, J.R., Prestegaard, K.L.,
0029 %   Richter, B.D., Sparks, R.E. and Stromberg, J.C., 1997. The natural flow
0030 %   regime. BioScience, 47(11), pp.769-784.
0031 %   Richter, B.D., Baumgartner, J.V., Powell, J. and Braun, D.P., 1996. A
0032 %   method for assessing hydrologic alteration within ecosystems.
0033 %
0034 %   Copyright (C) 2020
0035 %   This software is distributed under the GNU Public License Version 3.
0036 %   See <https://www.gnu.org/licenses/gpl-3.0.en.html> for details.
0037 
0038 % check input parameters
0039 if nargin < 2
0040     error('Not enough input arguments.')
0041 end
0042 
0043 ip = inputParser;
0044 ip.CaseSensitive = true;
0045 
0046 % required input arguments
0047 % Please input time series as a cell array of the following format:
0048 % {x_1; x_2; ...; x_n}, where each entry (1, 2, ..., n) corresponds to one
0049 % time series, e.g. from one catchment. For one catchment only, please
0050 % input {x}. Example: {Q_1; Q_2; ...; Q_n} for streamflow.
0051 addRequired(ip, 'Q_mat', @(Q_mat) iscell(Q_mat))
0052 addRequired(ip, 't_mat', @(t_mat) iscell(t_mat))
0053 
0054 parse(ip, Q_mat, t_mat)
0055 
0056 % calculate signatures
0057 
0058 % initialise arrays
0059 Q_mean = NaN(size(Q_mat,1),1);
0060 Q_mean_error_str = strings(size(Q_mat,1),1);
0061 Q5 = NaN(size(Q_mat,1),1);
0062 Q5_error_str = strings(size(Q_mat,1),1);
0063 Q95 = NaN(size(Q_mat,1),1);
0064 Q95_error_str = strings(size(Q_mat,1),1);
0065 Q_mean_monthly = NaN(size(Q_mat,1),12);
0066 Q_mean_monthly_error_str = strings(size(Q_mat,1),1);
0067 Q_7_day_min = NaN(size(Q_mat,1),1);
0068 Q_7_day_min_error_str = strings(size(Q_mat,1),1);
0069 BFI = NaN(size(Q_mat,1),1);
0070 BFI_error_str = strings(size(Q_mat,1),1);
0071 CoV = NaN(size(Q_mat,1),1);
0072 CoV_error_str = strings(size(Q_mat,1),1);
0073 x_Q_frequency = NaN(size(Q_mat,1),1);
0074 x_Q_frequency_error_str = strings(size(Q_mat,1),1);
0075 x_Q_duration = NaN(size(Q_mat,1),1);
0076 x_Q_duration_error_str = strings(size(Q_mat,1),1);
0077 HFD_mean = NaN(size(Q_mat,1),1);
0078 HFD_mean_error_str = strings(size(Q_mat,1),1);
0079 HFI_mean = NaN(size(Q_mat,1),1);
0080 HFI_mean_error_str = strings(size(Q_mat,1),1);
0081 AC1 = NaN(size(Q_mat,1),1);
0082 AC1_error_str = strings(size(Q_mat,1),1);
0083 FDC_slope = NaN(size(Q_mat,1),1);
0084 FDC_slope_error_str = strings(size(Q_mat,1),1);
0085 BaseflowRecessionK = NaN(size(Q_mat,1),1);
0086 BaseflowRecessionK_error_str = strings(size(Q_mat,1),1);
0087 
0088 % loop over all catchments
0089 for i = 1:size(Q_mat,1)
0090     
0091     [Q_mean(i),~,Q_mean_error_str(i)] = sig_Q_mean(Q_mat{i},t_mat{i});
0092     [Q5(i),~,Q5_error_str(i)] = sig_x_percentile(Q_mat{i},t_mat{i},[5]);
0093     [Q95(i),~,Q95_error_str(i)] = sig_x_percentile(Q_mat{i},t_mat{i},[95]);
0094     [Q_mean_monthly(i,:),~,Q_mean_monthly_error_str(i)] = sig_Q_mean_monthly(Q_mat{i},t_mat{i},[1:12]);
0095     [Q_7_day_min(i),~,Q_7_day_min_error_str(i)] = sig_Q_n_day_min(Q_mat{i},t_mat{i},7);
0096     [BFI(i),~,BFI_error_str(i)] = sig_BFI(Q_mat{i},t_mat{i});
0097     [CoV(i),~,CoV_error_str(i)] = sig_Q_CoV(Q_mat{i},t_mat{i});
0098     [x_Q_frequency(i),~,x_Q_frequency_error_str(i)] = sig_x_Q_frequency(Q_mat{i},t_mat{i},'low');
0099     [x_Q_duration(i),~,x_Q_duration_error_str(i)] = sig_x_Q_duration(Q_mat{i},t_mat{i},'low');
0100     [HFD_mean(i),~,HFD_mean_error_str(i)] = sig_HFD_mean(Q_mat{i},t_mat{i});
0101     [HFI_mean(i),~,HFI_mean_error_str(i)] = sig_HFI_mean(Q_mat{i},t_mat{i});
0102     [AC1(i),~,AC1_error_str(i)] = sig_Autocorrelation(Q_mat{i},t_mat{i},'lag',1);
0103     [FDC_slope(i),~,FDC_slope_error_str(i)] = sig_FDC_slope(Q_mat{i},t_mat{i});
0104     [BaseflowRecessionK(i),~,BaseflowRecessionK_error_str(i)] = sig_BaseflowRecessionK(Q_mat{i},t_mat{i},'eps',0.001*median(Q_mat{i},'omitnan'));
0105     
0106 end
0107 
0108 % add results to struct array
0109 results.Q_mean = Q_mean;
0110 results.Q_mean_error_str = Q_mean_error_str;
0111 results.Q5 = Q5;
0112 results.Q5_error_str = Q5_error_str;
0113 results.Q95 = Q95;
0114 results.Q95_error_str = Q95_error_str;
0115 results.Q_mean_monthly = Q_mean_monthly;
0116 results.Q_mean_monthly_error_str = Q_mean_monthly_error_str;
0117 results.Q_7_day_min = Q_7_day_min;
0118 results.Q_7_day_min_error_str = Q_7_day_min_error_str;
0119 results.BFI = BFI;
0120 results.BFI_error_str = BFI_error_str;
0121 results.CoV = CoV;
0122 results.CoV_error_str = CoV_error_str;
0123 results.x_Q_frequency = x_Q_frequency;
0124 results.x_Q_frequency_error_str = x_Q_frequency_error_str;
0125 results.x_Q_duration = x_Q_duration;
0126 results.x_Q_duration_error_str = x_Q_duration_error_str;
0127 results.HFD_mean = HFD_mean;
0128 results.HFD_mean_error_str = HFD_mean_error_str;
0129 results.HFI_mean = HFI_mean;
0130 results.HFI_mean_error_str = HFI_mean_error_str;
0131 results.AC1 = AC1;
0132 results.AC1_error_str = AC1_error_str;
0133 results.FDC_slope = FDC_slope;
0134 results.FDC_slope_error_str = FDC_slope_error_str;
0135 results.BaseflowRecessionK = BaseflowRecessionK;
0136 results.BaseflowRecessionK_error_str = BaseflowRecessionK_error_str;
0137 
0138 end

Generated on Wed 09-Dec-2020 16:59:11 by m2html © 2005