Bio.phenotype.phen_micro module

Classes to work with Phenotype Microarray data.

More information on the single plates can be found here: http://www.biolog.com/

Classes:
  • PlateRecord - Object that contain time course data on each well of the plate, as well as metadata (if any).

  • WellRecord - Object that contains the time course data of a single well

  • JsonWriter - Writer of PlateRecord objects in JSON format.

Functions:
  • JsonIterator - Incremental PM JSON parser, this is an iterator that returns PlateRecord objects.

  • CsvIterator - Incremental PM CSV parser, this is an iterator that returns PlateRecord objects.

  • _toOPM - Used internally by JsonWriter, converts PlateRecord objects in dictionaries ready to be serialized in JSON format.

class Bio.phenotype.phen_micro.PlateRecord(plateid, wells=None)

Bases: object

PlateRecord object for storing Phenotype Microarray plates data.

A PlateRecord stores all the wells of a particular phenotype Microarray plate, along with metadata (if any). The single wells can be accessed calling their id as an index or iterating on the PlateRecord:

>>> from Bio import phenotype
>>> plate = phenotype.read("phenotype/Plate.json", "pm-json")
>>> well = plate['A05']
>>> for well in plate:
...    print(well.id)
...
A01
...

The plate rows and columns can be queried with an indexing system similar to NumPy and other matrices:

>>> print(plate[1])
Plate ID: PM01
Well: 12
Rows: 1
Columns: 12
PlateRecord('WellRecord['B01'], WellRecord['B02'], WellRecord['B03'], ..., WellRecord['B12']')
>>> print(plate[:,1])
Plate ID: PM01
Well: 8
Rows: 8
Columns: 1
PlateRecord('WellRecord['A02'], WellRecord['B02'], WellRecord['C02'], ..., WellRecord['H02']')

Single WellRecord objects can be accessed using this indexing system:

>>> print(plate[1,2])
Plate ID: PM01
Well ID: B03
Time points: 384
Minum signal 0.00 at time 11.00
Maximum signal 76.25 at time 18.00
WellRecord('(0.0, 11.0), (0.25, 11.0), (0.5, 11.0), (0.75, 11.0), (1.0, 11.0), ..., (95.75, 11.0)')

The presence of a particular well can be inspected with the “in” keyword: >>> ‘A01’ in plate True

All the wells belonging to a “row” (identified by the first character of the well id) in the plate can be obtained:

>>> for well in plate.get_row('H'):
...     print(well.id)
...
H01
H02
H03
...

All the wells belonging to a “column” (identified by the number of the well) in the plate can be obtained:

>>> for well in plate.get_column(12):
...     print(well.id)
...
A12
B12
C12
...

Two PlateRecord objects can be compared: if all their wells are equal the two plates are considered equal:

>>> plate2 = phenotype.read("phenotype/Plate.json", "pm-json")
>>> plate == plate2
True

Two PlateRecord object can be summed up or subtracted from each other: the the signals of each well will be summed up or subtracted. The id of the left operand will be kept:

>>> plate3 = plate + plate2
>>> print(plate3.id)
PM01

Many Phenotype Microarray plate have a “negative control” well, which can be subtracted to all wells:

>>> subplate = plate.subtract_control()
__init__(plateid, wells=None)

Initialize the class.

__getitem__(index)

Access part of the plate.

Depending on the indices, you can get a WellRecord object (representing a single well of the plate), or another plate (representing some part or all of the original plate).

plate[wid] gives a WellRecord (if wid is a WellRecord id) plate[r,c] gives a WellRecord plate[r] gives a row as a PlateRecord plate[r,:] gives a row as a PlateRecord plate[:,c] gives a column as a PlateRecord

plate[:] and plate[:,:] give a copy of the plate

Anything else gives a subset of the original plate, e.g. plate[0:2] or plate[0:2,:] uses only row 0 and 1 plate[:,1:3] uses only columns 1 and 2 plate[0:2,1:3] uses only rows 0 & 1 and only cols 1 & 2

>>> from Bio import phenotype
>>> plate = phenotype.read("phenotype/Plate.json", "pm-json")

You can access a well of the plate, using its id.

>>> w = plate['A01']

You can access a row of the plate as a PlateRecord using an integer index:

>>> first_row = plate[0]
>>> print(first_row)
Plate ID: PM01
Well: 12
Rows: 1
Columns: 12
PlateRecord('WellRecord['A01'], WellRecord['A02'], WellRecord['A03'], ..., WellRecord['A12']')
>>> last_row = plate[-1]
>>> print(last_row)
Plate ID: PM01
Well: 12
Rows: 1
Columns: 12
PlateRecord('WellRecord['H01'], WellRecord['H02'], WellRecord['H03'], ..., WellRecord['H12']')

You can also access use python’s slice notation to sub-plates containing only some of the plate rows:

>>> sub_plate = plate[2:5]
>>> print(sub_plate)
Plate ID: PM01
Well: 36
Rows: 3
Columns: 12
PlateRecord('WellRecord['C01'], WellRecord['C02'], WellRecord['C03'], ..., WellRecord['E12']')

This includes support for a step, i.e. plate[start:end:step], which can be used to select every second row:

>>> sub_plate = plate[::2]

You can also use two indices to specify both rows and columns. Using simple integers gives you the single wells. e.g.

>>> w = plate[3, 4]
>>> print(w.id)
D05

To get a single column use this syntax:

>>> sub_plate = plate[:, 4]
>>> print(sub_plate)
Plate ID: PM01
Well: 8
Rows: 8
Columns: 1
PlateRecord('WellRecord['A05'], WellRecord['B05'], WellRecord['C05'], ..., WellRecord['H05']')

Or, to get part of a column,

>>> sub_plate = plate[1:3, 4]
>>> print(sub_plate)
Plate ID: PM01
Well: 2
Rows: 2
Columns: 1
PlateRecord(WellRecord['B05'], WellRecord['C05'])

However, in general you get a sub-plate,

>>> print(plate[1:5, 3:6])
Plate ID: PM01
Well: 12
Rows: 4
Columns: 3
PlateRecord('WellRecord['B04'], WellRecord['B05'], WellRecord['B06'], ..., WellRecord['E06']')

This should all seem familiar to anyone who has used the NumPy array or matrix objects.

__setitem__(key, value)
__delitem__(key)
__iter__()
__contains__(wellid)
__len__()

Return the number of wells in this plate.

__eq__(other)

Return self==value.

__add__(plate)

Add another PlateRecord object.

The wells in both plates must be the same

A new PlateRecord object is returned, having the same id as the left operand.

__sub__(plate)

Subtract another PlateRecord object.

The wells in both plates must be the same

A new PlateRecord object is returned, having the same id as the left operand.

get_row(row)

Get all the wells of a given row.

A row is identified with a letter (e.g. ‘A’)

get_column(column)

Get all the wells of a given column.

A column is identified with a number (e.g. ‘6’)

subtract_control(control='A01', wells=None)

Subtract a ‘control’ well from the other plates wells.

By default the control is subtracted to all wells, unless a list of well ID is provided

The control well should belong to the plate A new PlateRecord object is returned

__repr__()

Return a (truncated) representation of the plate for debugging.

__str__()

Return a human readable summary of the record (string).

The python built in function str works by calling the object’s __str__ method. e.g.

>>> from Bio import phenotype
>>> record = next(phenotype.parse("phenotype/Plates.csv", "pm-csv"))
>>> print(record)
Plate ID: PM01
Well: 96
Rows: 8
Columns: 12
PlateRecord('WellRecord['A01'], WellRecord['A02'], WellRecord['A03'], ..., WellRecord['H12']')

Note that long well lists are shown truncated.

__hash__ = None
class Bio.phenotype.phen_micro.WellRecord(wellid, plate=None, signals=None)

Bases: object

WellRecord stores all time course signals of a phenotype Microarray well.

The single time points and signals can be accessed iterating on the WellRecord or using lists indexes or slices:

>>> from Bio import phenotype
>>> plate = phenotype.read("phenotype/Plate.json", "pm-json")
>>> well = plate['A05']
>>> for time, signal in well:
...    print("Time: %f, Signal: %f" % (time, signal)) 
...
Time: 0.000000, Signal: 14.000000
Time: 0.250000, Signal: 13.000000
Time: 0.500000, Signal: 15.000000
Time: 0.750000, Signal: 15.000000
...
>>> well[1]
16.0
>>> well[1:5]
[16.0, 20.0, 18.0, 15.0]
>>> well[1:5:0.5]
[16.0, 19.0, 20.0, 18.0, 18.0, 18.0, 15.0, 18.0]

If a time point was not present in the input file but it’s between the minimum and maximum time point, the interpolated signal is returned, otherwise a nan value:

>>> well[1.3]
19.0
>>> well[1250]
nan

Two WellRecord objects can be compared: if their input time/signal pairs are exactly the same, the two records are considered equal:

>>> well2 = plate['H12']
>>> well == well2
False

Two WellRecord objects can be summed up or subtracted from each other: a new WellRecord object is returned, having the left operand id.

>>> well1 = plate['A05']
>>> well2 = well + well1
>>> print(well2.id)
A05

If SciPy is installed, a sigmoid function can be fitted to the PM curve, in order to extract some parameters; three sigmoid functions are available: * gompertz * logistic * richards The functions are described in Zwietering et al., 1990 (PMID: 16348228)

For example:

well.fit()
print(well.slope, well.model)
(61.853516785566917, 'logistic')

If not sigmoid function is specified, the first one that is successfully fitted is used. The user can also specify a specific function.

To specify gompertz:

well.fit('gompertz')
print(well.slope, well.model)
(127.94630059171354, 'gompertz')

If no function can be fitted, the parameters are left as None, except for the max, min, average_height and area.

__init__(wellid, plate=None, signals=None)

Initialize the class.

__setitem__(time, signal)

Assign a signal at a certain time point.

__getitem__(time)

Return a subset of signals or a single signal.

__iter__()
__eq__(other)

Return self==value.

__add__(well)

Add another WellRecord object.

A new WellRecord object is returned, having the same id as the left operand

__sub__(well)

Subtract another WellRecord object.

A new WellRecord object is returned, having the same id as the left operand

__len__()

Return the number of time points sampled.

__repr__()

Return a (truncated) representation of the signals for debugging.

__str__()

Return a human readable summary of the record (string).

The python built-in function str works by calling the object’s __str__ method. e.g.

>>> from Bio import phenotype
>>> plate = phenotype.read("phenotype/Plate.json", "pm-json")
>>> record = plate['A05']
>>> print(record)
Plate ID: PM01
Well ID: A05
Time points: 384
Minum signal 0.25 at time 13.00
Maximum signal 19.50 at time 23.00
WellRecord('(0.0, 14.0), (0.25, 13.0), (0.5, 15.0), (0.75, 15.0), (1.0, 16.0), ..., (95.75, 16.0)')

Note that long time spans are shown truncated.

get_raw()

Get a list of time/signal pairs.

get_times()

Get a list of the recorded time points.

get_signals()

Get a list of the recorded signals (ordered by collection time).

fit(function=('gompertz', 'logistic', 'richards'))

Fit a sigmoid function to this well and extract curve parameters.

If function is None or an empty tuple/list, then no fitting is done. Only the object’s .min, .max and .average_height are calculated.

By default the following fitting functions will be used in order:
  • gompertz

  • logistic

  • richards

The first function that is successfully fitted to the signals will be used to extract the curve parameters and update .area and .model. If no function can be fitted an exception is raised.

The function argument should be a tuple or list of any of these three function names as strings.

There is no return value.

__hash__ = None
Bio.phenotype.phen_micro.JsonIterator(handle)

Iterate over PM json records as PlateRecord objects.

Arguments:
  • handle - input file

Bio.phenotype.phen_micro.CsvIterator(handle)

Iterate over PM csv records as PlateRecord objects.

Arguments:
  • handle - input file

class Bio.phenotype.phen_micro.JsonWriter(plates)

Bases: object

Class to write PM Json format files.

__init__(plates)

Initialize the class.

write(handle)

Write this instance’s plates to a file handle.