Reference

fluxpart.api

fluxpart.api.fpread(saved_results)[source]
fluxpart.api.fvs_partition(file_or_dir, time_sorted=False, interval=None, hfd_format=None, hfd_options=None, meas_wue=None, wue_options=None, part_options=None, label=None, stdout=False, verbose=True)[source]

Partition CO2 & H2O fluxes into stomatal & nonstomatal components.

Provides a full implementation of the flux variance similarity partitioning algorithm; reads high frequency eddy covariance data; performs data transformations and data QA/QC; analyzes water vapor and carbon dioxide fluxes; and partitions the fluxes into stomatal (transpiration, photosynthesis) and nonstomatal (evaporation, respiration) components.

The following notation is used in variable naming and documentation to represent meteorological quantities:

u, v, w = wind velocities
q = water vapor mass concentration
c = carbon dioxide mass concentration
T = air temperature
P = total air pressure
Parameters:
  • file_or_dir (str or list of str) – Pathname for high frequency eddy covariance data file, data directory, or list of either. Paths can be relative or absolute, and can include unix-style wildcards and brace expansion. If pathname is a directory without any file pattern, then all files in the directory are read.

  • time_sorted (bool, optional) – Indicates whether the list of data files given by file_or_dir is sorted by time. Default is False.

  • interval (str, optional) – Time interval used to aggregate data and partition fluxes. Specified using the pandas string alias format (e.g., interval="30min"). Set to None to treat whole, individual data files as the analysis interval (e.g., if data files each hold 30 min of data). Default is None.

  • hfd_format (dict) – Dictionary of parameters specifying the high frequency data file format. See Other parameters below for an explanation of required and optional formatting parameters. Some pre-defined file formats can be specified by passing a string instead of a dictionary. Currently two formats are defined: “ec-TOB1” and “ec-TOA5”. These formats correspond to files commonly created when using Campbell Scientific eddy covariance data logger software. See Notes below for an explanation of these pre-defined formats. Default is “ec-TOB1”.

  • hfd_options (dict) – Dictionary of parameters specifying options for correcting high-frequency eddy covariance data and applying quality control measures. See Other parameters for a listing of options.

  • meas_wue (float or callable, optional) – Measured (or otherwise prescribed) leaf-level water use efficiency (kg CO2 / kg H2O). Note that by definition, meas_wue must be a negative value (< 0). If callable, it should take a datetime object as its sole argument and return a float value.

  • wue_options (dict) – Dictionary of parameters and options used to estimate water use efficiency if meas_wue is not provided. See Other parameters section for a description of valid fields. When passing wue_options, it is always required to provide entries “canopy_ht”, “meas_ht”, and “ppath”. Other entries are optional. The values for “canopy_ht” and “meas_ht” can be either a float or a callable. A callable must accept a datetime object as its sole argument and return a float value.

  • part_options (dict) – Dictionary of options for the fvs partitioning procedure. See Other parameters section for a listing of valid options.

  • label (optional) – Optional identifier to be appended to the results.

  • stdout (bool, optional) – If True (default), print to stdout information about the progress of the analysis.

  • verbose (bool, optional) – If True (default), print extra progress information to stdout.

  • hfd_format – High frequency data file format options. NOTE: When passing hfd_format, it is required at a minimum to provide values for “filetype” “cols”, and “time_col” (detailed below). “unit_convert” and “temper_unit” are also required if data are not in SI units.

  • hfd_format["filetype"] ({"csv", "tob1", "ghg", "pd.df"}) – “csv” = delimited text file; “tob1” = Campbell Scientific binary data table format; “ghg” = LI-COR raw data format; “pd.df” = pandas dataframe.

  • hfd_format["cols"] (7*(int,)) – 7-tuple of integers indicating the data column numbers that contain series data for (u, v, w, c, q, T, P), in that order. Uses 0-based indexing.

  • hfd_format["time_col"] (int or [int, int]) – The column index for the datetime column, or if a list, the indices for the date and time columns.

  • hfd_format["unit_convert"] (dict) – Dictionary of multiplication factors required to convert any u, v, w, c, q, or P data not in SI units to SI units (m/s, kg/m^3, Pa). (Note T is not in that list). The dictionary keys are the variable names. For example, if all data are in SI units except P and c, which are in kPa and mg/m^3, respectively, then set: hfd_options["unit_convert"] = {"P": 1e3, "c": 1e-6}, since it is necessary to multiply the kPa pressure data by 1e3 to obtain the SI pressure unit (Pa), and the mg/m^3 CO2 data by 1e-6 to obtain the SI concentration unit (kg/m^3). c and q data can be converted this way if they have units of mass or moles per volume. If instead the c and q data are given as mole fractions (e.g., mmol of q per mol dry air), then set the c and q converters to one of the following: “ppm_dry”, “ppt_dry”, “ppm_wet”, or “ppt_wet”. Use the “_dry” versions if the mole fraction data are w.r.t. to dry air (no water vapor) and “_wet” if they are w.r.t. moist air. For example, if q is in mmol q per mol dry air and c is umol c per mol dry air, then: {"P": 1e3, "q": "ppt_dry", "c": "ppm_dry"}.

  • hfd_format["temper_unit"] ({"K", "C"}) – Temperature data units. Default is “K” (Kelvin).

  • hfd_format["flags"] (2-tuple or list of 2-tuples) – Specifies that one or more data columns are used to flag bad data records. Each tuple is of the form (col, goodval), where col is an int specifying the column number containing the flag (0-based indexing), and goodval is the value of the flag that indicates a good data record.

  • hfd_format["to_datetime_kws"] (dict) – Dict of keyword arguments passed to pandas.to_datetime to read datafile dates and time. Generally needed only if a nonstandard format is used in the datafile.

  • ] (hfd_format[ other keys) – when hfd_format[“filetype”] is “csv” or “ghg”, all other key:value pairs in hfd_format are passed as keyword arguments to pandas.read_csv. Those keywords may be required to specify the details of the file formatting. Among the most commonly required are: “sep”, the str that is used to separate values or define column widths (default is sep=”,”); and “skiprows”, which will be needed if the file contains header rows. See pandas.read_csv for a full description of available format options.

  • hfd_options – Options for pre-processing high frequency data.

  • hfd_options["correct_external"] (bool, optional) – If True (default), the water vapor and carbon dioxide series data are corrected for external fluctuations associated with air temperature and vapor density according to [WPL80] and [DK07].

  • hfd_options["bounds"] (dict, optional) – Dictionary indicating any lower and upper bounds for valid data. Dictionary entries have the form varname: (float, float), where varname is one of “u”, “v”, “w”, “q”, “c”, “T”, or “P”, and the 2-tuple holds values for the lower and upper bounds: (lower, upper). Data records are rejected if a variable in the record is outside the prescribed bounds. Default is bounds = {"c": (0, np.inf), "q": (0, np.inf)} such that data records are rejected if c or q data are not positive values.

  • hfd_options["rd_tol"] (float, optional) – Relative tolerance for rejecting the datafile. Default is hfd_options[“rd_tol”] = 0.4. See hfd_options[“ad_tol”] for further explanation.

  • hfd_options["ad_tol"] (int, optional) – Absolute tolerance for rejecting the datafile. Default is ad_tol = 1024. If the datafile contains bad records (not readable, out-of-bounds, or flagged data), the partitioning analysis is performed using the longest stretch of consecutive good data records found, unless that stretch is too short, in which case the analysis is aborted. The criteria for judging “too short” can be specified in both relative and absolute terms: the datafile is rejected if the good stretch is a fraction of the total data that is less than rd_tol, and/or is less than ad_tol records long.

  • hfd_options["ustar_tol"] (float) – If the friction velocity (m/s) determined from the high frequency data is less than ustar_tol, the partitioning analysis is aborted due to insufficient turbulence. Defalult is hfd_options[“ustar_tol”] = 0.1 (m/s).

  • wue_options – Parameters for estimating water use efficiency. “canopy_ht”, “meas_ht”, and “ppath” are required keys, unless “heights” is given (in which case “canopy_ht” and “meas_ht” are ignored). See: water_use_efficiency() for more about water use efficiency estimation.

  • wue_options["canopy_ht"] (float or callable) – Vegetation canopy height (m). If callable must accept a date object as its sole argument and return a float value.

  • wue_options["meas_ht"] (float or callable) – Eddy covariance measurement height (m). If callable must accept a date object as its sole argument and return a float value.

  • wue_options["heights"] (str or callable, optional) – Alternative way to specify (canopy_ht, meas_ht). If str, is the name of a csv file with date strings in the first column, canopy heights in the second, and measurement heights in the third. If callable, accepts date as its sole argument and returns the tuple of heights.

  • wue_options["ppath"] ({"C3", "C4"}) – Photosynthetic pathway.

  • wue_options["ci_mod"] (str) – Valid values: “const_ratio”, “const_ppm”, “linear”, “sqrt”, “opt”. See: water_use_efficiency().

  • wue_options["ci_mod_param"] (float or (float, float)) – Paramter values to be used with ci_mod. See: water_use_efficiency().

  • wue_options["leaf_temper"] (float, str, or callable) – Canopy leaf temperature. The units (K or C) are whatever is indicated in hfd_format[“temper_unit”]. If not specified, the temperature is set to the air temperature. If a str, is the name of a csv file with datetime stamps in the first column and temperatures in the second. If callable, must accept a datetime as its sole argument and return a float value.

  • wue_options["leaf_temper_corr"] (float) – Offset adjustment applied to canopy temperature (C). Default is zero.

  • wue_options["diff_ratio"] (float, optional) – Ratio of molecular diffusivities for water vapor and CO2. Default is diff_ratio = 1.6.

  • part_options – Options for the fvs partitioning algorithm

  • part_options["adjust_fluxes"] (bool) – If True (default), the final partitioned fluxes are adjusted proportionally such that sum of the partitioned fluxes match exactly the total fluxes indicated in the original data.

  • part_options["daytime"] (2-tuple, str, or callable, optional) – A 2-tuple of python datetime.time objects or timestamps for (sunrise, sunset). If specified, all fluxes for time intervals ending before sunrise or starting after sunset will be assumed non-stomatal. If callable, function should take a datetime.date object as its sole argument and return a 2-tuple of time objects or stamps. If str, is the name of a csv file with date stamp in the first column, sunrise time in the second, and sunset time in the third.

Return type:

FluxpartResult

Notes

Three pre-defined hfd_formats are available.

“EC-TOA5”:

hfd_format = {
    "filetype": "csv",
    "skiprows": 4,
    "cols": (2, 3, 4, 5, 6, 7, 8),
    "time_col": 0,
    "temper_unit": "C",
    "unit_convert": {"q": 1e-3, "c": 1e-6, "P": 1e3},
    "na_values": "NAN",
}

“EC-TOB1”:

hfd_format = {
    "filetype": "tob1",
    "cols": (3, 4, 5, 6, 7, 8, 9),
    "temper_unit": "C",
    "unit_convert": {"q": 1e-3, "c": 1e-6, "P": 1e3},
}

“EC-GHG1”:

hfd_format = {
    "filetype": "ghg",
    "sep": "    ",
    "cols": (11, 12, 13, 7, 8, 9, 10),
    "time_col": [5, 6],
    "unit_convert": dict(q=1e-3 * MW.vapor, c=1e-3 * MW.co2, P=1e3),
    "temper_unit": "C",
    "skiprows": 8,
    "na_values": "NAN",
    "to_datetime_kws": {"format": "%Y-%m-%d %H:%M:%S:%f"},
}

fluxpart.constants

Physical constants

class fluxpart.constants.MolecularWeight(dryair, vapor, co2)[source]
dryair, vapor, co2

kg/mol

Type:

float

class fluxpart.constants.SpecificGasConstant(dryair, vapor, co2)[source]
dryair, vapor, co2

J/K/kg

Type:

float

class fluxpart.constants.SpecificHeatCapacity(dryair)[source]
dryair

J/K/kg

Type:

float

fluxpart.containers

Classes used to group data and results.

class fluxpart.containers.AllFluxes(Fq=nan, Fqt=nan, Fqe=nan, Fc=nan, Fcp=nan, Fcr=nan, temper_kelvin=nan, LE=False, LEt=False, LEe=False, Fq_mol=False, Fqt_mol=False, Fqe_mol=False, Fc_mol=False, Fcp_mol=False, Fcr_mol=False)[source]

Water vapor and CO2 fluxes.

Parameters:
  • Fq (float) – Total, transpiration, and evaporation H2O fluxes, kg/m^2/s.

  • Fqt (float) – Total, transpiration, and evaporation H2O fluxes, kg/m^2/s.

  • Fqe (float) – Total, transpiration, and evaporation H2O fluxes, kg/m^2/s.

  • Fc (float) – Total, photosynthesis, and respiration CO2 fluxes, kg/m^2/s.

  • Fcp (float) – Total, photosynthesis, and respiration CO2 fluxes, kg/m^2/s.

  • Fcr (float) – Total, photosynthesis, and respiration CO2 fluxes, kg/m^2/s.

  • temper_kelvin (float) – Temperature, K

common_units()[source]
common_units_labels()[source]
results_str()[source]
class fluxpart.containers.FVSPSolution(wqc_data=WQCData(var_q=nan, var_c=nan, corr_qc=nan, wq=nan, wc=nan), rootsoln=RootSoln(corr_cp_cr=nan, var_cp=nan, sig_cr=nan, co2soln_id=nan, valid_root=nan), wave_lvl=nan, valid_partition=nan, fvsp_mssg='')[source]

Result of FVS partitioning.

Parameters:

wave_lvl ((int, int)) – 2-tuple indicating the level of filtering applied (number of components removed from the series data). The second int is the maximum possible wavelet decompostion level given the length of the data. The first is the number of components remaining in the data. So when the first number is equal to the second, no components have been removed (no filtering applied). When the first number is 1, the maximum level of filtering was applied.

common_units()[source]
results_str()[source]
class fluxpart.containers.MassFluxes(Fq=nan, Fqt=nan, Fqe=nan, Fc=nan, Fcp=nan, Fcr=nan)[source]

H2O and CO2 mass flux components.

Parameters:
  • Fq (float) – Total, transpiration, and evaporation H2O fluxes, kg/m^2/s

  • Fqt (float) – Total, transpiration, and evaporation H2O fluxes, kg/m^2/s

  • Fqe (float) – Total, transpiration, and evaporation H2O fluxes, kg/m^2/s

  • Fc (float) – Total, photosynthesis, and respiration CO2 fluxes, kg/m^2/s

  • Fcp (float) – Total, photosynthesis, and respiration CO2 fluxes, kg/m^2/s

  • Fcr (float) – Total, photosynthesis, and respiration CO2 fluxes, kg/m^2/s

class fluxpart.containers.RootSoln(corr_cp_cr=nan, var_cp=nan, sig_cr=nan, co2soln_id=nan, valid_root=nan)[source]

Results from calculating the root (corr_cp_cr, var_cp).

corr_cp_cr

Correlation coefficient for CO2 concentrations connected with photosynthesis (cp) and respiration (cr).

Type:

float

var_cp

Variance of CO2 concentration connected with photosynthesis, (kg/m^3)^2.

Type:

float

sig_cr

Standard deviation of CO2 concentration connected with respiration, kg/m^3.

Type:

float

co2soln_id

Indicates the solution used for the quadratic Eq. 13b of [SS08]. Equal to 1 for the ‘+’ root, equal to 0 for the ‘-’ root.

Type:

{0 or 1}

valid_root

Indicates whether the obtained root (corr_cp_cr, var_cp) is physically plausible.

Type:

bool

root_mssg

Possibly informative message if isvalid = False.

Type:

str

common_units()[source]
common_units_labels()[source]
results_str(head=True)[source]
class fluxpart.containers.WQCData(var_q=nan, var_c=nan, corr_qc=nan, wq=nan, wc=nan)[source]

Summary stats for wind, water vapor, and CO2 data.

var_q, var_c

Variance of vapor (q) and CO2 (c) concentrations, (kg/m^3)^2.

Type:

float

corr_qc

Correlation coefficient for vapor and CO2 concentrations.

Type:

float

wq, wc

Mean vapor (wq) and CO2 (wc) fluxes, kg/m^2/s.

Type:

float

common_units()[source]
common_units_labels()[source]
results_str(head=True)[source]
class fluxpart.containers.WUE(wue=nan, inter_h2o=nan, inter_co2=nan, ambient_h2o=nan, ambient_co2=nan, vpd=nan, ci_mod=nan, ci_mod_param=nan, leaf_temper=nan, ppath=nan, meas_ht=nan, canopy_ht=nan, diff_ratio=nan)[source]

Summary of leaf-level water use efficiency calculation.

wue

Leaf-level water use efficiency, kg CO2 / kg H2O.

Type:

float

inter_h2o

Concentration of intercellular water vapor, kg/m^3

Type:

float

inter_co2

Concentrations of intercellular CO2, kg/m^3.

Type:

float

ambient_h2o

Concentrations of ambient atmospheric water vapor, kg/m^3.

Type:

float

ambient_co2

Concentration of ambient CO2, kg/m^3.

Type:

float

vpd

Atmospheric vapor pressure deficit, Pa.

Type:

float

ci_mod

The name of the model used to estimate the intercellular CO2 concentration.

Type:

str

ci_mod_param

Specific paramter values used with ci_mod.

Type:

float or 2-tuple of floats

leaf_temper
Type:

float

ppath

Photosynthetic pathway.

Type:

{‘C3’ or ‘C4’}

meas_ht, canopy_ht

Eddy covariance measument height and plant canopy height, m.

Type:

float

common_units()[source]
common_units_labels()[source]
results_str()[source]

fluxpart.fluxpart

exception fluxpart.fluxpart.Error[source]
class fluxpart.fluxpart.FVSResult(dataread=False, attempt_partition=False, partition_success=False, mssg=None, label=None, sunrise=None, sunset=None, fluxes=AllFluxes(Fq=nan, Fqt=nan, Fqe=nan, Fc=nan, Fcp=nan, Fcr=nan, temper_kelvin=nan, LE=nan, LEt=nan, LEe=nan, Fq_mol=nan, Fqt_mol=nan, Fqe_mol=nan, Fc_mol=nan, Fcp_mol=nan, Fcr_mol=nan), hfsummary=HFSummary(T=nan, P=nan, Pvap=nan, ustar=nan, wind_w=nan, var_w=nan, rho_vapor=nan, rho_co2=nan, var_vapor=nan, var_co2=nan, corr_q_c=nan, cov_w_q=nan, cov_w_c=nan, H=nan, rho_dryair=nan, rho_totair=nan, cov_w_T=nan, N=nan), wue=WUE(wue=nan, inter_h2o=nan, inter_co2=nan, ambient_h2o=nan, ambient_co2=nan, vpd=nan, ci_mod=nan, ci_mod_param=nan, leaf_temper=nan, ppath=nan, meas_ht=nan, canopy_ht=nan, diff_ratio=nan), fvsp_solution=FVSPSolution(wqc_data=WQCData(var_q=nan, var_c=nan, corr_qc=nan, wq=nan, wc=nan), rootsoln=RootSoln(corr_cp_cr=nan, var_cp=nan, sig_cr=nan, co2soln_id=nan, valid_root=nan), wave_lvl=nan, valid_partition=nan, fvsp_mssg=''))[source]

FVS partitioning result.

exception fluxpart.fluxpart.FluxpartError(message)[source]
class fluxpart.fluxpart.FluxpartResult(fp_results)[source]
istr(i)[source]

Return a string representation of the ith result

plot_co2(start=None, end=None, units='mass', components=(0, 1, 2), ax=None, **kws)[source]
plot_h2o(start=None, end=None, units='mass', components=(0, 1, 2), ax=None, **kws)[source]
save(filename)[source]
save_csv(filename)[source]
save_pickle(filename)[source]
fluxpart.fluxpart.flux_partition(*args, **kws)[source]
fluxpart.fluxpart.fvspart(file_or_dir, time_sorted=False, interval=None, hfd_format=None, hfd_options=None, meas_wue=None, wue_options=None, part_options=None, label=None, stdout=True, verbose=True)[source]

Partition CO2 & H2O fluxes into stomatal & nonstomatal components.

Provides a full implementation of the flux variance similarity partitioning algorithm [SS08]_[SAAS+18]_: reads high frequency eddy covariance data; performs data transformations and data QA/QC; analyzes water vapor and carbon dioxide fluxes; and partitions the fluxes into stomatal (transpiration, photosynthesis) and nonstomatal (evaporation, respiration) components.

The following notation is used in variable naming and documentation to represent meteorological quantities:

u, v, w = wind velocities
q = water vapor mass concentration
c = carbon dioxide mass concentration
T = air temperature
P = total air pressure
Parameters:

see (For parameters explanation) –

Return type:

FluxpartResult

fluxpart.hfdata

Read and model eddy covariance high-frequency time series data.

The following notation is used in variable naming to represent meteorological quantities (SI units):

u, v, w = wind velocities (m/s)
q = water vapor mass concentration (kg/m^3)
c = carbon dioxide mass concentration (kg/m^3)
T = air temperature (K)
P = total air pressure (Pa)
exception fluxpart.hfdata.Error[source]
class fluxpart.hfdata.HFData(hf_dataframe)[source]
Parameters:

hf_dataframe – High frequency eddy covariance dataframe. Must include columns for [“u”, “v”, “w”, “c”, “q”, “T”, “P”]. Normally, dataframe should have a datetime index. Dataframe may include also boolean mask columns.

cleanse(bounds=None, rd_tol=0.5, ad_tol=1024)[source]

Apply some data QC, remove bad data.

If problems are found (data not readable, out-of-bounds, or flagged), self.dataframe is modified to contain only the longest contiguous stretch of good data. An error is raised if the resulting data are too few. The criteria for judging ‘too few’ can be specified in both relative and absolute terms: the datafile is rejected if the good stretch is a fraction of the total data that is less than rd_tol, and/or is less than ad_tol records long.

Parameters:
  • bounds (dict, optional) – Dictionary specifying lower and upper bounds for legal data. Dict entries have the form varname: (lower, upper), where varname is one of ‘u’, ‘v’, ‘w’, ‘q’, ‘c’, ‘T’, or ‘P’, and the 2-tuple holds values for the lower and upper bounds. Default is None.

  • rd_tol (float, optional) – Relative tolerance for rejecting the datafile. Default is rd_tol = 0.5.

  • ad_tol (int, optional) – Absolute tolerance for rejecting the datafile. Default is ad_tol = 1024.

correct_external()[source]

Adjust q and c data series to correct for external effects.

Water vapor and carbon dioxide series data in the dataframe are corrected for external fluctuations associated with air temperature and vapor density. See: [WPL80] and [DK07].

summarize()[source]

Summarize high frequency dataframe statistics.

Return type:

HFSummary

truncate_pow2()[source]

Truncate dataframe length to largest possible power of 2.

exception fluxpart.hfdata.HFDataReadError(message)[source]
class fluxpart.hfdata.HFDataSource(files, filetype, cols, converters=None, time_col=None, flags=None, **kwargs)[source]

Reader for high-frequency eddy covariance data.

Parameters:
  • files (list or sequence of files) – Sorted sequence of data files

  • filetype ({'csv', 'tob1', 'ghg'}) – ‘csv’ = delimited text file; ‘tob1’ = Campbell Scientific binary format file; ‘ghg’ = LI-COR raw data format

  • cols (7*(int,)) – Column indices for (u, v, w, q, c, T, P) data, in that order. 0-based indexing.

  • converters (dict, required if datafile uses any non-SI units) – Dictionary of functions used to convert any non-SI data to SI units. Dict keys are ‘u’, ‘v’, ‘w’, ‘q’, ‘c’, ‘T’, or ‘P’. Functions take a single argument, e.g. converters={'P': lambda arg: 1e3 * arg}. If all data are SI units, set converters to None (default).

  • time_col (int, optional) – Datetime column for csv and ghg filetype. Default is None.

  • flags (2-tuple or list of 2-tuples, optional) – Indicate that one or more data columns are used to flag bad data records. Each tuple is of the form (col, goodval), where col is an int specifying the column number (0-based indexing), and goodval is the flag value indicating a good data record.

  • **kwargs – Passed to pandas.read_csv when filetype is csv or ghg. Should not include usecols or header keywords.

reader(interval, **kwargs)[source]

Consume data source and yield in chunks of the time interval.

Data chunks are returned in dataframe format.

Parameters:
  • interval (str) – Time interval used to chunk the data. Is specified using the pandas string alias format (e.g., interval="30min"). If set to -1, a single dataframe corresponding to a concatenation of all data files is returned on the first iteration. If set to None, dataframes corresponding to whole individual data files are returned with each iteration.

  • **kwargs – For csv and ghg filetypes, kwargs are passed to pandas.read_csv. Can be used to override or add to formatting specified in the initializer. Should not include usecols or header keywords. For tob1 filetype, kwargs is passed to numpy.fromfile. In this case, the only allowable kwarg is ‘count’, which can be used to limit the number of lines read.

class fluxpart.hfdata.HFSummary(T=nan, P=nan, Pvap=nan, ustar=nan, wind_w=nan, var_w=nan, rho_vapor=nan, rho_co2=nan, var_vapor=nan, var_co2=nan, corr_q_c=nan, cov_w_q=nan, cov_w_c=nan, H=nan, rho_dryair=nan, rho_totair=nan, cov_w_T=nan, N=nan)[source]

Summary of high frequency eddy covariance data.

Parameters:
  • T (float) – Mean air temperature, K.

  • P (float) – Mean total atmospheric (P) and vapor (Pvap) pressure, Pa.

  • Pvap (float) – Mean total atmospheric (P) and vapor (Pvap) pressure, Pa.

  • ustar (float) – Mean friction velocity, m/s.

  • wind_w (float) – Mean vertical wind velocity, m/s.

  • var_w (float) – Variance of vertical wind velocity, (m/s)^2.

  • rho_vapor (float) – Mean H2O vapor and CO2 concentrations, kg/m^3.

  • rho_co2 (float) – Mean H2O vapor and CO2 concentrations, kg/m^3.

  • var_vapor (float) – Variance of H2O vapor and CO2 concentrations, (kg/m^3)^2.

  • var_co2 (float) – Variance of H2O vapor and CO2 concentrations, (kg/m^3)^2.

  • corr_q_c (float) – Correlation coefficient for H2O and CO2 concentrations

  • cov_w_q (float) – Covariance of vertical wind velocity (w) with water vapor (q) and CO2 (c) mass densities, kg/m^2/s.

  • cov_w_c (float) – Covariance of vertical wind velocity (w) with water vapor (q) and CO2 (c) mass densities, kg/m^2/s.

  • H (float) – Sensible heat flux, W/m^2.

  • rho_dryair (float) – Dry and moist air densities, kg/m^3.

  • rho_totair (float) – Dry and moist air densities, kg/m^3.

  • cov_w_T (float) – Covariance of temperature and vertical wind velocity, K m/s.

  • N (int) – Length of data series.

common_units()[source]
common_units_labels()[source]
property fc_ov_fq
results_str()[source]
property sigc_ov_sigq
exception fluxpart.hfdata.TooFewDataError[source]

fluxpart.partition

TODO:

exception fluxpart.partition.Error[source]
exception fluxpart.partition.FVSError(message)[source]
fluxpart.partition.findroot(wqc_data, wue)[source]

Calculate (corr_cp_cr, var_cp).

Parameters:
  • wqc_data (namedtuple or equivalent namespace) – WQCData

  • wue (float) – Leaf-level water use efficiency, wue < 0, kg CO2 / kg H2O.

Returns:

RootSoln

Return type:

namedtuple

fluxpart.partition.flux_ratio(var_cp, corr_cp_cr, wqc_data, ftype, farg)[source]

Compute the nonstomatal:stomatal ratio of the H2O or CO2 flux.

The ratio (either wqe/wqt or wcr/wcp) is found by solving Eq. 13 of [SS08].

Parameters:
  • wqc_data (namedtuple or equivalent namespace) – WQCData

  • ftype ({'co2', 'h2o'}) – Specifies whether the flux is CO2 or H2O

  • farg (number) – If ftype = ‘co2’, then farg = {0 or 1} specifies the root of Eq. 13b to be used to calculate the CO2 flux ratio wcr/wcp: farg`=1 uses the ‘+’ solution, `farg`=0 uses the ‘-’ solution. If `ftype = ‘h2o’, then farg is a float equal to the water use efficiency (wue < 0), kg/kg.

Returns:

fratio – The requested flux ratio, wqe/wqt or wcr/wcp. Set to np.nan if solution is not real.

Return type:

float or np.nan

Notes

When solving for wqe/wqt, the ‘-’ solution of the quadratic Eq. 13a is not relevant because only the ‘+’ solution yields wqe/wqt > 0, which is required/assumed by the physical model in [SS08].

fluxpart.partition.fvspart_interval(wqc_data, wue, wipe_if_invalid=False)[source]

Partition H2O and CO2 fluxes using interval averaged data values.

Parameters:
  • wqc_data (WQCData) –

  • wue (float, kg CO2 / kg H2O) – Leaf-level water use efficiency, wue < 0

  • wipe_if_invalid (boolean) – If True, return default (nan) values for all mass fluxes if any calculated fluxes violate directional requirements. Default is False.

Return type:

FVSPSolution

fluxpart.partition.fvspart_progressive(w, q, c, wue, adjust_fluxes=True)[source]

FVS flux partitioning using high frequency eddy covariance data.

If a valid partitioning solution is not found for the passed series data, low-frequency (large-scale) components are progressively removed from the data until either a valid solution is found or the series decomposition is exhausted.

Parameters:
  • w (array_like) – 1D high frequency time series for vertical wind speed w (m/s), water vapor q (kg/m^3), and CO2 c (kg/m^3).

  • q (array_like) – 1D high frequency time series for vertical wind speed w (m/s), water vapor q (kg/m^3), and CO2 c (kg/m^3).

  • c (array_like) – 1D high frequency time series for vertical wind speed w (m/s), water vapor q (kg/m^3), and CO2 c (kg/m^3).

  • wue (float, wue < 0) – leaf-level water use efficiency, kg CO2 / kg H2O.

  • adjust_fluxes (bool, optional) – Indicates whether the obtained partitioned fluxes should be adjusted so that the totals match the original data. Default is True.

Return type:

FVSPSolution

Notes

If a valid partitioning is not found, the returned numersoln and wqc_data correspond to the final iteration attempted.

fluxpart.partition.fvspart_series(w, q, c, wue, wipe_if_invalid=False)[source]

FVS partition q and c fluxes using high frequency eddy cov data.

Parameters:
  • w (array_like) – 1D high frequency time series for vertical wind speed w (m/s), water vapor q (kg/m^3), and CO2 c (kg/m^3).

  • q (array_like) – 1D high frequency time series for vertical wind speed w (m/s), water vapor q (kg/m^3), and CO2 c (kg/m^3).

  • c (array_like) – 1D high frequency time series for vertical wind speed w (m/s), water vapor q (kg/m^3), and CO2 c (kg/m^3).

  • wue (float, wue < 0) – leaf-level water use efficiency, kg CO2 / kg H2O.

  • wipe_if_invalid (boolean) – If True, return default (nan) values for all mass fluxes if any calculated fluxes violate directional requirements. Default is False.

Return type:

FVSPSolution,

fluxpart.util

exception fluxpart.util.HFDataReadWarning[source]
fluxpart.util.cflux_mass_to_mol(massflux)[source]
fluxpart.util.chunked_df(dataframes, time_interval)[source]

Partition time-indexed dataframe sequence into time intervals.

fluxpart.util.dataframe_read_tob1(tobfile, count=-1)[source]
fluxpart.util.multifile_read_csv(files, *args, **kwargs)[source]

Buffered pd.read_csv of data split across multiple files.

fluxpart.util.multifile_read_ghg(files, *args, **kwargs)[source]

Buffered pd.read_csv of data split across multiple files.

fluxpart.util.multifile_read_tob1(tobfiles, count=-1)[source]

Buffered read of multiple tob1 files into dataframes.

fluxpart.util.ndarray_read_tob1(tobfile, count=-1)[source]

Read TOB1 data file into structured numpy array.

fluxpart.util.progressive_lowcut_series(series)[source]

Progressively remove low-frequency components of 1D series.

Yields sequence in which the low-frequency (large-scale) components of series are progressively removed. The sequence is obtained from reconstruction of a multilevel discrete haar wavelet decomposition of series.

N.B.: The length of series is assumed to be a power of 2 (does not check!)

Parameters:

series (array_like) – 1D data series with a length that is a power of 2

Yields:

lowcut_series (array) – Sequence of progressively lowcut filtered data series. The yielded series have the same length as series.

Notes

After an N level discrete wavelet decomposition, a data series S can be reconstructed in terms of wavelet ‘approximations’ (A) and ‘details’ (D):

S = A(N) + D(N) + D(N-1) ... D(2) + D(1)
  = A(N-1) + D(N-1) + ... D(2) + D(1)     [A(N-1) = A(N) + D(N)]
    ...
  = A(1) + D(1)                           [(A(1) = A(2) + D(2)]

where A(N) represents the ‘lowest’ level approximation. For the haar wavelet and data S having a length that is a power of 2, A(N) is equal to mean(S)]

The sequence returned by this function is:

S - A(N),
S - A(N-1),
S - A(N-2),
...,
S - A(1)

This sequence is computed by the equivalent:

S - A(N),
S - A(N) - D(N),
S - A(N) - D(N) - D(N-1),
...,
S - A(N) - D(N) - D(N-1) - ... - D(2),

i.e. the details are removed in sequence:

lowcut = S - A(N)
for j = N to 2
    lowcut -= D(j)
fluxpart.util.qflux_mass_to_heat(massflux, Tk)[source]
fluxpart.util.qflux_mass_to_mol(massflux)[source]
fluxpart.util.sat_vapor_press(t_kelvin)[source]
fluxpart.util.stats2(sarray, names=None)[source]

Calculate means and (co)variances for structured array data.

fluxpart.util.vapor_press_deficit(rho_vapor, t_kelvin)[source]
fluxpart.util.vapor_press_deficit_mass(rho_vapor, t_kelvin)[source]

fluxpart.wue

exception fluxpart.wue.Error[source]
exception fluxpart.wue.WUEError(message)[source]
fluxpart.wue.water_use_efficiency(hfs, ci_mod, ci_mod_param=None, leaf_temper=None, leaf_temper_corr=0, meas_ht=None, canopy_ht=None, ppath=None, diff_ratio=1.6)[source]

Estimate leaf-level water use efficiency.

Parameters:
  • hfs (HFSummary namedtuple) – Summary statistics for high frequency eddy covariance data interval. Possesses the following attributes (all floats): rho_vapor, mean vapor density (kg/m^3); rho_co2, mean carbon dioxide concentration (kg/m^3); T, mean air temperature (K); P, mean atmospheric pressure (Pa); cov_w_q, covariance of vert wind velocity and H2O density (kg/m^2/s); cov_w_q, covariance of vert wind velocity and CO2 density (kg/m^2/s); cov_w_T, covariance of wind and temperature (K m/s); ustar, friction velocity (m/s); rho_totair, moist air density (kg/m^3).

  • meas_ht (float or callable, optional) – Eddy covariance measurement height (m).

  • canopy_ht (float or callable, optional) – Vegetation canopy height (m).

  • ppath ({'C3', 'C4'}) – photosynthetic pathway

  • ci_mod ({'const_ratio', 'const_ppm', 'linear', 'sqrt', 'opt'}) – Specifies the model to be used to determine the leaf intercellular CO2 concentration. See Notes below for model descriptions.

  • ci_mod_param (float or 2-tuple of floats, optional) – Paramter values to be used with ci_mod. The number of parameters required depends on the model (see Notes). The default is ci_mod_param=None, in which case default values are used.

  • leaf_temper (float, optional) – Measured canopy tempeature in degrees K. If None (defalut), the leaf temperature is taken to be equal to the air temperature in hfs.

  • leaf_temper_corr (float, optional) – Optional adjustment to leaf temperature. The temperature used to calculate intercelluar vapor and CO2 concentrations is leaf_T + leaf_temp_corr, where leaf_T is leaf_temper if provided, and the air temperature in hfs otherwise. Default is leaf_temper_corr = 0.

  • diff_ratio (float, optional) – Ratio of molecular diffusivities for water vapor and CO2. Default is diff_ratio = 1.6.

Returns:

WUE

Return type:

namedtuple

Notes

Leaf-level water use efficiency is estimated as [CN98] (pg. 238):

wue = (1 / DR) * (ca - ci) / (qa - qi)

where:

DR = `diff_ratio`
ca, ci = ambient and intercellular CO2 concentration, resp.
qa, qi = ambient and intercellular H2O concentration, resp.

ca and qa are estimated from above-canopy tower measurements by extrapolating a logarithmic mean profile with stability corrections to the zero-plane displacement height [SS08].

qi corresponds to 100 percent relative humidity at leaf_temper.

To estimate ci, the following models are available:

‘const_ppm’

ci (kg/m^3) is determined from a specified constant ppm value, ci_ppm:

ci = f(ci_ppm; temperature, pressure)

Default parameter values for ci_ppm are 280 ppm when ppath=’C3’, and 130 ppm when ppath=’C4’. [CN98] (pg. 237).

‘const_ratio’

The ratio ci/ca is constant:

ci/ca = K.

Default parameter values are K=0.7 for C3 plants and K=0.44 for C4.

‘linear’

ci/ca is a linear function of vpd (the atmospheric vapor pressure deficit, which is calculated internally):

ci/ca = b - m * vpd

b is dimensionless with a value of ~1 while m has units of 1/Pa. The parameter pair (b, m) defaluts to (1, 1.6e-4) for C3 plants and (1, 2.7e-4) for C4. See [MG83].

‘sqrt’

ci/ca is a function of sqrt(vpd/ca) [KPO09]:

ci/ca = 1 - sqrt(1.6 * lambd *  vpd / ca)

The paramater lambd has units of kg-CO2 / m^3 / Pa, and defaults to 22e-9 for C3 plants. The sqrt model is not enabled for C4 plants.

‘opt’

ci/ca is determined based on an optimization model. The estimated value is given by Eq. 19 of [SSS19]. With this approach, no paramters (ci_mod_param) are specified. The model may not be suitable for C4 vegetation.