iotbx.managers

Overview

API documentation

iotbx.map_manager.add_tuples_int(t1, t2)
iotbx.map_manager.dummy_map_manager(crystal_symmetry, n_grid=12)

Make a map manager with crystal symmetry and unit sized map

iotbx.map_manager.get_indices_from_index(index=None, all=None)
iotbx.map_manager.get_sites_cart_from_index(indices_list=None, points=None, map_data=None, crystal_symmetry=None, all=None)

Get sites_cart from linear (1d) map indices. Supply either map_data or all to provide n_real crystal_symmetry is required Supply either 3D indices (i,j,k) or 1-D indices (points)

class iotbx.map_manager.map_manager(file_name=None, map_data=None, unit_cell_grid=None, unit_cell_crystal_symmetry=None, origin_shift_grid_units=None, ncs_object=None, wrapping=<libtbx.AutoType object>, experiment_type=<libtbx.AutoType object>, scattering_table=<libtbx.AutoType object>, resolution=<libtbx.AutoType object>, log=None)

map_manager, includes map_reader and write_ccp4_map

This class is intended to be the principal mechanism for reading and writing map information. It is intended to be used by the iotbx.data_manager for both of these purposes.

Use map_manager to read, write, and carry information about one map. Map_manager keeps track of the origin shifts and also the original full unit cell and cell dimensions. It writes out the map in the same place as it was read in.

Note on wrapping: Wrapping means that the map value outside the map boundaries can be obtained as the value inside the boundaries, (translated by some multiple of the unit cell translations.) Normally crystallographic maps can be wrapped and cryo EM maps cannot.

Wrapping should be specified on initialization if not read from a file. If read from a file, the value from the file labels is used if available, otherwise it is assumed to be wrapping = False unless specified (normal for a cryo-EM map. If not specified at all, it will need to be specified before a map_model_manager is created or the map_manager is written out.

Map_manager also keeps track of any changes in magnification. These are reflected in changes in unit_cell and crystal_symmetry cell dimensions and angles.

You normally create a new map_manager by initializing map_manager with a file name. Then you apply the shift_origin() method and the map is shifted to place the origin at (0, 0, 0) and the original origin is recorded as self.origin_shift_grid_units.

NOTE: do not set self.origin_shift_grid_units directly; instead use the method set_original_origin_and_gridding() or supply the value when initializing the map_manager.

You can also create a map_manager with a map_data object (3D flex.double() array) along with the meta-data below.

NOTE: MRC Maps may not represent the entire unit cell. Normally maps that

have an origin (corner with minimum i, j, k) that is not zero will be shifted at a later stage to have the origin at (0, 0, 0), along with any models and ncs objects (typically done with iotbx.map_and_model). To be able to write out a map in the same place as it was read in after shifting the origin and/or boxing the map, you need to keep track of 3 things. These are: 1. unit_cell_grid: grid representing one full unit cell as read in.

Saved in map_manager as self.unit_cell_grid

  1. unit_cell_crystal_symmetry: dimensions and space group of full unit cell

    Saved in map_manager as self._unit_cell_crystal_symmetry

  2. origin_shift_grid_units: the shift in grid units to apply to the working map to superimpose it on the original map. When you read the map in this is (0, 0, 0). If you shift the map origin from (i, j, k) to (0, 0, 0) then the origin_shift_grid_units is (i, j, k).

    Saved in map_manager as self.origin_shift_grid_units

Magnification (pixel size scaling) of a map: there is no general parameter describing magnification of an MRC map. Changes in scaling are recorded in map_manager as changes in the scaling matrix/translation that relates grid points in a map to real-space position.

Normal usage (NOTE: read/write should normally be done through data_manager):

Read in a map:

mm = map_manager(‘input_map.mrc’)

Summarize:

mm.show_summary()

Normally shift origin of map to (0, 0, 0) (you can do this here

or you can use iotbx.map_and_model to shift models and maps together):

mm.shift_origin()

Get the map_data (shifted if origin was shifted above):

map_data = mm.map_data()

Get the crystal_symmetry of the box of data that is present:

cs = mm.crystal_symmetry()

Get the crystal_symmetry of the whole unit cell (even if not present):

unit_cell_cs = mm.unit_cell_crystal_symmetry()

Write out the map in map_data() in original location:

mm.write_map(file_name = ‘output_map.ccp4’)

——– CONVENTIONS ————– See http://www.ccpem.ac.uk/mrc_format/mrc2014.php for MRC format See https://pypi.org/project/mrcfile/ for mrcfile library documentation

Same conventions as iotbx.ccp4_map

Default is to write maps with INTERNAL_STANDARD_ORDER of axes of [3, 2, 1]

corresponding to columns in Z, rows in Y, sections in X to match flex array layout. This can be modified by changing the values in output_axis_order.

Hard-wired to convert input maps of any order to

INTERNAL_STANDARD_ORDER = [3, 2, 1] before conversion to flex arrays This is not modifiable.

INTERNAL_STANDARD_ORDER = [3, 2, 1]

Standard limitations and associated message. These can be checked with: limitations = mrc.get_limitations()

which returns a group_args object with a list of limitations and a list of corresponding error messages, or None if there are none see phenix.show_map_info for example

These limitations can also be accessed with specific calls placed below:

For example: mrc.can_be_sharpened() returns False if “extract_unique” is present

Map labels that are not limitations can be accessed with:

additional_labels = mrc.get_additional_labels()

STANDARD_LIMITATIONS_DICT = {
“extract_unique”:

“This map is masked around unique region and not suitable for auto-sharpening.”,

“map_is_sharpened”:

“This map is auto-sharpened and not suitable for further auto-sharpening.”,

“map_is_density_modified”: “This map has been density modified.”,

}

NOTES ON ORDER OF AXES

Phenix standard order is 3, 2, 1 (columns Z, rows Y, sections in X).

Convert everything to this order.

This is the order that allows direct conversion of a numpy 3D array

with axis order (mapc, mapr, maps) to a flex array.

For reverse = True, supply order that converts flex array to numpy 3D array

with order (mapc, mapr, maps)

Note that this does not mean input or output maps have to be in this order.

It just means that before conversion of numpy to flex or vice-versa the array has to be in this order.

Note that MRC standard order for input/ouput is 1, 2, 3.

NOTE: numpy arrays indexed from 0 so this is equivalent to

order of 2, 1, 0 in the numpy array

NOTE: MRC format allows data axes to be swapped using the header

mapc mapr and maps fields. However the mrcfile library does not attempt to swap the axes and assigns the columns to X, rows to Y and sections to Z. The data array is indexed C-style, so data values can be accessed using mrc.data[z][y][x].

NOTE: normal expectation is that phenix will read/write with the

order 3, 2, 1. This means X-sections (index = 3), Y rows (index = 2), Z columns (index = 1). This correxponds to

mapc (columns) = 3 or Z mapr (rows) = 2 or Y maps (sections) = 1 or X

In the numpy array (2, 1, 0 instead of 3, 2, 1):

To transpose, specify i0, i1, i2 where:

i0 = 2 means input axis 0 becomes output axis 2 NOTE: axes are 0, 1, 2 etc, not 1, 2, 3

Examples:

np.transpose(a, (0, 1, 2)) does nothing np.transpose(a, (1, 2, 0)) output axis 0 is input axis 1

We want output axes to always be 2, 1, 0 and input axes for numpy array are

(mapc-1, mapr-1, maps-1):

For example, in typical phenix usage, the transposition is:

i_mapc = 3 i_mapc_np = 2 i_mapr = 2 i_mapr_np = 1 i_maps = 1 i_maps_np = 0

——– END CONVENTIONS ————–

absolute_center_cart(use_assumed_end=False, place_on_grid_point=False, use_unit_cell_grid=False)

Return center of map (absolute position) in Cartesian coordinates A little tricky because for example the map goes from 0 to nx-1, not nx

If use_assumed_end, go to nx

Also map could start at non-zero origin If place_on_grid_point then guess the end by whether the center ends

on a grid point

If use_unit_cell_grid just find center of full unit cell

add_label(label=None, verbose=False)

Add a label (up to 80-character string) to write to output map. Default is to specify the program name and date

add_limitation(limitation=None)

Add a limitation from STANDARD_LIMITATIONS_DICT Supply the key (such as “map_is_sharpened”)

apply_mask(set_outside_to_mean_inside=False)

Replace map_data with masked version based on current mask Just uses method in create_mask_around_atoms

apply_spectral_scaling(d_min=None, d_max=None, n_bins=100)
as_full_size_map()

Create a full-size map with the current map inside it, padded by zero

A little tricky because the starting map is going to have its origin at (0, 0, 0) but the map we are creating will have that point at self.origin_shift_grid_units.

First use box.with_bounds to create map from -self.origin_shift_grid_units

to -self.origin_shift_grid_units+self.unit_cell_grid-(1, 1, 1). Then

shift that map to place origin at (0, 0, 0)

If the map is full size already, return the map as is If the map is bigger than full size stop as this is not suitable

as_map_model_manager()

Return a map_model_manager

binary_filter(threshold=0.5)

Apply a binary filter to the map (value at pixel i,j,k=1 if average of all 27 pixels within 1 of this one is > threshold, otherwise 0) Changes and overwrites contents of this map_manager.

cart_to_grid_units(xyz)
cc_to_other_map_manager(other_map_manager)
create_mask_around_atoms(model, mask_atoms_atom_radius=None, invert_mask=None)

Use cctbx.maptbx.mask.create_mask_around_atoms to create a mask around atoms in model

Does not apply the mask (use apply_mask_to_map etc for that)

mask_atoms_atom_radius default is max(3, resolution)

invert_mask makes outside 1 and inside 0

create_mask_around_density(resolution=None, molecular_mass=None, sequence=None, solvent_content=None)
Use cctbx.maptbx.mask.create_mask_around_density to create a

mask automatically

Does not apply the mask (use apply_mask_to_map etc for that)

Parameters are:
resolutionresolution of map, taken from self.resolution() if not

specified

molecular_mass: optional mass (Da) of object in density sequence: optional sequence of object in density solvent_content : optional solvent_content of map

create_mask_around_edges(boundary_radius=None)

Use cctbx.maptbx.mask.create_mask_around_edges to create a mask around edges of map. Does not make a soft mask. For a soft mask, follow with soft_mask(boundary_radius =boundary_radius) The radius is to define the boundary around the map.

Does not apply the mask (use apply_mask_to_map etc for that)

create_mask_with_map_data(map_data)

Set mask to be map_data

Does not apply the mask (use apply_mask_to_map etc for that)

Uses cctbx.maptbx.mask.create_mask_with_mask_data to do it

Requires origin to be zero of both self and new mask

customized_copy(map_data=None, origin_shift_grid_units=None, use_deep_copy_for_map_data=True, crystal_symmetry_space_group_number=None, wrapping=None)

Return a customized deep_copy of this map_manager, replacing map_data with supplied map_data.

The map_data and any _created_mask will be deep_copied before using them unless use_deep_copy_for_map_data = False

Normally this customized_copy is applied with a map_manager that has already shifted the origin to (0, 0, 0) with shift_origin.

Normally the new map_data will have the same dimensions of the current map_data. Normally new map_data will also have origin at (0, 0, 0).

NOTE: It is permissible for map_data to have different bounds or origin than the current self.map_data. In this case you must specify a new value of origin_shift_grid_units corresponding to this new map_data. This new origin_shift_grid_units specifies the original position in the full unit cell grid of the most-negative corner grid point of the new map_data. The new map_manager will still have the same unit cell dimensions and grid as the original.

NOTE: It is permissible to get a customized copy before shifting the origin. Applying with non-zero origin requires that:

self.origin_shift_grid_units == (0, 0, 0) origin_shift_grid_units = (0, 0, 0) map_data.all() (size in each direction) of current and new maps

are the same.

origins of current and new maps are the same

NOTE: wrapping is normally copied from original map, but if new map is not full size then wrapping is always set to False.

If crystal_symmetry_space_group_number is specified, use it

deep_copy()

Return a deep copy of this map_manager object Uses customized_copy to deepcopy everything including map_data

Origin does not have to be at (0, 0, 0) to apply

delete_mask()
density_at_sites_cart(sites_cart)
Return flex.double list of density values corresponding to sites (cartesian

coordinates in A)

experiment_type()
external_origin_as_grid_units(as_inverse=False)
external_origin_is_compatible_with_gridding()
find_map_symmetry(include_helical_symmetry=False, symmetry_center=None, min_ncs_cc=None, symmetry=None, ncs_object=None, check_crystal_symmetry=True, only_proceed_if_crystal_symmetry=False)

Use run_get_symmetry_from_map tool in segment_and_split_map to find map symmetry and save it as an mmtbx.ncs.ncs.ncs object

Here map symmetry is the reconstruction symmetry used to generate the map. Normally it is essentially perfect symmetry and normally the principal axes are aligned with x,y,z and normally the center is at the original center of the map.

Sets self._warning_message if failure, sets self._ncs_object and

self._ncs_cc if success

This procedure may fail if the above assumptions do not hold. Optional center of map can be supplied, and minimum NCS correlation can also be supplied

Requires that map_manager is already shifted to place origin at (0, 0, 0)

Assumes that center of symmetry is at (1/2, 1/2, 1/2) in the full map

It is optional to include search for helical symmetry. Reason is that this is much slower than other symmetries.

symmetry (symbol such as c1, O, D7) can be supplied and search will be limited to that symmetry

ncs_object can be supplied in which case it is just checked

If check_crystal_symmetry, try to narrow down possibilities by looking for space-group symmetry first

If only_proceed_if_crystal_symmetry, skip looking if nothing comes up

with check_crystal_symmetry

find_n_highest_grid_points_as_sites_cart(n=0, n_tolerance=0, max_tries=100)

Return the n highest grid points in the map as sites_cart

fourier_coefficients_as_map_manager(map_coeffs)
Convert Fourier coefficients into to a real-space map with gridding

matching this existing map_manager. Returns a map_manager object.

Requires that this map_manager has origin at (0, 0, 0) (i.e., shift_origin() has been applied if necessary)

NOTE: Assumes that the map_coeffs are in the same frame of reference as this map_manager (i.e., similar to those that would be written out using map_as_fourier_coefficients).

gaussian_filter(smoothing_radius)

Gaussian blur the map in map_manager with given smoothing radius. Changes and overwrites contents of this map_manager.

get_boxes_to_tile_map(target_for_boxes=24, box_cushion=3, get_unique_set_for_boxes=None, dist_min=None, do_not_go_over_target=None, target_xyz_center_list=None)

Return a group_args object with a list of lower_bounds and upper_bounds corresponding to a set of boxes that tiles the part of the map that is present. The boxes may not be the same size but will tile to exactly cover the existing part of the map. Approximately target_for_boxes will be returned (may be fewer or greater) Also return boxes with cushion of box_cushion If get_unique_set_for_boxes is set, try to use map symmetry to identify

duplicates and set ncs_object

If target_xyz_center_list is set, use these points as centers but try

to use standard box size.

get_density_along_line(start_site=None, end_site=None, n_along_line=10, include_ends=True)

Return group_args object with density values and coordinates along a line segment from start_site to end_site (cartesian coordinates in A) with n_along_line sampling points. Optionally include/exclude ends.

get_mask_as_map_manager()
get_n_real_for_grid_spacing(grid_spacing=None)
grid_units_to_cart(grid_units)

Convert grid units to cartesian coordinates

initialize_map_data(map_value=0)

Set all values of map_data to map_value

invert_hand()

Invert the hand of this map by swapping the order of all sections in Z

is_compatible_model(model, require_match_unit_cell_crystal_symmetry=True, absolute_angle_tolerance=0.01, absolute_length_tolerance=0.01, shift_tol=0.001)

Model is compatible with this map_manager if it is not specified as being different.

They are different if:
  1. original and current symmetries are present and different from each

other and do not match

  1. model current symmetry does not match map original or current

  2. model has a shift_cart (shift applied) different than map shift_cart

NOTE: a True result does not mean that the model crystal_symmetry matches the map crystal_symmetry. It does mean that it is reasonable to set the model crystal_symmetry to match the map ones.

If require_match_unit_cell_crystal_symmetry is True, then they are different if anything is different. If it is False, allow original

symmetries do not have to match

is_compatible_ncs_object(ncs_object, tol=0.001)

ncs_object is compatible with this map_manager if shift_cart is the same as map

is_consistent_with_wrapping(relative_sd_tol=0.01)

Report if this map looks like it is a crystallographic map and can be wrapped

If it is not full size…no wrapping If origin is not at zero…no wrapping If it is not periodic, no wrapping If very small or resolution_factor for map is close to 0.5…cannot tell If has all zeroes (or some other constant on edges) … no wrapping

relative_sd_tol defines how close to constant values at edges must be to qualify as “constant”

Returns True, False, or None (unsure)

is_dummy_map_manager()

Is this a dummy map manager

is_full_size()

Report if map is full unit cell

is_mask()

Is this a mask

is_similar(other=None, absolute_angle_tolerance=0.01, absolute_length_tolerance=0.01)
map_as_fourier_coefficients(d_min=None, d_max=None, box=True, resolution_factor=0.3333333333333333)

Convert a map to Fourier coefficients to a resolution of d_min, if d_min is provided, otherwise box full of map coefficients will be created.

Filter results with low resolution of d_max if provided

NOTE: Fourier coefficients are relative the working origin (not original origin). A map calculated from the Fourier coefficients will superimpose on the working (current map) without origin shifts.

This method and fourier_coefficients_as_map_manager interconvert map_data and map_coefficients without changing origin. Both are intended for use with map_data that has an origin at (0, 0, 0).

The map coefficients are always in space group P1.

map_map_cc(other_map_manager)

Return simple map correlation to other map_manager

minimum_resolution(set_minimum_resolution=True)

Get minimum resolution. If set previously, use that value

ncs_cc()
ncs_object()
origin_is_zero()
randomize(d_min=None, low_resolution_fourier_noise_fraction=0.01, high_resolution_fourier_noise_fraction=2, low_resolution_real_space_noise_fraction=0, high_resolution_real_space_noise_fraction=0, low_resolution_noise_cutoff=None, random_seed=None)

Randomize a map.

Unique aspect of this noise generation is that it can be specified whether the noise is local in real space (every point in a map gets a random value before Fourier filtering), or local in Fourier space (every Fourier coefficient gets a complex random offset). Also the relative contribution of each type of noise vs resolution can be controlled.

Parameters:

d_min: high-resolution limit in Fourier transformations

low_resolution_fourier_noise_fraction (float, 0): Low-res Fourier noise high_resolution_fourier_noise_fraction (float, 0): High-res Fourier noise low_resolution_real_space_noise_fraction(float, 0): Low-res

real-space noise

high_resolution_real_space_noise_fraction (float, 0): High-res

real-space noise

low_resolution_noise_cutoff (float, None): Low resolution where noise

starts to be added

remove_labels()

Remove all labels

resample_on_different_grid(n_real=None, target_grid_spacing=None)

Resample the map on a grid of n_real and return new map_manager If an ncs_object is present, set its shift_cart too

Allows map to be boxed; if so returns new boxed map of similar size and location. If boxed map has cubic gridding, keep that after resampling.

resolution(force=False, method='d99', set_resolution=True)

Get nominal resolution Return existing if present unless force is True choices:

d9: resolution correlated at 0.9 with original d99: resolution correlated at 0.99 with original d999: resolution correlated at 0.999 with original d_min: d_min_from_map

resolution_filter(d_min=None, d_max=None)

High- or low-pass filter the map in map_manager. Changes and overwrites contents of this map_manager. Remove all components with resolution < d_min or > d_max Either d_min or d_max or both can be None. To make a low_pass filter with cutoff at 3 A, set d_min=3 To make a high_pass filter with cutoff at 2 A, set d_max=2

scattering_table()
set_crystal_symmetry_to_p1(space_group_number=1)

Change the working crystal symmetry to P1 This changes map in place Do a deep_copy first if you do not want it changed (Actually can set space group number to anything so you can set it back)

set_experiment_type(experiment_type)

Set the experiment type xray,neutron, or cryo_em If scattering_table is not defined, it is guessed from experiment_type

set_is_mask(value=True)

define if this is a mask

set_log(log=<colorama.ansitowin32.StreamWrapper object>)

Set output log file

set_map_data(map_data=None)

Replace self.data with map_data. The two maps must have same gridding

NOTE: This uses selections to copy all the data in map_data into self.data. The map_data object is not associated with self.data, the data is simply copied. Also as self.data is modified in place, any objects that currently are just pointers to self.data are affected.

set_mean_zero_sd_one()

Function to normalize the map If the map has a constant value, do nothing

set_model_symmetries_and_shift_cart_to_match_map(model)

Set the model original and working crystal_symmetry to match map.

Overwrites any information in model on symmetry and shift_cart Modifies model in place

NOTE: This does not shift the coordinates in model. It is used to fix crystal symmetry and set shift_cart, not to actually shift a model. For shifting a model, use:

model.shift_model_and_set_crystal_symmetry(shift_cart=shift_cart)

set_ncs_object(ncs_object)

set the ncs object for this map_manager. Incoming ncs_object must

be compatible (shift_cart values must match or be defined).

Incoming ncs_object is deep_copied.

set_ncs_object_shift_cart_to_match_map(ncs_object)

Set the ncs_object shift_cart to match map

Overwrites any information in ncs_object on shift_cart Modifies ncs_object in place. Does not return anything

Do not use this to try to shift the ncs object. That is done in the ncs object itself with ncs_object.coordinate_offset(shift_cart), which returns a new ncs object

You can use ncs_object =
self.shift_ncs_object_to_match_map_and_return_new_ncs_object(ncs_object)

to shift the ncs object and set its shift_cart and get a new copy.

set_original_origin_and_gridding(original_origin=None, gridding=None)

Specify the location in the full unit cell grid where the origin of the map that is present is to be placed to match its original position. This is referred to here as the “original” origin, as opposed to the current origin of this map.

Note that this method does not actually shift the origin of the working map. It just defines where that origin is going to be placed when restoring the map to its original position.

Also optionally redefine the definition of the unit cell, keeping the grid spacing the same.

This allows redefining the location of the map that is present within the full unit cell. It also allows redefining the unit cell itself. Only use this to create a new partial map in a defined location.

Previous definition of the location of the map that is present is discarded.

Fundamental parameters set:

self.origin_shift_grid_units: shift to place origin in original location self._unit_cell_crystal_symmetry: dimensions of full unit cell self.unit_cell_grid: grid units of full unit cell

At end, recheck wrapping

set_output_external_origin(value)
set_program_name(program_name=None)

Set name of program doing work on this map_manager for output (string)

set_resolution(resolution)

Set the nominal resolution of map

set_scattering_table(scattering_table)

Set the scattering table (type of scattering) electron: cryo_em n_gaussian x-ray (standard) wk1995: x-ray (alternative) it1992: x-ray (alternative) neutron: neutron scattering

set_unit_cell_crystal_symmetry(crystal_symmetry)

Specify the dimensions and space group of unit cell. This also changes the crystal_symmetry of the part that is present and the grid spacing. Also resets crystal_symmetry of ncs object

Purpose is to redefine the dimensions of the map without changing values of the map. Normally used to correct the dimensions of a map where something was defined incorrectly.

Does not change self.unit_cell_grid or self.origin_shift_grid_units

Does change the result of self.shift_cart(), which is based on

self.origin_shift_grid_units and self.unit_cell_grid and self.crystal_symmetry

Fundamental parameters set:

self._unit_cell_crystal_symmetry: dimensions of full unit cell self._crystal_symmetry: dimensions of part of unit cell that is present

set_wrapping(wrapping_value)

Set wrapping to be wrapping_value

shift_aware_rt(from_obj=None, to_obj=None, working_rt_info=None, absolute_rt_info=None)

Returns shift_aware_rt object

Uses rt_info objects (group_args with members of r, t).

Simplifies keeping track of rotation/translation between two

objects that each may have an offset from absolute coordinates.

absolute rt is rotation/translation when everything is in original,

absolute cartesian coordinates.

working_rt is rotation/translation of anything in “from_obj” object

to anything in “to_obj” object using working coordinates in each.

Usage: shift_aware_rt = self.shift_aware_rt(absolute_rt_info = rt_info) shift_aware_rt = self.shift_aware_rt(working_rt_info = rt_info,

from_obj=from_obj, to_obj = to_obj)

apply RT using working coordinates in objects sites_cart_to_obj = shift_aware_rt.apply_rt(sites_cart_from_obj,

from_obj=from_obj, to_obj=to_obj)

apply RT absolute coordinates sites_cart_to = shift_aware_rt.apply_rt(sites_cart_from)

shift_cart()

Return the shift_cart of this map from its original location.

(the negative of the origin shift ) in cartesian coordinates

shift_model_to_match_map(model)

Move the model to match this map. Note difference from set_model_symmetries_and_shift_cart_to_match_map

which sets model symmetry and shift_cart but does not move the model

shift_ncs_object_to_match_map_and_return_new_ncs_object(ncs_object)

Move the ncs_object to match this map. Also sets ncs_object shift_cart Returns new copy of ncs_obect and does not affect the original

Note difference from set_ncs_object_shift_cart_to_match_map which

sets the shift_cart but does not move the object

shift_origin(desired_origin=(0, 0, 0), apply_external_origin_if_present=True)
Shift the origin of the map to desired_origin

(normally desired_origin is (0, 0, 0), so just update origin_shift_grid_units

Origin is the value of self.map_data().origin() origin_shift_grid_units is the shift to apply to self.map_data() to

superimpose it on the original map.

If you shift the origin by (i, j, k) then add -(i, j, k) to

the current origin_shift_grid_units.

If current origin is at (a, b, c) and

desired origin = (d, e, f) and existing origin_shift_grid_units is (g, h, i):

the shift to make is (d, e, f) - (a, b, c)

the new value of origin_shift_grid_units will be:

(g, h, i)+(a, b, c)-(d, e, f) or new origin_shift_grid_units is: (g, h, i)- shift

the new origin of map_data will be (d, e, f)

Note on external origin: If input map had ORIGIN specified,

so that the value of self.external_origin is not (0,0,0) and not None, and apply_external_origin_if_present is set, then: determine if self.external_origin is on a grid point and if so, convert

and use negative of it as origin. Then set self.external_origin to zero.

Does not apply if origin is already not (0,0,0).

shift_origin_to_match_original()

Shift origin by self.origin_shift_grid_units to place origin in its original location

soft_mask(soft_mask_radius=None)

Make mask a soft mask. Just uses method in cctbx.maptbx.mask Use resolution for soft_mask radius if not specified

trace_atoms_in_map(dist_min, n_atoms)

Utility to find positions where n_atoms atoms separated by dist_min can be placed in density in this map

warning_message()
wrapping()

Report if map can be wrapped

write_map(file_name)

Simple version of write

file_name is output file name map_data is map_data object with 3D values for map. If not supplied,

use self.map_data()

Normally call with file_name (file to be written) Output labels are generated from existing self.labels, self.program_name, and self.limitations

If self.output_external_origin is specified, write that value to file

iotbx.map_manager.remove_site_with_most_neighbors(sites_cart)
iotbx.map_manager.select_n_in_biggest_cluster(sites_cart, dist_min=None, n=None, dist_min_ratio=1.0, dist_min_ratio_min=0.5, minimize_density_of_points=None)

Select n of sites_cart, taking those near biggest cluster if possible If minimize_density_of_points, remove those with the most neighbors

class iotbx.map_manager.shift_aware_rt(from_obj=None, to_obj=None, working_rt_info=None, absolute_rt_info=None)

Class to simplify keeping track of rotation/translation between two objects that each may have an offset from absolute coordinates.

Basic idea: absolute rt is rotation/translation when everything is in original, absolute cartesian coordinates.

working_rt is rotation/translation of anything in “from_obj” object to anything

in “to_obj” object using working coordinates in each.

The from_obj and to objects must have a shift_cart method

absolute_rt_info()
apply_rt(site_cart=None, sites_cart=None, from_obj=None, to_obj=None)

Apply absolute rt if from and to not specified. Apply relative if specified

get_absolute_rt_info(working_rt_info=None, from_obj=None, to_obj=None)

working_rt_info describes how to map from_xyz -> to_xyz in local coordinates from_xyz is shifted from absolute by from.shift_cart() to_xyz is shifted from absolute by to.shift_cart()

We have:

r from_xyz + t = to_xyz in working frame of reference

We want to describe how to map:

(from_xyz - from.shift_cart()) -> (to_xyz - to.shift_cart())

where r is going to be the same and T will be different than t

r ((from_xyz - from.shift_cart()) + T = (to_xyz - to.shift_cart()) T = (to_xyz - to.shift_cart() - r from_xyz + r from.shift_cart()

but: to_xyz - r from_xyz = t

T = t - to.shift_cart() + r from.shift_cart()

Note reverse:

t = T + to.shift_cart() - r from.shift_cart()

inverse()
is_similar(other_shift_aware_rt_info, tol=0.001)
working_rt_info(from_obj=None, to_obj=None)

Get rt in working frame of reference

iotbx.map_manager.subtract_tuples_int(t1, t2)
iotbx.map_model_manager.apply_aniso_b_cart_to_f_array_info(f_array_info, b_iso, d_min, aniso_b_cart)
iotbx.map_model_manager.apply_ncs_to_dv_results(direction_vectors=None, xyz=None, scaling_group_info=None, ncs_object=None)
iotbx.map_model_manager.convert_tlso_group_info_to_lists(tlso_group_info)
iotbx.map_model_manager.create_fine_spacing_array(unit_cell, cell_ratio=10)
iotbx.map_model_manager.create_map_manager_with_value_list(n_real=None, crystal_symmetry=None, value_list=None, sites_cart_list=None, target_spacing=None, max_iterations=None, default_value=None)

Create a map_manager with values set with a set of sites_cart and values Use nearest available value for each grid point, done iteratively

with radii in shells of target_spacing/2 and up to max_iterations shells

If default_value is set, use that for all empty locations after max_iterations

iotbx.map_model_manager.cutoff_values(inside=True)
iotbx.map_model_manager.get_average_scale_factors(scale_factor_info)
iotbx.map_model_manager.get_center_of_box_frac(lower_bounds=None, upper_bounds=None, n_real=None, crystal_symmetry=None)

get center of this box

iotbx.map_model_manager.get_map_coeffs_as_fp_phi(map_coeffs, d_min=None, n_bins=None)

Get map_coeffs as fp and phi. also set up binner if n_bins is not None

iotbx.map_model_manager.get_map_counts(map_data, crystal_symmetry=None)
iotbx.map_model_manager.get_map_histograms(data, n_slots=20, data_1=None, data_2=None)
iotbx.map_model_manager.get_normalization_data_for_unit_binning(f_array)

Get normalizations for each reflection in a scheme for interpolating a top-hat function over bins

iotbx.map_model_manager.get_pointer_to_old_dv_id_dict(working_dv_list=None, dv_list=None, very_similar=0.95, allow_multiple_use=True)

For each member of working_dv_list, identify best match to member of dv_list. Only use each dv_list member once unless allow_multiple_use. ID by abs(dot product) allow_multiple_use is for matching any to dv_list, False is for # rearranging only

iotbx.map_model_manager.get_selection_inside_box(lower_bounds=None, upper_bounds=None, n_real=None, model=None, crystal_symmetry=None)

get selection for all the atoms inside this box

iotbx.map_model_manager.get_selections_and_boxes_to_split_model(model_id=None, map_model_manager=None, selection_method='by_chain', selection_list=None, skip_waters=False, skip_hetero=False, target_for_boxes=24, box_cushion=3, select_final_boxes_based_on_model=None, skip_empty_boxes=True, mask_around_unselected_atoms=None, mask_all_maps_around_edges=None, mask_radius=3, masked_value=-10, get_unique_set_for_boxes=True, mask_id=None, exclude_points_outside_density=None, minimum_boxes_inside_density=25, log=<colorama.ansitowin32.StreamWrapper object>)

Split up model into pieces using selection_method Choices are [‘by_chain’, ‘by_segment’,’all’, ‘boxes’] by_chain: each chain is a selection by_segment: each unbroken part of a chain is a selection boxes: map is split into target_for_boxes boxes, all atoms in

each box selected requires map_model_manager to be present

Skip waters or hetero atoms in selections if specified If select_final_boxes_based_on_model and selection_method == ‘boxes’ then

make the final boxes just go around the selected parts of the model and not tile the map.

If skip_empty_boxes then skip anything with no model. if get_unique_set_for_boxes then get a unique set for ‘boxes’ method If mask_id is set and exclude_points_outside_density , skip boxes

outside of mask for boxes method

If exclude_points_outside_density,

try to add boxes inside density (basically add the proportional number of boxes but put them definitely inside the density instead of evenly spaced.

iotbx.map_model_manager.get_selections_for_segments(model, skip_waters=True, skip_hetero=True, no_water_or_het_with_and=None, skip_n_residues_on_ends=None, minimum_length=1, return_as_group_args=False)

‘ Generate selections corresponding to each segment (chain or part of a chain that is separate from remainder of chain)

Use skip_waters and skip_hetero to specify whether to include them If skip_n_residues_on_ends is set, skip residues within

skip_n_residues_on_ends of an end

iotbx.map_model_manager.get_selections_from_boxes(box_info=None, model=None, overall_selection=None, skip_empty_boxes=None, mask_map_manager=None, exclude_points_outside_density=None)
Generate a list of selections that covers all the atoms in model,

grouped by the boxes defined in box_info

iotbx.map_model_manager.get_skip_waters_and_hetero_lines(skip_waters=True, skip_hetero=True)
iotbx.map_model_manager.get_split_maps_and_models(map_model_manager=None, box_info=None, first_to_use=None, last_to_use=None)

Apply selections and boxing in box_info to generate a set of small map_model_managers

if mask_around_unselected_atoms is set, then mask within each box

around all the atoms that are not selected (including waters/hetero) with a mask_radius of mask_radius and set the value inside the mask to

masked_value

mask_around_unselected_atoms = mask_around_unselected_atoms, mask_radius = mask_radius, masked_value = masked_value,

iotbx.map_model_manager.get_tlso_group_info_from_model(model, nproc=1, log=<colorama.ansitowin32.StreamWrapper object>)

Extract tlso_group_info from aniso records in model

iotbx.map_model_manager.get_tlso_resid(T, L, S, cm, u_cart, xyz)
iotbx.map_model_manager.get_weights_for_unit_binning(f_array, i_bin)

get weighting for each reflection in a scheme for interpolating a top-hat function over bins

iotbx.map_model_manager.is_inside_mask(mask_map_manager, site_frac=None, inside=True)
class iotbx.map_model_manager.map_model_manager(model=None, map_manager=None, map_manager_1=None, map_manager_2=None, extra_model_list=None, extra_model_id_list=None, extra_map_manager_list=None, extra_map_manager_id_list=None, ncs_object=None, ignore_symmetry_conflicts=None, wrapping=None, absolute_angle_tolerance=0.01, absolute_length_tolerance=0.01, log=None, stop_file=None, make_cell_slightly_different_in_abc=False, name='map_model_manager', verbose=False)

Class for shifting origin of map(s) and model to (0, 0, 0) and keeping track of the shifts.

Typical use: mam = map_model_manager(

model = model, map_manager = map_manager, ncs_object = ncs_object)

mam.box_all_maps_around_model_and_shift_origin(

box_cushion=3)

shifted_model = mam.model() # at (0, 0, 0), knows about shifts shifted_map_manager = mam.map_manager() # also at (0, 0, 0) knows shifts shifted_ncs_object = mam.ncs_object() # also at (0, 0, 0) and knows shifts

Optional after boxing: apply soft mask to map (requires soft_mask_radius)

The maps allowed are:
map_dict has four special ids with interpretations:

map_manager: full map map_manager_1, map_manager_2: half-maps 1 and 2 map_manager_mask: a mask as a map_manager

All other ids are any strings and are assumed to correspond to other maps

Note: It is permissible to call with no map_manager, but supplying

both map_manager_1 and map_manager_2. In this case, the working map_manager will be the average of map_manager_1 and map_manager_2. This will be created the first time map_manager is referenced.

Note: mam.map_manager() contains mam.ncs_object(), so it is not necessary to keep both.

Note: model objects may contain internal ncs objects. These are separate from those in the map_managers and separate from ncs_object in the call to map_model_manager.

The ncs_object describes the NCS of the map and is a property of the map. It must be shared by all maps

The model ncs objects describe the NCS of the individual model. They can differ between models and between models and maps.

Note: set wrapping of all maps to match map_manager if they differ. Set all to be wrapping if it is set

If a model has no crystal_symmetry or not unit_cell, it receives the

crystal symmetry of the map_manager and the shift_cart of the map_manager The same result is obtained if ignore_symmetry_conflicts is set.

add_crystal_symmetry_to_model_if_necessary(model, map_manager=None)

Take any model and add crystal symmetry if it is missing Changes model in place

Parameters: model

Also add unit_cell_crystal_symmetry if missing

add_map_from_fourier_coefficients(map_coeffs, map_id='map_from_fourier_coefficients')

Create map_manager from map_coeffs and add it to maps with map_id The map_coeffs must refer to a map with origin at (0, 0, 0) such as is produced by map_as_fourier_coefficients.

add_map_manager_by_id(map_manager, map_id, overwrite=True, force=False)

Add a new map_manager Must be similar to existing Overwrites any existing with the same id unless overwrite = False Is a mask if is_mask is set

add_model_by_id(model, model_id, overwrite=True)

Add a new model Must be similar to existing map_managers Overwrites any existing with the same id unless overwrite = False If model is None, removes model unless overwrite is False

add_to_info(item_name=None, item=None)
apply_mask_to_map(map_id, mask_id='mask', set_outside_to_mean_inside=False)

Apply the mask in ‘mask’ to map specified by map_id

Optionally set the value outside the mask equal to the mean inside,

changing smoothly from actual values inside the mask to the constant value outside (otherwise outside everything is set to zero)

Optionally use any mask specified by mask_id

NOTE: Does not change the gridding or shift_cart of the map

apply_mask_to_maps(map_ids=None, mask_id='mask', set_outside_to_mean_inside=False)

Apply the mask in ‘mask’ to maps specified by map_ids. If map_ids is None apply to all

Optionally set the value outside the mask equal to the mean inside,

changing smoothly from actual values inside the mask to the constant value outside (otherwise outside everything is set to zero)

Optionally use any mask specified by mask_id

NOTE: Does not change the gridding or shift_cart of the maps

as_map_model_manager()

Return this object (allows using .as_map_model_manager() on both map_model_manager objects and others including box.around_model() etc.

as_map_model_manager_with_resampled_maps(sampling_ratio=2)

Return a new map_model_manager with maps sampled more finely Parameter: sampling_ratio must be an integer Creates new maps, keeps same models Sampling ratio must be greater than 1

as_match_map_model_ncs()

Return this object as a match_map_model_ncs

Includes only the map_manager and model and ncs object, ignores all other maps and models (match_map_model_ncs takes only one of each).

box_all_maps_around_density_and_shift_origin(box_cushion=5.0, threshold=0.05, map_id='map_manager', get_half_height_width=True, model_can_be_outside_bounds=None, soft_mask_radius=None, soft_mask_around_edges=None, boundary_to_smoothing_ratio=2.0, stay_inside_current_map=None, use_cubic_boxing=False, require_match_unit_cell_crystal_symmetry=False, extract_box=False)

Box all maps around the density in map_id map (default is map_manager) shift origin of maps, model

If extract_box=True: Creates new object with deep_copies. Otherwise: replaces existing map_managers and shifts model in place

Replaces existing map_managers and shifts model in place

NOTE: This changes the gridding and shift_cart of the maps and model

Can be used in map_model_manager to work with boxed maps and model or in map_model_manager to re-box all maps and model

Does not require a model, but a model can be supplied. If model is supplied, it is possible that the model will be outside the density after boxing. To avoid this, use box_all_maps_around_model_and_shift_origin instead.

The box_cushion defines how far away from the nearest density the new box boundaries will be placed

The threshold defines how much (relative to maximum in map) above mean value of map near edges is significant and should count as density.

if soft_mask_around_edges, makes a bigger box and makes a soft mask around

the edges. Use this option if you are going to calculate a FT of the map or otherwise manipulate it in reciprocal space.

If use_cubic_box, make a cubic box (in grid units). If also

stay_inside_current_map is set, keep the cubic box inside current map

If require_match_unit_cell_crystal_symmetry is False, do not require unit_cell crystal symmetry to match.

box_all_maps_around_mask_and_shift_origin(box_cushion=5.0, mask_id='mask', model_can_be_outside_bounds=None, soft_mask_radius=None, soft_mask_around_edges=None, boundary_to_smoothing_ratio=2.0, stay_inside_current_map=True, use_cubic_boxing=False, require_match_unit_cell_crystal_symmetry=False, extract_box=False)

Box all maps around specified mask, shift origin of maps, model Replaces existing map_managers and shifts model in place

If extract_box=True: Creates new object with deep_copies. Otherwise: replaces existing map_managers and shifts model in place

NOTE: This changes the gridding and shift_cart of the maps and model

Requires a mask

The box_cushion defines how far away from the edge of the mask the new box boundaries will be placed

if soft_mask_around_edges, makes a bigger box and makes a soft mask around

the edges. Use this option if you are going to calculate a FT of the map or otherwise manipulate it in reciprocal space.

If use_cubic_box, make a cubic box (in grid units). If also

stay_inside_current_map is set, keep the cubic box inside current map

If require_match_unit_cell_crystal_symmetry is False, do not require unit_cell crystal symmetry to match.

box_all_maps_around_model_and_shift_origin(selection_string=None, selection=None, box_cushion=5.0, select_unique_by_ncs=False, model_can_be_outside_bounds=None, stay_inside_current_map=None, soft_mask_radius=None, soft_mask_around_edges=None, boundary_to_smoothing_ratio=2.0, use_cubic_boxing=False, require_match_unit_cell_crystal_symmetry=None, extract_box=False)

Box all maps around the model, shift origin of maps, model If extract_box=True: Creates new object with deep_copies. Otherwise: replaces existing map_managers and shifts model in place

NOTE: This changes the gridding and shift_cart of the maps and model

Can be used in map_model_manager to work with boxed maps and model or in map_model_manager to re-box all maps and model

Requires a model

The box_cushion defines how far away from the nearest atoms the new box boundaries will be placed

The selection_string defines what part of the model to keep (‘ALL’ is

default)

If selection is specified, use instead of selection_string

If select_unique_by_ncs is set, select the unique part of the model automatically. Any selection in selection_string or selection

will not be applied.

if soft_mask_around_edges, makes a bigger box and makes a soft mask around

the edges. Use this option if you are going to calculate a FT of the map or otherwise manipulate it in reciprocal space. Do not use this option if you are going to mask around atoms, density, mask or anything else afterwards as you should apply a mask only once.

If use_cubic_box, make a cubic box (in grid units). If also

stay_inside_current_map is set, keep the cubic box inside current map

If require_match_unit_cell_crystal_symmetry is False, do not require unit_cell crystal symmetry to match.

box_all_maps_around_unique_and_shift_origin(resolution=None, solvent_content=None, sequence=None, molecular_mass=None, soft_mask=True, chain_type='PROTEIN', box_cushion=5, target_ncs_au_model=None, regions_to_keep=None, keep_low_density=True, symmetry=None, mask_expand_ratio=1, soft_mask_radius=None, soft_mask_around_edges=None, boundary_to_smoothing_ratio=2.0, keep_this_region_only=None, residues_per_region=None, use_symmetry_in_extract_unique=True, stay_inside_current_map=True, use_cubic_boxing=False, extract_box=False, require_match_unit_cell_crystal_symmetry=False)

Box all maps using bounds obtained with around_unique, shift origin of maps, model, and mask around unique region

If extract_box=True: Creates new object with deep_copies. Otherwise: replaces existing map_managers and shifts model in place

Replaces existing map_managers and shifts model in place

NOTE: This changes the gridding and shift_cart of the maps and model and masks the map

Normally supply just sequence; resolution will be taken from map_manager resolution if present. other options match all possible ways that segment_and_split_map can estimate solvent_content

Must supply one of (sequence, solvent_content, molecular_mass)

Symmetry is optional symmetry (i.e., D7 or C1). Used as alternative to ncs_object supplied in map_manager

Use_symmetry can be set to False to ignore symmetry found in ncs_object.

NCS object is still kept and shifted however.

if soft_mask_around_edges, makes a bigger box and makes a soft mask around

the edges. Use this option if you are going to calculate a FT of the map or otherwise manipulate it in reciprocal space.

If use_cubic_box, make a cubic box (in grid units). If also

stay_inside_current_map is set, keep the cubic box inside current map

Additional parameters:
mask_expand_ratio: allows increasing masking radius beyond default at

final stage of masking

solvent_content: fraction of cell not occupied by macromolecule. Can

be None in which case it is estimated from map

sequence: one-letter code of sequence of unique part of molecule chain_type: PROTEIN or RNA or DNA. Used with sequence to estimate

molecular_mass

molecular_mass: Molecular mass (Da) of entire molecule used to

estimate solvent_content

target_ncs_au_model: model marking center of location to choose as

unique

box_cushion: buffer around unique region to be boxed soft_mask: use soft mask keep_low_density: keep low density regions regions_to_keep: Allows choosing just highest-density contiguous

region (regions_to_keep=1) or a few

residues_per_region: Try to segment with this many residues per region keep_this_region_only: Keep just this region (first one is 0 not 1) require_match_unit_cell_crystal_symmetry: require model unit_cell_

crystal_symmetry to match when extracting.

box_all_maps_with_bounds_and_shift_origin(lower_bounds, upper_bounds, model_can_be_outside_bounds=None, soft_mask_radius=None, soft_mask_around_edges=None, boundary_to_smoothing_ratio=2.0, stay_inside_current_map=True, use_cubic_boxing=False, require_match_unit_cell_crystal_symmetry=False, extract_box=False)

Box all maps using specified bounds, shift origin of maps, model Replaces existing map_managers and shifts model in place

If extract_box=True: Creates new object with deep_copies. Otherwise: replaces existing map_managers and shifts model in place

NOTE: This changes the gridding and shift_cart of the maps and model Also changes space group to p1

Can be used in map_model_manager to work with boxed maps and model or in map_model_manager to re-box all maps and model

The lower_bounds and upper_bounds define the region to be boxed. These bounds are relative to the current map with origin at (0, 0, 0).

if soft_mask_around_edges, still uses the same bounds, but makes

a soft mask around the edges. Use this option if you are going to calculate a FT of the map or otherwise manipulate it in reciprocal space. Do not use this option if you are going to mask around atoms, density, mask or anything else afterwards as you should apply a mask only once.

If use_cubic_box, make a cubic box (in grid units). If also

stay_inside_current_map is set, keep the cubic box inside current map

If require_match_unit_cell_crystal_symmetry is False, do not require unit_cell crystal symmetry to match.

check_stop_file()
choose_best_set(dd, max_dist=None)
counts()
create_mask_around_atoms(model=None, mask_atoms_atom_radius=3, soft_mask=False, soft_mask_radius=None, invert_mask=None, mask_id='mask')

Generate mask based on model. Does not apply the mask to anything. Normally follow with apply_mask_to_map or apply_mask_to_maps

Optional: radius around atoms for masking Optional: soft mask (default = True)

Radius will be soft_mask_radius (default radius is self.resolution() or resolution calculated

from gridding)

If soft mask is set, mask_atoms_atom_radius increased by

soft_mask_radius

Optional: invert_mask: keep outside atoms instead of inside NOTE: if space group is not p1 and wrapping is True then

atoms are expanded to P1 before calculating mask

Generates new entry in map_manager dictionary with id of mask_id (default=’mask’) replacing any existing entry with that id

create_mask_around_density(resolution=None, solvent_content=None, soft_mask=True, soft_mask_radius=None, mask_id='mask', map_id='map_manager')

Generate mask based on density in map_manager (map_id defines it). Does not apply the mask to anything. Normally follow with apply_mask_to_map or apply_mask_to_maps

Optional: supply working resolution Optional: supply approximate solvent fraction

Optional: soft mask (default = True)

Radius will be soft_mask_radius (default radius is resolution calculated from gridding)

Generates new entry in map_manager dictionary with id of mask_id (default=’mask’) replacing any existing entry with that id

create_mask_around_edges(soft_mask_radius=None, boundary_to_smoothing_ratio=2.0, boundary_radius=None, mask_id='mask')

Generate new mask map_manager with soft mask around edges of mask Does not apply the mask to anything. Normally follow with apply_mask_to_map or apply_mask_to_maps

Optional: radius around edge for masking
(default radius is self.resolution() or resolution calculated

from gridding)

Generates new entry in map_manager dictionary with id of mask_id (default=’mask’) replacing any existing entry with that id

create_masked_copies_of_maps(map_id_list=None, mask_id='mask')

Create masked copies of all maps identified by map_id_list (default is all) Return list of map_id for masked versions

create_spherical_mask(soft_mask=True, mask_center_cart=None, mask_radius=None, soft_mask_radius=None, boundary_radius=None, boundary_to_smoothing_ratio=2.0, mask_id='mask', minimum_mask_radius_ratio=0.25)

Generate spherical mask with radius mask_radius around the cartesian point mask_center_cart. The value of

mask_center_cart is relative to the

shifted (origin at (0,0,0) ) position of the map. Does not apply the mask to anything. Normally follow with apply_mask_to_map or apply_mask_to_maps Default: calculate spherical soft mask centered at center of map

soft_mask radius default is resolution() boundary between mask and closest edge is

soft_mask_radius * boundary_to_smoothing_ratio

Optional: not soft mask. Same as soft mask in dimensions but the

soft_mask will not be applied

Optional: mask_center_cart (default is center of map)

Generates new entry in map_manager dictionary with id of mask_id (default=’mask’) replacing any existing entry with that id

crystal_symmetry()

Get the working crystal_symmetry

customized_copy(model_dict=None, map_dict=None, name=None)

Produce a copy of this map_model object, replacing nothing, maps or models, or both

deep_copy()

Return a deep_copy of this map_manager Use customized copy with default map_dict and model_dict (from self)

density_at_model_sites(map_id='map_manager', model_id='model', selection_string=None, model=None)

Return density at sites in model

duplicate_id(dd)
duplicate_map_manager(map_id='map_manager', new_map_id='new_map_manager')

Duplicate (deep_copy) map_manager Overwrites any existing with the new id

expand_mask(buffer_radius=5, resolution=None, soft_mask=True, soft_mask_radius=None, mask_id='mask')
experiment_type()
external_sharpen(map_id='map_manager', map_id_external_map='external_map', map_id_to_be_scaled_list=None, map_id_scaled_list=None, exclude_points_outside_density=None, minimum_boxes_inside_density=None, resolution=None, d_min=None, k_sol=None, b_sol=None, n_bins=None, n_boxes=None, core_box_size=None, box_cushion=None, smoothing_radius=None, local_sharpen=None, anisotropic_sharpen=None, expected_ssqr_list=None, expected_ssqr_list_rms=None, tlso_group_info=None, get_tls_from_u=None, overall_sharpen_before_and_after_local=False, get_scale_as_aniso_u=None, use_dv_weighting=None, n_direction_vectors=None, run_analyze_anisotropy=True, sharpen_all_maps=False, nproc=None)
Scale map_id with scale factors identified from map_id vs

map_id_external_map

Changes the working map_manager

resolution is nominal resolution of map d_min is minimum resolution to use in calculation of Fourier coefficients

extract_all_maps_around_density(box_cushion=5.0, threshold=0.05, get_half_height_width=True, model_can_be_outside_bounds=None, boundary_to_smoothing_ratio=2.0, soft_mask_around_edges=None, soft_mask_radius=None, stay_inside_current_map=True, require_match_unit_cell_crystal_symmetry=None, use_cubic_boxing=False, map_id='map_manager')

Runs box_all_maps_around_density_and_shift_origin with extract_box=True

extract_all_maps_around_mask(box_cushion=5.0, model_can_be_outside_bounds=None, boundary_to_smoothing_ratio=2.0, soft_mask_around_edges=None, soft_mask_radius=None, stay_inside_current_map=True, require_match_unit_cell_crystal_symmetry=None, use_cubic_boxing=False, mask_id='mask')

Runs box_all_maps_around_mask_and_shift_origin with extract_box=True

extract_all_maps_around_model(selection_string=None, selection=None, select_unique_by_ncs=False, model_can_be_outside_bounds=None, stay_inside_current_map=None, box_cushion=5.0, boundary_to_smoothing_ratio=2.0, soft_mask_around_edges=None, soft_mask_radius=None, require_match_unit_cell_crystal_symmetry=None, use_cubic_boxing=False)

Runs box_all_maps_around_model_and_shift_origin with extract_box=True

extract_all_maps_around_unique(resolution=None, solvent_content=None, sequence=None, molecular_mass=None, soft_mask=True, chain_type='PROTEIN', box_cushion=5, target_ncs_au_model=None, regions_to_keep=None, require_match_unit_cell_crystal_symmetry=None, keep_low_density=True, symmetry=None, boundary_to_smoothing_ratio=2.0, soft_mask_around_edges=None, keep_this_region_only=None, residues_per_region=None, soft_mask_radius=None, mask_expand_ratio=1, stay_inside_current_map=True, use_cubic_boxing=False, use_symmetry_in_extract_unique=True)

Runs box_all_maps_around_unique_and_shift_origin with extract_box=True

extract_all_maps_with_bounds(lower_bounds, upper_bounds, boundary_to_smoothing_ratio=2.0, soft_mask_around_edges=None, soft_mask_radius=None, stay_inside_current_map=True, use_cubic_boxing=False, require_match_unit_cell_crystal_symmetry=None, model_can_be_outside_bounds=None)

Runs box_all_maps_with_bounds_and_shift_origin with extract_box=True

find_k_sol_b_sol(model=None, d_min=None, model_map_id=None, comparison_map_id=None, n_bins=5)

Routine to guess k_sol and b_sol by low-resolution Fc calculation

generate_map(d_min=None, origin_shift_grid_units=None, file_name=None, model=None, n_residues=None, b_iso=30, k_sol=None, b_sol=None, box_cushion=5, scattering_table=None, fractional_error=0.0, gridding=None, wrapping=False, map_id=None, f_obs_array=None, resolution_factor=None)

Simple interface to cctbx.development.generate_map allowing only a small subset of keywords. Useful for quick generation of models, map coefficients, and maps

For full functionality use cctbx.development.generate_model, cctbx.development.generate_map_coeffs, and cctbx.development.generate_map

Summary:

If no map_manager is present, use supplied or existing model to

generate map_manager and model.

If map_manager is present and is not a dummy map_manager,

use supplied or existing model as model and create new entry in this this map_model_manager with name map_id. If map_id is None, use ‘model_map’

If no existing or supplied model, use default model from library, box with box_cushion around it and choose n_residues to include (default=10).

Parameters:

model (model.manager object, None): model to use (as is) file_name (path , None): file containing coordinates to use (instead

of default model)

n_residues (int, 10): Number of residues to include (from default

model or file_name)

b_iso (float, 30): B-value (ADP) to use for all atoms if model

is not supplied

box_cushion (float, 5): Buffer (A) around model d_min (float, 3): high_resolution limit (A) gridding (tuple (nx, ny, nz), None): Gridding of map (optional) origin_shift_grid_units (tuple (ix, iy, iz), None): Move location of

origin of resulting map to (ix, iy, iz) before writing out

wrapping: Defines if map is to be specified as wrapped resolution_factor: Defines ratio of resolution to gridding if gridding is not set scattering_table (choice, ‘electron’): choice of scattering table

All choices: wk1995 it1992 n_gaussian neutron electron

fractional_error: resolution-dependent fractional error, ranging from

zero at low resolution to fractional_error at d_min. Can be more than 1.

map_id: ID of map_manager to be created with model-map information (only

applies if there is an existing map_manager)

get_any_map_manager()

Return any map manager

get_cb_resseq_to_skip(cb_resseq_list, far_away=None, far_away_n=None, very_far_away=None, very_far_away_n=None)

Identify numbers in resseq_list that are way different than all others by at least max_gap

get_counts_and_histograms()
get_diffs_for_matching_target_and_model(matching_info=None, max_dist=None, ca_only=None, reverse=False, allow_reverse=False)
get_info(item_name=None)
get_map_data_by_id(map_id)

Get map_data from a map_manager with the name map_id

get_map_manager_by_id(map_id)

Get a map_manager with the name map_id If map_id is ‘map_manager’ specifically return self.map_manager() so that it will create a map_manager from map_manager_1 and map_manager_2 if map_manager is not present

get_map_model_manager_with_selected(map_id_list=None, model_id_list=None, deep_copy=False)
get_map_rmse(map_id='map_manager', mask_id='mask', fc_map_id='fc_for_model_comparison', fofc_map_id='fofc_for_model_comparison', model=None, d_min=None, r_free=None, r_work=None, match_overall_b=True)

Estimate map rms error by comparison with model in region containing the model. Suitable for cases where model is not refined against this map. If refined against this map, correct for degrees of freedom by estimating overfitting from ratio of Rfree/Rwork. If match_overall_b is True, adjust average B of model to match fall-off of the Fourier coefficients representing the map Return group_args object with:

map_rmse = fofc_info_inside.sd # map rms error map_rmse_to_sd_ratio = ratio_error_to_input_map # rmse / sd of map map_sd = map_info.sd # SD of the map (rms after setting mean to zero)

get_model_by_id(model_id)

Get a model with the name model_id

get_model_from_other(other, other_model_id='model')
Take a model with id other_model_id from other_map_model_manager with any

boxing and origin shifts allowed, and put it in the same reference frame as the current model. Used to build up a model from pieces that were worked on in separate boxes.

Changes model from other in place

Parameters: other: Other map_model_manager containing a model

get_ncs_from_map(use_existing=True, include_helical_symmetry=False, symmetry_center=None, min_ncs_cc=None, symmetry=None, ncs_object=None)

Use existing ncs object in map if present or find ncs from map Sets ncs_object in self.map_manager() Sets self._ncs_cc which can be retrieved with self.ncs_cc()

get_ncs_from_model()

Return model NCS as ncs_spec object if available Does not set anything. If you want to save it use:

self.set_ncs_object(self.get_ncs_from_model()) This will set the ncs object in the map_manager (if present)

get_rms_f_list(map_id='map_manager', d_min=None, n_bins=None, resolution=None)

Return list of rms amplitude by bins

get_selection_string_from_chain_dict(chain_dict=None, max_gap=None, minimum_length=None, return_as_group_args_list=False)

Return a selection string for the segments represented in chain_dict, allowing gaps of up to max-gap (fill them in) Require minimum_length if set

half_map_sharpen(map_id='map_manager', map_id_1='map_manager_1', map_id_2='map_manager_2', map_id_scaled_list=None, map_id_to_be_scaled_list=None, exclude_points_outside_density=None, minimum_boxes_inside_density=None, resolution=None, d_min=None, k_sol=None, b_sol=None, n_bins=None, n_boxes=None, core_box_size=None, box_cushion=None, smoothing_radius=None, rmsd=None, local_sharpen=None, anisotropic_sharpen=None, minimum_low_res_cc=None, get_scale_as_aniso_u=None, use_dv_weighting=None, n_direction_vectors=None, run_analyze_anisotropy=True, spectral_scaling=True, expected_rms_fc_list=None, expected_ssqr_list=None, expected_ssqr_list_rms=None, tlso_group_info=None, get_tls_from_u=None, model_id_for_rms_fc=None, replace_aniso_with_tls_equiv=None, max_abs_b=None, nproc=None, optimize_b_eff=None, equalize_power=None, overall_sharpen_before_and_after_local=False, get_tls_info_only=None, coordinate_shift_to_apply_before_tlso=None, sharpen_all_maps=False, remove_overall_anisotropy=True)

Scale map_id with scale factors identified from map_id_1 vs map_id_2 Changes the working map_manager unless map_id_scaled_list is set.

max_abs_b applies if get_scale_as_aniso_u and anisotropic_sharpen and

local_sharpen are set. It limits range of anisotropic B. Default is 100 at 4 A, proportional to resolution squared

resolution is nominal resolution of map d_min is minimum resolution to use in calculation of Fourier coefficients

hierarchy()
histograms()
info()
initialize_maps(map_value=0)

Set values of all maps to map_value Used to set up an empty set of maps for filling in from boxes

local_fsc(map_id='map_manager', map_id_1='map_manager_1', map_id_2='map_manager_2', model_id=None, map_id_to_be_scaled_list=None, map_id_scaled_list=None, mask_id=None, exclude_points_outside_density=None, minimum_boxes_inside_density=None, resolution=None, d_min=None, k_sol=None, b_sol=None, max_resolution_ratio=None, min_bin_width=20, n_bins=None, fsc_cutoff=0.143, n_boxes=None, core_box_size=None, box_cushion=None, rmsd=None, smoothing_radius=None, nproc=None, is_model_based=None, optimize_b_eff=None, equalize_power=None, is_external_based=None, return_scale_factors=False, direction_vectors=None, minimum_low_res_cc=None, get_scale_as_aniso_u=None, use_dv_weighting=None, n_direction_vectors=None, run_analyze_anisotropy=None, spectral_scaling=None, expected_rms_fc_list=None, expected_ssqr_list=None, expected_ssqr_list_rms=None, tlso_group_info=None, get_tls_from_u=None, model_id_for_rms_fc=None, replace_aniso_with_tls_equiv=None, max_abs_b=None, get_tls_info_only=None, coordinate_shift_to_apply_before_tlso=None, sharpen_all_maps=None, n_bins_default=2000, b_iso=None, aniso_b_cart=None)

Calculates local Fourier Shell Correlations to estimate local resolution Creates map with smoothed local resolution

Optionally estimates scale factors vs resolution at each point in map to apply to yield a locally-scaled map (return_scale_factors = True).

If direction_vector is specified, weight scale factor calculation by dot product of reflection directions with direction_vector

map_as_fourier_coefficients(d_min=None, d_max=None, map_id='map_manager')

Return Miller array to resolution specified based on map with id map_id

Note that the map_manager is always zero-based (origin at (0,0,0)). The Fourier coefficients represent the map in this location at (0, 0, 0)

map_data()
map_data_1()
map_data_2()
map_data_list()
map_dict()

Get the dictionary of all maps and masks as map_manager objects

map_id_list()

Get all the names (ids) for all map_managers that are present

map_info(map_id='map_manager', sigma_cutoff=1, quiet=False)

Summarizes info about map and returns group_args

map_manager()

Get the map_manager

If not present, calculate it from map_manager_1 and map_manager_2 and set it.

map_manager_1()

Get half_map 1 as a map_manager object

map_manager_2()

Get half_map 2 as a map_manager object

map_manager_mask()

Get the mask as a map_manager object

map_managers()

Get all the map_managers as a list

map_map_cc(map_id='map_manager_1', other_map_id='map_manager_2', mask_id=None, mask_cutoff=0.5, resolution=None)
map_map_fsc(map_id_1='map_manager_1', map_id_2='map_manager_2', resolution=None, mask_id=None, mask_cutoff=0.5, min_bin_width=20, n_bins=2000, fsc_cutoff=0.143)

Return the map-map FSC for these two maps, optionally masked with mask_id Returns fsc object which contains d_min which is d_min where fsc

drops to fsc_cutoff, and sub-object fsc with arrays d, d_inv and fsc which are the FSC curve

map_model_cc(resolution=None, map_id='map_manager', model_id='model', selection_string=None, model=None, use_b_zero=True)
mask_all_maps_around_atoms(model=None, mask_atoms_atom_radius=3, set_outside_to_mean_inside=False, soft_mask=False, soft_mask_radius=None, skip_n_residues_on_ends=None, invert_mask=None, mask_id='mask')
mask_all_maps_around_density(solvent_content=None, soft_mask=True, soft_mask_radius=None, mask_id='mask', map_id='map_manager')

Apply a soft mask around density. Mask calculated using map_id and written to mask_id . Overwrites values in maps Default is to use ‘mask’ as the mask id

NOTE: Does not change the gridding or shift_cart of the maps and model

mask_all_maps_around_edges(soft_mask_radius=None, boundary_to_smoothing_ratio=2.0, mask_id='mask')

Apply a soft mask around edges of all maps. Overwrites values in maps Use ‘mask’ as the mask id

NOTE: Does not change the gridding or shift_cart of the maps and model

mask_info(mask_id='mask', cutoff=0.5, quiet=False)

Summarizes info about mask and returns group_args

match_cb_to_ca(ca_sites=None, cb_as_list=None, ca_residue_names=None, cb_residue_names=None, max_dist=2.0)

Identify cb sites that match ca sites If ca_residue_names and cb_residue names, require that residue names match

match_model_b_to_map(map_id='map_manager', model=None, d_min=None)
merge_split_maps_and_models(model_id=None, box_info=None, replace_coordinates=True, replace_u_aniso=False, allow_changes_in_hierarchy=False, output_model_id=None)
Replaces coordinates in working model with those from the

map_model_managers in box_info. The box_info object should come from running split_up_map_and_model in this instance of the map_model_manager.

If allow_changes_in_hierarchy is set, create a new working model where the hierarchy has N “models”, one from each box. This allows changing the hierarchy structure. This will create one model in a new hierarchy for each box, numbered by box number. The new hierarchy will be placed in a model with id output_model_id (default is model_id, replacing existing model specified by model_id; usually this is just ‘model’, the default model.)

minimum_resolution()

Return d_min, normally minimum available but if set, return value of d_min

model()

Get the model

model_building(nproc=None, soft_zero_boundary_mask=True, soft_zero_boundary_mask_radius=None, model_id='model', normalize=True)

Return this object as a local_model_building object The model-building object has pointers to model and map_manager, not

copies

resolution is resolution for Fourier coefficients is_xray_map is True for x-ray map nproc is number of processors to use

If no map_manager is present, resolution, experiment type are not required

model_dict()

Get the dictionary of all models

model_from_hierarchy(hierarchy, return_as_model=False, model_id='model_from_hierarchy')

Convenience method to convert a hierarchy into a model, where the model has symmetry and shift cart matching this manager

model_from_sites_cart(sites_cart, model_id='model_from_sites', return_model=None, **kw)

Convenience method to use the from_sites_cart method of model manager to create a model.

Model is saved with the id of model_id (default ‘model_from_sites’)

Unique feature of this method is that it sets up the crystal_symmetry,

unit_cell_crystal_symmetry, and shift_cart in the model to match this map_model_manager.

Note: sites_cart are assumed to be in the same coordinate frame as this

map_model_manager.

Any keywords used in from_sites_cart are allowed (i.e.,

atom_name and scatterer or atom_name_list and scatterer_list occ or occ_list b_iso or b_iso_list resname or resname_list resseq_list

model_from_text(text, return_as_model=False, model_id='model_from_text')

Convenience method to convert txt into a model, where the model has symmetry and shift cart matching this manager

model_id_list()

Get all the names (ids) for all models

model_sharpen(map_id='map_manager', model_id='model', map_id_scaled_list=None, map_id_to_be_scaled_list=None, exclude_points_outside_density=True, minimum_boxes_inside_density=True, resolution=None, d_min=None, k_sol=None, b_sol=None, find_k_sol_b_sol=True, d_min_for_k_sol_b_sol=6.0, n_bins=None, n_boxes=None, core_box_size=None, box_cushion=None, smoothing_radius=None, rmsd=None, local_sharpen=None, anisotropic_sharpen=None, minimum_low_res_cc=0.2, get_scale_as_aniso_u=None, use_dv_weighting=None, n_direction_vectors=None, run_analyze_anisotropy=True, spectral_scaling=True, expected_rms_fc_list=None, expected_ssqr_list=None, expected_ssqr_list_rms=None, tlso_group_info=None, get_tls_from_u=None, find_tls_from_model=None, model_id_for_rms_fc=None, replace_aniso_with_tls_equiv=None, max_abs_b=None, nproc=None, optimize_b_eff=None, equalize_power=None, map_id_model_map='model_map_for_scaling', optimize_with_model=None, overall_sharpen_before_and_after_local=False, mask_around_model=True, get_tls_info_only=None, coordinate_shift_to_apply_before_tlso=None, sharpen_all_maps=False, remove_overall_anisotropy=True, save_model_map=False)

Scale map_id with scale factors identified from map_id vs model Changes the working map_manager unless map_id_scaled is set.

max_abs_b applies if get_scale_as_aniso_u and anisotropic_sharpen and

local_sharpen are set. It limits range of anisotropic B. Default is 100 at 4 A, proportional to resolution squared

resolution is nominal resolution of map d_min is minimum resolution to use in calculation of Fourier coefficients

models()

Get all the models as a list

ncs_cc()
ncs_object()
nproc()
propagate_model_from_other(other, model_id='model', other_model_id='model')

Import a model from other with get_model_from_other (other_model_id), then set coordinates of corresponding atoms in model_id

The model in other must have been extracted from the model in this object or one just like it with select_unique_by_ncs=True, and no atoms can have been added or removed.

remove_anisotropy(d_min=None, map_coeffs=None, aniso_b_cart=None, map_id='map_manager', map_ids=None, remove_from_all_maps=False, model_map_ids_to_leave_as_is=None, b_iso=None)
Remove anisotropy from map, optionally remove anisotropy specified by

aniso_b_cart and b_iso

remove_map_manager_by_id(map_id='extra')

Remove this map manager Note: you cannot remove ‘map_manager’ … you can only replace it

remove_model_by_id(model_id='extra')

Remove this model

remove_model_outside_map(model=None, boundary=3, return_as_new_model=False)

Remove all the atoms in the model that are well outside the map (more than boundary). Boundary can be negative (remove inside box near edges)

remove_origin_shift_and_unit_cell_crystal_symmetry()

Set the origin of the current map to be (0,0,0) and reset the unit_cell_crystal_symmetry to match the current crystal_symmetry.

Do not use this method if you want to use the map_model_manager or any of its maps and models in the usual way. This is a method to create a set of maps and model that have no reference to their original locations.

Typical use:

box_mmm = mmm.extract_all_maps_around_model() box_mmm.remove_origin_shift_and_unit_cell_crystal_symmetry() model_file = data_manager.write_model_file(box_mmm.model(), model_file) data_manager.write_real_map_file(box_mmm.map_manager(), map_file)

Now model_file and map_file have no origin shift and the model CRYST1 matches the crystal_symmetry and map_file unit_cell_crystal_symmetry of map_file

resolution(use_fsc_if_no_resolution_available_and_maps_available=True, map_id_1='map_manager_1', map_id_2='map_manager_2', fsc_cutoff=0.143)
resolution_filter(d_min=None, d_max=None, map_id='map_manager', new_map_id='map_manager')

Resolution-filter a map with range of d_min to d_max and place in new map (can be the same)

Typically used along with duplicate_map_manager to create a new map and filter it:

rm.duplicate_map_manager(map_id=’map_manager’,

new_map_id=’resolution_filtered’)

rm.resolution_filter(map_id = ‘resolution_filtered’,)

rmsd_of_matching_residues(target_model_id='model', matching_model_id=None, target_model=None, matching_model=None, max_dist=None, minimum_length=None, chain_type=None, atom_name=None, element=None, ignore_element=None, ca_only=True, matching_info_list=None, allow_reverse=None, quiet=True)

Get rmsd of all or ca(P) atoms in matching residues

scattering_table()
select_matching_segments(target_model_id='model', matching_model_id=None, target_model=None, matching_model=None, chain_type=None, atom_name=None, element=None, ignore_element=False, max_dist=None, max_dist_extra=None, minimum_length=None, max_gap=5, far_away=7, far_away_n=2, very_far_away=20, very_far_away_n=1000, one_to_one=False, residue_names_must_match=False, minimum_match_length=2, shift_without_deep_copy=False, quiet=True)
sequence(model=None, model_id='model', selection_string=None)

Return sequence of model

set_experiment_type(experiment_type)

Set nominal experiment_type

set_info(info)
set_log(log=<colorama.ansitowin32.StreamWrapper object>)

Set output log file

set_map_id_lists(kw)
set_map_manager(map_manager)

Overwrites existing map_manager with id ‘map_manager’

set_minimum_resolution(d_min)

Set minimum resolution used in calculations

set_model(model, overwrite=True)

Overwrites existing model with id ‘model’ Allows setting model to None if overwrite is True

set_model_symmetries_and_shift_cart_to_match_map(model)

Method to set the crystal_symmetry, unit_cell_crystal_symmetry and shift_cart of any model to match those of this map_model_manager

This method changes the model in place.

Assumes that the sites in this model are already in the same coordinate system as this map_model_manager

This is the preferred method of making sure that a model has the same symmetry and shifts as a map_model_manager if the model coordinates already match the map_model_manager.

If you have a model that is relative to some other coordinate system, first make sure that model has a complete set of symmetry and shifts for that coordinate system (any model coming from a map_model_manager will have these). Then use either self.shift_any_model_to_match(model) to shift that model directly to match this map_model_manager, or use self.add_model_by_id(model_id=model_id, model=model) which will import that model into this manager and shift it to match in the process.

set_multiprocessing(nproc=None, multiprocessing=None, queue_run_command=None)

Set multiprocessing parameters

set_name(name=None)

Set name

set_ncs_object(ncs_object)

Set the ncs object of map_manager

set_original_origin_grid_units(original_origin_grid_units=None)
Reset (redefine) the original origin of the maps and models (apply an

origin shift in effect).

Procedure is: calculate shift_cart and set origin_shift_grid_units and

shift_cart everywhere

set_resolution(resolution)

Set nominal resolution. Pass along to any map managers

set_scattering_table(scattering_table)
Set nominal scattering_table. Overrides anything in map_managers

electron: cryo_em n_gaussian x-ray (standard) wk1995: x-ray (alternative) it1992: x-ray (alternative) neutron: neutron scattering

set_stop_file(file_name=None)

Define file name that means “STOP”

set_up_map_dict(map_manager=None, map_manager_1=None, map_manager_2=None, extra_map_manager_list=None, extra_map_manager_id_list=None)
map_dict has four special ids with interpretations:

map_manager: full map map_manager_1, map_manager_2: half-maps 1 and 2 map_manager_mask: a mask in a map_manager

All other ids are any strings and are assumed to correspond to other maps map_manager must be present

set_up_model_dict(model=None, extra_model_list=None, extra_model_id_list=None)
map_dict has one special id with interpretation:

model: standard model

All other ids are any strings and are assumed to correspond to other models.

set_verbose(verbose=None)

Set verbose

shift_any_model_to_match(model, map_manager=None, set_unit_cell_crystal_symmetry=False)

Take any model and shift it to match the working shift_cart Also sets crystal_symmetry. Changes model in place

Parameters: model

set_unit_cell_crystal_symmetry: optionally set this as well

shift_aware_rt(from_obj=None, to_obj=None, working_rt_info=None, absolute_rt_info=None)

Returns shift_aware_rt object

Uses rt_info objects (group_args with members of r, t).

Simplifies keeping track of rotation/translation between two

objects that each may have an offset from absolute coordinates.

absolute rt is rotation/translation when everything is in original,

absolute cartesian coordinates.

working_rt is rotation/translation of anything in “from_obj” object

to anything in “to_obj” object using working coordinates in each.

Usage: shift_aware_rt = self.shift_aware_rt(absolute_rt_info = rt_info) shift_aware_rt = self.shift_aware_rt(working_rt_info = rt_info,

from_obj=from_obj, to_obj = to_obj)

apply RT using working coordinates in objects sites_cart_to_obj = shift_aware_rt.apply_rt(sites_cart_from_obj,

from_obj=from_obj, to_obj=to_obj)

apply RT absolute coordinates sites_cart_to = shift_aware_rt.apply_rt(sites_cart_from)

shift_aware_rt_to_superpose_other(other, selection_string=None)
Identify rotation/translation to transform model from other on to

model in this object.

Optionally apply selection_string to both models before doing the

mapping

shift_cart()

get the shift_cart (shift since original location)

shift_origin_to_match_original()

Shift the origin of all maps and models to match original

show_summary(log=<colorama.ansitowin32.StreamWrapper object>)
split_up_map_and_model_by_boxes(model_id='model', skip_waters=False, skip_hetero=False, write_files=False, target_for_boxes=24, select_final_boxes_based_on_model=True, box_cushion=3, mask_around_unselected_atoms=None, mask_radius=3, masked_value=-10, skip_empty_boxes=True, apply_box_info=True, mask_id=None, exclude_points_outside_density=None, minimum_boxes_inside_density=None)

Split up the map, creating boxes that time the entire map.

Try to get about target_for_boxes boxes. Do not go over this target

If select_final_boxes_based_on_model then make the final boxes just go

around the selected parts of the model with cushion defined by box_cushion and not tile the map.

Otherwise select atoms inside the boxes and afterwards expand the boxes

with box_cushion

If skip_empty_boxes then skip boxes with no model.

Note that this procedure just selects by atom so you can get a single atom

in a box

Returns a group_args object containing list of the map_model_manager

objects and a list of the selection objects that define which atoms from the working model are in each object.

Normally do work on each map_model_manager to create a new model with

the same atoms, then use merge_split_maps_and_models() to replace coordinates in the original model with those from all the component models.

Optionally carry out the step box_info = get_split_maps_and_models(…)

separately with the keyword apply_box_info=False

skip_waters and skip_hetero define whether waters and hetero atoms are

ignored

split_up_map_and_model_by_chain(model_id='model', skip_waters=False, skip_hetero=False, box_cushion=3, mask_around_unselected_atoms=None, mask_radius=3, masked_value=-10, write_files=False, apply_box_info=True)

Split up the map, boxing around each chain in the model.

Returns a group_args object containing list of the map_model_manager

objects and a list of the selection objects that define which atoms from the working model are in each object.

Normally do work on each map_model_manager to create a new model with

the same atoms, then use merge_split_maps_and_models() to replace coordinates in the original model with those from all the component models.

Optionally carry out the step box_info = get_split_maps_and_models(…)

separately with the keyword apply_box_info=False

skip_waters and skip_hetero define whether waters and hetero atoms are

ignored

box_cushion is the padding around the model atoms when creating boxes

split_up_map_and_model_by_ncs_groups(model_id='model', box_cushion=3, mask_around_unselected_atoms=None, mask_radius=3, masked_value=-10, write_files=False, apply_box_info=True)

Split up the map, boxing around atoms selected with each ncs group in ncs_groups obtained from supplied model

Returns a group_args object containing list of the map_model_manager

objects and a list of the selection objects that define which atoms from the working model are in each object.

Normally do work on each map_model_manager to create a new model with

the same atoms, then use merge_split_maps_and_models() to replace coordinates in the original model with those from all the component models.

Optionally carry out the step box_info = get_split_maps_and_models(…)

separately with the keyword apply_box_info=False

box_cushion is the padding around the model atoms when creating boxes

split_up_map_and_model_by_segment(model_id='model', skip_waters=False, skip_hetero=False, box_cushion=3, mask_around_unselected_atoms=None, mask_radius=3, masked_value=-10, write_files=False, apply_box_info=True)
Split up the map, boxing around each segment (each unbroken part of

each chain) in the model

Returns a group_args object containing list of the map_model_manager

objects and a list of the selection objects that define which atoms from the working model are in each object.

Normally do work on each map_model_manager to create a new model with

the same atoms, then use merge_split_maps_and_models() to replace coordinates in the original model with those from all the component models.

Optionally carry out the step box_info = get_split_maps_and_models(…)

separately with the keyword apply_box_info=False

skip_waters and skip_hetero define whether waters and hetero atoms are

ignored

box_cushion is the padding around the model atoms when creating boxes

split_up_map_and_model_by_supplied_selections(selection_list, model_id='model', box_cushion=3, mask_around_unselected_atoms=None, mask_radius=3, masked_value=-10, write_files=False, apply_box_info=True)
Split up the map, boxing around atoms selected with each selection in

selection_list Note: a selection can be obtained with:

self.model().selection(selection_string)

Returns a group_args object containing list of the map_model_manager

objects and a list of the selection objects that define which atoms from the working model are in each object.

Normally do work on each map_model_manager to create a new model with

the same atoms, then use merge_split_maps_and_models() to replace coordinates in the original model with those from all the component models.

Optionally carry out the step box_info = get_split_maps_and_models(…)

separately with the keyword apply_box_info=False

box_cushion is the padding around the model atoms when creating boxes

superposed_map_manager_from_other(other, working_rt_info=None, absolute_rt_info=None, shift_aware_rt_info=None, selection_string=None)
Identify rotation/translation to transform model from other on to

model in this object.

Optionally apply selection_string to both models before doing the

mapping

Then extract map from other to cover map in this object, Fill in with zero where undefined if wrapping is False.

Allow specification of working_rt (applies to working coordinates in

other and self), or absolute_rt_info (applies to absolute, original coordinates)

tls_from_map(map_id_1=None, map_id_2=None, map_id=None, model_id=None, mask_id=None, tls_by_chain=True, apply_tls_to_model=True, iterations=1, skip_waters=True, skip_hetero=True, coordinate_shift_to_apply_before_tlso=None, core_box_size_ratio=None, box_cushion_ratio=None, exclude_points_outside_density=True, minimum_boxes_inside_density=True, d_min=None, **kw)
unit_cell_crystal_symmetry()

Get the unit_cell_crystal_symmetry (full or original symmetry)

warning_message()
write_map(file_name, map_id='map_manager')
write_model(file_name, model_id=None, model=None, data_manager=None, format=None)
xray_structure()
class iotbx.map_model_manager.match_map_model_ncs(log=None, ignore_symmetry_conflicts=None, absolute_angle_tolerance=0.01, absolute_length_tolerance=0.01)

Use: Container to hold map, model, ncs object and check consistency and shift origin

Normal usage:

Initialize empty, then read in or add a group of model.manager, map_manager, and ncs objects

Read in the models, maps, ncs objects

Shift origin to (0, 0, 0) and save position of this (0, 0, 0) point in the

original coordinate system so that everything can be written out superimposed on the original locations. This is origin_shift_grid_units in grid units

NOTE: modifies the model, map_manager, and ncs objects. Call with deep_copy() of these if originals need to be preserved.

Input models, maps, and ncs_object must all match in crystal_symmetry, original (unit_cell) crystal_symmetry, and shift_cart for maps)

If map_manager contains an ncs_object and an ncs_object is supplied, the map_manager receives the supplied ncs_object

absolute_angle_tolerance and absolute_length_tolerance are tolerances for crystal_symmetry.is_similar_symmetry()

add_map_manager(map_manager)
add_model(model, set_model_log_to_null=True)
add_ncs_object(ncs_object)
as_map_model_manager()

Return map_model_manager object with contents of this class (not a deepcopy)

check_model_and_set_to_match_map_if_necessary()
crystal_symmetry()
deep_copy()
get_coordinate_shift(reverse=False)
map_manager()

Return the map_manager

model()
ncs_object()
read_map(file_name)
read_model(file_name)
read_ncs_file(file_name)
set_log(log=<colorama.ansitowin32.StreamWrapper object>)

Set output log file

set_original_origin_and_gridding(original_origin=None, gridding=None)

Use map_manager to reset (redefine) the original origin and gridding of the map. You can supply just the original origin in grid units, or just the gridding of the full unit_cell map, or both.

Update shift_cart for model and ncs object if present.

shift_model_to_match_original_map(model=None)
shift_model_to_match_working_map(model=None, reverse=False, coordinate_shift=None, new_shift_cart=None, final_crystal_symmetry=None, final_unit_cell_crystal_symmetry=None)

Shift a model based on the coordinate shift for the working map.

Also match the crystal_symmetry and unit_cell_crystal_symmetry

of the model to the map, unless specified as final_crystal_symmetry and final_unit_cell_crystal_symmetry.

Optionally specify the shift to apply (coordinate shift) and the new value of the shift recorded in the model (new_shift_cart)

shift_ncs_to_match_original_map(ncs_object=None)
shift_ncs_to_match_working_map(ncs_object=None, reverse=False, coordinate_shift=None, new_shift_cart=None)

Shift an ncs object to match the working map (based on self._map_manager.origin_shift_grid_units)

The working map is the current map in its current location. Typically origin is at (0,0,0).

This shifts an ncs object (typically is in its original location) to match this working map.

If the ncs object was already shifted (as reflected in its shift_cart()) it will receive the appropriate additional shift to match current map.

If coordinate_shift is specified, it is the target final coordinate shift instead of the shift_cart() for the working map.

shift_origin(desired_origin=(0, 0, 0))
show_summary()
unit_cell_crystal_symmetry()
write_map(file_name=None)
write_model(file_name=None, data_manager=None, format=None)
iotbx.map_model_manager.residue_group_is_linked_to_previous(rg, previous_rg)
class iotbx.map_model_manager.run_anisotropic_scaling_as_class(map_model_manager=None, direction_vectors=None, scale_factor_info=None, setup_info=None)
class iotbx.map_model_manager.run_fsc_as_class(map_model_manager=None, run_list=None, box_info=None)
iotbx.map_model_manager.set_nearby_empty_values(map_manager, set_values_map_manager, xyz_list, radius, value)
Set values within radii of xyz_list points to value if not already

set