API Documentation

Core modules

Model

A model with thermodynamic properties

class thermo_flux.core.model.ThermoModel(model, pH: ~typing.Dict[str, ~pint.Quantity] = None, I: ~typing.Dict[str, ~pint.Quantity] = None, pMg: ~typing.Dict[str, ~pint.Quantity] = None, T: ~typing.Dict[str, ~pint.Quantity] = None, phi: ~typing.Dict[str, ~pint.Quantity] = None, gdiss_lim: ~pint.Quantity = None, dfG0_cov_sqrt=None, update_thermo_info=False, rmse_inf=<Quantity(100000.0, 'kilojoule / mole')>, cc=None, lc=None, m=None, max_drG=None, split_biomass=True, update_biomass_dfG0=False, add_charge_exchange=True, phi_dict: ~typing.Dict[str, ~pint.Quantity] = None, compartment_parents: ~typing.Dict[str, str] = None, **kwargs)

Bases: Model

Thermodynamic model object : thermodynamic extension of a COBRA model.

This class extends the base Model cobra class adds attributes/methods required for thermodynamic calculations :
  • Compartmental pH, ionic strength, pMg, temperature, and membrane potential differences.

  • containers for the proton, charge, and magnesium metabolites of the model.

  • Optional biomass splitting, biomass ΔfG’° estimation, and adding charge-exchange reactions.

The main method is ‘’update_thermo_info’’ which calls equilibrator to calculate the thermodynamic properties of the metabolites and reactions (see drg_tools module).

property I
property T
add_TFBA_variables(m=None, conds=[''], error_type='covariance', qnorm='sep_norm', alpha=0.95, epsilon=0.5, nullspace=None, gdiss_constraint=False, sigmac_limit=12.3, split_v=False, big_M=False)
add_charge_exchange()

function to add charge transport and exchange reactions to the model. Skips compartment if there is already a proton transport reaction

add_proton_exchange()

function to add charge transport and exchange reactions to the model. Skips compartment if there is already a proton transport reaction

property cc
property charge_dict
property compartment_parents
class compartment_parents_dict(model)

Bases: dict

A dict subclass to represent compartment parents with custom string representation.

property dfG0_cov_sqrt

Gets the square root of the covariance of the standard formation energies

property gdiss_lim
get_charge_transporters()

Get the charge transporters in the model

get_compounds(search=False, update_annotations=False, check_consisitency=True)
get_proton_transporters()

Get the proton transporters in the model

property inner_compartments
property lc
property m
property max_drG
property mg_dict
property pH
property pMg
property phi
property phi_dict
property proton_dict
regression(conds, flux_data, metabolite_data, volume_data, conc_fit=True, flux_fit=True, drG_fit=True, resnorm=1, qm_resnorm=2, error_type='covariance', conc_units=None, extracellular=None)
property rmse_inf
solution()

return a thermo optimization solution object for the model

split_biomass(biomass_rxn=None)

Automatically expand the biomass reaction to include transport of biomass out of the cell. This ensures the gibbs energy balance is maintained as biomass contains protons that have been transported into the cell, the energy of those protons leaving the cell must also be accounted for.

update_biomass_dfG0()

update the dfG0 of the biomass reaction based on the formula of biomass_c if no formula is defined it will be automatically calcualted from the biomass_c reaction

update_parameters()
update_thermo_info(fit_unknown_dfG0=False, search=False, round_dp=False, report=False)

Compute thermodynamic parameters for all metabolites and all reactions.

Steps: - Identify compounds and initialize charge, proton, and magnesium lookup tables. - Estimate metabolite ΔfG’° (mean) and the square root of ΔfG’° covariance

using drg_tools.calc_dfG0prime. Store ΔfG’°, ΔfG’ (pH-corrected), and standard errors.

  • Compute reaction ΔrG° and ΔrG’° means, and ΔrGm′ (including physicochemical

    corrections: proton, Mg, ionic strength, membrane potential).

  • Using the stoichiometric matrix S, form the square root of drG0 covariancestandard_dgr_Q = Sᵀ @ dfG0_cov_sqrt
    • Zero out numerically small elements in standard_dgr_Q and drop all-zero degrees of freedom.

  • For transport processes (2 compartment reactions), compute transport contributions to ΔrG0prime using drg_tools.calc_drGtransport.

  • Stores standard_dgr_Q as _drG0_cov_sqrt for further use in the gurobi model.

Returns: - If report=True, returns a thermodynamic summary report.

thermo_flux.core.model.copy(self) ThermoModel

Modified from COBRApy: Provide a partial ‘deepcopy’ of the Model.

All the Metabolite, Gene, and Reaction objects are created anew but in a faster fashion than deepcopy.

Returns

cobra.Model: new model copy

thermo_flux.core.model.custom_formatwarning(message, category, filename, lineno, line=None)

Reactions

A reaction with thermodynamic properties

class thermo_flux.core.reaction.ThermoReaction(reaction, model=None, drG0: float = 0, ignore_snd: bool = False)

Bases: Reaction

Thermodynamic reaction object: thermodynamic extension of a COBRA reaction.

This class adds attributes/methods required for thermodynamic calculations: - Standard, transformed, and physicochemically corrected Gibbs energies

(ΔrG°, ΔrG′°, ΔrGm′).

  • Transport processes- related reaction energy terms (total transport, proton, and charge components).

  • Flags for ignoring second-law constraints in boundary, water transport, or biomass reactions.

  • methods to balance reactions are called from the drg_tools module (net_elements,transported_c_h)

Initialization: - Copies non-structural attributes from the COBRA reaction. - Wraps participating metabolites in ThermoMetabolite objects. - Assigns transport- and second-law-related flags based on reaction content.

property balanced

Defines if a reaction has already been charge and proton balanced

check_atom_balance(round_dp=False, electrons=False, pMg=<Quantity(14, 'dimensionless')>)

Check the mass balance of a reaction

Uses the eQuilibrator compound atom bag to check the mass balance. If this is not available the COBRA model formula is used. Note magnesium balance is ignored by default.

Parameters

rxncobra.Reaction

Reaction to check

round_dpint

Number of decimal places to round to

electronsbool

Include electrons in the balance

pMgfloat

pMg of the reaction

Returns

net_elementsdict

a dict of {element: amount} for unbalanced elements.

This should be empty for balanced reactions.

property drG

Gibbs free energy of the reaction. drG0’ + RTlnC

property drG0

Standard Gibbs energy of a reaction before transformation.

property drG0prime

Transformed Gibbs energy of a reaction.

property drG_SE

standard error on the gibbs free energy estimate. Can be set explicitly or calcualted from the covariance

property drG_c_transport

Charge transport component of the transport component of the gibbs energy of a transporter reaction

property drG_h_transport

Proton transport component of the transport component of the gibbs energy of a transporter reaction

property drGmprime

Transformed Gibbs energy of a reaction.

property drGtransform

Gibbs energy to be added to the standard gibbs energy of a reaction to get the transformed gibbs free energy

property drGtransport

Transport componenent of gibbs energy of a transporter reaction.

property ignore_snd

Flag to ignore the second law constraint for a reaction

property inner_compartment

Return the inner compartment of a transport reaction

net_elements(balance_mg=False, round_dp=False, rxn_already_balanced=True)
split_reaction()

Split a multiple compartment reaction into subunits of 2 compartments or less

transported_c_h(round_dp=False)
property transported_charge

Free or additional charge transported in a transport reaction. Represents a generic positive ion

property transported_h

Free or additional protons transported in a transport reaction

property transported_mets

explicitly define metabolite transported in a reaction. This is useful for transporters that include chemical transformation during transport. Leave empty for automatic calculation of transported metabolites

Metabolites

A metabolite with thermodynamic properties

class thermo_flux.core.metabolite.ThermoMetabolite(metabolite, model=None, upper_bound: ~pint.Quantity = <Quantity(10, 'millimolar')>, lower_bound: ~pint.Quantity = <Quantity(0.1, 'micromolar')>, concentration: ~pint.Quantity = <Quantity(1, 'molar')>, accession: str = None, dfG0: ~pint.Quantity = None, dfG0prime: ~pint.Quantity = None, compound=None, dfG_SE: ~pint.Quantity = None, redox: bool = False, biomass: bool = False, unknown: bool = False, ignore_conc: bool = False)

Bases: Metabolite

Thermodynamic metabolite object : metabolite subclass that adds thermodynamic state and concentration bounds.

This class extends the base Metabolite and adds attributes/methods required for thermodynamic calculations: - standard Gibbs energies (dfG0, dfG0prime) - concentration bounds (lower_bound and upper_bound) - equilibrator compound (compound), used for chemical species calculations

Some thermodynamic exceptions applied at initialization: - proton and water concentration are fixed at 1 M - redox and biomass attributes enforce .ignore_conc

property accession
average_charge_protons(pH=None, pMg=None, ionic_strength=None, temperature=None, accuracy=0.1, round_dp=False, cobra_formula=False)
property biomass
check_consistency(ignore_H=True)

Check the consistency of a metabolite between the original model definiton and the metabolite identified in equilibrator database for thermodynamic analysis

property compound
property concentration
property dfG0
property dfG0prime
property dfG_SE
property dfGprime
property ignore_conc
property lower_bound
major_microspecies(pH=None, pMg=None, ionic_strength=None, temperature=None)
property redox
property unknown
property upper_bound

Solver modules : gurobi

thermo_flux.solver.gurobi.add_TFBA_variables(tmodel, m, conds=[''], error_type='linear', qnorm=1, alpha=0.95, epsilon=0.5, nullspace=None, gdiss_constraint=False, sigmac_limit=12.3, split_v=False, big_M=False)

Generates gurobi model from a thermodynamic model. The gurobi model contains :

  • Flux variables, Gibbs energy variables (ΔrG°, ΔrG, ΔrG_conc), direction binaries, and (optionally) split forward/backward fluxes.

  • Concentration log-variables with bounds from metabolite constraints.

  • Mixed-integer constraints linking reaction direction, flux sign, and ΔrG.

  • Mass-balance constraints for each condition.

  • Error models for uncertainty in drG0 using covariance (recommended) or box errors.

  • Optional Gibbs dissipation rate constraint (Niebel et al., 2019).

  • Objective is defined using the original FBA objective.

Returns:
  • Updated Gurobi model with TFBA variables and constraints.

  • Dictionary of created variable blocks.Contains constraints for FBA optimization

thermo_flux.solver.gurobi.calc_conc_bounds(tmodel, conds, metabolite_data=None, extracellular_data=None, volume_data=None, extracellular='e', CI=2.6, conc_units=None)

Calcualte concentration bounds for individual metabolites based on data from whole cell metabolomics and extracellular measurements.

Parameters

tmodel: ThermoModel

Thermodynamic model

conds: list

List of conditions to calculate bounds for (must match data)

metabolite_data: pd.DataFrame

Dataframe of metabolite concentrations (mean and sd) for each condition

extracellular_data: pd.DataFrame

Dataframe of extracellular concentrations (lo and up) for each condition (to match legacy GAMS data)

volume_data: pd.DataFrame

Dataframe of volume fractions for each compartment for each condition

extracellular: str

Name of extracellular compartment (default = ‘e’)

CI: float

Zvalue for Confidence interval for data (default = 2.6 to match GAMS legacy ~99.9% CI)

conc_units: str

Units of concentration data (default = mM)

Returns

data_bounds_conc_df: pd.DataFrame

Dataframe of bounds for each metabolite for each condition

thermo_flux.solver.gurobi.compute_IIS(tmodel)

Creates an Irreducible Infeasible Subsystem (IIS) for the gurobi model, solves it and prints the IIS constraints and bounds. Gurobi documentation : an IIS is a subset of the constraints and variable bounds in the infeasible model with the following properties :

  • It is still infeasible.

  • If a single constraint or bound is removed, the subsystem becomes feasible.

https://docs.gurobi.com/projects/optimizer/en/current/concepts/logging/iis.html Could be used on any gurobi model (except multiscenario models).

thermo_flux.solver.gurobi.drG_bounds(tmodel, concentration=True, drG_error=True, alpha=0.95, condition_index=0)

calculate the drg bounds for a given model

modified from equilibrator-api examples https://equilibrator.readthedocs.io/en/latest/equilibrator_examples.html

thermo_flux.solver.gurobi.gdiss_model(tmodel)

Print the individual gibbs energy balance for different components of drG. This can be used to identify why the Gibbs enery balance might not be met

thermo_flux.solver.gurobi.gdiss_var(tmodel, var, verbose=False)

Return the gibbs energy dissiaption of a specific variable

thermo_flux.solver.gurobi.get_solution(tmodel)
thermo_flux.solver.gurobi.model_start(tmodel, sol_file, ignore_vars=[], fix_vars=[], fix='start')

Set the model start point from a saved gurobi solution fix_vars: the variables that you want to import the starting points from fix: [‘start’, hint, bounds], Note that for some solutions numerical issue make an imported start point ifeasible. Therefore it is beneficial to only fix some start points and allow the remining unknwon variables to be calcualted

thermo_flux.solver.gurobi.multi_scenario_sol(tmodel, var)

Return the solution for a variable from a multiple scenario optimization

thermo_flux.solver.gurobi.read_bounds(bounds_file)

Read in bounds text file from hpc optimisation and returns a dataframe of the bounds. Parameters ———-

bounds_file: str

Path to bounds file with format variable_name: [lb, ub]

Returns

drG_bounds: pandas.DataFrame

DataFrame of bounds for each drG

thermo_flux.solver.gurobi.regression(tmodel, m, mvars, conds, flux_data, metabolite_data, volume_data, conc_fit=True, flux_fit=True, drG_fit=True, resnorm=1, qm_resnorm=2, error_type='covariance', conc_units=None, extracellular=None)

Add variables and constraints for regression to data. If conc_fit is True add constraints to regress metabolite concentration using total_cell_conc. Modifies the objective of the gurobi model : minimize sum of residuals.

Parameters

tmodel: ThermoModel

ThermoModel object to be used for regression

m: gurobipy.Model

gurobipy model object to add variables and constraints to

mvars: dict

dictionary of variables to be updated with new variables

conds: list

list of conditions to be used for regression

flux_data: pandas.DataFrame

dataframe of flux data to be regressed to

metabolite_data: pandas.DataFrame

dataframe of metabolite data to be regressed to

volume_data: pandas.DataFrame

dataframe of volume data to be used for total cell concentration constraints

conc_fit: bool

if True add constraints for metabolite concentration regression

flux_fit: bool

if True add constraints for flux regression

drG_fit: bool

if True add constraints for drG regression

resnorm: int

1 for linear residual formulation (sum of absolute differences), 2 for quadratic residual formulation (sum of squared differences)

error_type: str

‘covariance’ for covariance based error, ‘linear’ for linear error

conc_units: str

units of concentration data to be regressed to

Returns

m: gurobipy.Model

gurobipy model object with added variables and constraints

mvars: dict

dictionary of variables now updated with new variables

thermo_flux.solver.gurobi.regression_legacy(tmodel, m, mvars, conds, flux_data, metabolite_data, volume_data, conc_fit=True, flux_fit=True, drG_fit=True, resnorm=1, error_type='linear')

legacy function to set up regession optimization to replicated GAMS results included incorect assigning of cytosolic volume to entire cell

thermo_flux.solver.gurobi.total_cell_conc(tmodel, conds=[''], metabolites=[], volume_data=None, extracellular=None, GAMS_style=False)

Add total cellular concentration constraints (whole-cell) to a TFBA model, used to set bounds from metabolome data or for metabolite concentration regression. Metabolome data indeed provide whole-cell concentrations of metabolites whereas in the model we have compartmental concentrations. To define whole-cell concentrations, we add :

  • Variables for whole-cell concentrations for the metabolites in argument (typically measured ones) and a list of conditions.

  • Variables for exponential of the corresponding compartmental concentrations (linked with ln c = exp (new_var), where new var is called met_conc_dist)

  • Constraints linking compartment concentrations to total cell concentrations, scaled by compartment volumes : C_wholecell=Sum_compartment(C_compartment*vol_compartment)

Returns:
  • Updated Gurobi model with added concentration variables and constraints.

  • Dictionary of model variables including ‘cell_conc’ and ‘met_conc_dist’.

thermo_flux.solver.gurobi.variability_analysis(m, vars=[])

Set up a multiscenario optimisation problem to perform variability analysis on a given variable. The model will have 2 scenarios for each variable, one for the lower bound and one for the upper bound.

Parameters

tmodel: ThermoModel

ThermoModel object to be used for variability analysis

vars: list

list of variables to perform variability analysis on

Returns

m: gurobipy.Model

gurobipy model object with added variables and constraints

thermo_flux.solver.gurobi.variability_results(m)

gets results from a local variability analysis

Parameters

m: gurobipy.Model

gurobipy model object

Returns

obj_val: dict

dictionary of objective values (actual best incumbent) for each variable

obj_bound: dict

dictionary of objective bounds (best known bound) for each variable

optimal_bounds: dict

dictionary of optimal bounds for each variable

MIPGaps: dict

dictionary of MIPGaps (MIPGap = abs((ObjBound-ObjVal)/ObjVal)) for each variable

thermo_flux.solver.gurobi.variable_scan(tmodel, scan_range, var)

Scan a variable over a range of values and return the solution for each value

Tools modules : drg_tools

thermo_flux.tools.drg_tools.add_transporter_charge_varaint(reaction: Any, charge_state: int, round_dp: bool | int = False) Any

Add a transporter variant with a specific net charge transport. Protons are added to either side of the equation to represent different charge states.

Parameters

reactionThermoReaction

Reaction to add variants of.

charge_stateint

Desired net charge to be transported.

round_dpUnion[bool, int], optional

Decimal places for rounding (default is False).

Returns

ThermoReaction

ThermoReaction with the desired net charge transport.

thermo_flux.tools.drg_tools.add_transporter_varaints(reaction: Any, add_charge_neutral: bool = True, balance_charge: bool = False, round_dp: bool | int = False) Set[Any]

Add all transporter variants of a reaction. Variants are added for all species of transported metabolites with an abundance of >10% in the inner compartment. Additional transporters to represent a charge neutral transporter can also be added.

Parameters

reactionThermoReaction

Reaction to add transporter variants of.

add_charge_neutralbool, optional

If True, a charge neutral variant will be added (default is True).

balance_chargebool, optional

If True, the reaction will be balanced for charge (default is False).

round_dpUnion[bool, int], optional

Decimal places for rounding (default is False).

Returns

Set[ThermoReaction]

Set of ThermoReaction objects representing the transporter variants.

thermo_flux.tools.drg_tools.calc_average_charge_protons(compound: Any, pH: Any, pMg: Any, ionic_strength: Any, temperature: Any, accuracy: float = 0.1, round_dp: bool | int = False, cobra_formula: bool = False) Tuple[float | None, float | None, float | None, DataFrame | None]

Calculate the average charge, protons, and magnesiums for a compound under specific conditions. Modified from Elad Noor.

Parameters

compoundThermoMetabolite or eQ_compound

The compound to analyze.

pHQuantity

The pH of the environment.

pMgQuantity

The pMg of the environment.

ionic_strengthQuantity

The ionic strength of the environment.

temperatureQuantity

The temperature of the environment.

accuracyfloat, optional

The threshold for microspecies abundance to be considered (default is 0.1).

round_dpUnion[bool, int], optional

The number of decimal places to round to, or False to skip rounding (default is False).

cobra_formulabool, optional

Whether to force using the COBRA model formula (default is False).

Returns

Tuple[Optional[float], Optional[float], Optional[float], Optional[pd.DataFrame]]

A tuple containing: - Average charge. - Average number of protons. - Average number of magnesium atoms. - DataFrame of microspecies distribution.

thermo_flux.tools.drg_tools.calc_biomass_formula(biomass_rxn: Any) Tuple[str, Dict[str, float]]

Calculate biomass elemental composition from the biomass equation.

Parameters

biomass_rxnThermoReaction

The biomass reaction.

Returns

Tuple[str, Dict[str, float]]

A tuple containing: - The biomass formula string. - The biomass atom bag (dictionary of elements).

thermo_flux.tools.drg_tools.calc_dfG0(tmodel: Any, fit_unknown_dfG0: bool = False) Tuple[Any, Any, Any]

Calculate the standard Gibbs energy of formation (dfG0) for metabolites in the model.

Parameters

tmodelThermoModel

The thermodynamic model containing metabolites.

fit_unknown_dfG0bool, optional

Whether to fit unknown dfG0 values (default is False).

Returns

Tuple[Quantity, Quantity, Quantity]

A tuple containing: - Mean dfG0 values. - Square root of the covariance matrix of dfG0. - Basis for unknown metabolites.

thermo_flux.tools.drg_tools.calc_dfG0prime(tmodel: Any, fit_unknown_dfG0: bool = False) Tuple[Any, Any, Any, Any]

Calculate the standard transformed Gibbs energy of formation (dfG0’) for metabolites.

Parameters

tmodelThermoModel

The thermodynamic model.

fit_unknown_dfG0bool, optional

Whether to fit unknown dfG0 values (default is False).

Returns

Tuple[Quantity, Quantity, Quantity, Quantity]

A tuple containing: - Mean dfG0 values. - Mean dfG0’ values. - Square root of the covariance matrix of dfG0. - Basis for unknown metabolites.

thermo_flux.tools.drg_tools.calc_dfG_transform(met: Any) Any

Calculate the Legendre transform for a metabolite based on model compartment information.

Parameters

metThermoMetabolite

The metabolite to calculate the transform for.

Returns

Quantity

The calculated Gibbs energy transform.

thermo_flux.tools.drg_tools.calc_drG0(S: ndarray, dfG0: Any) Any

Calculate the standard Gibbs energy of a reaction matrix(drG0).

Parameters

Snp.ndarray

Stoichiometric matrix.

dfG0Quantity

Standard Gibbs energy of formation.

Returns

Quantity

Standard Gibbs energy of reaction.

thermo_flux.tools.drg_tools.calc_drGtransport(reaction: Any, round_dp: bool | int = False, rxn_already_balanced: bool = True) Tuple[Any, Any, Any]

Calculate the Gibbs energy of transport (drGtransport). Note this will balance the reaction if rxn_already_balanced is False as the reaction must be balanced for accurate calculation of drGtransport.

Parameters

reactionThermoReaction

The transport reaction.

round_dpUnion[bool, int], optional

Decimal places for rounding or stoichiometry of protons when automatically balancing the reaction (default is False).

rxn_already_balancedbool, optional

Whether the reaction is already balanced (default is True).

Returns

Tuple[Quantity, Quantity, Quantity]

A tuple containing: - Total Gibbs energy of transport. - Gibbs energy of proton transport. - Gibbs energy of electrostatic potential.

thermo_flux.tools.drg_tools.calc_model_drG0(tmodel: Any) Any

Calculate the standard Gibbs energy of reaction (drG0) for all reactions in the model.

Parameters

tmodelThermoModel

The thermodynamic model.

Returns

Quantity

Standard Gibbs energy of reaction for all reactions.

thermo_flux.tools.drg_tools.calc_model_drG0prime(tmodel: Any) Any

Calculate the standard transformed Gibbs energy of reaction (drG0’) for all reactions in the model.

Parameters

tmodelThermoModel

The thermodynamic model.

Returns

Quantity

Standard transformed Gibbs energy of reaction for all reactions.

thermo_flux.tools.drg_tools.calc_phys_correction(tmodel: Any) Any

Calculate the physiological concentration correction for Gibbs energy assuming metabolites are at 1 mM not 1 M.

Parameters

tmodelThermoModel

The thermodynamic model.

Returns

Quantity

The physiological concentration correction term.

thermo_flux.tools.drg_tools.calc_transported_mets(reaction: Any) Dict[Any, float]

Calculate metabolites transported in a transport reaction.

Parameters

reactionThermoReaction

The transport reaction.

Returns

Dict[ThermoMetabolite, float]

Dictionary of transported metabolites and their stoichiometry.

thermo_flux.tools.drg_tools.calculate_biomass_dfG0(biomass: Any) Any

Calculate the formation energy of biomass from the biomass formula.

Parameters

biomassThermoMetabolite

The biomass metabolite object.

Returns

Quantity

Formation energy of biomass. Units are defined as kJ/mol for compatibility with other reactions, but actual units correspond to J/gDW due to conversion in biomass equation.

thermo_flux.tools.drg_tools.charge_dict(tmodel: Any) Dict[str, Any]

Identify charge metabolites in the model for balancing.

Parameters

tmodelThermoModel

The thermodynamic model.

Returns

Dict[str, ThermoMetabolite]

Dictionary mapping compartment IDs to charge metabolites.

thermo_flux.tools.drg_tools.comp_split(reaction: Any, compartment: str) Dict[Any, float]

Separate reaction metabolites into those belonging to a specific compartment.

Parameters

reactionThermoReaction

The reaction to split.

compartmentstr

The compartment ID to filter by.

Returns

Dict[ThermoMetabolite, float]

Dictionary of metabolites in the specified compartment and their stoichiometry.

thermo_flux.tools.drg_tools.custom_formatwarning(message, category, filename, lineno, line=None)
thermo_flux.tools.drg_tools.dfGbm(formula: Dict[str, int | float] = {}, units: str = 'kJ/mol', Mr_bio: float | None = None) Tuple[Any, Any, float, ndarray]

Calculate the formation energy of biomass or macromolecules based on their empirical formula. Modified from Saadat et. al Entropy 2020, 22(3), 277. https://doi.org/10.3390/e22030277 https://gitlab.com/qtb-hhu/thermodynamics-in-genome-scale-models/-/blob/master/EnergyOfFormationBiomass.py?ref_type=heads

Parameters

formulaDict[str, Union[int, float]], optional

Empirical formula of macromolecule (default is empty dict).

unitsstr, optional

Units of the output (default is ‘kJ/mol’).

Mr_biofloat, optional

Molecular weight of biomass in units carbon mol/gDW. Default is None and will be automatically calculated from the elemental composition.

Returns

Tuple[Quantity, Quantity, float, np.ndarray]

A tuple containing: - Gibbs energy of formation (Gf). - Gibbs energy of combustion (Gc). - Degree of reduction (y). - Stoichiometry of combustion reaction.

thermo_flux.tools.drg_tools.formula_dict_to_string(formula: Dict[str, int | float]) str

Convert a formula dictionary to a string representation.

Parameters

formulaDict[str, Union[int, float]]

Dictionary where keys are element symbols and values are counts.

Returns

str

String representation of the formula (e.g., “C6H12O6”).

thermo_flux.tools.drg_tools.generate_combinations(dictionary: Dict[Any, List[Any]]) List[Dict[Any, Any]]

Generate all combinations of values from a dictionary of lists.

Parameters

dictionaryDict[Any, List[Any]]

Dictionary where values are lists of options.

Returns

List[Dict[Any, Any]]

List of dictionaries representing all combinations.

thermo_flux.tools.drg_tools.get_compound(met: str | Compound | Any) Compound | None

Retrieve an equilibrator compound from an identifier string, ThermoMetabolite, or eQ_compound.

Parameters

metUnion[str, eQ_compound, ThermoMetabolite]

The input identifier or object.

Returns

Optional[eQ_compound]

The corresponding equilibrator compound object, or None if not found.

thermo_flux.tools.drg_tools.get_suitable_ids(met: Any, search: bool = False, update_annotations: bool = False) Tuple[Compound | None, Dict, str | None, str | None, bool]

Find suitable identifiers for a metabolite in the equilibrator database.

Parameters

metThermoMetabolite

The metabolite to search for.

searchbool, optional

Whether to perform a broader search using common names (default is False).

update_annotationsbool, optional

Whether to update the metabolite’s annotations with found identifiers (default is False). Warning update_annotations can be slow to extract all annotations from eQuilibrator

Returns

Tuple[Optional[eQ_compound], Dict, Optional[str], Optional[str], bool]

A tuple containing: - The equilibrator compound object (or None). - A dictionary of annotations. - The chemical formula (or None). - The InChI string (or None). - A boolean indicating if a search was performed.

thermo_flux.tools.drg_tools.leading_zeros(decimal: float) bool | int

Calculate the number of leading zeros in a decimal.

Parameters

decimalfloat

The decimal number.

Returns

Union[bool, int]

The number of leading zeros, or False if the number is too small.

thermo_flux.tools.drg_tools.major_microspecies(met: Any, pH: Any, pMg: Any, ionic_strength_M: Any, T_in_K: Any, cobra_formula: bool = False) Tuple[float | None, float | None, float | None]

Calculate the major microspecies of a compound. Useful for transporter calculations.

Parameters

metThermoMetabolite or eQ_compound

The metabolite to analyze.

pHQuantity

The pH of the environment.

pMgQuantity

The pMg of the environment.

ionic_strength_MQuantity

The ionic strength in Molar.

T_in_KQuantity

The temperature in Kelvin.

cobra_formulabool, optional

Whether to force using the COBRA model formula (default is False).

Returns

Tuple[Optional[float], Optional[float], Optional[float]]

A tuple containing: - Charge of the major microspecies. - Number of protons in the major microspecies. - Number of magnesium atoms in the major microspecies.

thermo_flux.tools.drg_tools.mg_dict(tmodel: Any) Dict[str, Any]

Identify magnesium metabolites in the model for balancing.

Parameters

tmodelThermoModel

The thermodynamic model.

Returns

Dict[str, ThermoMetabolite]

Dictionary mapping compartment IDs to magnesium metabolites.

thermo_flux.tools.drg_tools.net_elements(reaction: Any, balance_mg: bool = True, round_dp: bool | int = False, rxn_already_balanced: bool = True) Tuple[Dict[Any, float], Dict[str, float]]

Calculate the net protons, charge, and magnesium of a reaction.

Parameters

reactionThermoReaction

The reaction to analyze.

balance_mgbool, optional

Whether to balance magnesium (default is True).

round_dpUnion[bool, int], optional

Decimal places for rounding (default is False).

rxn_already_balancedbool, optional

Whether the reaction is already balanced (default is True).

Returns

Tuple[Dict[ThermoMetabolite, float], Dict[str, float]]

A tuple containing: - Dictionary of net elements (protons, charge, Mg) to add/remove. - Dictionary of transported free protons.

thermo_flux.tools.drg_tools.new_eq_metabolite(met: Any) Compound | None

Add a new metabolite to the local equilibrator database if it does not already exist.

Parameters

metThermoMetabolite

The metabolite object to add.

Returns

Optional[eQ_compound]

The equilibrator compound object if found or created, else None.

thermo_flux.tools.drg_tools.new_reaction_name(reaction: Any, charge_states: List[int]) str

Generate a new reaction name for a transporter variant with a specific charge state.

Parameters

reactionThermoReaction

Reaction to add variants of.

charge_statesList[int]

List of charge states of the transported metabolites.

Returns

str

Name of the new reaction.

thermo_flux.tools.drg_tools.pka_graph(metabolite: Any, pMg: Any | None = None, ionic_strength: Any | None = None, temperature: Any | None = None, accuracy: float = 0, round_dp: bool | int = False) DataFrame

Return a dataframe of the charge distribution of a metabolite at different pHs.

Parameters

metaboliteThermoMetabolite

The metabolite object.

pMgfloat or Quantity, optional

The pMg value (default is None, uses model default).

ionic_strengthQuantity, optional

The ionic strength (default is None, uses model default).

temperatureQuantity, optional

The temperature (default is None, uses model default).

accuracyfloat, optional

The accuracy threshold for microspecies (default is 0).

round_dpUnion[bool, int], optional

Decimal places for rounding (default is False).

Returns

pd.DataFrame

Dataframe with pH as index and charge as columns.

thermo_flux.tools.drg_tools.proton_dict(tmodel: Any) Dict[str, Any]

Identify protons in the model so they can be added to reactions for balancing.

Parameters

tmodelThermoModel

The thermodynamic model.

Returns

Dict[str, ThermoMetabolite]

Dictionary mapping compartment IDs to proton metabolites.

thermo_flux.tools.drg_tools.reaction_balance(reaction: Any, balance_charge: bool = True, balance_mg: bool = True, round_dp: bool | int = False, rxn_already_balanced: bool = False) None

Balance a reaction for protons, charge, and magnesium.

Parameters

reactionThermoReaction

The reaction to balance.

balance_chargebool, optional

Whether to balance charge (default is True).

balance_mgbool, optional

Whether to balance magnesium (default is True).

round_dpUnion[bool, int], optional

Decimal places for rounding (default is False).

rxn_already_balancedbool, optional

Whether the reaction is already balanced (default is False).

thermo_flux.tools.drg_tools.round_and_normalize(numbers: List[float] | ndarray, round_dp: int = 2) List[float]

Round a list of numbers to a specified decimal place and normalize them so they sum to 1.

Parameters

numbersUnion[List[float], np.ndarray]

The numbers to round and normalize.

round_dpint, optional

The number of decimal places to round to (default is 2).

Returns

List[float]

The rounded and normalized numbers.

thermo_flux.tools.drg_tools.transported_c_h(reaction: Any, round_dp: bool | int = False, verbose: bool = False, rxn_already_balanced: bool = True) Tuple[float, float, float, float, str, str, bool]

Calculate the transported protons and charge for a reaction.

Parameters

reactionThermoReaction

The reaction to analyze.

round_dpUnion[bool, int], optional

Decimal places for rounding (default is False).

verbosebool, optional

Whether to print verbose output (default is False).

rxn_already_balancedbool, optional

Whether the reaction is already balanced (default is True).

Returns

Tuple[float, float, float, float, str, str, bool]

A tuple containing: - Net protons in inner compartment. - Net charge in inner compartment. - Net protons in outer compartment. - Net charge in outer compartment. - Inner compartment ID. - Outer compartment ID. - Whether the reaction is balanced.