Package specific_fns :: Module n_state_model :: Class N_state_model
[hide private]
[frames] | no frames]

Class N_state_model

source code


Class containing functions for the N-state model.

Instance Methods [hide private]
 
__init__(self)
Initialise the class by placing API_common methods into the API.
source code
numpy array
_assemble_param_vector(self, sim_index=None)
Assemble all the parameters of the model into a single array.
source code
numpy rank-2 array
_assemble_scaling_matrix(self, data_types=None, scaling=True)
Create and return the scaling matrix.
source code
list of str
_base_data_types(self)
Determine all the base data types.
source code
float
_calc_ave_dist(self, atom1, atom2, exp=1)
Calculate the average distances.
source code
 
_CoM(self, pivot_point=None, centre=None)
Centre of mass analysis.
source code
 
_cone_pdb(self, cone_type=None, scale=1.0, file=None, dir=None, force=False)
Create a PDB file containing a geometric object representing the various cone models.
source code
 
_disassemble_param_vector(self, param_vector=None, data_types=None, sim_index=None)
Disassemble the parameter vector and place the values into the relevant variables.
source code
 
_elim_no_prob(self)
Remove all structures or states which have no probability.
source code
tuple of len 2 of a numpy rank-2, size NxM matrix and numpy rank-1, size N array
_linear_constraints(self, data_types=None, scaling_matrix=None)
Function for setting up the linear constraint matrices A and b.
source code
 
_minimise_bc_data(self, model)
Extract and unpack the back calculated data.
source code
numpy rank-3 array, numpy rank-1 array.
_minimise_setup_atomic_pos(self)
Set up the atomic position data structures for optimisation using PCSs and PREs as base data sets.
source code
tuple of (numpy rank-2 array, numpy rank-2 array, numpy rank-2 array, numpy rank-1 array, numpy rank-1 array)
_minimise_setup_pcs(self, sim_index=None)
Set up the data structures for optimisation using PCSs as base data sets.
source code
tuple of (numpy rank-2 array, numpy rank-2 array, numpy rank-2 array)
_minimise_setup_rdcs(self, sim_index=None)
Set up the data structures for optimisation using RDCs as base data sets.
source code
tuple of (list, numpy rank-1 array, numpy rank-1 array, numpy rank-1 array)
_minimise_setup_tensors(self, sim_index=None)
Set up the data structures for optimisation using alignment tensors as base data sets.
source code
numpy rank-1 array.
_minimise_setup_fixed_tensors(self)
Set up the data structures for the fixed alignment tensors.
source code
int
_num_data_points(self)
Determine the number of data points used in the model.
source code
 
_number_of_states(self, N=None)
Set the number of states in the N-state model.
source code
str
_param_model_index(self, param)
Return the N-state model index for the given parameter.
source code
int
_param_num(self)
Determine the number of parameters in the model.
source code
 
_ref_domain(self, ref=None)
Set the reference domain for the '2-domain' N-state model.
source code
 
_select_model(self, model=None)
Select the N-state model type.
source code
 
_target_fn_setup(self, sim_index=None, scaling=True)
Initialise the target function for optimisation or direct calculation.
source code
(int, AlignTensorData instance)
_tensor_loop(self, red=False)
Generator method for looping over the full or reduced tensors.
source code
 
_update_model(self)
Update the model parameters as necessary.
source code
list of str
base_data_loop(self)
Loop over the base data of the spins - RDCs, PCSs, and NOESY data.
source code
 
calculate(self, spin_id=None, verbosity=1, sim_index=None)
Calculation function.
source code
list of floats
create_mc_data(self, data_id=None)
Create the Monte Carlo data by back-calculation.
source code
float
default_value(self, param)
The default N-state model parameter values.
source code
 
grid_search(self, lower=None, upper=None, inc=None, constraints=False, verbosity=0, sim_index=None)
The grid search function.
source code
bool
is_spin_param(self, name)
Determine whether the given parameter is spin specific.
source code
list of float
map_bounds(self, param, spin_id=None)
Create bounds for the OpenDX mapping function.
source code
 
minimise(self, min_algor=None, min_options=None, func_tol=None, grad_tol=None, max_iterations=None, constraints=False, scaling=True, verbosity=0, sim_index=None, lower=None, upper=None, inc=None)
Minimisation function.
source code
tuple of (int, int, float)
model_statistics(self, model_info=None, spin_id=None, global_stats=None)
Return the k, n, and chi2 model statistics.
source code
str
return_data_name(self, param)
Return a unique identifying string for the N-state model parameter.
source code
list of floats
return_error(self, data_id=None)
Create and return the spin specific Monte Carlo Ri error structure.
source code
str
return_grace_string(self, param)
Return the Grace string representation of the parameter.
source code
str
return_units(self, param)
Return a string representing the parameters units.
source code
 
set_error(self, model_info, index, error)
Set the parameter errors.
source code
 
set_param_values(self, param=None, value=None, spin_id=None, force=True)
Set the N-state model parameter values.
source code
 
sim_init_values(self)
Initialise the Monte Carlo parameter values.
source code
 
sim_pack_data(self, data_id, sim_data)
Pack the Monte Carlo simulation data.
source code
list of float
sim_return_param(self, model_info, index)
Return the array of simulation parameter values.
source code

Inherited from api_base.API_base: back_calc_ri, bmrb_read, bmrb_write, data_init, data_names, data_type, deselect, duplicate_data, eliminate, get_param_names, get_param_values, has_errors, model_desc, model_loop, model_type, molmol_macro, num_instances, overfit_deselect, pymol_macro, read_columnar_results, return_conversion_factor, return_data, return_data_desc, return_value, set_selected_sim, set_update, sim_return_chi2, sim_return_selected, skip_function, test_grid_ops

Inherited from object: __delattr__, __format__, __getattribute__, __hash__, __new__, __reduce__, __reduce_ex__, __repr__, __setattr__, __sizeof__, __str__, __subclasshook__

Class Variables [hide private]
  default_value_doc = Desc_container("N-state model default valu...
  return_data_name_doc = Desc_container("N-state model data type...
  _table = uf_tables.add_table(label= "table: N-state data type ...
  set_doc = Desc_container("N-state model set details")

Inherited from api_base.API_base: eliminate_doc

Properties [hide private]

Inherited from object: __class__

Method Details [hide private]

__init__(self)
(Constructor)

source code 

Initialise the class by placing API_common methods into the API.

Overrides: object.__init__

_assemble_param_vector(self, sim_index=None)

source code 

Assemble all the parameters of the model into a single array.

Parameters:
  • sim_index (None or int) - The index of the simulation to optimise. This should be None if normal optimisation is desired.
Returns: numpy array
The parameter vector used for optimisation.

_assemble_scaling_matrix(self, data_types=None, scaling=True)

source code 

Create and return the scaling matrix.

Parameters:
  • data_types (list of str) - The base data types used in the optimisation. This list can contain the elements 'rdc', 'pcs' or 'tensor'.
  • scaling (bool) - If False, then the identity matrix will be returned.
Returns: numpy rank-2 array
The square and diagonal scaling matrix.

_base_data_types(self)

source code 

Determine all the base data types.

The base data types can include:

   - 'rdc', residual dipolar couplings.
   - 'pcs', pseudo-contact shifts.
   - 'noesy', NOE restraints.
   - 'tensor', alignment tensors.
Returns: list of str
A list of all the base data types.

_calc_ave_dist(self, atom1, atom2, exp=1)

source code 

Calculate the average distances.

The formula used is:

             _N_
         / 1 \                  \ 1/exp
   <r> = | -  > |p1i - p2i|^exp |
         \ N /__                /
              i

where i are the members of the ensemble, N is the total number of structural models, and p1 and p2 at the two atom positions.

Parameters:
  • atom1 (str) - The atom identification string of the first atom.
  • atom2 (str) - The atom identification string of the second atom.
  • exp (int) - The exponent used for the averaging, e.g. 1 for linear averaging and -6 for r^-6 NOE averaging.
Returns: float
The average distance between the two atoms.

_CoM(self, pivot_point=None, centre=None)

source code 

Centre of mass analysis.

This function does an analysis of the centre of mass (CoM) of the N states. This includes calculating the order parameter associated with the pivot-CoM vector, and the associated cone of motions. The pivot_point argument must be supplied. If centre is None, then the CoM will be calculated from the selected parts of the loaded structure. Otherwise it will be set to the centre arg.

Parameters:
  • pivot_point (list of float of length 3) - The pivot point in the structural file(s).
  • centre (list of float of length 3) - The optional centre of mass vector.

_cone_pdb(self, cone_type=None, scale=1.0, file=None, dir=None, force=False)

source code 

Create a PDB file containing a geometric object representing the various cone models.

Currently the only cone types supported are 'diff in cone' and 'diff on cone'.

Parameters:
  • cone_type (str) - The type of cone model to represent.
  • scale (float) - The size of the geometric object is eqaul to the average pivot-CoM vector length multiplied by this scaling factor.
  • file (str) - The name of the PDB file to create.
  • dir (str) - The name of the directory to place the PDB file into.
  • force (int) - Flag which if set to True will cause any pre-existing file to be overwritten.

_disassemble_param_vector(self, param_vector=None, data_types=None, sim_index=None)

source code 

Disassemble the parameter vector and place the values into the relevant variables.

For the 2-domain N-state model, the parameters are stored in the probability and Euler angle data structures. For the population N-state model, only the probabilities are stored. If RDCs are present and alignment tensors are optimised, then these are stored as well.

Parameters:
  • data_types (list of str) - The base data types used in the optimisation. This list can contain the elements 'rdc', 'pcs' or 'tensor'.
  • param_vector (numpy array) - The parameter vector returned from optimisation.
  • sim_index (None or int) - The index of the simulation to optimise. This should be None if normal optimisation is desired.

_linear_constraints(self, data_types=None, scaling_matrix=None)

source code 

Function for setting up the linear constraint matrices A and b.

Standard notation

The N-state model constraints are:

   0 <= pc <= 1,

where p is the probability and c corresponds to state c.

Matrix notation

In the notation A.x >= b, where A is an matrix of coefficients, x is an array of parameter values, and b is a vector of scalars, these inequality constraints are:

   | 1  0  0 |                   |    0    |
   |         |                   |         |
   |-1  0  0 |                   |   -1    |
   |         |                   |         |
   | 0  1  0 |                   |    0    |
   |         |     |  p0  |      |         |
   | 0 -1  0 |     |      |      |   -1    |
   |         |  .  |  p1  |  >=  |         |
   | 0  0  1 |     |      |      |    0    |
   |         |     |  p2  |      |         |
   | 0  0 -1 |                   |   -1    |
   |         |                   |         |
   |-1 -1 -1 |                   |   -1    |
   |         |                   |         |
   | 1  1  1 |                   |    0    |

This example is for a 4-state model, the last probability pn is not included as this parameter does not exist (because the sum of pc is equal to 1). The Euler angle parameters have been excluded here but will be included in the returned A and b objects. These parameters simply add columns of zero to the A matrix and have no effect on b. The last two rows correspond to the inequality:

   0 <= pN <= 1.

As:

           N-1
                       pN = 1 - >  pc,
           /__
           c=1

then:

   -p1 - p2 - ... - p(N-1) >= -1,

    p1 + p2 + ... + p(N-1) >= 0.
Parameters:
  • data_types (list of str) - The base data types used in the optimisation. This list can contain the elements 'rdc', 'pcs' or 'tensor'.
  • scaling_matrix (numpy rank-2 square matrix) - The diagonal scaling matrix.
Returns: tuple of len 2 of a numpy rank-2, size NxM matrix and numpy rank-1, size N array
The matrices A and b.

_minimise_bc_data(self, model)

source code 

Extract and unpack the back calculated data.

Parameters:
  • model (class instance) - The instantiated class containing the target function.

_minimise_setup_atomic_pos(self)

source code 

Set up the atomic position data structures for optimisation using PCSs and PREs as base data sets.

Returns: numpy rank-3 array, numpy rank-1 array.
The atomic positions (the first index is the spins, the second is the structures, and the third is the atomic coordinates) and the paramagnetic centre.

_minimise_setup_pcs(self, sim_index=None)

source code 

Set up the data structures for optimisation using PCSs as base data sets.

Parameters:
  • sim_index (None or int) - The index of the simulation to optimise. This should be None if normal optimisation is desired.
Returns: tuple of (numpy rank-2 array, numpy rank-2 array, numpy rank-2 array, numpy rank-1 array, numpy rank-1 array)
The assembled data structures for using PCSs as the base data for optimisation. These include:
  • the PCS values.
  • the unit vectors connecting the paramagnetic centre (the electron spin) to
  • the PCS weight.
  • the nuclear spin.
  • the pseudocontact shift constants.

_minimise_setup_rdcs(self, sim_index=None)

source code 

Set up the data structures for optimisation using RDCs as base data sets.

Parameters:
  • sim_index (None or int) - The index of the simulation to optimise. This should be None if normal optimisation is desired.
Returns: tuple of (numpy rank-2 array, numpy rank-2 array, numpy rank-2 array)
The assembled data structures for using RDCs as the base data for optimisation. These include:
  • rdc, the RDC values.
  • rdc_err, the RDC errors.
  • rdc_weight, the RDC weights.
  • vectors, the heteronucleus to proton vectors.
  • rdc_const, the dipolar constants.

_minimise_setup_tensors(self, sim_index=None)

source code 

Set up the data structures for optimisation using alignment tensors as base data sets.

Parameters:
  • sim_index (None or int) - The index of the simulation to optimise. This should be None if normal optimisation is desired.
Returns: tuple of (list, numpy rank-1 array, numpy rank-1 array, numpy rank-1 array)
The assembled data structures for using alignment tensors as the base data for optimisation. These include:
  • full_tensors, the data of the full alignment tensors.
  • red_tensor_elem, the tensors as concatenated rank-1 5D arrays.
  • red_tensor_err, the tensor errors as concatenated rank-1 5D arrays.
  • full_in_ref_frame, flags specifying if the tensor in the reference frame is the full or reduced tensor.

_minimise_setup_fixed_tensors(self)

source code 

Set up the data structures for the fixed alignment tensors.

Returns: numpy rank-1 array.
The assembled data structures for the fixed alignment tensors.

_num_data_points(self)

source code 

Determine the number of data points used in the model.

Returns: int
The number, n, of data points in the model.

_number_of_states(self, N=None)

source code 

Set the number of states in the N-state model.

Parameters:
  • N (int) - The number of states.

_param_model_index(self, param)

source code 

Return the N-state model index for the given parameter.

Parameters:
  • param (str) - The N-state model parameter.
Returns: str
The N-state model index.

_param_num(self)

source code 

Determine the number of parameters in the model.

Returns: int
The number of model parameters.

_ref_domain(self, ref=None)

source code 

Set the reference domain for the '2-domain' N-state model.

Parameters:
  • ref (str) - The reference domain.

_select_model(self, model=None)

source code 

Select the N-state model type.

Parameters:
  • model (str) - The N-state model type. Can be one of '2-domain', 'population', or 'fixed'.

_target_fn_setup(self, sim_index=None, scaling=True)

source code 

Initialise the target function for optimisation or direct calculation.

Parameters:
  • sim_index (None or int) - The index of the simulation to optimise. This should be None if normal optimisation is desired.
  • scaling (bool) - If True, diagonal scaling is enabled during optimisation to allow the problem to be better conditioned.

_tensor_loop(self, red=False)

source code 

Generator method for looping over the full or reduced tensors.

Parameters:
  • red (bool) - A flag which if True causes the reduced tensors to be returned, and if False the full tensors are returned.
Returns: (int, AlignTensorData instance)
The tensor index and the tensor.

base_data_loop(self)

source code 

Loop over the base data of the spins - RDCs, PCSs, and NOESY data.

This loop iterates for each data point (RDC, PCS, NOESY) for each spin, returning the identification information.

Returns: list of str
A list of the spin ID string, the data type ('rdc', 'pcs', 'noesy'), and the alignment ID if required.
Overrides: api_base.API_base.base_data_loop

calculate(self, spin_id=None, verbosity=1, sim_index=None)

source code 

Calculation function.

Currently this function simply calculates the NOESY flat-bottom quadratic energy potential, if NOE restraints are available.

Parameters:
  • spin_id (None or str) - The spin identification string (unused).
  • verbosity (int) - The amount of information to print. The higher the value, the greater the verbosity.
  • sim_index (None) - The MC simulation index (unused).
Overrides: api_base.API_base.calculate

create_mc_data(self, data_id=None)

source code 

Create the Monte Carlo data by back-calculation.

Parameters:
  • data_id (str) - The list of spin ID, data type, and alignment ID, as yielded by the base_data_loop() generator method.
Returns: list of floats
The Monte Carlo Ri data.
Overrides: api_base.API_base.create_mc_data

default_value(self, param)

source code 

The default N-state model parameter values.

Parameters:
  • param (str) - The N-state model parameter.
Returns: float
The default value.
Overrides: api_base.API_base.default_value

grid_search(self, lower=None, upper=None, inc=None, constraints=False, verbosity=0, sim_index=None)

source code 

The grid search function.

Parameters:
  • lower (array of numbers) - The lower bounds of the grid search which must be equal to the number of parameters in the model.
  • upper (array of numbers) - The upper bounds of the grid search which must be equal to the number of parameters in the model.
  • inc (array of int) - The increments for each dimension of the space for the grid search. The number of elements in the array must equal to the number of parameters in the model.
  • constraints (bool) - If True, constraints are applied during the grid search (elinating parts of the grid). If False, no constraints are used.
  • verbosity (int) - A flag specifying the amount of information to print. The higher the value, the greater the verbosity.
Overrides: api_base.API_base.grid_search

is_spin_param(self, name)

source code 

Determine whether the given parameter is spin specific.

Parameters:
  • name (str) - The name of the parameter.
Returns: bool
False
Overrides: api_base.API_base.is_spin_param

map_bounds(self, param, spin_id=None)

source code 

Create bounds for the OpenDX mapping function.

Parameters:
  • param (str) - The name of the parameter to return the lower and upper bounds of.
  • spin_id (None) - The spin identification string (unused).
Returns: list of float
The upper and lower bounds of the parameter.
Overrides: api_base.API_base.map_bounds

minimise(self, min_algor=None, min_options=None, func_tol=None, grad_tol=None, max_iterations=None, constraints=False, scaling=True, verbosity=0, sim_index=None, lower=None, upper=None, inc=None)

source code 

Minimisation function.

Parameters:
  • min_algor (str) - The minimisation algorithm to use.
  • min_options (array of str) - An array of options to be used by the minimisation algorithm.
  • func_tol (None or float) - The function tolerance which, when reached, terminates optimisation. Setting this to None turns of the check.
  • grad_tol (None or float) - The gradient tolerance which, when reached, terminates optimisation. Setting this to None turns of the check.
  • max_iterations (int) - The maximum number of iterations for the algorithm.
  • constraints (bool) - If True, constraints are used during optimisation.
  • scaling (bool) - If True, diagonal scaling is enabled during optimisation to allow the problem to be better conditioned.
  • verbosity (int) - A flag specifying the amount of information to print. The higher the value, the greater the verbosity.
  • sim_index (None or int) - The index of the simulation to optimise. This should be None if normal optimisation is desired.
  • lower (array of numbers) - The lower bounds of the grid search which must be equal to the number of parameters in the model. This optional argument is only used when doing a grid search.
  • upper (array of numbers) - The upper bounds of the grid search which must be equal to the number of parameters in the model. This optional argument is only used when doing a grid search.
  • inc (array of int) - The increments for each dimension of the space for the grid search. The number of elements in the array must equal to the number of parameters in the model. This argument is only used when doing a grid search.
Overrides: api_base.API_base.minimise

model_statistics(self, model_info=None, spin_id=None, global_stats=None)

source code 

Return the k, n, and chi2 model statistics.

k - number of parameters. n - number of data points. chi2 - the chi-squared value.

Parameters:
  • model_info (None) - The data returned from model_loop() (unused).
  • spin_id (None or str) - The spin identification string. This is ignored in the N-state model.
  • global_stats (None or bool) - A parameter which determines if global or local statistics are returned. For the N-state model, this argument is ignored.
Returns: tuple of (int, int, float)
The optimisation statistics, in tuple format, of the number of parameters (k), the number of data points (n), and the chi-squared value (chi2).
Overrides: api_base.API_base.model_statistics

return_data_name(self, param)

source code 

Return a unique identifying string for the N-state model parameter.

Parameters:
  • param (str) - The N-state model parameter.
Returns: str
The unique parameter identifying string.
Overrides: api_base.API_base.return_data_name

return_error(self, data_id=None)

source code 

Create and return the spin specific Monte Carlo Ri error structure.

Parameters:
  • data_id (str) - The list of spin ID, data type, and alignment ID, as yielded by the base_data_loop() generator method.
Returns: list of floats
The Monte Carlo simulation data errors.
Overrides: api_base.API_base.return_error

return_grace_string(self, param)

source code 

Return the Grace string representation of the parameter.

This is used for axis labelling.

Parameters:
  • param (str) - The specific analysis parameter.
Returns: str
The Grace string representation of the parameter.
Overrides: api_base.API_base.return_grace_string

return_units(self, param)

source code 

Return a string representing the parameters units.

Parameters:
  • param (str) - The name of the parameter to return the units string for.
Returns: str
The parameter units string.
Overrides: api_base.API_base.return_units

set_error(self, model_info, index, error)

source code 

Set the parameter errors.

Parameters:
  • model_info (int) - The global model index originating from model_loop().
  • index (int) - The index of the parameter to set the errors for.
  • error (float) - The error value.
Overrides: api_base.API_base.set_error

set_param_values(self, param=None, value=None, spin_id=None, force=True)

source code 

Set the N-state model parameter values.

Parameters:
  • param (list of str) - The parameter name list.
  • value (list) - The parameter value list.
  • spin_id (None) - The spin identification string (unused).
  • force (bool) - A flag which if True will cause current values to be overwritten. If False, a RelaxError will raised if the parameter value is already set.
Overrides: api_base.API_base.set_param_values

sim_init_values(self)

source code 

Initialise the Monte Carlo parameter values.

Overrides: api_base.API_base.sim_init_values

sim_pack_data(self, data_id, sim_data)

source code 

Pack the Monte Carlo simulation data.

Parameters:
  • data_id (list of str) - The list of spin ID, data type, and alignment ID, as yielded by the base_data_loop() generator method.
  • sim_data (list of float) - The Monte Carlo simulation data.
Overrides: api_base.API_base.sim_pack_data

sim_return_param(self, model_info, index)

source code 

Return the array of simulation parameter values.

Parameters:
  • model_info (int) - The global model index originating from model_loop().
  • index (int) - The index of the parameter to return the array of values for.
Returns: list of float
The array of simulation parameter values.
Overrides: api_base.API_base.sim_return_param

Class Variable Details [hide private]

default_value_doc

Value:
Desc_container("N-state model default values")

return_data_name_doc

Value:
Desc_container("N-state model data type string matching patterns")

_table

Value:
uf_tables.add_table(label= "table: N-state data type patterns", captio\
n= "N-state model data type string matching patterns.")