jLM.CMEPostProcessing#

Post-processing for the Well-stirred CME model

Functions

autocurve(graph[, attribute, default])

Calculates curvature values for each of the edges in the graph to make sure that multiple edges are shown properly on a graph plot.

closeLMFile(f)

Close a Lattice Microbes File

color_name_to_rgb(color[, palette])

Converts a color given in one of the supported color formats to R-G-B values.

color_name_to_rgba(color[, palette])

Converts a color given in one of the supported color formats to R-G-B-A values.

community_to_membership(merges, nodes, steps)

compare_communities(comm1, comm2[, method, ...])

Compares two community structures using various distance measures.

convex_hull(vs[, coords])

Calculates the convex hull of a given point set.

disjoint_union(graphs)

Graph disjoint union.

getAvgVarTrace(f, specie[, doublingTime])

Get the average and variance of the specie trace over time

getHistogram(f, species)

Get the histogram for a specie

getNumTrajectories(f)

Get the number of trajectories (replicates) in the file

getOccupancyKymograph(f[, species, replicate])

Compute the specie density(occupancy) among a slice of the simulation domain as a function over time for the given direciton

getPhaseSpace(f, species, replicate)

Get the nD phase space associated with the traces of species.

getSpecieTrace(f, specie[, replicate, ...])

Extract data for a particular specie for the specified replicate(s)

getSpecies(f)

Exract the species names

getTimesteps(f)

Extract the timestep times

getTrajectory(f, replicate, specie)

Get trajectory data for a specific replicate and species

get_include()

Returns the folder that contains the C API headers of the Python interface of igraph.

hsl_to_rgb(h, s, l)

Converts a color given by its HSL coordinates (hue, saturation, lightness) to RGB coordinates.

hsla_to_rgba(h, s, l[, alpha])

Converts a color given by its HSLA coordinates (hue, saturation, lightness, alpha) to RGBA coordinates.

hsv_to_rgb(h, s, v)

Converts a color given by its HSV coordinates (hue, saturation, value) to RGB coordinates.

hsva_to_rgba(h, s, v[, alpha])

Converts a color given by its HSVA coordinates (hue, saturation, value, alpha) to RGB coordinates.

intersection(graphs[, byname, keep_all_vertices])

Graph intersection.

is_bigraphical(degrees1, degrees2[, multiple])

Returns whether two sequences of integers can be the degree sequences of a bipartite graph.

is_degree_sequence(out_deg[, in_deg])

Deprecated since 0.9 in favour of L{is_graphical()}.

is_graphical(out_deg[, in_deg, loops, multiple])

Returns whether a list of degrees can be a degree sequence of some graph, with or without multiple and loop edges, depending on the allowed edge types in the remaining arguments.

is_graphical_degree_sequence(out_deg[, in_deg])

Deprecated since 0.9 in favour of L{is_graphical()}.

load(filename, *args, **kwds)

Loads a graph from the given filename.

mean(xs)

Returns the mean of an iterable.

median(xs[, sort])

Returns the median of an unsorted or sorted numeric vector.

openLMFile(filename)

Open a Lattice Microbes File for reading

percentile(xs[, p, sort])

Returns the pth percentile of an unsorted or sorted numeric vector.

plot(obj[, target, bbox])

Plots the given object to the given target.

plotAvgVar(f[, species, filename])

Plot a specific species average over time and variance

plotAvgVarFromFile(filename, species, ...)

Plot species from an output file

plotCMEReactionNetwork(sim, filename[, ...])

Plot the reaction scheme as a network

plotOccupancyKymograph(f[, species, ...])

Plot a specific specie density(occupancy) among a slice of the simulation domain as a function of a direction over time

plotPhaseSpace(f[, species, replicate, ...])

Plot the 2D or 3D phase space associated with the given species over the replicates indicated

plotTrace(f[, species, replicate, filename])

Plot a specific species trace

plotTraceFromFile(filename, species, ...)

Plot species from an output file

power_law_fit(data[, xmin, method, p_precision])

Fitting a power-law distribution to empirical data

quantile(xs[, q, sort])

Returns the qth quantile of an unsorted or sorted numeric vector.

read(filename, *args, **kwds)

Loads a graph from the given filename.

rescale(values[, out_range, in_range, ...])

Rescales a list of numbers into a given range.

rgb_to_hsl(r, g, b)

Converts a color given by its RGB coordinates to HSL coordinates (hue, saturation, lightness).

rgb_to_hsv(r, g, b)

Converts a color given by its RGB coordinates to HSV coordinates (hue, saturation, value).

rgba_to_hsla(r, g, b[, alpha])

Converts a color given by its RGBA coordinates to HSLA coordinates (hue, saturation, lightness, alpha).

rgba_to_hsva(r, g, b[, alpha])

Converts a color given by its RGBA coordinates to HSVA coordinates (hue, saturation, value, alpha).

save(graph, filename, *args, **kwds)

Saves a graph to the given file.

setLMLogConsole([level])

Set the logger to write to the console as the code is working

setLMLogFile(filename[, level])

Set up file handler to print log to file

setLMLoggerLevel(level)

Set the level of the logger for the application

set_progress_handler(handler)

Sets the handler to be called when igraph is performing a long operation. @param handler: the progress handler function. It must accept two arguments, the first is the message informing the user about what igraph is doing right now, the second is the actual progress information (a percentage).

set_random_number_generator(generator)

Sets the random number generator used by igraph. @param generator: the generator to be used. It must be a Python object with at least three attributes: C{random}, C{randint} and C{gauss}. Each of them must be callable and their signature and behaviour must be identical to C{random.random}, C{random.randint} and C{random.gauss}. Optionally, the object can provide a function named C{getrandbits} with a signature identical to C{randpm.getrandbits} that provides a given number of random bits on demand. By default, igraph uses the C{random} module for random number generation, but you can supply your alternative implementation here. If the given generator is C{None}, igraph reverts to the default PCG32 generator implemented in the C layer, which might be slightly faster than calling back to Python for random numbers, but you cannot set its seed or save its state.

set_status_handler(handler)

Sets the handler to be called when igraph tries to display a status message.

showAvgVar(f, species, **kwargs)

Show a specific species average over time and variance

showAvgVarFromFile(filename, species, **kwargs)

Show species from an output file :param filename: The name of an HDF5 output file generated by LatticeMicrobes :param species: A list of species to plot :param kwargs: Additional arguments to be passed on to matplotlib.plot.

showOccupancyKymograph(f[, species, replicate])

Show a specific specie density(occupancy) among a slice of the simulation domain as a function of a direction over time

showTrace(f, species, replicate, **kwargs)

Show a specific species trace

showTraceFromFile(filename, species, ...)

Show species trace from a particular replicate

split_join_distance(comm1, comm2[, remove_none])

Calculates the split-join distance between two community structures.

summary(obj[, stream])

Prints a summary of object o to a given stream

umap_compute_weights(graph, dist)

Compute undirected UMAP weights from directed distance graph.

union(graphs[, byname])

Graph union.

write(graph, filename, *args, **kwds)

Saves a graph to the given file.

Classes

ARPACKOptions

Class representing the parameters of the ARPACK module.

AdvancedGradientPalette(colors[, indices, n])

Advanced gradient that consists of more than two base colors.

BFSIter()

igraph BFS iterator object

BoundingBox(*args)

Class representing a bounding box (a rectangular area) that encloses some objects.

CairoGraphDrawer(context[, bbox, ...])

Class implementing the default visualisation of a graph.

ClusterColoringPalette(n)

A palette suitable for coloring vertices when plotting a clustering.

Clustering(membership[, params])

Class representing a clustering of an arbitrary ordered set.

CohesiveBlocks(graph[, blocks, cohesion, parent])

The cohesive block structure of a graph.

Configuration([filename])

Class representing igraph configuration details.

Cover(clusters[, n])

Class representing a cover of an arbitrary ordered set.

Cut(graph[, value, cut, partition, partition2])

A cut of a given graph.

DFSIter()

igraph DFS iterator object

DefaultGraphDrawer

alias of CairoGraphDrawer

Dendrogram(merges)

The hierarchical clustering (dendrogram) of some dataset.

DyadCensus([iterable])

Dyad census of a graph.

Edge

Class representing a single edge in a graph.

EdgeSeq

Class representing a sequence of edges in the graph.

FittedPowerLaw(continuous, alpha, xmin, L, D, p)

Result of fitting a power-law to a vector of samples

Flow(graph, value, flow, cut, partition)

A flow of a given graph.

GradientPalette(color1, color2[, n])

Base class for gradient palettes

Graph(*args, **kwds)

Generic graph.

GraphBase

Low-level representation of a graph.

GraphSummary(graph[, verbosity, width, ...])

Summary representation of a graph.

Histogram([bin_width, data])

Generic histogram class for real numbers

Layout([coords, dim])

Represents the layout of a graph.

Matching(graph, matching[, types])

A matching of vertices in a graph.

MatplotlibGraphDrawer(ax[, ...])

Graph drawer that uses a pyplot.Axes as context

Matrix([data])

Simple matrix data type.

Palette(n)

Base class of color palettes.

Plot

alias of CairoPlot

Point(x, y)

Class representing a point on the 2D plane.

PrecalculatedPalette(items)

A palette that returns colors from a pre-calculated list of colors

RainbowPalette([n, s, v, start, end, alpha])

A palette that varies the hue of the colors along a scale.

Rectangle(*args)

Class representing a rectangle.

RunningMean([items, n, mean, sd])

Running mean calculator.

TriadCensus([iterable])

Triad census of a graph.

UniqueIdGenerator([id_generator, initial])

A dictionary-like class that can be used to assign unique IDs to names (say, vertex names).

Vertex

Class representing a single vertex in a graph.

VertexClustering(graph[, membership, ...])

The clustering of the vertex set of a graph.

VertexCover(graph[, clusters])

The cover of the vertex set of a graph.

VertexDendrogram(graph, merges[, ...])

The dendrogram resulting from the hierarchical clustering of the vertex set of a graph.

VertexSeq

Class representing a sequence of vertices in the graph.

Exceptions

InternalError

jLM.CMEPostProcessing.closeLMFile(f)[source]#

Close a Lattice Microbes File

Parameters:

f – A previously opened lattice microbes file

jLM.CMEPostProcessing.getAvgVarTrace(f, specie, doublingTime=None)[source]#

Get the average and variance of the specie trace over time

Parameters:
  • f – The HDF5 file handle to extract from or the name of a file to open

  • specie – The specie to extract

  • doublingTime – An optional doubling time parameter that will normalize average each trace by time in the cell cycle assuming exponentially growing cell (2log2 * x/2^(t/DT)), effectively normalizing against cell size. Default: no averaging is performed

Returns:

avg, var, time

jLM.CMEPostProcessing.getHistogram(f, species)[source]#

Get the histogram for a specie

Parameters:
  • f – The HDF5 file handle to extract from or the name of a file to open

  • specie – An array of specie/s to extract. returns a nD histogram of the data.

Returns:

bins, edges (len(bins)+1)

jLM.CMEPostProcessing.getNumTrajectories(f)[source]#

Get the number of trajectories (replicates) in the file

Parameters:

f – The HDF5 file handle to extract from or the name of a file to open

Returns:

Number of trajectories/replicates in the file

Return type:

int

jLM.CMEPostProcessing.getOccupancyKymograph(f, species=None, replicate=1)[source]#

Compute the specie density(occupancy) among a slice of the simulation domain as a function over time for the given direciton

Parameters:
  • filename – Name of file to extract data from

  • specie – A particular specie to plot density for

  • replicate – The replicate to show trace for

jLM.CMEPostProcessing.getPhaseSpace(f, species, replicate)[source]#

Get the nD phase space associated with the traces of species.

If a single replicate is specified, a single trace will be returned, otherwise a 2/3D density matrix will be returned.

Parameters:
  • f – The HDF5 file handle to extract from or the name of a file to open

  • specie – An iterable of 2 more specie names

  • replicate – The replicate for which to extract the phase space

Returns:

A numpy array of a trace from a single replicate (specie1, specie2, …)

jLM.CMEPostProcessing.getSpecieTrace(f, specie, replicate=None, doublingTime=None)[source]#

Extract data for a particular specie for the specified replicate(s)

Parameters:
  • f – The HDF5 file handle to extract from or the name of a file to open

  • specie – The specie to extract data

  • replicate – The number of the replicate to extract from. If None (default), returns data for all replicates

  • doublingTime – An optional doubling time parameter that will normalize average each trace by time in the cell cycle assuming exponentially growing cell (2log2 * x/2^(t/DT)), effectively normalizing against cell size. Default: no averaging is performed

Returns:

The species time trace in a numpy array If replicate is None: A dictionary with replicate numbers as keys and time traces as values

Return type:

If replicate is specified

jLM.CMEPostProcessing.getSpecies(f)[source]#

Exract the species names

Parameters:

f – The HDF5 file handle to extract from or the name of a file to open

jLM.CMEPostProcessing.getTimesteps(f)[source]#

Extract the timestep times

Parameters:

f – The HDF5 file handle to extract from or the name of a file to open

Returns:

The timestep times in a numpy array

jLM.CMEPostProcessing.getTrajectory(f, replicate, specie)[source]#

Get trajectory data for a specific replicate and species

Parameters:
  • f – The HDF5 file handle to extract from or the name of a file to open

  • replicate – The replicate number (0-indexed for compatibility with your usage)

  • specie – The species name to extract

Returns:

The species trajectory data

Return type:

numpy.ndarray

jLM.CMEPostProcessing.openLMFile(filename)[source]#

Open a Lattice Microbes File for reading

Parameters:

filename – Name of the file

Returns:

a handle to the file

jLM.CMEPostProcessing.plotAvgVar(f, species=None, filename=None, **kwargs)[source]#

Plot a specific species average over time and variance

Parameters:
  • f – An h5py object handle

  • species – A list of species to show

  • filename – A filename to print to (default None, gives the same behavior as showAvgVar(…))

  • kwargs

    Additional arguments to be passed on to matplotlib.plot.

    These are any arguments that are valid with matplotlib.plot

Returns:

A handle to the figure object created which allows customization of plot attributes

jLM.CMEPostProcessing.plotAvgVarFromFile(filename, species, outfile, **kwargs)[source]#

Plot species from an output file

Parameters:
  • filename – The name of an HDF5 output file generated by LatticeMicrobes

  • species – A list of species to plot

  • kwargs – Additional arguments to be passed on to matplotlib.plot. These are any arguments that are valid with matplotlib.plot

  • outfile – A filename to plot to

jLM.CMEPostProcessing.plotCMEReactionNetwork(sim, filename, withRxnNodes=False, collapseReversible=False)[source]#

Plot the reaction scheme as a network

Parameters:
  • sim – An CMESimulation object

  • filename – A file to which to output

  • withRxnNodes – Plot the graph with reactions as nodes (default: false)

  • collapseReversible – Collapsee reversible reactions into one node (default: false)

jLM.CMEPostProcessing.plotOccupancyKymograph(f, species=None, replicate=1, filename=None)[source]#

Plot a specific specie density(occupancy) among a slice of the simulation domain as a function of a direction over time

Parameters:
  • f – The name of the LM file

  • specie – A particular specie to plot density for

  • filename – A filename to print to (default None, gives the same behavior as showAvgVar(…))

  • replicate – The replicate to show trace for the kymograph from the file

jLM.CMEPostProcessing.plotPhaseSpace(f, species=None, replicate=1, withHistogram=False)[source]#

Plot the 2D or 3D phase space associated with the given species over the replicates indicated

Parameters:
  • f – An h5py object handle

  • specie – An iterable of 2 or 3 specie names indicating whether to plot in 2D or 3D space

  • replicate – The replicate to show trace for

  • withHistogram – If set to true, a heatmap of the phase space over all replicates will be plotted in the background. (NOTE: this only works in 2D.)

jLM.CMEPostProcessing.plotTrace(f, species=None, replicate=1, filename=None, **kwargs)[source]#

Plot a specific species trace

Parameters:
  • f – An h5py object handle

  • species – A specie name or a list of species to show; can be a single string or a iterable list of species

  • filename – A filename to print to (default None, gives the same behavior as showAvgVar(…))

  • replicate – The replicate to show trace for. Can be an integer or an iterable list of replicte numbers

  • kwargs

    Additional arguments to be passed on to matplotlib.plot.

    These are any arguments that are valid with matplotlib.plot

Returns:

A handle to the figure object created which allows customization of plot attributes

jLM.CMEPostProcessing.plotTraceFromFile(filename, species, replicate, outfile, **kwargs)[source]#

Plot species from an output file

Parameters:
  • filename – The patch to an HDF5 output file generated by LatticeMicrobes

  • species – A list of species to plot

  • replicate – The replicate to show trace for

  • kwargs – Additional arguments to be passed on to matplotlib.plot. These are any arguments that are valid with matplotlib.plot

  • outfile – A filename to plot to

jLM.CMEPostProcessing.showAvgVar(f, species, **kwargs)[source]#

Show a specific species average over time and variance

Parameters:
  • f – An h5py object handle

  • species – A list of species to show

  • kwargs – Additional arguments to be passed on to matplotlib.plot. These are any arguments that are valid with matplotlib.plot

jLM.CMEPostProcessing.showAvgVarFromFile(filename, species, **kwargs)[source]#

Show species from an output file :param filename: The name of an HDF5 output file generated by LatticeMicrobes :param species: A list of species to plot :param kwargs: Additional arguments to be passed on to matplotlib.plot. These are any arguments that are valid with matplotlib.plot

jLM.CMEPostProcessing.showOccupancyKymograph(f, species=None, replicate=1)[source]#

Show a specific specie density(occupancy) among a slice of the simulation domain as a function of a direction over time

Parameters:
  • f – An h5py object handle

  • specie – A particular specie to plot density for

  • replicate – The replicate to show trace for

jLM.CMEPostProcessing.showTrace(f, species, replicate, **kwargs)[source]#

Show a specific species trace

Parameters:
  • f – An h5py object handle

  • species – A list of species to show

  • replicate – The replicate to show trace for

  • kwargs – Additional arguments to be passed on to matplotlib.plot. These are any arguments that are valid with matplotlib.plot

jLM.CMEPostProcessing.showTraceFromFile(filename, species, replicate, **kwargs)[source]#

Show species trace from a particular replicate

Parameters:
  • filename – The patch to an HDF5 output file generated by LatticeMicrobes

  • species – A list of species to plot

  • replicate – The replicate to show trace for

  • kwargs – Additional arguments to be passed on to matplotlib.plot. These are any arguments that are valid with matplotlib.plot