mfapy package¶
Submodules¶
mfapy.carbonsource module¶
carbonsource.py:CarbonSource class in mfapy
The module includes:
CarbonSource class
Todo
Subtraction of natural isotope
- class mfapy.carbonsource.CarbonSource(carbon_sources, emu_dict)¶
Bases:
object
Class for carbon source information.
Instance of this class is generated in the MetabolicModel instance
- generate_carbonsource_MDV(carbonsource=[], correction='no')¶
Generator of MDVs of all EMUs of each carbon source.
This function is called in the mfapy.metabolicmodel.MetabolicModel
- Parameters:
carbonsource (array) – List of carbon source metabolite (optional)
correction (str) – (yes/no) Correction of isotopomer distribution considering natural 13C occurence
- Returns:
Dictionary of mdv data of carbon source.
- Return type:
dict
Examples
>>> cs.generate_carbonsource_MDV()
History:
Correcton of IDV by natural 13C method was improved.
labeled by 2 13C and 3 13C was taken into consideration.
- generate_dict()¶
Generator of a dictionary of MDVs of all EMUs of carbon sources
- Parameters:
required. (Not) –
- Returns:
Dictionary of MDVs of all EMUs of carbon sources
- Return type:
dict
Examples
>>> mdvs = cs.generate_dict()
- set_all_isotopomers(compound, list, correction='no')¶
Setter of IDV data by list of all mass isotopomer distribution
- Parameters:
compound (str) – Name of carbon source.
list (array) – list of mass isotopomer distribution (from 000, 001, 010, to 111)
correction (str) – (yes/no) Correction of isotopomer distribution considering natural 13C occurence
- Returns:
True/False
- Return type:
Booleans
Examples
>>> cs.set_all_isotopomers('AcCoA', [0.3, 0.3, 0.3, 0.1])
- set_carbonsources(filename, correction='no', format='text', output='normal')¶
Setter of isotopomer data of multiple carbon sourses from text file.
- Parameters:
filename (str) –
filename of MDV data with following format:
Name Isotopomer Ratio Asp #0000 0.5 Asp #1111 0.5 AcCoA #00 0.5 AcCoA #11 0.5
correction (str) – (yes/no) Correction of isotopomer distribution considering natural 13C occurence
format (str) – “text” (defalut) or “csv”
output (str) – “normal” (defalut) or “debug”
- Returns:
True/False
- Return type:
Boolean
Examples
>>> cs2.set_isotopomers_from_file('Example_1_carbonsource2.txt', correction = "yes")
- set_each_isotopomer(compound, dict, correction='no')¶
Setter of IDV data of selected mass isotopomers
- Parameters:
compound (str) – Name of carbon source
dict (dict) –
Dictionary of mass isotopomer and its relative abundance:
{'#111': 0.5, '#001': 0.5}
correction (str) – (yes/no) Correction of isotopomer distribution considering natural 13C occurence
- Returns:
True/False
- Return type:
Boolean
Examples
>>> cs.set_each_isotopomers('AcCoA', {'#11':0.5, '#10':0.25}, correction = 'yes')
- set_labeled_compounds(compound, dict, correction='no')¶
Setter of IDV data by distribution of selected mass isotopomers
Following symbols are available:
"[13C]CO2": "#1", "[12C]CO2": "#0", "[13C]THF": "#1", "[12C]THF": "#0", "[1-13C]glucose": "#100000", "[2-13C]glucose": "#010000", "[1,2-13C]glucose": "#110000", "[U-13C]glucose": "#111111", "[1-13C]glutamine": "#10000", "[2-13C]glutamine": "#01000", "[5-13C]glutamine": "#00001", "[U-13C]glutamine": "#11111",
- Parameters:
compound (str) – Name of carbon source.
dict (dict) – Dictionary of mass isotopomer and its relative abundance {‘#111’: 0.5, ‘#001’: 0.5}
correction (str) – (yes/no) Correction of isotopomer distribution considering natural 13C occurence
- Returns:
True/False
- Return type:
Boolean
Examples
>>> cs.set_labeled_compounds('Glc',{'[1_13C]glucose': 0.5, '[U_13C]glucose':0.5}, correction = 'yes')
- show()¶
Method to display contents of CarbonSource instance
- Parameters:
required. (Not) –
- Returns:
Nothing.
Examples
>>> cs.show() Name: Asp Carbon number: 4 #0000 1.0 Name: AcCoA Carbon number: 2 #00 0.5 #10 0.25 #11 0.25
- mfapy.carbonsource.main()¶
mfapy.mdv module¶
mdv.py:MdvData class in mfapy
The module includes:
MdvData class
MdvTimeCourseData class
Todo
I/O
- mfapy.mdv.INV_correcting(formula)¶
- class mfapy.mdv.MdvData(target_fragments)¶
Bases:
object
Class of mass distribution vector (MDV) information. A instance of the MdvData class could be generate from an instance of MetabolicModel class using model.generate_mdv_templete()
- add_gaussian_noise(stdev, iteration=1, method='absolute', normalize='on')¶
Addition of gausian noise to MDV data
- Parameters:
stdev (float) – noise level (standard deviation of normal distribution)
iteration (int) – number of sampleing (experimental)
method (str) –
method to add gausian noise
”relative” intentsity = (randn() * stdev + 1) * original_intensity
”absolute” (default) intentsity = randn() * stdev) + original_intensity
normalize (str) – on (default)/off sum of ratio is set to 1.0
- Returns:
True/False
- Return type:
boolean
Examples
>>> mdv.add_gaussian_noise(0.01)
- add_natural_isotope()¶
Addition of natural isotope effect to fragments
Molecular formula information in “target fragment” of model file is used.
- Parameters:
required. (Not) –
- Returns:
Nothing.
Example
>>> mdv.add_natural_isotope()
- check(output='debug')¶
Function to check missing values in the MDV data before MDV comparison
- Parameters:
output (str) – Show details when “debug” mode
- Returns:
True/False
- Return type:
boolean
Examples
>>> mdv.check()
- comparison_with_another_mdv(mdv)¶
Comparison of two MDV data
- Parameters:
mdv – instanse of another MdvData class for comparison
- Returns:
Nothing
Examples
>>> mdv.compare_mdv(mdv_fitted)
- correct_natural_isotope(mode='normal')¶
Subtraction of natural isotope effect from fragments
Molecular formula information of “target fragment” in the model definition file is used.
- Parameters:
mode (str) – “normal”(defalut) and “correction” (forced removal of negatives value and sum(ratio) is set at 1.0)
- Returns:
Nothing
Examples
>>> mdv.correct_natural_isotope()
- generate_observed_mdv()¶
Generator of MDV related inforamtion in list format.
This function is onely used in model.set_experiment()
- Parameters:
required. (Nor) –
- Returns:
id_array (array) List of ids
ratio_array (array) List of ratio data
std_array (array) List of std data
use_array (array) List of use data
observed_fragments (array) List of observed fragments
data (array) List of other information
Examples
>>> id_array, ratio_array, std_array, use_array, observed_fragments, data = mdv.generate_observed_mdv()
- get_data(fragment, number)¶
Getter of [number] th isotopomer in the MDV of [fragment].
- Parameters:
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
- Returns:
ratio (float) Relative abundance of mass spectra data
std (float) Standard deviation of measured MDV data
use (str) Data to be used for MDV comparison. ‘use’ or ‘no’
Examples
>>> ratio, std, use = mdv.get_data('Leu85', 2)
- get_fragment_mdv(fragment)¶
Getter of MDV data of [fragment].
- Parameters:
fragment (str) – name of target_fragment
- Returns:
mdv_pattern (array) Array of MDV of fragment
Examples
>>> mdv_pattern = mdv.get_fragment_mdv('Leu85')
- get_fragments_for_mdv_calculation()¶
Getter of fragment list whose MDVs are calculated in calmdv fucntion
- Parameters:
required. (Not) –
- Returns:
List of fragment names calculated in calmdv fucntion
- Return type:
array
Examples
>>> list_of_fragment = mdv.get_fragments_for_mdv_calculation()
- get_number_of_measurement()¶
Getter of number of measurement of the MDV data set.
- Parameters:
required. (Nor) –
- Returns:
Number of measurement
- Return type:
int
Examples
>>> number_of_measurement = mdv.get_number_of_measurement()
- has_data(fragment, number)¶
Checker whether the MdvData instance has data of [number] th isotopomer in the MDV of [fragment].
- Parameters:
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
- Returns:
TrueFalse
- Return type:
Boolean
Examples
>>> if mdv.has_data('Leu85', 2): brabra
- load(filename, format='text', output='normal')¶
Method to load MDV data from text/csv file
- Parameters:
filename (str) –
filename of MDV data with the format:
Name Spectrum Select MDV Std Ala57 m0 1 0.566990778 0.000774686 Ala57 m1 1 0.148623963 0.000774686 Ala57 m2 1 0.039467636 0.000774686 Ala57 m3 1 0.244917622 0.000774686
format (str) –
file format:
’csv’ : CSV.
’text’ (defalut) : tab-deliminated text.
output – “normal” (defalut) or “debug”
- Returns:
True/False
- Return type:
Boolean
Examples
>>> mdv.load('filename', format = 'text',output = "normal")
- save(filename, format='text')¶
Method to save MDV data in text/csv file
- Parameters:
filename (str) –
filename of MDV data with the format:
Name Spectrum Select MDV Std Ala57 m0 1 0.566990778 0.000774686 Ala57 m1 1 0.148623963 0.000774686 Ala57 m2 1 0.039467636 0.000774686 Ala57 m3 1 0.244917622 0.000774686
format (str) –
file format:
’csv’ : CSV.
’text’ (defalut) : tab-deliminated text.
- Returns:
True/False
- Return type:
Boolean
Examples
>>> mdv.save('filename', format = "csv")
- set_data(fragment, number, ratio, std, use)¶
Setter of [number] th isotopomer in the MDV of [fragment].
- Parameters:
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
ratio (float) – Relative abundance of mass spectra data
std (float) – Standard deviation of measured MDV data
use (str) – Data to be used for MDV comparison. ‘use’ or ‘no’
- Returns:
TrueFalse
- Return type:
Boolean
Examples
>>> mdv.set_data('Leu85', 2, 0.763, 0.0054, 'use')
- set_mdv_for_comparison(fragment, number)¶
Setter of mass data to be used for MDV comparison
- Parameters:
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
- Returns:
Nothing.
Examples
>>> mdv.set_mdv_for_comparison(fragment, number)
- set_mdvs_for_comparison(threshold_ratio, threshold_std=1.0)¶
Setters of multiple data to be used for MDV comparison
- Parameters:
threshold_ratios (float) – MDV data whose ratio (0.0-1.0) is above threshold_ratios are used for MDV comparison
threshold_std (float) – MDV data whose std is below threshold_std are used for MDV comparison
- Returns:
True/False
- Return type:
boolean
Examples
>>> mdv.set_mdvs_for_comparison(0.05, 0.005)
- set_observed_fragments(fragments)¶
Setter of observed fragment.
When new MdvData instance is generated ‘observed_fragment’ == ‘target_fragment’.
- Parameters:
fragments (Array) – List of fragment names calculated in calmdv fucntion
- Returns:
Nothing.
Examples
>>> mdv.set_observed_fragments(['Phe_85', 'Ala_57'])
- set_std(value, method='absolute')¶
Setter of standard deviation level
- Parameters:
value (float) – level of standard deviation.
method (str) –
Method to calculate stdev levels.
’relative’: Levels are set by stdev = [value] * signal intensity.
’absolute’ (detault) : Levels are set by stdev = [value].
- Returns:
Nothing.
Examples
>>> mdv.set_std(0.05, method = 'absolute')
- set_unused_mdv_for_comparison(fragment, number)¶
Setter of mass data to be ignored for MDV comparison
- Parameters:
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
- Returns:
Nothing.
Examples
>>> mdv.set_unused_mdv_for_comparison(fragment, number)
- class mfapy.mdv.MdvTimeCourseData¶
Bases:
object
- add_gaussian_noise(stdev, iteration, method='absolute', normalize='on')¶
Addition of gausian noise to MDV data
- Parameters:
stdev (float) – noise level (standard deviation of normal distribution)
iteration (int) – number of sampleing (experimental)
method (str) –
method to add gausian noise
”relative” intentsity = (randn() * stdev + 1) * original_intensity
”absolute” (default) intentsity = randn() * stdev) + original_intensity
normalize (str) – on (default)/off sum of ratio is set to 1.0
- Returns:
True/False
- Return type:
boolean
Examples
>>> mdv.add_gaussian_noise(0.01)
- add_time_point(time, mdv_instance)¶
Addition of new time point
- Parameters:
time (float) – time point
mdv_instance (MdvData) – MdvData instance
- Returns:
nothing
Examples
>>> mdv.add_time_point(3.0, mdv3)
- generate_observed_mdv()¶
Generator of MDV related inforamtion in list format.
This function is onely used in model.set_experiment()
- Parameters:
required. (Nor) –
- Returns:
id_array (array) Array of list of ids at each time point
ratio_array (array) Array of list of ratio data at each time point
std_array (array) Array of list of std data at each time point
use_array (array) Array of list of use data at each time point
observed_fragments (array) Array of list of observed fragments at each time point
data (array) Array of list of other information at each time point
Examples
>>> id_array, ratio_array, std_array, use_array, observed_fragments, data = mdv.generate_observed_mdv()
- get_data(time, fragment, number)¶
Getter of [number] th isotopomer in the MDV of [fragment] at [time].
- Parameters:
time (float) – time point
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
- Returns:
ratio (float) Relative abundance of mass spectra data
std (float) Standard deviation of measured MDV data
use (str) Data to be used for MDV comparison. ‘use’ or ‘no’
Examples
>>> ratio, std, use = mdv.get_data('Leu85', 2)
- get_number_of_measurement()¶
Getter of number of measurement of the MDV data set.
- Parameters:
required. (Nor) –
- Returns:
Number of measurement
- Return type:
int
Examples
>>> number_of_measurement = mdv.get_number_of_measurement()
- History
Revised at 27/5/2018
- get_timecourse_data(fragment, number)¶
Getter of time course MDV data
Examples
>>> tc = mdv.get_timecourse_data('Leu85', 2)
- Parameters:
fragment (name of target_fragment) –
number (Number of isotope) –
Reterns –
---------- –
tc (Array of time course of Mass isotopomer) –
- get_timepoints()¶
Getter of Time points
- Parameters:
Reterns –
---------- –
timepoints (list of time points) –
Examples
>>> timepoints = mdv.get_timepoints() >>> timepoints >>> [0., 1., 5,]
- has_data(time, fragment, number)¶
Checker whether the MdvData instance has data of [number] th isotopomer in the MDV of [fragment] at [time].
- Parameters:
time (float) – time point
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
- Returns:
TrueFalse
- Return type:
Boolean
Examples
>>> if mdv.has_data(0, 'Leu85', 2): brabra
- save(filename, format='text')¶
Method to save MDV data in text/csv file
Files of each timepoint are generated with names filename+”timepoint”.txt
- Parameters:
filename (str) –
filename of MDV data with the format:
Name Spectrum Select MDV Std Ala57 m0 1 0.566990778 0.000774686 Ala57 m1 1 0.148623963 0.000774686 Ala57 m2 1 0.039467636 0.000774686 Ala57 m3 1 0.244917622 0.000774686
format (str) –
file format:
’csv’ : CSV.
’text’ (defalut) : tab-deliminated text.
- Returns:
True/False
- Return type:
Boolean
Examples
>>> mdv.save('filename', format = "csv")
- set_data(time, fragment, number, ratio, std, use)¶
Setter of [number] th isotopomer in the MDV of [fragment] at [time].
- Parameters:
time (float) – time point
fragment (str) – name of target_fragment
number (int) – [number] th isotopomer
ratio (float) – Relative abundance of mass spectra data
std (float) – Standard deviation of measured MDV data
use (str) – Data to be used for MDV comparison. ‘use’ or ‘no’
- Returns:
TrueFalse
- Return type:
Boolean
Examples
>>> mdv.set_data(0, 'Leu85', 2, 0.763, 0.0054, 'use')
- set_mdvs_for_comparison(threshold_ratio, threshold_std=1.0)¶
Setters of multiple data to be used for MDV comparison
- Parameters:
threshold_ratios (float) – MDV data whose ratio (0.0-1.0) is above threshold_ratios are used for MDV comparison
threshold_std (float) – MDV data whose std is below threshold_std are used for MDV comparison
- Returns:
True/False
- Return type:
boolean
Examples
>>> mdv.set_mdvs_for_comparison(0.05, 0.005)
- set_observed_fragments(fragments)¶
Setter of observed fragment.
When new MdvData instance is generated ‘observed_fragment’ == ‘target_fragment’.
- Parameters:
fragments (Array) – List of fragment names calculated in calmdv fucntion
- Returns:
Nothing.
Examples
>>> mdv.set_observed_fragments(['Phe_85', 'Ala_57'])
- set_std(value, method='absolute')¶
Set standard deviation levels from mass sepctral intensity
- Parameters:
method (Method to calculate stdev levels:) – ‘relative’ (detault): Levels are set by stdev = [value] * signal intensity ‘abusolute’: Levels are set by stdev = [value]
Examples
>>> mdv.set_std(0.05, method = 'absolute')
- mfapy.mdv.binomial(n, k)¶
- mfapy.mdv.count_atom_number(formula)¶
- mfapy.mdv.transition_matrix(formula)¶
mfapy.metabolicmodel module¶
metabolicmodel.py:MetabolicModel class in mfapy
This is a core module of mfapy.
The module includes:
MetabolicModel class
Todo
Rewriting model construction part.
Check callbacklevel.
- class mfapy.metabolicmodel.MetabolicModel(reactions, reversible, metabolites, target_fragments, mode='normal')¶
Bases:
object
Class of MetabolicModel
An instances of this class has fuctions to analyze mass spectrometry based 13C MFA and INST MFA data.
An metabolicModel instance generates MDVData and CarbonSouece classes
- Returns:
instance of metabolic model
- Return type:
instance
Examples
>>> reactions, reversible, metabolites, target_fragments = mfapy.mfapyio.load_metabolic_model("example_1_toymodel_model.txt", format = "text") >>> model = mfapy.metabolicmodel.MetabolicModel(reactions, reversible, metabolites, target_fragments)
Attributes (incomplete):
self.reactions = reactions self.metabolites = metabolites self.reversible = reversible self.target_fragments = target_fragments self.symmetry = {} self.carbon_source = {} self.experiments = {}
- add_perturbation(flux)¶
New state dict is generarate from a given state dict by an addition of small random perturbation to an independent flux
- Parameters:
state_dict (dict) – Dictionary of metabolic state data.
- Returns:
check (boolean)
perturbed_state (dict)Dictionary of metabolic state data.
flux_vector (array) List of metabolic flux, metabolite concentration data.
reaction_name_list (array) List of ids
Examples
>>> check, perturbed_state, flux_vector, reaction_name_list = model.add_perturbation(state_dict)
- History
Newly created at 21/1/2021
- calc_idv(flux, carbon_sources)¶
Calc IDV from given flux and carbon_sources.
Isotopomer distribution vector (IDV) was used to calc MDV at time = 0 in the INST mode.
Examples
>>> idv = model.calc_idv(flux, carbon_sources)
- Parameters:
flux (dict) – Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
carbon_sources (instance) – Instance of carbon source object
- Returns:
array of idv
- Return type:
array
See also
set_experiment() generate_mdv()
- calc_rss(fluxes, mode='flux')¶
Getter to calc residual sum of square (RSS) level between estimated MDVs of ‘flux’ and measured MDVs in ‘experiment(s)’
- Parameters:
flux (dict) – Dictionary of metabolic flux distribution or independent flux vector
mode (str) –
“independent”:Calculation of RSS from independent flux vector
”flux”:Calculation of RSS from state dictionary
- Returns:
Residual sum of square (RSS)
- Return type:
float
Examples
>>> rss = model.calc_rss(flux_opt) >>> print rss 7564.123
- History:
Call of calmdv funcion is removed.12/9/2014
- check_independents(independents)¶
Checker method to test whether a given vector of independent fluxes is able to produce feasible status (flux) vector within lower and upper boundaries.
- Parameters:
independents (array) – list of independent fluxes generated by self.get_independents
- Returns:
check (Boolean) True means that a given vector of independent fluxes is able to produce feasible status (flux) vector within lower and upper boundaries.
flux_vector (array) List of metabolic flux, metabolite concentration data.
reaction_name (array) List of ids
Examples
>>> check, flux_vector, reaction_name = self.check_independents(independents)
History:
Newly created at 21/1/2021
- clear_experiment()¶
Clear all ‘experiment set(s)’ in the metabolic model.
- Parameters:
required (Not) –
Examples
>>> model.clear_experiment()
See also
model.set_experiment()
- History:
Newly developed at 4/9/2014
- fitting_flux(method='SLSQP', flux=[], output='result')¶
Method for fitting a metabolic model to multiple experiments
Used parameters in self.configuration().
‘iteration_max’: Maximal number of interation in the optimizers. Example: model.set_configuration(iteration_max = 1000)
‘number_of_repeat’: Number of repeated execution by ‘deep’ functions. Example: model.set_configuration(number_of_repeat = 2)
‘ncpus’: Number of CPUs used in parallel processing. Example: model.set_configuration(ncpus = 16)
- Parameters:
method (str) –
Algorism used for optimization.
’deep’: Repeated execution of SLSQP & LN_PRAXIS
’SLSQP’: Sequential Least SQuares Programming algorithm implemented in scipy
’COBYLA’: Constrained Optimization By Linear Approximation implemented by scipy
”LN_BOBYQA”:The BOBYQA algorithm for bound constrained optimization without derivatives by nlopt
”LN_SBPLX”: Sbplx (based on Subplex, (a variant of Nelder-Mead that uses Nelder-Mead on a sequence of subspaces) by nlopt
”LN_NEWUOA”: The NEWUOA software for unconstrained optimization without derivatives by nlopt
”LN_COBYLA”: Constrained Optimization By Linear Approximation implemented by nlopt
”LN_PRAXIS”: Gradient-free local optimization via the “principal-axis method” implemented by nlopt
”LN_NELDERMEAD”: the original Nelder-Mead simplex algorithm by nlopt
”GN_ESCH”: Evolutionary Algorithm for global optimization by nlopt
”GN_CRS2_LM”: Controlled random searchimplemented for global optimizatoin by nlopt
”GN_DIRECT_L”: DIviding RECTangles algorithm for global optimization by nlopt
”GN_IRES”: Improved Stochastic Ranking Evolution Strategy by nlopt
”GN_AGS”: AGS NLP solver by nlopt
flux (dict) –
Initial flux distribution generated by self.generate_initial_states():
When flux is a dict of one state, single fitting trial is executed.
When flux is a array of multiple dists of states, muptiple fitting trial is exected by using the parallel processing.
output (str) –
Output method.
’result’ (default): Fitting problem is solved inside of the function.
’for_parallel’: Fitting problems is NOT solved inside of the function. In this mode, a tupple of parameters for fit_r_mdv_pyopt() functions are generated for parallel processinng.
- Returns:
stop condition
[number_of_initial_states == 1 ] Str of stop condition.
[number_of_initial_states > 1 ] Array of str of stop condition.
flux_opt (array or dict): Initial flux distribution generated by self.generate_initial_states()
[number_of_initial_states == 1 ] dict of optimzed metabolic state.
[number_of_initial_states > 1 ] Array of dict of optimzed metabolic state. Sorted by ascending order of their RSS levels
RSS_bestfit (array or float) : Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
[number_of_initial_states == 1 ] Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
[number_of_initial_states > 1 ] Array of dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
- Return type:
state (str or array)
Examples:
# # Single fitting trial # >>> state, flux_initial = model.generate_initial_states(10, 1) >>> state, RSS_bestfit, flux_opt = model.fitting_flux(method = 'deep', flux = flux_initial) >>> results= [("deep", flux_opt)] >>> model.show_results(results, pool_size = "off") # # Multipe fitting trials by using parallel python # >>> state, flux_initial = model.generate_initial_states(50, 4) >>> state, RSS_bestfit, flux_opt_slsqp = model.fitting_flux(method = "SLSQP", flux = flux_initial) >>> results= [("deep", flux_opt[0])] >>> model.show_results(results, pool_size = "off") # # Generating parameters for parallel processing # >>> parameters = model.fitting_flux(method = "SLSQP", flux = flux_initial, output = "for_parallel")
See also
generate_initial_states()
- generate_calmdv(mode='normal')¶
Constructer of the python function to calculate MDV of target fragments from a metabolic flux distribution.
This function is performed in reconstruct()
- Parameters:
mode (str) – optimize EMU network reduction (experimental)
- Returns:
string, cython_string (str) Texts of python code describing calmdv and diffmdv functions for python 3 and Cython (Experimental)
Examples
>>> calmdv_text, ccalmdv_text = self.generate_calmdv(mode)
- generate_carbon_source_template()¶
Generator of a templete of CarbonSourse instance. Labeling information of carbon sources will be added to the instance.
Carbon source compounds are derived from the metabolic model (//Metabolites)
- Parameters:
required. (Not) –
Examples
>>> carbon_source_idv = model.generate_carbon_source_templete()
See also
CarbonSource.py
- generate_carbon_source_templete()¶
- generate_ci_template(targets='normal')¶
Generator of a template dictionary to keep confidence interval search results.
- The templatedisctionary generated by this method is used for:
Setting reactions, metabolites, and reversible reactions for CI searching in self.search_ci method.
Store all record of CI searching by self.search_ci method.
- Parameters:
targets (str) –
‘normal’ Free reversible and inreversible reactions
’independent’ Independent reactions
’all’ all reactions
’without_reversible_reactions’ inreversible reactions
- Returns:
template dictionary
- Return type:
dict
Examples
>>> ci = model.generate_ci_templete(targets = [("reaction","r29_g6pdh"),("reaction","r27_pc"),("reaction","r28_mae")]) >>> ci = model.generate_ci_templete(targets = 'normal')
See also
search_ci.
- generate_ci_templete(targets='normal')¶
- generate_initial_states(iterations=100, initial_states=1, method='normal', template=[])¶
Generator of random initial metabolic states for model fitting
Initial metabolic states are randomly generated for “iterations” times from which better states with lower RSS (number_of_initial_states) were selected.
Used parameters in self.configuration().
‘initial_search_iteration_max’: Maximal number of interation in the optimizers. Example: model.set_configuration(initial_search_iteration_max = 1000)
‘ncpus’: Number of CPUs used in parallel processing. Example: model.set_configuration(ncpus = 16)
- Parameters:
method (str) – computational method * ‘normal’ : Generation of flux distributions by single thread. * ‘parallel’ : Parallel generation of flux distributions using parallel processing.
iterations (int) – Number of trials to generate initial metabolic state. A trial sometimes failed due to a bad optimization state.
initial_states (int) – Number of initial metabolic state(s) with lower RSS to be generated.
template (array) – Array of templete metabolic state(s) dictionaries. When template is available, initial metabolic state(s) similiar to the templete metabolic state are generated. This function is used in the grid search fucntion.
- Returns:
Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
[number_of_initial_states > 1 ] Array of dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
[number_of_initial_states == 1 ] Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
- state (str or array)stop condition
[number_of_initial_states > 1 ] Str of stop condition.
[number_of_initial_states == 1 ] Array of Str of stop condition.
- Return type:
flux (array or dict)
Examples
>>> state, dict = model.generate_initial_states(50,1) >>> state, dict = model.generate_initial_states(50, 2, template = flux)
See also
Fitting
- generate_mdv(flux, carbon_sources, timepoint=[], startidv=[])¶
Generator of a MdvData instance including MDV data generated from given flux and carbon sources.
- Parameters:
flux (dict) – Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
carbon_sources (instance) – Instance of CarbonSource
timepoint (array) – For INST-13C MFA. An array of time points of measurement in MDV
startidv (array) – array of idv as starting isotope distribution for INST
- Returns:
instance* MdvData instance
Examples
>>> mdv = model.generate_mdv(flux, mdv_carbon_sources)
See also
generate_carbonsource_MDV
- History:
140912 calc_MDV_from_flux instead of calmdv is used.
- generate_state(template=[])¶
Generator of a random metabolic flux distribution.
Used parameters in self.configuration().
‘initial_search_iteration_max’: Maximal number of interation in the optimizers. Example: model.set_configuration(initial_search_iteration_max = 1000)
- Parameters:
template (dict) – Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels. When template is available, metabolic state most similar to the template is generated. The function is used in the grid search.
- Returns:
flux (dict) Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
independent_flux (array) List of independent flux
Examples
>>> flux, independent_flux = model.generate_intial_flux_distribution()
- generate_state_dict(tmp_r)¶
Generator of a state dict from a given vector.
- Parameters:
tmp_r (array) – tmp_r = numpy.dot(matrixinv, Rm_intial)
- Returns:
Dictionary of metabolic state inclucing metabolic flux and metabolite pool size levels.
- Return type:
dict
Examples
>>> state_dict = model.generate_state_dict(tmp_r)
- get_constrain(group, id)¶
Getter of constrains of metabolic reaction.
- Parameters:
group (str) – group of parameters (“reaction”, “reversible”, “metabolite”)
id (str) – id in model definition
- Returns:
type (str) type of constrain “fitting”, “fixed”, “free”, “pseudo”
value (float) Level of metabolic flux
stdev (float) Standard deviation of metabolic flux level used for ‘fitting’ type reactions.
Examples
>>> type, value, stdev = model.get_constrain("reaction","pgi")
See also
set_constrain
- get_degree_of_freedom()¶
Getter of degree of freedom of the model
Examples
>>> degree_of_freedom = model.get_degree_of_freedom()
- Returns:
degree of freedom.
- Return type:
int
- get_independents(flux)¶
Getter of independent status (flux) vector, its lower and upper boundaries and a list of independent flux ids of given state (flux) dict.
Examples
>>> independent_vector, lb, ub, independent_list = model.get_independents(state_dict)
- Parameters:
flux (dict) – Dictionary of metabolic state data.
- Returns:
independent_vector (array) list of independent flux values
lb (array) list of independent flux values
ub (array) list of independent flux values
independent_list (array) list of independent fluxes, metabolites, and reversible reactions.
- History:
Newly created at 21/1/2021
- get_number_of_independent_measurements()¶
Getter of number of independent measurements of the model
- Returns:
number of independent measurements.
- Return type:
int
- History:
29/9/2014 Modified to consider “fitting” flux information.
Example
>>> number_of_independent_measurements = model.get_number_of_independent_measurements()
- get_thres_confidence_interval(flux, alpha=0.05, dist='F_dist')¶
Getter to obtain a threshold level of RSS for searching confidence interval
- Parameters:
flux (dict) – Dictionary of best fitted flux distribution.
alpha (float) – confidence level. For 95% confience interval, please set alpha = 0.05
dist (str) – Probability distribution for test * ‘F-dist’ (default): Using F-distribution * ‘chai-dist’:Using Chai-square distribution
- Returns:
thres (float) a threshold level of RSS for searching confidence interval
number_of_measurements (int) number of measurements
degree_of_freedom (int) degree of freedom
Examples
>>> thres, number_of_reactions, degree_of_freedom = model.get_thres_confidence_interval(RSS, alpha = 0.05, dist = 'F_dist')
- goodness_of_fit(flux, alpha=0.05)¶
Getter to calculate goodness-of-fit of a given flux distribution
- Parameters:
flux (dict) – Dictionary of metabolic flux distribution
alpha (float) – Confidence level. Default is 0.05.
- Returns:
pvalue (float) p-value determined by chi-squere test.
rss_thres (float) Threshold RSS at the alpha levels.
Examples
>>> pvalue, rss_thres= model.goodness_of_fit(flux_opt_fitted, alpha = 0.05)
- load_mdv_data(filename, format='text', output='normal')¶
Load MDV data from the text file. This function generate new instance of mfapy.mdv class from a instance of mfapy.model.
- Parameters:
filename (str) –
filename of MDV data with following format:
Name Spectrum Select MDV Std Ala57 m0 1 0.566990778 0.000774686 Ala57 m1 1 0.148623963 0.000774686 Ala57 m2 1 0.039467636 0.000774686 Ala57 m3 1 0.244917622 0.000774686
format (str) – “text” (defalut) or “csv”
output (str) – “normal” (defalut) or “debug”
- Returns:
True/False
- Return type:
Boolean
Examples
>>> mdv = load_mdv_data('filename')
- load_states(filename, format='text')¶
Load a text/csv file with ‘reacton type’ information to generate new states dict.
- Parameters:
filename (str) –
filename of flux data with following format:
type Id type flux_calue flux_std lb ub reaction v1 fixed 100 1 0.001 100 reaction v2 free 100 1 0.001 100 reaction v3 free 50 1 0.001 100 metabolites AcCoA fixed 1 1 0.001 100 metabolites OAC fixed 1 1 0.001 100 metabolites OACs fixed 1 1 0.001 100 metabolites OACx fixed 1 1 0.001 100 reversible FUM free 1 1 0.001 100
format – ‘csv’ CSV or ‘text’tab-deliminated text.
- Returns:
Dictionary of states dict.
- Return type:
dict
Examples
>>> states = model.load_states("test.csv", format = 'csv')
- modelcheck()¶
Method to check medel definition data.
- Parameters:
nothing –
- Returns:
Nothing.
Examples
>>> model.datacheck()
See also
reconstruct
- History:
Newly developed for version 059
- posterior_distribution(state_dict, number_of_steps=1000000)¶
Generator of posterior distribution by the Metropolis-Hastings algorism
- Parameters:
state_dict (dist) – Dictionary of initial metabolic state
number_of_steps (int) – Length of Markov-chain
- Returns
array: 2d array. Markov-chain of metabolic state vector
Examples
>>> record = model.posterior_distribution(perturbed_state, number_of_steps = numberofsteps)
History:
Newly created at 21/1/2021
- reconstruct(mode='no')¶
Method to reconstruct the metabolic model to generate new stoichiometry matrix and a function to calculate MDV data.
- Parameters:
mode (str) – experimental
- Returns:
calmdv_text (dict) The dictionary has two keys. “calmdv” and “diffmdv” are instance of functions for 13C-MFA and INST-13C-MFA to generate MDV from a given metabolic flux, respectively.
Examples
>>> model.reconstruct()
See also
update()
- save_states(dict, filename, format='text')¶
Save state dict to a text/csv file with ‘type’ information.
- Parameters:
filename (str) – filename of MDV data with the format.
dict (dict) – dictinary of metabolic state (flux)
format – ‘csv’ CSV or ‘text’tab-deliminated text.
- Returns:
True/False
- Return type:
Boolean
Examples
>>> model.save_states(flux_dict, "test.csv", format = 'csv')
- search_ci(ci, flux, method='grid', alpha=0.05, dist='F_dist', outputthres=2.0)¶
Method to estimate confidence intervals using information in a ‘ci’ dictionary.
Used parameters in self.configuration().
‘grid_search_iterations’ (default = 1): MNumber of trials for model fitting at each point in grid search. Example: model.set_configuration(grid_search_iterations = 1)
‘number_of_repeat’(default = 50): Number of initial metabolic flux disributions generated for each trial in grid search.Best initial metabolic flux disribution with smallest RSS is employed.
‘number_of_repeat’: Number of repeated execution by ‘deep’ functions. Example: model.set_configuration(number_of_repeat = 2)
‘ncpus’: Number of CPUs used in parallel processing. Example: model.set_configuration(ncpus = 16)
- Argss:
ci (dict): dictionary for storeing confidence interval results. ci is generated by self.generate_ci_templete()
flux (dict): dict of metabolic state (best fitted)
method (str): ‘grid’ (Grid search method) is available
alpha (float): (default 0.05) Confidence interval level. 0.05 means a 95% confidence interval level
- dist (str):
‘F-dist’ (default): Using F-distribution
‘chai-dist’:Using Chai-square distribution
outputthres (float): The grid search results larger than outputthres * threshold rss were removed from output
- Returns:
dictionary of confidence interval results.
- Return type:
dict
Examples
>>> ci = model.generate_ci_templete(targets = 'independent') >>> ci = model.search_confidence_interval_parallel(ci, flux, method = 'grid')
- set_boundary(group, id, lb, ub)¶
Setter of a lower and upper boundaries of metablic states
- Parameters:
group (str) – metabolite, reaction, reversible
id (str) – id of state parameter
lb (float) – lower boundary
ub (float) – upper boundary
- Returns:
True/False.
- Return type:
Booleans
Examples
>>> model.set_boundary("reaction", 'HEX1', 0.001, 100.0)
See also
set_constrain
- set_configuration(*args, **kwargs)¶
Setter of a configuration of model parameters.
- Parameters:
*args (list) – “parameter_name = number”
- Returns:
Nothing.
Parameters:
'callbacklevel': default = 0, 'callbacklevel' determines a frequency level of callbacks from metabolic model. # 0: No call back # 1: Error report # 2: Simple pregress report for grid search # 3: Detailed pregress report for grid search # 4: Simple pregress report for flux fitting # 5: Detailed pregress report for flux fitting # 7: Detailed report for model constrution # 8: Debug 'iteration_max': default = 1000, Maximal number of interations (steps) in each fitting task. 'default_reaction_lb': default = 0.001, 'default_reaction_ub': default = 5000.0, 'default_metabolite_lb': default = 0.001, 'default_metabolite_ub': default = 500.0, 'default_reversible_lb': default = -300, 'default_reversible_ub': default = 300.0, Default lower and upper boundaries used only when there is no data in the model definition file. 'grid_search_iterations': default = 1, Number of trials for model fitting at each point in grid search. 'initial_search_repeats_in_grid_search': default = 50, Number of initial metabolic flux disributions generated for each trial in grid search. Best initial metabolic flux disribution with smallest RSS is employed. 'initial_search_iteration_max': default = 1000, Maximal number of interations (steps) allowed in each task to find feasible initial metabolic flux distribution. 'add_naturalisotope_in_calmdv':default = "no" For INST-13C MFA. Natural isotope is added to MDVs of all intermediate at time = 0. 'number_of_repeat':default = 3, For "deep" function of self.fitting flux. "Global" => "Local" optimization methods are repeated for n times. 'ncpus':default = 3, Number of CPUs used in parallel processing 'odesolver': default = "scipy" # or "sundials" Only "scipy" is avaliable due to memory leak problem of sundials
Examples
>>> model.set_configuration(callbacklevel = 1) >>> model.set_configuration(iteration_max = 1000)
- set_constrain(group, id, type, value=0.1, stdev=1.0)¶
Setter of a metabolic constrains.
Please perform update() after model modification.
- Parameters:
group (str) – group of parameters (“reaction”, “reversible”, “metabolite”)
id (str) – id in model definition
type (str) – type of constrain “fitting”, “fixed”, “free”, “pseudo”: * “fitting” : For the calclation of RSS by using “value” and stdev * “fixed” : Fixed to “value”. Not used for the calculation of RSS. * “free” : Free between lower and upper boundary. Not used for the calculation of RSS. * “pseudo” : Pseudo reaction is a special “free” reaction to consider G-index. The reaction is ignored in a stoichiometry matrix. However, the reaction is considered in the MDV calculation.
value (float) – Level of metabolic flux
stdev (float) – Standard deviation of metabolic flux level used for ‘fitting’ type reactions.
- Returns:
True/False.
- Return type:
Booleans
Examples
>>> model.set_fixed_reaction("reaction", 'pgi', "fixed", 100.0,1)
See also
get_constrain
- set_constraints_from_state_dict(dict)¶
Reaction types, value, stdev, lb, and ub are set from the state dict data at once
- Parameters:
dict (dict) – Dictionary of flux. Data in ‘value’, ‘stdev’, ‘lb’, ‘ub’ and ‘type’ fields are used.
Examples
>>> model.set_constraints_from_state_dict(flux_dict)
- set_experiment(name, mdv, carbon_sources, startidv=[])¶
Setter of an ‘experiment’ to metabolic model.
Here an ‘experiment’ indicated a set of a carbon source labeling pattern and a measured MDV data.
Parallel labeling experiment can be performed by setting multiple sets of ‘experiment’.
- Parameters:
name (str) – name of the experiment (unique)
mdv (instance) – MDVdata instance including measured MDV data of target fragments.
carbon_sources (instance) – Instance of carbon source object
startidv (array) – array of idv as starting isotope distribution for INST
Examples
>>> model.set_experiment('ex1', mdv1, cs1) >>> model.set_experiment('ex2', mdv2, cs2, startidv)
- show_flux_balance(results, metabolite, filename='', format='csv')¶
List of metabolic flux levels of reactions related to given metabolite.
- Parameters:
input – Tuple of (“name”, dic of metabolic statet) [(‘name1’, state1),(‘name2’, state2),(‘name3’, state3),]
metabolite – ID of metabolite of interest
filename – Results are saved when file name is not “”.
format – “csv” or “text”
Examples
>>> model.show_metabolite_balance(results, "Fum") >>> model.show_metabolite_balance(results, "Fum", filename = "metabolite_balance.csv", format = "csv")
- show_results(input, flux='on', rss='on', mdv='on', pool_size='off', checkrss='off', filename='', format='csv')¶
Method to output metabolic state data.
- Parameters:
input (tuple) – Tuple of (“name”, dic of metabolic statet) [(‘name1’, state1),(‘name2’, state2),(‘name3’, state3),]
flux (str) – show flux data “on”/”off”
rss (check) – show rss data “on”/”off”
mdv (str) – show mdv data “on”/”off”
pool_size (str) – show metabolite data “on”/”off”
rss – show each RSS value “on”/”off”
filename (str) – Results are saved when file name is not “”/
format (str) – “csv” or “text”
Examples
>>> model.show_results(results, flux = "on", rss = "on", mdv = "on", pool_size = "off", checkrss = "on")# Output to console >>> model.show_results(results, filename = "show_results.csv", format = "csv")# Output to .csv file
- show_results_in_map(flux, mapfilename, outputfilename, formattext='.1f')¶
Method to project metabolic state data into metabolic file (.GML).
- Parameters:
flux (dict) – dict of metabolic state
mapfilename (str) – name of blank GML file
outputfilename (str) – name of uutput GML file
format (str) – format string of “format” method of string
Examples
>>> model.show_results_in_map(flux_opt1, "Example_2_cancer_map_new.gml", "Example_2_cancer_map_mapped.gml") >>> model.show_results_in_map(flux_opt1, "Example_2_cancer_map_new.gml", "Example_2_cancer_map_mapped.gml", formattext = ".2f")
- update()¶
Method to generate stoichiometry matrix.
A metabolic model have to be updated when any types (fixed, free, fitting, pseudo) of reactions, metabolites, and reversible reactions are modified.
- Parameters:
required. (Not) –
Examples
>>> model.update()
See also
reconstruct
mfapy.mfapyio module¶
mfapyio.py:I/O funtions for mfapy
This module includes functions to read model description files.
Following functions are available:
load_metabolic_model_reactions
load_metabolic_model_metabolites
load_metabolic_model_reversibles
load_metabolic_model_fragments
load_metabolic_model
Example
>>> reactions, reversible, metabolites, target_fragments = mfapy.mfapyio.load_metabolic_model("Explanation_1_13CMFA_toymodel_model.txt", format = "text")
Todo
Support other file formats
- mfapy.mfapyio.load_metabolic_model(filename, format='text', mode='normal')[source]¶
Function to load metabolic model information from a text file
CAUTION: Now this function has no error checking.
- Parameters:
filename (str) – metabolic model file
format (str) – “text” (defalut, tab deliminated) or “csv”
mode (str) – “normal” (defalut) or “debug” (to show loaded metabolic file data)
- Returns:
reactions (dict), Dictionary describing metabolite reactions
reversible (dict), Dictionary for defining reversible reactions
metabolites (dict), Dictionary including metabolite information
target_fragments (dict), Dictionary of target_fragments
Examples
>>> reactions, reversible, metabolites, target_fragments = load_metabolic_model("filename.txt') # The obtaind data (dictionaries) are directly used for generation of new Metabolic Model object >>> model = MetabolicModel(reactions, reversible, metabolites, target_fragments)
- mfapy.mfapyio.load_metabolic_model_fragments(filename, format='text', mode='normal')[source]¶
Function to load mass fragment information from a metabolic model file with following format.
Examples:
//Target_fragments Glue gcms Glu_1:2:3:4:5 use C5H10N2O3 Gluee gcms Glu_1:2:3+Glu_4:5 use C5H10N2O3
- Parameters:
filename (str) – metabolic model file
format (str) – “text” (defalut, tab deliminated) or “csv”
mode (str) – “normal” (defalut) or “debug” (to show loaded metabolic file data)
- Returns:
Dictionary of target_fragments
- Return type:
dict
Examples
>>> target_fragments = load_metabolic_model_fragments('filename.txt') >>> print target_fragments {'Thr302': {'atommap': 'Thr_12', 'use': 'no', 'type': 'gcms', 'order': 35, 'number': 3}, ...}
- mfapy.mfapyio.load_metabolic_model_metabolites(filename, format='text', mode='normal')[source]¶
Function to load Metabolite information from a file with following format.
Examples:
//Metabolites CO2ex 1 no no excreted (kegg:C00011) 0.0 300 AcCoA 2 no carbonsource no (kegg:C00024) 0.0 300 OAC 4 no no no (kegg:C00036) 0.0 300 Cit 6 no no no (kegg:C00158) 0.0 300 AKG 5 no no no (kegg:C00026) 0.0 300 Suc 4 symmetry no no (kegg:C00042) 0.0 300 Fum 4 symmetry no no (kegg:C00122) 0.0 300 Glu 5 no no no (kegg:C00025) 0.0 300 Asp 4 no carbonsource no (kegg:C00049) 0.0 300
- Parameters:
filename (str) – metabolic model file
format (str) – “text” (defalut, tab deliminated) or “csv”
mode (str) – “normal” (defalut) or “debug” (to show loaded metabolic file data)
- Returns:
Dictionary including metabolite information
- Return type:
dict
Examples
>>> metabolites = load_metabolic_model_metabolites('filename.txt') >>> print metabolites {'GLUEX': {'excreted': 'no', 'carbonsource': 'carbonsource', 'C_number': 5, 'symmetry': 'no', 'order': 4}, ...}
- mfapy.mfapyio.load_metabolic_model_reactions(filename, format='text', mode='normal')[source]¶
Function to load metabolic reaction information from a text or CSV file with following format.
Examples:
//Reactions v1 AcCoA + OAC --> Cit AcCoA + OAC --> Cit AB + CDEF --> FEDBAC (kegg:R00351) 0.0 300 v2 Cit --> AKG + CO2ex Cit --> AKG + CO2ex ABCDEF --> ABCDE + F (kegg:R00709) 0.0 300 v3 AKG --> Glu AKG --> Glu ABCDE --> ABCDE (kegg:R00243) 0.0 300 v4 AKG --> Suc + CO2ex AKG --> Suc + CO2ex ABCDE --> BCDE + A (kegg:R01197) 0.0 300 v5 Suc --> Fum Suc --> Fum ABCD --> ABCD (kegg:R02164) 0.0 300 v6 Fum --> OAC Fum --> OAC ABCD --> ABCD (kegg:R01082) 0.0 300 v7 OAC --> Fum OAC --> Fum ABCD --> ABCD (kegg:R01082) 0.0 300 v8 Asp --> OAC Asp --> OAC ABCD --> ABCD (kegg:R00355) 0.0 300
- Parameters:
filename (str) – metabolic model file name
format (str) – “text” (defalut, tab deliminated) or “csv”
mode (str) – “normal” (defalut) or “debug” (to show loaded metabolic file data).
- Returns:
Dictionary describing metabolite reactions
- Return type:
dict
Examples
>>> reactions = load_metabolic_model_reaction('filename.txt') >>> print reactions {'r4': {'atommap': 'AB+CDEF-->FEDBAC', 'reaction': 'ACCOA+OAA-->ICI', 'use': 'use', 'lb': 0.1, 'flux_value': 0.0, 'reversible': 'no', 'flux_var': 1.0, 'type': 'free', 'order': 3, 'ub': 300.0},...}
- mfapy.mfapyio.load_metabolic_model_reversibles(filename, format='text', mode='normal')[source]¶
Function to load definitions of reversible reactions from a metabolic model file with following format.
Examples:
//Reversible_reactions FUM v6 v7 (kegg:R01082) 0.0 300
- Parameters:
filename (str) – metabolic model file
format (str) – “text” (defalut, tab deliminated) or “csv”
mode (str) – “normal” (defalut) or “debug” (to show loaded metabolic file data)
- Returns:
Dictionary for defining reversible reactions
- Return type:
dict
Examples
>>> reversible = load_metabolic_model_reversibles('filename.txt') >>> print reversible {'MDH': {'flux_value': 0.0, 'reverse': 'r27', 'flux_var': 1.0, 'forward': 'r26', 'type': 'free', 'order': 6}, ...}
mfapy.optimize module¶
optimize.py:low level optimizer functions used in mfapy.
These functions were separated from model instance for the parallel execution.
Todo
Cleaning-up and support other optimizers
- mfapy.optimize.calc_MDV_from_flux(tmp_r, target_fragments, mdv_carbon_sources, func, timepoint=[], y0temp=[])¶
Low level function to calculate mdv vector and mdv hash from metabolic flux and carbon source MDV using calmdv.
This funcition is called from mfapy.metabolicmodel.show_results.
- Parameters:
tmp_r (array) – list of metabolic state data (tmp_r = numpy.dot(matrixinv, Rm_temp)
target_fragments (array) – list of targed mdvs for MDV calclation, model.target_fragments.keys()
mdv_carbon_sources (dict) – dict of mdv_carbon_sources in model.experiments[ex_id][‘mdv_carbon_sources’]
func (dict) – Dict of functions for MDV calclation in model.func
timepoint (array) – For INST mode only. timepoints for MDV comparison in model.experiments[ex_id][‘timepoint’] When the length of timepoint array >= 1, INST mode is used.
y0temp (dict) – Start IDV state for INST mode
- Returns:
mdv (array) list of MDV data
mdv_hash (dict) dict of MDV data
- INST-MFA mode:
mdv (array) array of mdv at each time point
mdv_hash (array) array of mdv_hash at each time point
- Return type:
13C-MFA mode
Example
>>> mdv_exp, mdv_hash = optimize.calc_MDV_from_flux(tmp_r, target_fragments_temp, mdv_carbon_sources_temp, self.func)
See also
mfapy.metabolicmodel.show_results
- mfapy.optimize.calc_MDV_residue(x, *args, **kwargs)¶
Low level function for residual sum of square calculation from mfapy.metabolicmodel.MetaboliModel.calc_rss
- Parameters:
x (array) – vector of independent flux.
*args (array) – list of parameters.
**kwargs (dict) – dic of parameters.
- Returns:
RSS + Penalty score (When out side of the lower and upper boundaries)
- Return type:
float
See also
fit_r_mdv_scipy
- mfapy.optimize.calc_MDV_residue_nlopt(x, grad, kwargs)¶
Low level function for residual sum of square calculation for model fitting using nlopt.nlopt
- Parameters:
x (array) – vector of independent flux.
*args (array) – list of parameters.
- Returns:
RSS + Penalty score (When out side of the lower and upper boundaries)
- Return type:
float
See also
fit_r_mdv_scipy
- mfapy.optimize.calc_MDV_residue_scipy(x, *args)¶
Low level function for residual sum of square calculation for model fitting using scipy.optimize.minimize
- Parameters:
x (array) – vector of independent flux.
*args (array) – list of parameters.
- Returns:
RSS + Penalty score (When out side of the lower and upper boundaries)
- Return type:
float
See also
fit_r_mdv_scipy
- mfapy.optimize.calc_protrude_nlopt(independent_flux, grad, kwargs)¶
Objective function used in initializing_Rm_fitting (nlpot)
Calc penalty score of metabolic state out side of the feasible space.
- Parameters:
independent_flux (array) – vector of independent flux
grad – not used
*args (array) – list of parameters.
- Returns:
Penalty score
- Return type:
float
See also
initializing_Rm_fitting
- mfapy.optimize.calc_protrude_scipy(independent_flux, *args)¶
Objective function used in initializing_Rm_fitting (SLSQP)
This function calculates penalty score of metabolic state out side of the feasible space.
- Parameters:
independent_flux (array) – vector of independent flux
*args (list) – list of parameters.
- Returns:
Penalty score
- Return type:
float
See also
initializing_Rm_fitting
- mfapy.optimize.fit_r_mdv_deep(configure, experiments, numbers, vectors, matrixinv, func, flux)¶
Low level function for model fitting by iterative fittings.
1st iteration: GN_CRS2_LM (global optimizer)
2n th iterations: SLSQP (local)
2n + 1 th iterations: LN_SBPLX (local)
This combination is empirically best
- Parameters:
configures (dict) – “model.configures” including various configulatoins of the model.
experiments (dict) – “model.experiments” including experiments defined in the model.
numbers (dict) – “model.numbers” including various numbers of the model.
vectors (dict) – “model.vector” including various vectors of the model.
matrixinv (2d array) – “model.matrixinv” is a inversed matrix for the flux calculation.
func (dict) – Dict of functions for MDV calclation in model.func
flux (dict) – Dictionary of initial metabolic state.
- Returns:
state (str) finishing condition
kai (float) Residual sum of square of fitted metabolic state
opt_flux (array) list of fitted metabolix state
Rm_ind_sol (array) list of fitted independent flux
Example
>>> state, kai, opt_flux, Rm_ind_sol = optimize.fit_r_deep(configure, self.experiments, numbers, vectors, self.matrixinv, self.func, flux)
See also
optimize.fit_r_nlopt optimize.fit_r_scipy
- mfapy.optimize.fit_r_mdv_nlopt(configure, experiments, numbers, vectors, matrixinv, func, flux, method='LN_PRAXIS')¶
Low level function for model fitting using nlopt.opt
- Parameters:
configures (dict) – “model.configures” including various configulatoins of the model.
experiments (dict) – “model.experiments” including experiments defined in the model.
numbers (dict) – “model.numbers” including various numbers of the model.
vectors (dict) – “model.vector” including various vectors of the model.
matrixinv (2d array) – “model.matrixinv” is a inversed matrix for the flux calculation.
func (dict) – Dict of functions for MDV calclation in model.func
flux (dict) – Dictionary of initial metabolic state.
method (str) – “LN_COBYLA”, “LN_BOBYQA”, “LN_NEWUOA”, “LN_PRAXIS”, “LN_SBPLX”, “LN_NELDERMEAD”, “GN_DIRECT_L”, “GN_CRS2_LM”,”GN_ESCH”
- Returns:
state (str) finishing condition
kai (float) Residual sum of square of fitted metabolic state
opt_flux (array) list of fitted metabolix state
Rm_ind_sol (array) list of fitted independent flux
Example
>>> state, kai, opt_flux, Rm_ind_sol = optimize.fit_r_mdv_nlopt(configure, self.experiments, numbers, vectors, self.matrixinv, self.func, flux, method = "LN_PRAXIS")
See also
calc_MDV_residue_nlopt
- mfapy.optimize.fit_r_mdv_scipy(configure, experiments, numbers, vectors, matrixinv, func, flux, method='SLSQP')¶
Low level function for model fitting using scipy.optimize.minimize
- Parameters:
configures (dict) – “model.configures” including various configulatoins of the model.
experiments (dict) – “model.experiments” including experiments defined in the model.
numbers (dict) – “model.numbers” including various numbers of the model.
vectors (dict) – “model.vector” including various vectors of the model.
matrixinv (2d array) – “model.matrixinv” is a inversed matrix for the flux calculation.
func (dict) – Dict of functions for MDV calclation in model.func
flux (dict) – Dictionary of initial metabolic state.
method (str) – “SLSQP” and “COBYLA” are available
- Returns:
state (str) finishing condition
kai (float) Residual sum of square of fitted metabolic state
opt_flux (array) list of fitted metabolix state
Rm_ind_sol (array) list of fitted independent flux
Example
>>> state, kai, opt_flux, Rm_ind_sol = optimize.fit_r_mdv_scipy(configure, self.experiments, numbers, vectors, self.matrixinv, self.func, flux, method = "SLSQP")
See also
calc_MDV_residue_scipy
- mfapy.optimize.initializing_Rm_fitting(numbers, vectors, matrixinv, template, initial_search_iteration_max, method='fitting')¶
Funcition to generate randomized initial flux dixtribution using scipy.optimize.minimize SLSQP
- Parameters:
numbers (dict) – “model.numbers” including various number related data of the model.
vectors (dict) – “model.vector” including various vector related data of the model.
matrixinv (numpy 2d array) – “model.matrixinv” is a inversed matrix of stoichiometry matrix for flux calculation.
template (dict) – Dictionary of metabolic state. When template is available, metabolic state most similar to the template is generated. The function is used in the grid search.
initial_search_iteration_max (int) – “configure[“initial_search_iteration_max”]”. Maximal number of interations (steps) allowed in each task to find feasible initial metabolic flux distribution.
method (str) – “fitting” is only available.
- Returns:
tmp_r (list) list of metabolic state data (tmp_r = numpy.dot(matrixinv, Rm_temp)
Rm_temp (list) metabolic state vector
Rm_ind (list) independent flux vector
state (str) State of finishing condition “Failed”/”Determined”
Examples
>>> tmp_r, Rm_temp, Rm_ind, state = optimize.initializing_Rm_fitting(numbers, vectors, matrixinv, template ,initial_search_iteration_max)
See also
calc_protrude_scipy