6.1.2.2. equilibrator_api.component_contribution

A wrapper for the GibbeEnergyPredictor in component-contribution.

6.1.2.2.1. Module Contents

equilibrator_api.component_contribution.find_most_abundant_ms(cpd: equilibrator_cache.Compound, p_h: equilibrator_api.Q_, p_mg: equilibrator_api.Q_, ionic_strength: equilibrator_api.Q_, temperature: equilibrator_api.Q_) equilibrator_cache.CompoundMicrospecies[source]

Find the most abundant microspecies based on transformed energies.

equilibrator_api.component_contribution.predict_protons_and_charge(rxn: equilibrator_api.phased_reaction.PhasedReaction, p_h: equilibrator_api.Q_, p_mg: equilibrator_api.Q_, ionic_strength: equilibrator_api.Q_, temperature: equilibrator_api.Q_) Tuple[float, float, float][source]

Find the #protons and charge of a transport half-reaction.

class equilibrator_api.component_contribution.ComponentContribution(rmse_inf: equilibrator_api.Q_ = default_rmse_inf, ccache: equilibrator_cache.CompoundCache | None = None, predictor: component_contribution.predict.GibbsEnergyPredictor | None = None)[source]

Bases: object

A wrapper class for GibbsEnergyPredictor.

Also holds default conditions for compounds in the different phases.

property p_h: equilibrator_api.Q_[source]

Get the pH.

property p_mg: equilibrator_api.Q_[source]

Get the pMg.

property ionic_strength: equilibrator_api.Q_[source]

Get the ionic strength.

property temperature: equilibrator_api.Q_[source]

Get the temperature.

property RT: equilibrator_api.Q_[source]

Get the value of RT.

static legacy() ComponentContribution[source]

Initialize a ComponentContribution object with the legacy version.

The legacy version is intended for compatibility with older versions of equilibrator api (0.2.x - 0.3.1). Starting from 0.3.2, there is a significant change in the predictions caused by an improved Mg2+ concentration model.

Return type:

A ComponentContribution object

static initialize_custom_version(rmse_inf: equilibrator_api.Q_ = default_rmse_inf, ccache_settings: component_contribution.ZenodoSettings = DEFAULT_COMPOUND_CACHE_SETTINGS, cc_params_settings: component_contribution.ZenodoSettings = DEFAULT_CC_PARAMS_SETTINGS) ComponentContribution[source]

Initialize ComponentContribution object with custom Zenodo versions.

Parameters:
  • rmse_inf (Quantity, optional) – Set the factor by which to multiply the error covariance matrix for reactions outside the span of Component Contribution. (Default value: 1e-5 kJ/mol)

  • settings (ZenodoSettings) – The doi, filename and md5 of the

Return type:

A ComponentContribution object

get_compound(compound_id: str) equilibrator_cache.Compound | None[source]

Get a Compound using the DB namespace and its accession.

Returns:

cpd

Return type:

Compound

get_compound_by_inchi(inchi: str) equilibrator_cache.Compound | None[source]

Get a Compound using InChI.

Returns:

cpd

Return type:

Compound

search_compound_by_inchi_key(inchi_key: str) List[equilibrator_cache.Compound][source]

Get a Compound using InChI.

Returns:

cpd

Return type:

Compound

search_compound(query: str) None | equilibrator_cache.Compound[source]

Try to find the compound that matches the name best.

Parameters:

query (str) – an (approximate) compound name

Returns:

cpd – the best match

Return type:

Compound

parse_reaction_formula(formula: str) equilibrator_api.phased_reaction.PhasedReaction[source]

Parse reaction text using exact match.

Parameters:

formula (str) – a string containing the reaction formula

Returns:

rxn

Return type:

PhasedReaction

search_reaction(formula: str) equilibrator_api.phased_reaction.PhasedReaction[source]

Search a reaction written using compound names (approximately).

Parameters:

formula (str) – a string containing the reaction formula

Returns:

rxn

Return type:

PhasedReaction

balance_by_oxidation(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.phased_reaction.PhasedReaction[source]

Convert an unbalanced reaction into an oxidation reaction.

By adding H2O, O2, Pi, CO2, and NH4+ to both sides.

get_oxidation_reaction(compound: equilibrator_cache.Compound) equilibrator_api.phased_reaction.PhasedReaction[source]

Generate an oxidation Reaction for a single compound.

Generate a Reaction object which represents the oxidation reaction of this compound using O2. If there are N atoms, the product must be NH3 (and not N2) to represent biological processes. Other atoms other than C, N, H, and O will raise an exception.

standard_dg_formation(compound: equilibrator_cache.Compound) Tuple[float | None, numpy.ndarray | None][source]

Get the (mu, sigma) predictions of a compound’s formation energy.

Parameters:

compound (Compound) – a Compound object

Returns:

  • mu (float) – the mean of the standard formation Gibbs energy estimate

  • sigma_fin (array) – a vector representing the square root of the covariance matrix (uncertainty)

  • sigma_inf (array) – a vector representing the infinite-uncertainty eigenvalues of the covariance matrix

standard_dg(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the chemical reaction energies of a reaction.

Returns:

standard_dg – the dG0 in kJ/mol and standard error. To calculate the 95% confidence interval, use the range -1.96 to 1.96 times this value

Return type:

ureg.Measurement

standard_dg_multi(reactions: List[equilibrator_api.phased_reaction.PhasedReaction], uncertainty_representation: str = 'cov') Tuple[numpy.ndarray, numpy.ndarray][source]

Calculate the chemical reaction energies of a list of reactions.

Using the major microspecies of each of the reactants.

Parameters:
  • reactions (List[PhasedReaction]) – a list of PhasedReaction objects to estimate

  • uncertainty_representation (str) – which representation to use for the uncertainties. cov would return a full covariance matrix. precision would return the precision matrix (i.e. the inverse of the covariance matrix). sqrt would return a sqaure root of the covariance, based on the uncertainty vectors. fullrank would return a full-rank square root of the covariance which is a compressed form of the sqrt result. (Default value: cov)

Returns:

  • standard_dg (Quantity) – the estimated standard reaction Gibbs energies based on the the major microspecies

  • dg_uncertainty (Quantity) – the uncertainty matrix (in either ‘cov’, ‘sqrt’ or ‘fullrank’ format)

standard_dg_prime(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the transformed reaction energies of a reaction.

Returns:

standard_dg – the dG0_prime in kJ/mol and standard error. To calculate the 95% confidence interval, use the range -1.96 to 1.96 times this value

Return type:

ureg.Measurement

dg_prime(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the dG’0 of a single reaction.

Returns:

dg – the dG_prime in kJ/mol and standard error. To calculate the 95% confidence interval, use the range -1.96 to 1.96 times this value

Return type:

ureg.Measurement

standard_dg_prime_multi(reactions: List[equilibrator_api.phased_reaction.PhasedReaction], uncertainty_representation: str = 'cov', minimize_norm: bool = False) Tuple[equilibrator_api.Q_, equilibrator_api.Q_][source]

Calculate the transformed reaction energies of a list of reactions.

Parameters:
  • reactions (List[PhasedReaction]) – a list of PhasedReaction objects to estimate

  • uncertainty_representation (str) – which representation to use for the uncertainties. cov would return a full covariance matrix. precision would return the precision matrix (i.e. the inverse of the covariance matrix). sqrt would return a sqaure root of the covariance, based on the uncertainty vectors. fullrank would return a full-rank square root of the covariance which is a compressed form of the sqrt result. (Default value: cov)

  • minimize_norm (bool) – if True, use an orthogonal projection to minimize the norm2 of the result vector (keeping it within the finite-uncertainty sub-space, i.e. only moving along eigenvectors with infinite uncertainty).

Returns:

  • standard_dg_prime (Quantity) – the CC estimation of the reactions’ standard transformed energies

  • dg_uncertainty (Quantity) – the uncertainty co-variance matrix (in either ‘cov’, ‘sqrt’ or ‘fullrank’ format)

physiological_dg_prime(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the dG’m of a single reaction.

Assume all aqueous reactants are at 1 mM, gas reactants at 1 mbar and the rest at their standard concentration.

Returns:

  • standard_dg_primes (ndarray) – a 1D NumPy array containing the CC estimates for the reactions’ physiological dG’

  • dg_sigma (ndarray) – the second is a 2D numpy array containing the covariance matrix of the standard errors of the estimates. one can use the eigenvectors of the matrix to define a confidence high-dimensional space, or use dg_sigma as the covariance of a Gaussian used for sampling (where ‘standard_dg_primes’ is the mean of that Gaussian).

dgf_prime_sensitivity_to_p_h(compound: equilibrator_cache.Compound) equilibrator_api.ureg.Quantity[source]

Calculate the sensitivity of the chemical formation energy to pH.

Returns:

The derivative of \(\Delta G_f\) with respect to pH, in kJ/mol.

Return type:

Quantity

dg_prime_sensitivity_to_p_h(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Quantity[source]

Calculate the sensitivity of the chemical reaction energy to pH.

Returns:

The derivative of \(\Delta G_r\) with respect to pH, in kJ/mol.

Return type:

Quantity

ln_reversibility_index(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the reversibility index (ln Gamma) of a single reaction.

Returns:

ln_RI – the reversibility index (in natural log scale).

Return type:

ureg.Measurement

standard_e_prime(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the E’0 of a single half-reaction.

Returns:

  • standard_e_prime (ureg.Measurement)

  • the estimated standard electrostatic potential of reaction and

  • E0_uncertainty is the standard deviation of estimation. Multiply it

  • by 1.96 to get a 95% confidence interval (which is the value shown on

  • eQuilibrator).

physiological_e_prime(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the E’0 of a single half-reaction.

Returns:

  • physiological_e_prime (ureg.Measurement)

  • the estimated physiological electrostatic potential of reaction and

  • E0_uncertainty is the standard deviation of estimation. Multiply it

  • by 1.96 to get a 95% confidence interval (which is the value shown on

  • eQuilibrator).

e_prime(reaction: equilibrator_api.phased_reaction.PhasedReaction) equilibrator_api.ureg.Measurement[source]

Calculate the E’0 of a single half-reaction.

Returns:

  • e_prime (ureg.Measurement)

  • the estimated electrostatic potential of reaction and

  • E0_uncertainty is the standard deviation of estimation. Multiply it

  • by 1.96 to get a 95% confidence interval (which is the value shown on

  • eQuilibrator).

dg_analysis(reaction: equilibrator_api.phased_reaction.PhasedReaction) Iterable[Dict[str, object]][source]

Get the analysis of the component contribution estimation process.

Return type:

the analysis results as a list of dictionaries

is_using_group_contribution(reaction: equilibrator_api.phased_reaction.PhasedReaction) bool[source]

Check whether group contribution is needed to get this reactions’ dG.

Return type:

true iff group contribution is needed

multicompartmental_standard_dg_prime(reaction_inner: equilibrator_api.phased_reaction.PhasedReaction, reaction_outer: equilibrator_api.phased_reaction.PhasedReaction, e_potential_difference: equilibrator_api.Q_, p_h_outer: equilibrator_api.Q_, ionic_strength_outer: equilibrator_api.Q_, p_mg_outer: equilibrator_api.Q_ = default_physiological_p_mg, tolerance: float = 0.0) equilibrator_api.ureg.Measurement[source]

Calculate the transformed energies of a multi-compartmental reaction.

Based on the equations from Harandsdottir et al. 2012 (https://doi.org/10.1016/j.bpj.2012.02.032)

Parameters:
  • reaction_inner (PhasedReaction) – the inner compartment half-reaction

  • reaction_outer (PhasedReaction) – the outer compartment half-reaction

  • e_potential_difference (Quantity) – the difference in electro-static potential between the outer and inner compartments

  • p_h_outer (Quantity) – the pH in the outside compartment

  • ionic_strength_outer (Quantity) – the ionic strength outside

  • p_mg_outer (Quantity (optional)) – the pMg in the outside compartment

  • tolerance (Float (optional)) – tolerance for identifying inbalance between inner and outer reactions (default = 0)

Returns:

standard_dg_prime – the transport reaction Gibbs free energy change

Return type:

Measurement

static parse_formula_side(s: str) Dict[str, float][source]

Parse one side of the reaction formula.

static parse_formula(formula: str) Dict[str, float][source]

Parse a two-sided reaction formula.

static create_stoichiometric_matrix_from_reaction_formulas(formulas: Iterable[str]) pandas.DataFrame[source]

Build a stoichiometric matrix.

Parameters:

formulas (Iterable[str]) – String representations of the reactions.

Returns:

The stoichiometric matrix as a DataFrame whose indexes are the compound IDs and its columns are the reaction IDs (in the same order as the input).

Return type:

DataFrame

create_stoichiometric_matrix_from_reaction_objects(reactions: Iterable[equilibrator_api.phased_reaction.PhasedReaction]) pandas.DataFrame[source]

Build a stoichiometric matrix.

Parameters:

reactions (Iterable[PhasedReaction]) – The collection of reactions to build a stoichiometric matrix from.

Returns:

The stoichiometric matrix as a DataFrame whose indexes are the compounds and its columns are the reactions (in the same order as the input).

Return type:

DataFrame