Generate a 3D models

pytadbit.imp.imp_modelling.generate_3d_models(zscores, resolution, start=1, n_models=5000, n_keep=1000, close_bins=1, n_cpus=1, keep_all=False, verbose=0, outfile=None, config={'maxdist': 600, 'upfreq': 0.3, 'reference': 'victor corces dataset 2013', 'kforce': 5, 'lowfreq': -0.7, 'scale': 0.005}, values=None)[source]

This function generates three-dimensional models starting from Hi-C data. The final analysis will be performed on the n_keep top models.

Parameters:
  • zscores – the dictionary of the Z-score values calculated from the Hi-C pairwise interactions
  • resolution – number of nucleotides per Hi-C bin. This will be the number of nucleotides in each model’s particle
  • n_models (5000) – number of models to generate
  • n_keep (1000) – number of models used in the final analysis (usually the top 20% of the generated models). The models are ranked according to their objective function value (the lower the better)
  • keep_all (False) – whether or not to keep the discarded models (if True, models will be stored under StructuralModels.bad_models)
  • close_bins (1) – number of particles away (i.e. the bin number difference) a particle pair must be in order to be considered as neighbors (e.g. 1 means consecutive particles)
  • n_cpus – number of CPUs to use
  • verbose (False) – if set to True, information about the distance, force and Z-score between particles will be printed
  • values (None) – the normalized Hi-C data in a list of lists (equivalent to a square matrix)
  • config (CONFIG[‘dmel_01’]) –

    a dictionary containing the standard parameters used to generate the models. The dictionary should contain the keys kforce, lowrdist, maxdist, upfreq and lowfreq. Examples can be seen by doing:

    from pytadbit.imp.CONFIG import CONFIG
    
    where CONFIG is a dictionary of dictionaries to be passed to this function:

    ::

    CONFIG = {
     'dmel_01': {
         # Paramaters for the Hi-C dataset from:
         'reference' : 'victor corces dataset 2013',
    
         # Force applied to the restraints inferred to neighbor particles
         'kforce'    : 5,
    
         # Minimum distance between two non-bonded particles
         'lowrdist'  : 100,
    
         # Maximum experimental contact distance
         'maxdist'   : 600, # OPTIMIZATION: 500-1200
    
         # Maximum threshold used to decide which experimental values have to be
         # included in the computation of restraints. Z-score values greater than upfreq
         # and less than lowfreq will be included, while all the others will be rejected
         'upfreq'    : 0.3, # OPTIMIZATION: min/max Z-score
    
         # Minimum thresholds used to decide which experimental values have to be
         # included in the computation of restraints. Z-score values bigger than upfreq
         # and less that lowfreq will be include, whereas all the others will be rejected
         'lowfreq'   : -0.7 # OPTIMIZATION: min/max Z-score
    
         # Space occupied by a nucleotide (nm)
         'scale'     : 0.005
    
         }
     }
Returns:

a TheeDeeModels object

IMPmodel class

class pytadbit.imp.impmodel.IMPmodel[source]

A container for the IMP modeling results. The container is a dictionary with the following keys:

  • log_objfun: The list of IMP objective function values

  • objfun: The final objective function value of the corresponding model

    (from log_objfun). This value will be used to rank all the generated models

  • rand_init: The random number generator feed (needed for model

    reproducibility)

  • x, y, z: The spatial 3D coordinates of each particles. Each coordinate is

    represented as a list

objective_function(log=False, smooth=True, axe=None, savefig=None)[source]

This function plots the objective function value per each Monte-Carlo step.

Parameters:
  • log (False) – log plot
  • smooth (True) – curve smoothing
view_model(tool='chimera', savefig=None, cmd=None)[source]

Visualize a selected model in the three dimensions.

Parameters:
  • model_num – model to visualize
  • tool (‘chimera’) – path to the external tool used to visualize the model
  • savefig (None) – path to a file where to save the image OR movie generated (depending on the extension; accepted formats are png, mov and webm). if set to None, the image or movie will be shown using the default GUI.
  • cmd (None) –

    list of commands to be passed to the viewer. The chimera list is:

    focus
    set bg_color white
    windowsize 800 600
    bonddisplay never #0
    represent wire
    shape tube #0 radius 5 bandLength 100 segmentSubdivisions 1 followBonds on
    clip yon -500
    ~label
    set subdivision 1
    set depth_cue
    set dc_color black
    set dc_start 0.5
    set dc_end 1
    scale 0.8

    Followed by the movie command to record movies:

    movie record supersample 1
    turn y 3 120
    wait 120
    movie stop
    movie encode output SAVEFIG

    Or the copy command for images:

    copy file SAVEFIG png

    Passing as the following list as ‘cmd’ parameter:

    cmd = ['focus', 'set bg_color white', 'windowsize 800 600', 'bonddisplay never #0', 'shape tube #0 radius 10 bandLength 200 segmentSubdivisions 100 followBonds on', 'clip yon -500', '~label', 'set subdivision 1', 'set depth_cue', 'set dc_color black', 'set dc_start 0.5', 'set dc_end 1', 'scale 0.8']
    

    will return the default image (other commands can be passed to modified the final image/movie).

write_cmm(directory, color=<function color_residues at 0x38945f0>, rndname=True, model_num=None)[source]

Save a model in the cmm format, read by Chimera (http://www.cgl.ucsf.edu/chimera).

Note: If none of model_num, models or cluster parameter are set, ALL the models will be written.

Parameters:
  • directory – location where the file will be written (note: the name of the file will be model_1.cmm if model number is 1)
  • model_num (None) – the number of the model to save
  • rndname (True) – If True, file names will be formatted as: model.RND.cmm, where RND is the random number feed used by IMP to generate the corresponding model. If False, the format will be: model_NUM_RND.cmm where NUM is the rank of the model in terms of objective function value
  • color (color_residues) – either a coloring function like pytadbit.imp.imp_model.color_residues() or a list of (r, g, b) tuples (as long as the number of particles)
write_xyz(directory, model_num=None, get_path=False, rndname=True)[source]

Writes a xyz file containing the 3D coordinates of each particle in the model.

Note: If none of model_num, models or cluster parameter are set, ALL the models will be written.

Parameters:
  • directory – location where the file will be written (note: the file name will be model.1.xyz, if the model number is 1)
  • model_num (None) – the number of the model to save
  • rndname (True) – If True, file names will be formatted as: model.RND.xyz, where RND is the random number feed used by IMP to generate the corresponding model. If False, the format will be: model_NUM_RND.xyz where NUM is the rank of the model in terms of objective function value
  • get_path (False) – whether to return, or not, the full path where the file has been written

IMPoptimizer class

class pytadbit.imp.impoptimizer.IMPoptimizer(experiment, start, end, n_models=500, cutoff=300, n_keep=100, close_bins=1)[source]

This class optimizes a set of paramaters (scale, maxdist, lowfreq and upfreq) in order to maximize the correlation between the models generated by IMP and the input data.

Parameters:
  • experiment – an instance of the class pytadbit.experiment.Experiment
  • start – first bin to model (bin number)
  • end – last bin to model (bin number)
  • n_models (5000) – number of models to generate
  • n_keep (1000) – number of models used in the final analysis (usually the top 20% of the generated models). The models are ranked according to their objective function value (the lower the better)
  • close_bins (1) – number of particles away (i.e. the bin number difference) a particle pair must be in order to be considered as neighbors (e.g. 1 means consecutive particles)
load_from_file(f_name)[source]

Loads the optimized parameters from a file generated with the function: pytadbit.imp.impoptimizer.IMPoptimizer.write_result. This function does not overwrite the parameters that were already loaded or calculated.

Parameters:f_name – file name with the absolute path
plot_2d(axes=('scale', 'maxdist', 'upfreq', 'lowfreq'), show_best=0, skip=None)[source]

A grid of heatmaps representing the result of the optimization.

Parameters:
  • axes (‘scale’,’maxdist’,’upfreq’,’lowfreq’) – list of axes to be represented in the plot. The order will define which parameter will be placed on the x, y, z or w axe.
  • show_best (0) – number of best correlation values to highlight in the plot
  • skip (None) – if passed (as a dictionary), fix a given axe, e.g.: {‘scale’: 0.001, ‘maxdist’: 500}
plot_3d(axes=('scale', 'maxdist', 'upfreq', 'lowfreq'))[source]

A grid of heatmaps representing the result of the optimization.

Parameters:axes (‘scale’,’maxdist’,’upfreq’,’lowfreq’) – tuple of axes to be represented in the plot. The order will define which parameter will be placed on the x, y, z or w axe.

This function calculates the correlation between the models generated by IMP and the input data for the four main IMP parameters (scale, maxdist, lowfreq and upfreq) in the given ranges of values.

Parameters:
  • n_cpus – number of CPUs to use
  • lowfreq_range ((-1,0,0.1)) – range of lowfreq values to be optimized. The last value of the input tuple is the incremental step for the lowfreq values
  • upfreq_range ((0,1,0.1)) – range of upfreq values to be optimized. The last value of the input tuple is the incremental step for the upfreq values
  • maxdist_range ((400,1400,100)) – upper and lower bounds used to search for the optimal maximum experimental distance. The last value of the input tuple is the incremental step for maxdist values
  • scale_range ([0.01]) – upper and lower bounds used to search for the optimal scale parameter (nm per nucleotide). The last value of the input tuple is the incremental step for scale parameter values
  • verbose (True) – print the results to the standard output
write_result(f_name)[source]

This function writes a log file of all the values tested for each parameter, and the resulting correlation value.

This file can be used to load or merge data a posteriori using the function pytadbit.imp.impoptimizer.IMPoptimizer.load_from_file

Parameters:f_name – file name with the absolute path

StructuralModels class

class pytadbit.imp.structuralmodels.StructuralModels(nloci, models, bad_models, resolution, original_data=None, zscores=None, clusters=None, config=None)[source]

This function generates three-dimensional models starting from Hi-C data. The final analysis will be performed on the n_keep top models.

Parameters:
  • nloci – number of particles in the selected region
  • models – a dictionary containing the generated pytadbit.imp.impmodel.IMPmodel to be used as ‘best models’
  • bad_models – a dictionary of pytadbit.imp.impmodel.IMPmodel, these model will not be used, just stored in case the set of ‘best models’ needs to be extended later-on ( pytadbit.imp.structuralmodels.StructuralModels.define_best_models() ).
  • resolution – number of nucleotides per Hi-C bin. This will be the number of nucleotides in each model’s particle.
  • original_data (None) – a list of list (equivalent to a square matrix) of the normalized Hi-C data, used to build this list of models.
  • clusters (None) – a dictionary of type pytadbit.imp.structuralmodels.ClusterOfModels
  • config (None) – a dictionary containing the parameter to be used for the generation of three dimensional models.
average_model(models=None, cluster=None, verbose=False)[source]
Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • verbose (False) – prints the distance of each model to average model (in stderr)
Returns:

the average model of a given group of models (a new and ARTIFICIAL model)

centroid_model(models=None, cluster=None, verbose=False)[source]
Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • verbose (False) – prints the distance of each model to average model (in stderr)
Returns:

the centroid model of a given group of models (the most model representative)

cluster_analysis_dendrogram(n_best_clusters=None, color=False, axe=None, savefig=None)[source]

Representation of the clustering results. The length of the leaves if proportional to the final objective function value of each model. The branch widths are proportional to the number of models in a given cluster (or group of clusters, if it is an internal branch).

Parameters:
  • n_best_clusters (None) – number of clusters to represent (by default all clusters will be shown)
  • color (False) – color the dendrogram based on the significance of the clustering (basically it depends of the internal branch lengths)
cluster_models(fact=0.75, dcutoff=200, var='score', method='mcl', mcl_bin='mcl', tmp_file=None, verbose=True)[source]

This function performs a clustering analysis of the generated models based on structural comparison. The result will be stored in StructuralModels.clusters

Parameters:
  • fact (0.75) – factor to define the percentage of equivalent positions to be considered in the clustering
  • dcutoff (200) – distance threshold (nm) to determine if two particles are in contact
  • var (‘score’) –

    value to return, which can be either (i) ‘drmsd’ (symmetry independent: mirrors will show no differences), or (ii) ‘score’, defined as:

                          dRMSD[i] / max(dRMSD)
    score[i] = eqvs[i] * -----------------------
                           RMSD[i] / max(RMSD)

    where eqvs[i] is the number of equivalent position for the ith pairwise model comparison

  • method (‘mcl’) – clustering method to use, which can be either ‘mcl’ or ‘ward’. The last one uses scipy implementation, and is NOT RECOMMENDED.
  • mcl_bin (‘mcl’) – path to the mcl executable file, in case of the ‘mcl is not in the PATH’ warning message
  • tmp_file (None) – path to a temporary file created during the clustering computation. Default will be created in /tmp/ folder
  • verbose (True) – same as print StructuralModels.clusters
contact_map(models=None, cluster=None, cutoff=150, axe=None, savefig=None, savedata=None)[source]

Plots a contact map representing the frequency of interaction (defined by a distance cutoff) between two particles.

Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • cutoff (150) – distance cutoff (nm) to define whether two particles are in contact or not
  • axe (None) – a matplotlib.axes.Axes object to define the plot appearance
  • savefig (None) – path to a file where to save the image generated; if None, the image will be shown using matplotlib GUI
  • savedata (None) – path to a file where to save the contact map data generated, in three columns format (particle1, particle2, percentage of models where these two particles are in contact)
correlate_with_real_data(models=None, cluster=None, cutoff=200, plot=False, axe=None, savefig=None)[source]
Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • cutoff (200) – distance cutoff (nm) to define whether two particles are in contact or not
  • savefig (None) – path to a file where to save the image generated; if None, the image will be shown using matplotlib GUI
  • plot (False) – to display the plot
  • axe (None) – a matplotlib.axes.Axes object to define the plot appearance
Returns:

correlation coefficient rho, between the two matrices. A rho value greater than 0.7 indicates a very good correlation

define_best_models(nbest)[source]

Defines the number of top models (based on the objective function) to keep. If keep_all is set to True in pytadbit.imp.imp_model.generate_3d_models() or in pytadbit.experiment.Experiment.model_region(), then the full set of models (n_models parameter) will be used, otherwise only the n_keep models will be available.

Parameters:nbest – number of top models to keep (usually 20% of the generated models).
density_plot(models=None, cluster=None, steps=(1, 2, 3, 4, 5), error=False, axe=None, savefig=None, savedata=None)[source]

Plots the number of nucleotides per nm of chromatin vs the modeled region bins.

Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • 2, 3, 4, 5) steps ((1,) – how many particles to group for the estimation. By default 5 curves are drawn
  • error (False) – represent the error of the estimates
  • savedata (None) – path to a file where to save the density data generated (1 column per step + 1 for particle number).
dihedral_angle(parta, partb, partc, partd, models=None, cluster=None)[source]
fetch_model_by_rand_init(rand_init, all_models=False)[source]

Models are stored according to their objective function value (first best), but in order to reproduce a model, we need its initial random number. This method helps to fetch the model corresponding to a given initial random number stored under StructuralModels.models[N][‘rand_init’].

Parameters:
  • rand_init – the wanted rand_init number.
  • all_models (False) – whether or not to use ‘bad’ models
Returns:

index of 3d model

get_contact_matrix(models=None, cluster=None, cutoff=150)[source]

Returns a matrix with the number of interactions observed below a given cutoff distance.

Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • cutoff (150) – distance cutoff (nm) to define whether two particles are in contact or not
Returns:

matrix frequency of interaction

measure_angle_3_particles(parta, partb, partc, models=None, cluster=None, radian=False, all_angles=False)[source]

Given three particles A, B and C, the angle g (angle ACB, shown below):

      A
     /|
    /i|
  c/  |
  /   |
 /    |
B )g  |b
 \    |
  \   |
  a\  |
    \h|
     \|
      C

is given by the theorem of Al-Kashi:

\[b^2 = a^2 + c^2 - 2ac\cos(g)\]
Parameters:
  • parta – A particle number
  • partb – A particle number
  • partc – A particle number
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • radian (False) – if True, return value in radians (in degrees otherwise)
Returns:

an angle, either in degrees or radians. If all_angles is true returns a list of the angle g, h, i (see picture above)

median_3d_dist(part1, part2, models=None, cluster=None, plot=True, median=True, axe=None, savefig=None)[source]

Computes the median distance between two particles over a set of models

Parameters:
  • part1 – number corresponding to the first particle
  • part2 – number corresponding to the second particle
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • plot (True) – if True, display a histogram and a box-plot of the distribution of the calculated distances. If False, return either the full list of the calculated distances or their median value
  • median (True) – return either the full list of the calculated distances (False) or their median value (True), when ‘plot’ is set to False
Returns:

if ‘plot’ is False, return either the full list of the calculated distances or their median value distances, either the list of distances.

model_consistency(cutoffs=(50, 100, 150, 200), models=None, cluster=None, axe=None, savefig=None, savedata=None)[source]

Plots the particle consistency, over a given set of models, vs the modeled region bins. The consistency is a measure of the variability (or stability) of the modeled region (the higher the consistency value, the higher stability).

Parameters:
  • cutoffs ((50,100,150,200)) – list of distance cutoffs (nm) used to compute the consistency. Two particle are considered consistent if their distance is less than the given cutoff
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • tmp_path (‘/tmp/tmp_cons’) – location of the input files for TM-score program
  • tmsc (‘’) – path to the TMscore_consistency script (assumed to be installed by default)
  • savedata (None) – path to a file where to save the consistency data generated (1 column per cutoff + 1 for particle number).
objective_function_model(model, log=False, smooth=True, axe=None, savefig=None)[source]

This function plots the objective function value per each Monte-Carlo step

Parameters:
  • model – the number of the model to plot
  • log (False) – log plot
  • smooth (True) – curve smoothing
particle_coordinates(part, models=None, cluster=None)[source]
save_models(path_f)[source]

Saves all the models in pickle format (python object written to disk).

Parameters:path_f – path where to save the pickle file
view_model(models=None, cluster=None, tool='chimera', savefig=None, cmd=None)[source]

Visualize a selected model in the three dimensions.

Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • tool (‘chimera’) – path to the external tool used to visualize the model
  • savefig (None) – path to a file where to save the image OR movie generated (depending on the extension; accepted formats are png, mov and webm). if set to None, the image or movie will be shown using the default GUI.
  • cmd (None) –

    list of commands to be passed to the viewer. The chimera list is:

    focus
    set bg_color white
    windowsize 800 600
    bonddisplay never #0
    represent wire
    shape tube #0 radius 5 bandLength 100 segmentSubdivisions 1 followBonds on
    clip yon -500
    ~label
    set subdivision 1
    set depth_cue
    set dc_color black
    set dc_start 0.5
    set dc_end 1
    scale 0.8

    Followed by the movie command to record movies:

    movie record supersample 1
    turn y 3 120
    wait 120
    movie stop
    movie encode output SAVEFIG

    Or the copy command for images:

    copy file SAVEFIG png

    Passing as the following list as ‘cmd’ parameter:

    cmd = ['focus', 'set bg_color white', 'windowsize 800 600', 'bonddisplay never #0', 'shape tube #0 radius 10 bandLength 200 segmentSubdivisions 100 followBonds on', 'clip yon -500', '~label', 'set subdivision 1', 'set depth_cue', 'set dc_color black', 'set dc_start 0.5', 'set dc_end 1', 'scale 0.8']
    

    will return the default image (other commands can be passed to modified the final image/movie).

walking_angle(models=None, cluster=None, steps=(1, 3), signed=True, savefig=None, savedata=None, axe=None)[source]

Plots the angle between successive loci in a given model or set of models. In order to limit the noise of the measure angle is calculated between 3 loci, between each are two other loci. E.g. in the scheme bellow, angle are calculated between loci A, D and G.

Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • 3) steps ((1,) – how many particles to group for the estimation. By default 2 curves are drawn
  • signed (True) – whether to compute the sign of the angle according to a normal plane, or not.
  • savefig (None) – path to a file where to save the image generated; if None, the image will be shown using matplotlib GUI
  • savedata (None) – path to a file where to save the angle data generated (1 column per step + 1 for particle number).
                      C..........D
                   ...            ...
                ...                 ...
             ...                       ...
  A..........B                            .E
 ..                                        .
.                                          .
                                           .                 .
                                           .                .
                                           F...............G
walking_dihedral(models=None, cluster=None, steps=(1, 3), plot=True, savefig=None, axe=None)[source]

Plots the dihedral angle between successive plans. A plan is formed by 3 successive loci.

Parameters:
  • models (None) – if None (default) the contact map will be computed using all the models. A list of numbers corresponding to a given set of models can be passed
  • cluster (None) – compute the contact map only for the models in the cluster number ‘cluster’
  • 3) steps ((1,) – how many particles to group for the estimation. By default 2 curves are drawn
  • signed (True) – whether to compute the sign of the angle according to a normal plane, or not.
  • savefig (None) – path to a file where to save the image generated; if None, the image will be shown using matplotlib GUI
  • savedata (None) – path to a file where to save the angle data generated (1 column per step + 1 for particle number).
                      C..........D
                   ...            ...
                ...                 ...
             ...                       ...
  A..........B                            .E
 ..                                        .
.                                          .
                                           .                 .
                                           .                .
                                           F...............G
write_cmm(directory, model_num=None, models=None, cluster=None, color=<function color_residues at 0x38945f0>, rndname=True)[source]

Save a model in the cmm format, read by Chimera (http://www.cgl.ucsf.edu/chimera).

Note: If none of model_num, models or cluster parameter are set, ALL the models will be written.

Parameters:
  • directory – location where the file will be written (note: the name of the file will be model_1.cmm if model number is 1)
  • model_num (None) – the number of the model to save
  • models (None) – a list of numbers corresponding to a given set of models to save
  • cluster (None) – save the models in the cluster number ‘cluster’
  • rndname (True) – If True, file names will be formatted as: model.RND.cmm, where RND is the random number feed used by IMP to generate the corresponding model. If False, the format will be: model_NUM_RND.cmm where NUM is the rank of the model in terms of objective function value
  • color (color_residues) – either a coloring function like pytadbit.imp.imp_model.color_residues() or a list of (r, g, b) tuples (as long as the number of particles)
write_xyz(directory, model_num=None, models=None, cluster=None, get_path=False, rndname=True)[source]

Writes a xyz file containing the 3D coordinates of each particle in the model.

Note: If none of model_num, models or cluster parameter are set, ALL the models will be written.

Parameters:
  • directory – location where the file will be written (note: the file name will be model.1.xyz, if the model number is 1)
  • model_num (None) – the number of the model to save
  • models (None) – a list of numbers corresponding to a given set of models to be written
  • cluster (None) – save the models in the cluster number ‘cluster’
  • rndname (True) – If True, file names will be formatted as: model.RND.xyz, where RND is the random number feed used by IMP to generate the corresponding model. If False, the format will be: model_NUM_RND.xyz where NUM is the rank of the model in terms of objective function value
  • get_path (False) – whether to return, or not, the full path where the file has been written
zscore_plot(axe=None, savefig=None)[source]

Generate 3 plots. Two heatmaps of the Z-scores used for modeling, one of which is binary showing in red Z-scores higher than upper cut-off; and in blue Z-scores lower than lower cut-off. Last plot is an histogram of the distribution of Z-scores, showing selected regions.

Parameters:
  • axe (None) – a matplotlib.axes.Axes object to define the plot appearance
  • savefig (None) – path to a file where to save the image generated; if None, the image will be shown using matplotlib GUI