Module Index
RetentionParameterEstimator.Program
RetentionParameterEstimator.check_measurement
RetentionParameterEstimator.elution_temperature
RetentionParameterEstimator.estimate_parameters
RetentionParameterEstimator.estimate_start_parameter
RetentionParameterEstimator.estimate_start_parameter_mean_elu_temp
RetentionParameterEstimator.estimate_start_parameter_single_ramp
RetentionParameterEstimator.load_chromatograms
RetentionParameterEstimator.load_chromatograms
RetentionParameterEstimator.loss
RetentionParameterEstimator.loss
RetentionParameterEstimator.method_m1
RetentionParameterEstimator.method_m2
RetentionParameterEstimator.method_m3
RetentionParameterEstimator.opt_Kcentric
RetentionParameterEstimator.opt_dKcentric
RetentionParameterEstimator.optimize_Kcentric
RetentionParameterEstimator.optimize_dKcentric
RetentionParameterEstimator.prog_list
RetentionParameterEstimator.reference_holdup_time
RetentionParameterEstimator.separate_error_columns
RetentionParameterEstimator.stderror
RetentionParameterEstimator.substance_list
RetentionParameterEstimator.tR_calc
RetentionParameterEstimator.tR_τR_calc
Detailed API
RetentionParameterEstimator.Program
— MethodProgram(TP, FpinP, L; pout="vacuum", time_unit="min")
Construct the structure Program
with conventional formulation (see function conventional_program
in GasChromatographySimulator.jl
) of programs for the case without a thermal gradient.
Arguments
TP
: conventional formulation of a temperature program.FpinP
: conventional formulation of a Flow (in m³/s) resp. inlet pressure (in Pa(a)) program.L
: Length of the capillary measured in m (meter).pout
: Outlet pressure, "vacuum" (default), "atmosphere" or the outlet pressure in Pa(a).time_unit
: unit of time in the programs, "min"` (default) times are measured in minutes, "s" times are measured in seconds.
The argument L
is used to construct the temperature interpolation T_itp(x,t)
.
Examples
julia> Program((40.0, 1.0, 5.0, 280.0, 2.0, 20.0, 320.0, 2.0),
(400000.0, 10.0, 5000.0, 500000.0, 20.0),
10.0)
RetentionParameterEstimator.check_measurement
— Methodcheck_measurement(meas, col_input; min_th=0.1, loss_th=1.0, se_col=true)
Similar to method_m1
(method_m1
) estimate the three retention parameters $T_{char}$, $θ_{char}$ and $ΔC_p$ including standard errors, see stderror
. In addition, if the found optimized minima is above a threshold min_th
, it is flagged and the squared differences of single measured retention times and calculated retention times above another threshold loss_th
are recorded.
Arguments
meas
... Tuple with the loaded measurement data, seeload_chromatograms
.col_input
... Named tuple withcol_input.L
the column length in m andcol_input.d
the column diameter in mm. If this parameter is not gicen, than these parameters are taken frommeas
.se_col=true
... Iftrue
the standard errors (from the Hessian matrix, seestderror
) of the estimated parameters are added as separate columns to the result dataframe. Iffalse
the standard errors are added to the values asMasurement
type.
Output
check
... Boolean.true
if all values are below the thresholds,false
if not.msg
... String. Description ofcheck
df_flag
... Dataframe containing name of the flagged measurements, solutes and the corresponding mesured and calculated retention times.index_flag
... Indices of flagged resultsres
... Dataframe with the optimized parameters and the found minima.Telu_max
... The maximum of elution temperatures every solute experiences in the measured programs.
RetentionParameterEstimator.elution_temperature
— Methodelution_temperature(tRs, prog)
Calculate the elution temperatures from retention times tRs
and GC programs prog
.
RetentionParameterEstimator.estimate_parameters
— Methodestimate_parameters()
Calculate the estimates for the K-centric parameters and (optional) the column diameter.
Arguments
chrom
... Tuple of the loaded chromatogram, seeload_chromatograms
Options
method=NewtonTrustRegion()
... used optimization methodopt=std_opt
... general options,std_opt = GasChromatographySimulator.Options(abstol=1e-8, reltol=1e-5, ng=true, odesys=false)
maxiters=10000
... maximum number of iterations for every single optimizationmaxtime=600.0
... maximum time for every single optimizationmode="dKcentric"
... mode of the estimation. Possible options:- "Kcentric_single" ... optimization for the three K-centric retention parameters separatly for every solute
- "Kcentric" ... optimization for the three K-centric retention parameters together for all solutes
- "dKcentric_single" ... optimization for the column diameter and the three K-centric retention parameters separatly for every solute
- "dKcentric" ... optimization for the column diameter and the three K-centric retention parameters together for all solutes
metric="squared"
... used metric for the loss function ("squared" or "abs")
Output
df
... DataFrame with the columnsName
(solute names),d
(estimated column diameter, optional),Tchar
(estimated Tchar),θchar
(estimated θchar),ΔCp
(estimated ΔCp) andmin
(value of the loss function at the found optima)sol
... Array ofSciMLBase.OptimizationSolution
with the results of the optimization with some additional informations.
RetentionParameterEstimator.estimate_start_parameter
— Methodestimate_start_parameter(tRs::DataFrame, col, prog; time_unit="min", control="Pressure")
Estimation of initial parameters for Tchar
, θchar
and ΔCp
based on the elution temperatures calculated from the retention times tR
and GC programs prog
for column col
. The initial value of Tchar
is estimated from the elution temperatures of the measurements. Based on this estimated Tchar
estimates for the initial values of θchar
and ΔCp
are calculated as $\theta_{char,init} = 22 \left(\frac{T_{char,init}}{T_{st}}\right)^{0.7} \left(1000\frac{d_f}{d}\right)^{0.09} °C$ and $\Delta C_p = (-52 + 0.34 T_{char,init}) \mathrm{J mol^{-1} K^{-1}}$
Output
Tchar_est
... estimate for initial guess of the characteristic temperatureθchar_est
... estimate for initial guess of θcharΔCp_est
... estimate for initial guess of ΔCpTelu_max
... the maximum of the calculated elution temperatures of the solutes
RetentionParameterEstimator.estimate_start_parameter_mean_elu_temp
— Methodestimate_start_parameter_mean_elu_temp(tRs::DataFrame, col, prog; time_unit="min", control="Pressure")
Estimation of initial parameters for Tchar
, θchar
and ΔCp
based on the elution temperatures calculated from the retention times tR
and GC programs prog
for column col
. This function is used, if the temperature program is not a single ramp heating program. The elution temperatures of all measurements are calculated and the mean value of the elution temperatures is used as the initial characteristic temperature of a substance. Based on this estimated Tchar
estimates for the initial values of θchar
and ΔCp
are calculated as $\theta_{char,init} = 22 \left(\frac{T_{char,init}}{T_{st}}\right)^{0.7} \left(1000\frac{d_f}{d}\right)^{0.09} °C$ and $\Delta C_p = (-52 + 0.34 T_{char,init}) \mathrm{J mol^{-1} K^{-1}}$
Output
Tchar_est
... estimate for initial guess of the characteristic temperatureθchar_est
... estimate for initial guess of θcharΔCp_est
... estimate for initial guess of ΔCpTelu_max
... the maximum of the calculated elution temperatures of the solutes
RetentionParameterEstimator.estimate_start_parameter_single_ramp
— Methodestimate_start_parameter_single_ramp(tRs::DataFrame, col, prog; time_unit="min", control="Pressure")
Estimation of initial parameters for Tchar
, θchar
and ΔCp
based on the elution temperatures calculated from the retention times tR
and GC programs prog
for column col
. For this function it is assumed, that single ramp heating programs are used. The elution temperatures of all measurements are calculated and than interpolated over the heating rates. For a dimensionless heating rate of 0.6 the elution temperature and the characteristic temperature of a substance are nearly equal. Based on this estimated Tchar
estimates for the initial values of θchar
and ΔCp
are calculated as $\theta_{char,init} = 22 \left(\frac{T_{char,init}}{T_{st}}\right)^{0.7} \left(1000\frac{d_f}{d}\right)^{0.09} °C$ and $\Delta C_p = (-180 + 0.63 T_{char,init}) \mathrm{J mol^{-1} K^{-1}}$
Output
Tchar_est
... estimate for initial guess of the characteristic temperatureθchar_est
... estimate for initial guess of θcharΔCp_est
... estimate for initial guess of ΔCpTelu_max
... the maximum of the calculated elution temperatures of the solutes
RetentionParameterEstimator.load_chromatograms
— Methodload_chromatograms(file; delim=";")
Loading of the chromatographic data (column information, GC program information, retention time information, see also "Structure of input data") from a file.
Arguments
file
... path to the file.delim=";"
... Delimiter for the imported .csv file.
Output
A tuple of the following quantities:
col
... settings of the column asGasChromatographySimulator.Column
structure.prog
... Array of the GC programs asGasChromatographySimulator.Program
structure.tRs
... DataFrame of the retention times.solute_names
... Vector of the solute names.pout
... outlet pressure (detector pressure), "vacuum" or "atmospheric".time_unit
... unit of time scale used in the retention times and GC programs, "min" or "s".
RetentionParameterEstimator.load_chromatograms
— Methodload_chromatograms(file::Dict{Any, Any}; path=joinpath(dirname(pwd()), "data", "exp_pro"))
Loading of the chromatographic data (column information, GC program information, retention time information, see also "Structure of input data") from a file selected by the FilePicker in a Pluto notebook.
Arguments
file
... file dictionary from the FilePicker.path
... if the temperature programs are defined by measured temperatures over time, define the path to these files.
Output
A tuple of the following quantities:
col
... settings of the column asGasChromatographySimulator.Column
structure.prog
... Array of the GC programs asGasChromatographySimulator.Program
structure.tRs
... DataFrame of the retention times.solute_names
... Vector of the solute names.pout
... outlet pressure (detector pressure), "vacuum" or "atmospheric".time_unit
... unit of time scale used in the retention times and GC programs, "min" or "s".
RetentionParameterEstimator.loss
— Methodloss(tR, Tchar, θchar, ΔCp, L, d, prog, gas; opt=std_opt, metric="squared")
Loss function as sum of squares of the residuals between the measured and calculated retention times.
Arguments
tR
... mxn-array of the measured retention times in seconds.Tchar
... n-array of characteristic temperatures in K.θchar
... n-array of characteristic constants in °C.ΔCp
... n-array of the change of adiabatic heat capacity in J mol^-1 K^-1.L
... number of the length of the column in m.d
... number of the diameters of the column in m.prog
... m-array of structure GasChromatographySimulator.Programs containing the definition of the GC-programs.gas
... string of name of the mobile phase gas.
Output
The output is a tuple of the following quantites:
sum((tR.-tRcalc).^2)
... sum of the squared residuals over m GC-programs and n solutes.
RetentionParameterEstimator.loss
— Methodloss(tR::Array{T, 1}, Tchar, θchar, ΔCp, substance_list::Array{String, 1}, L, d, prog::Array{GasChromatographySimulator.Program, 1}, gas; opt=std_opt, metric="squared") where T<:Number
Calculates the loss between the measured retention times tR
and the calculated retention times tRcalc
for a list (vector) of substances.
Arguments
tR::Array{T, 1}
: Array of measured retention times.Tchar
: Array of characteristic temperatures for each substance.θchar
: Array of characteristic parameters for each substance.ΔCp
: Array of heat capacity changes for each substance.substance_list::Array{String, 1}
: List of substances.L
: Length of the column.d
: Diameter of the column.prog::Array{GasChromatographySimulator.Program, 1}
: Array of programs for gas chromatography.gas
: The gas used.opt=std_opt
: (Optional) Options for calculating retention times.metric="squared"
: (Optional) The loss metric to use. Can be "abs" for absolute loss or "squared" for squared loss.
Return Value
l
: The calculated loss value.
Errors
- Throws an error if the lengths of
tR
,substance_list
, andprog
do not match.
RetentionParameterEstimator.method_m1
— Methodmethod_m1(meas, col_input; se_col=true)
Estimation of the three retention parameters $T_{char}$, $θ_{char}$ and $ΔC_p$ including standard errors, see stderror
.
Arguments
meas
... Tuple with the loaded measurement data, seeload_chromatograms
.col_input
... Named tuple withcol_input.L
the column length in m andcol_input.d
the column diameter in mm.se_col=true
... Iftrue
the standard errors (from the Hessian matrix, seestderror
) of the estimated parameters are added as separate columns to the result dataframe. Iffalse
the standard errors are added to the values asMasurement
type.
Output
res
... Dataframe with the optimized parameters and the found minima.Telu_max
... The maximum of elution temperatures every solute experiences in the measured programs.
RetentionParameterEstimator.method_m2
— Methodmethod_m2(meas; se_col=true)
Estimation of the column diameter $d$ and three retention parameters $T_{char}$, $θ_{char}$ and $Δ C_p$ including standard errors, see stderror
. In a first run all four parameters are estimated for every substance separatly, resulting in different optimized column diameters. The mean value of the column diameter is used for a second optimization using this mean diameter and optimize the remainig thre retention parameters $T_{char}$, $θ_{char}$ and $Δ C_p$.
Arguments
meas
... Tuple with the loaded measurement data, seeload_chromatograms
.se_col=true
... Iftrue
the standard errors (from the Hessian matrix, seestderror
) of the estimated parameters are added as separate columns to the result dataframe. Iffalse
the standard errors are added to the values asMasurement
type.
Output
res
... Dataframe with the optimized parameters and the found minima.Telu_max
... The maximum of elution temperatures every solute experiences in the measured programs.
RetentionParameterEstimator.method_m3
— Methodmethod_m3(meas; se_col=true)
Estimation of the column diameter $d$ and three retention parameters $T_{char}$, $θ_{char}$ and $Δ C_p$ including standard errors, see stderror
. Brute-force method, where all parameters (3n+1
for n
substances) are estimate in one optimization.
ATTENTION: This method can take long time to finish. The more substances, the longer it takes.
Arguments
meas
... Tuple with the loaded measurement data, seeload_chromatograms
.se_col=true
... Iftrue
the standard errors (from the Hessian matrix, seestderror
) of the estimated parameters are added as separate columns to the result dataframe. Iffalse
the standard errors are added to the values asMasurement
type.
Output
res
... Dataframe with the optimized parameters and the found minima.Telu_max
... The maximum of elution temperatures every solute experiences in the measured programs.
RetentionParameterEstimator.opt_Kcentric
— Methodopt_Kcentric(x, p)
Function used for optimization of the loss-function in regards to the three K-centric parameters.
Arguments
x
... 3n-vector of the three K-centric parameters of n solutes. Elements 1:n are Tchar, n+1:2n are θchar and 2n+1:3n are ΔCp values.p
... vector containing the fixed parameters:tR = p[1]
... vector of the measured retention times in seconds.substance_list = p[2]
... vector of the names of the solutes related totR
andprog
.L = p[3]
... number of the length of the column in m.d = p[4]
... number of the diameters of the column in m.prog = p[5]
... vector of structure GasChromatographySimulator.Programs containing the definition of the GC-programs.opt = p[6]
... struture GasChromatographySimulator.Options containing the settings for options for the simulation.gas = p[7]
... string of name of the mobile phase gas.metric = p[8]
... string of the metric used for the loss function (squared
orabs
).
Output
sum((tR.-tRcalc).^2)
... sum of the squared residuals over m GC-programs and n solutes.
RetentionParameterEstimator.opt_dKcentric
— Methodopt_dKcentric(x, p)
Function used for optimization of the loss-function in regards to the colum diameter and the three K-centric parameters.
Arguments
x
... 1+3n-vector of column diameter and the three K-centric parameters of n solutes. Element 1 is the diameter, 2:n+1 are Tchar, n+2:2n+1 are θchar and 2n+2:3n+1 are ΔCp values.p
... vector containing the fixed parameters:tR = p[1]
... vector of the measured retention times in seconds.substance_list = p[2]
... vector of the names of the solutes related totR
andprog
.L = p[3]
... number of the length of the column in m.prog = p[4]
... vector of structure GasChromatographySimulator.Programs containing the definition of the GC-programs.opt = p[5]
... struture GasChromatographySimulator.Options containing the settings for options for the simulation.gas = p[6]
... string of name of the mobile phase gas.metric = p[7]
... string of the metric used for the loss function (squared
orabs
).
Output
sum((tR.-tRcalc).^2)
... sum of the squared residuals over m GC-programs and n solutes.
RetentionParameterEstimator.optimize_Kcentric
— Methodoptimize_Kcentric(tR::Vector{T}, substance_list, col, prog, Tchar_e::Vector{T}, θchar_e::Vector{T}, ΔCp_e::Vector{T}; method=RetentionParameterEstimator.NewtonTrustRegion(), opt=RetentionParameterEstimator.std_opt, maxiters=10000, maxtime=600.0, metric="squared") where T<:Number
Optimize the K-centric retention parameters for a given set of substances, measured retention times and used programs.
Arguments
tR::Vector{T}
: Vector of retention times.substance_list
: List of substances to be optimized.col
: Column characteristics including lengthL
, diameterd
, and gas typegas
.prog
: Program conditions.Tchar_e::Vector{T}
: Initial estimates for characteristic temperatures.θchar_e::Vector{T}
: Initial estimates for characteristic phase ratios.ΔCp_e::Vector{T}
: Initial estimates for heat capacity changes.method
: Optimization method to be used (default:RetentionParameterEstimator.NewtonTrustRegion()
).opt
: Optimization options (default:RetentionParameterEstimator.std_opt
).maxiters
: Maximum number of iterations (default: 10000).maxtime
: Maximum time allowed for optimization in seconds (default: 600.0).metric
: Metric to be used for optimization (default: "squared").
Returns
opt_sol
: The solution of the optimization problem.
RetentionParameterEstimator.optimize_dKcentric
— Methodoptimize_dKcentric(tR::Vector{T}, substance_list, col, prog, Tchar_e::Vector{T}, θchar_e::Vector{T}, ΔCp_e::Vector{T}; method=RetentionParameterEstimator.NewtonTrustRegion(), opt=RetentionParameterEstimator.std_opt, maxiters=10000, maxtime=600.0, metric="squared") where T<:Number
Optimize the column diameter and K-centric retention parameters for a given set of substances, measured retention times and used programs.
Arguments
tR::Vector{T}
: Vector of retention times.substance_list
: List of substances to be optimized.col
: Column characteristics including lengthL
and gas typegas
.prog
: Program conditions.Tchar_e::Vector{T}
: Initial estimates for characteristic temperatures.θchar_e::Vector{T}
: Initial estimates for characteristic phase ratios.ΔCp_e::Vector{T}
: Initial estimates for heat capacity changes.method
: Optimization method to be used (default:RetentionParameterEstimator.NewtonTrustRegion()
).opt
: Optimization options (default:RetentionParameterEstimator.std_opt
).maxiters
: Maximum number of iterations (default: 10000).maxtime
: Maximum time allowed for optimization in seconds (default: 600.0).metric
: Metric to be used for optimization (default: "squared").
Returns
opt_sol
: The solution of the optimization problem.
RetentionParameterEstimator.prog_list
— Methodprog_list(prog, tRs)
Generate a vector of program conditions, excluding rows with missing values in tRs
.
Arguments
prog
: A program condition to be repeated.tRs
: A 2D array where each row corresponds to a set of retention times and columns to different substances.
Returns
- A vector of program conditions corresponding to non-missing rows in
tRs
.
RetentionParameterEstimator.reference_holdup_time
— Methodreference_holdup_time(prog, L, d, gas; control="Pressure")
Calculate the reference holdup time for the GC program prog
for a column with length L
and diameter d
and gas
as mobile phase. The reference holdup time is the holdup time at the reference temperature 150°C.
RetentionParameterEstimator.separate_error_columns
— Methodseparate_error_columns(res)
If the result dataframe res
contains columns with Measurements
typed values, these columns will be split in two. The first column contains the value and uses the original column name. The second column contains the uncertainty and the name of the column is the original name with an added "_uncertainty". If the column is not of the Measurements
type, it will be copied as is for the new dataframe.
RetentionParameterEstimator.stderror
— Methodstderror(meas, res)
Calculation of the standard error of the found optimized parameters using the hessian matrix at the optima.
Attention
The used loss-function is hard coded in the function opt_Kcentric
and has to be changed if another loss-function is used.
Arguments
meas
... Tuple with the loaded measurement data, seeload_chromatograms
.res
... Dataframe with the result of the optimization, seeestimate_parameters
.
Optional parameters
col_input
... Named tuple withcol_input.L
the column length in m andcol_input.d
the column diameter in mm. If this parameter is not gicen, than these parameters are taken frommeas
.
Output
stderrors
... Dataframe with the standard errors of the optimized parameters.Hessian
... The hessian matrix at the found optims.
RetentionParameterEstimator.substance_list
— Methodsubstance_list(substance_names, tRs)
Generate a vector of substance names, excluding rows with missing values in tRs
.
Arguments
substance_names
: A vector of substance names.tRs
: A 2D array where each row corresponds to a set of retention times and columns to different substances.
Returns
- A vector of substance names corresponding to non-missing rows in
tRs
.
RetentionParameterEstimator.tR_calc
— MethodtR_calc(Tchar, θchar, ΔCp, L, d, prog, gas; opt=std_opt)
Calculates the retention time tR for a solute with the K-centric parameters Tchar
θchar
and ΔCp
for a column with length L
, internal diameter d
, the (conventional) program prog
, options opt
and mobile phase gas
. For this calculation only the ODE for the migration of a solute in a GC column is solved, using the function GasChromatographySimulator.solving_migration
.
RetentionParameterEstimator.tR_τR_calc
— MethodtR_τR_calc(Tchar, θchar, ΔCp, L, d, df, prog, Cag, t₀, τ₀, gas; opt=std_opt)
Calculates the retention time tR for a solute with the K-centric parameters Tchar
θchar
and ΔCp
for a column with length L
, internal diameter d
, film thickness df
, the (conventional) program prog
, diffusirivity coefficient Cag
, start time t₀
, initial peak width τ₀
, options opt
and mobile phase gas
. For this calculation the ODE-system for the migration of a solute and peak width in a GC column is solved, using the function GasChromatographySimulator.solving_odesystem_r
. The result is a tuple of retention time tR
and peak width τR
.