TEMPy.maps.em_map

class TEMPy.maps.em_map.Map(fullMap, origin, apix, filename, header=[], ext_header=[])[source]

A class representing information from a density map file.

Parameters:
  • fullMap – 3D Numpy array containing the EM density

  • origin – Origin of the EM map, in [x, y, z] format.

  • apix – Pixel size (Å) of the EM map.

  • filename – Filename of the parsed map.

  • header – Python array, containing the header information, ordered based on the mrc header format.

  • ext_header

    Information present in the extended header, as per the mrc header format , if present.

Returns:

Map Instance.

Variables:
  • Map.fullMap – 3D Numpy array containing the EM density

  • Map.origin – Origin of the EM map, in [x, y, z] format.

  • Map.apix – Pixel size (Å) of the EM map.

  • Map.filename – Filename of the parsed map.

  • Map.header

    Python array, containing the header information, ordered based on the mrc header format.

  • Map.ext_header

    Information present in the extended header, as per the mrc header format , if present.

box_size()[source]
Returns:

Size of the map array, in Numpy (Z, Y, X) format.

centre()[source]
Returns:

Vector instance of the centre of the map in Angstroms.

change_origin(x_origin, y_origin, z_origin)[source]

Change the origin of the map.

Parameters:
  • x_origin – New origin co-ordinate.

  • y_origin – New origin co-ordinate.

  • z_origin – New origin co-ordinate.

Returns:

Map instance with new origin

copy()[source]
Returns:

copy of the Map Instance.

fourier_transform()[source]

Apply Fourier transform on the density map.

The zero frequency component is in the centre of the spectrum.

Returns:

Fourier Transformed Map instance.

getMap()[source]
Returns:

Numpy array containing the map density data.

get_com()[source]
Returns:

the centre of mass as a Vector instance

get_cropped_box(centroid, new_box_size)[source]

Returns a new Map object with specified shape centred an arbitrary point.

Note: centroid and new_box_size should follow numpy indexing rules, i.e. follow [z, y, x] convention.

Parameters:
  • centroid – List, tuple or np.ndarray with shape (3,) defining the central point in the Map to crop around. Must be (i, j, k) coordinates, i.e. an index into the Map.fullMap ndarray.

  • new_box_size – List, tuple or np.ndarray with shape (3,) defining the shape for the new Map object

Returns:

Map Instance for the cropped out density.

get_cropped_box_around_atom(atom, new_box_size)[source]

Returns a new Map object with specified shape centred on the CA of a specific atom.

Parameters:
  • atom – An Atom (or subclass) instance.

  • new_box_size – List, tuple or np.ndarray with shape (3,) defining the shape for the new Map object

Returns:

Map Instance for the cropped out density.

get_line_between_points(point1, point2)[source]

Return an array of float values representing a line of density values between two points on the map.

Parameters:
  • point1 – Vector instances defining one of the line to extract density values.

  • point2 – Vector instances defining one of the line to extract density values.

Returns:

Array of density values along the defined line.

Return type:

Numpy array

get_normal_vector(x_pos, y_pos, z_pos)[source]

Calculate the normal vector at the point specified.

Point calculated using 3SOM algorithm used by Ceulemans H. & Russell R.B. (2004).

Parameters:
  • x_pos – Voxel in map on which to calculate normal vector.

  • y_pos – Voxel in map on which to calculate normal vector.

  • z_pos – Voxel in map on which to calculate normal vector.

Returns:

The Normal vector at the point specified.

Return type:

Vector instance

get_point_map(min_thr, percentage=0.2)[source]

Calculates the amount of points to use for the Normal vector score and Chamfer distance score.

Minimum number of returned points is 100.

Parameters:
  • min_thr – Minimum threshold, recommended to get value by running the get_primary_boundary on target map.

  • percentage – Percentage of the protein volume.

Returns:

Number of points for Chamfer and Normal Vector scores.

get_pos(minDens, maxDens)[source]

Find the index for all voxels in the Map whose density values fall within the specified thresholds.

Parameters:
  • minDens – Minimum density threshold.

  • maxDens – Maximum density threshold.

Returns:

An array of 3-tuples (indices of the voxels in x,y,z format)

get_primary_boundary(molWeight, low, high, vol_factor=1.21)[source]

Estimates a threshold value to create an envelope that encapsulates a volume expected for a protein with a certain molecular size.

The estimation uses a recursive algorithm, so some parameters can be set as approximate guesses.

NOTE: when used to calculated the NV score this value should be calculated for the experimental map.

Parameters:
  • molWeight – molecular weight of all molecules in model. Can be calculated using get_prot_mass_from_atoms.

  • low – Minimum threshold to consider. Can be a crude first guess e.g. Map.min <TEMPy.maps.em_map.Map.min()

  • high – Minimum threshold to consider. Can be a crude first guess e.g. Map.min <TEMPy.maps.em_map.Map.min()

  • vol_factor – Approximate value for cubic Angstroms per Dalton for globular proteins. Value used in Chimera (Petterson et al, 2004), taken from Harpaz 1994, is 1.21. Other recommended volume factor are 1.5 (1.1-1.9) cubic Angstroms per Dalton in EMAN Volume/mass conversions assume a density of 1.35 g/ml (0.81 Da/A3) (~1.23A3/Da)

Returns:

Primary boundary density value

Return type:

float

get_second_boundary(primary_boundary, noOfPoints, low, high, err_percent=1)[source]

Calculate the second bound density value. For a given boundary value, it calculates the second bound density value such that a specified number of points whose density values fall between the defined boundaries. Uses recursive algorithm.

Parameters:
  • *primary_boundary* – primary threshold, normally given by get_primary_boundary method based on protein molecular weight.

  • *noOfPoints* – Number of points to use in the normal vector score - try first with 10% (maybe 5%better) of the number of points in the map ( round((self.map_size())*0.1)

  • *low – minimum and maximum values between which the threshold will be taken. low should be equal to the value returned by the get_primary_boundary() method and high is the maximum density values in map.

  • high* – minimum and maximum values between which the threshold will be taken. low should be equal to the value returned by the get_primary_boundary() method and high is the maximum density values in map.

  • *err_percent* – default value of 1. Allowed to find a secondary boundary that includes a 1% error.

Returns:

secondary boundary density value

get_significant_points()[source]

Retrieve all points with a density greater than one standard deviation above the mean.

Returns:

An array of 4-tuple (indices of the voxels in z, y, x format and density value)

Return type:

Numpy array

get_vectors()[source]

Retrieve all non-zero density points as Vector instances.

Returns:

An array of 2-tuple ( Vector Instance with respect to origin, and density value)

Return type:

Numpy array

laplace_filtered()[source]

Laplacian filter Map.

Returns:

Laplacian filtered Map instance

makeKDTree(minDens, maxDens)[source]

Make a k-dimensional tree for points in the Map within specified density thresholds.

The KD-Tree can be used for quick nearest neighbour look-ups.

Parameters:
  • minDens – Minimum density value to include in k-dimensional tree.

  • maxDens – Maximum density value to include in k-dimensional tree.

Returns:

Scipy KDTree containing all relevant points.

make_bin_map(cutoff)[source]

Make a binarised version of Map, where values above cutoff = 1 and below cutoff = 0.

Parameters:

cutoff – Cutoff density value

Returns:

Binarised Map instance

map_rotate_by_axis_angle(x, y, z, angle, CoM, rad=False)[source]

Rotate Map instance around pivot point.

Parameters:
  • angle – Angle (in radians if rad == True, else in degrees) of rotation

  • x – Axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.

  • y – Axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.

  • z – Axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.

  • CoM – Point around which map will be rotated.

Returns:

Rotated Map instance

map_size()[source]
Returns:

Number of voxels in the fullMap array.

matrix_transform(mat, x=0, y=0, z=0)[source]

“ Apply affine transform to the map.

Parameters:
  • mat – Affine 3x3 transformation matrix

  • shape – New box dimensions

  • x – Translation in angstroms in respective plane.

  • y – Translation in angstroms in respective plane.

  • z – Translation in angstroms in respective plane.

Returns:

Transformed Map instance

max()[source]
Returns:

Maximum density value of the map.

mean()[source]
Returns:

Mean density value of map.

median()[source]
Returns:

Median density value of map.

min()[source]
Returns:

Minimum density value of the map.

normalise()[source]

Returns: 0-mean normalised Map instance.

normalise_to_1_minus1(in_place=False)[source]
old_write_to_MRC_file(mrcfilename, imod=False)[source]

Write out Map instance as an MRC file

Old implementation that uses Python write(filename, ‘wb’) to write binary data. Extended header not saved using this implementation.

Parameters:

mrcfilename – Filename for the output mrc file

origin_change_maps(MapRef)[source]

Translate Map such that the origin is moved to the origin of a reference Map.

Parameters:

MapRef – Reference Map

Returns:

Translated Map instance

pad_map(nx, ny, nz)[source]

Pad a map (in place) by specified number of voxels along each dimension.

Parameters:
  • nx – Number of voxels to pad in each dimension.

  • ny – Number of voxels to pad in each dimension.

  • nz – Number of voxels to pad in each dimension.

Returns:

Padded Map instance

pixel_centre()[source]
Returns:

Vector instance of the centre of the map in pixels.

represent_normal_vectors(min_threshold, max_threshold)[source]

Create a Structure instance representing normal vectors of density points specified.

Parameters:
  • min_threshold – Minimum/maximum values to include in normal vector representation.

  • max_threshold – Minimum/maximum values to include in normal vector representation.

Returns:

Structure Instance

resample_by_apix(new_apix)[source]

Resample the map to a new pixel size.

Parameters:

new_apix – New Angstrom per pixel sampling

Returns:

Resampled Map instance

resample_by_box_size(new_box_size)[source]

Resample the map based on new box size.

Parameters:

new_box_size – New box dimensions in (Z, Y, X) format

Returns:

Resampled Map instance

resize_map(new_size)[source]

Resize Map instance by cropping/padding (with zeros) the box.

Parameters:

new_size – 3-tuple (x,y,z) giving the box size.

Returns:

Map instance with new box size.

rotate_by_axis_angle(x, y, z, angle, CoM, rad=False)[source]

Rotate map around pivot, given by CoM, using Euler angles.

Parameters:
  • x – Axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.

  • y – Axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.

  • z – Axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.

  • angle – Angle (in radians if rad == True, else in degrees) of rotation.

  • CoM – Centre of mass around which map will be rotated.

Returns:

Rotated Map instance

rotate_by_euler(x, y, z, CoM, rad=False)[source]

Rotate map around pivot, given by CoM, using Euler angles.

Parameters:
  • x – Euler angles (in radians if rad == True, else in degrees) used to rotate map.

  • y – Euler angles (in radians if rad == True, else in degrees) used to rotate map.

  • z – Euler angles (in radians if rad == True, else in degrees) used to rotate map.

  • CoM – Centre of mass around which map will be rotated.

  • *x – translation in angstroms.

  • y – translation in angstroms.

  • z* – translation in angstroms.

Returns:

Rotated Map instance

rotate_by_matrix(mat, CoM)[source]

Rotate map around pivot, given by CoM, using a rotation matrix

Parameters:
  • mat – 3x3 matrix used to rotate map.

  • CoM – Rotation pivot, usually the centre of mass around which map will be rotated.

Returns:

Rotated Map instance

scale_map(scaling)[source]

Scaling Map by scaling factor

Returns:

new Map instance

shift_origin(x_shift, y_shift, z_shift)[source]

Shift the Map origin.

Parameters:
  • x_shift – Shift (in voxels) applied to origin in respective plane.

  • y_shift – Shift (in voxels) applied to origin in respective plane.

  • z_shift – Shift (in voxels) applied to origin in respective plane.

Returns:

Map instance with new origin.

std()[source]
Returns:

Standard deviation of density values in map.

threshold_map(minDens, maxDens)[source]

Create a Map instance containing only density values between the specified min and max values.

Parameters:
  • minDens – Minimum density threshold

  • maxDens – Maximum density threshold

Returns:

Thresholded Map instance

translate(x, y, z)[source]

Translate the map.

Parameters:
  • x – Translation (in A) applied to the map in respective plane.

  • y – Translation (in A) applied to the map in respective plane.

  • z – Translation (in A) applied to the map in respective plane.

Returns:

Shifted Map instance

update_header(binarise=False)[source]

Updates all entries in Map.header.

Parameters:

binarse – If True, returns binary version of map header data. Used in old MRC writing function.

Returns:

Binary version of Map.header if binarise == True

vectorise_point(x, y, z)[source]

Get a Vector instance for a specific voxel coordinate (index) of the map.

Parameters:
  • x – Co-ordinate of the density point to be vectorised.

  • y – Co-ordinate of the density point to be vectorised.

  • z – Co-ordinate of the density point to be vectorised.

Returns:

Vector instance for the point, in angstrom units, relative to the origin.

write_to_MRC_file(mrcfilename, overwrite=True)[source]

Write out Map instance as an MRC file

Uses mrcfile library for file writing.

Parameters:

mrcfilename – Filename for the output mrc file

x_origin()[source]
Returns:

x-coordinate of the origin.

x_size()[source]
Returns:

Size of the map array in x direction.

y_origin()[source]
Returns:

y-coordinate of the origin.

y_size()[source]
Returns:

Size of the map array in y direction.

z_origin()[source]
Returns:

z-coordinate of the origin.

z_size()[source]
Returns:

Size of the map array in z direction.