mmtbx.maps package

Submodules

mmtbx.maps.composite_omit_map module

Methods and utilities for generating composite omit maps. The actual end-user application, phenix.composite_omit_map, is part of the phenix sources.

mmtbx.maps.composite_omit_map.combine_maps(map_arrays, omit_groups, background_map_coeffs, resolution_factor, flatten_background=False, sigma_scaling=False, control_map=False)

For each box, FFT the corresponding omit map coefficients, extract the omit regions, and copy them to the combined map, using Marat’s asymmetric map class to handle symmetry expansion internally.

mmtbx.maps.composite_omit_map.create_omit_regions(xray_structure, selection=None, fraction_omit=0.05, optimize_binning=True, box_cushion_radius=2.5, even_boxing=False, asu_buffer_thickness=1e-05, log=None)

Divide the asymmetric unit into boxes of atoms to omit. This will include a cushion around the actual omit region. Although the step size is uniform, by default boxes will be grouped as needed to make the distribution roughly even and keep the total number of omit regions to a minimum. If a region does not contain any atoms, it will be skipped and filled in with the background map.

Returns a list of omit_regions objects, which specify the selection of atoms to omit along with the boxes of interest (in fractional coordinates).

mmtbx.maps.composite_omit_map.grid_max(grid_coords, n_real)
mmtbx.maps.composite_omit_map.join_selections(sel1, sel2)
class mmtbx.maps.composite_omit_map.omit_box(**kwds)

Bases: slots_getstate_setstate_default_initializer

Defines a region in fractional coordinates containing a selection of atoms to be omitted.

frac_max
frac_min
property n_atoms
selection
serial
show(out=<colorama.ansitowin32.StreamWrapper object>, prefix='')
class mmtbx.maps.composite_omit_map.omit_region_results(**kwds)

Bases: slots_getstate_setstate_default_initializer

map_coeffs_list
r_free
r_work
serial
class mmtbx.maps.composite_omit_map.omit_regions(serial, selection=None)

Bases: slots_getstate_setstate

Groups together multiple omit_box objects (without any implicit spatial relationship); this is used to reduce the total number of bins of omitted atoms and to make the fraction omitted per bin approximately similar.

add_box(box)
boxes
combine_with(other)
property n_atoms
property n_boxes
selection
serial
show(out=<colorama.ansitowin32.StreamWrapper object>, prefix='')
class mmtbx.maps.composite_omit_map.run(fmodel, crystal_gridding, box_size_as_fraction=0.1, max_boxes=2000, neutral_volume_box_cushion_width=2, full_resolution_map=True, log=<colorama.ansitowin32.StreamWrapper object>)

Bases: object

Composite residual OMIT map. Details:

Afonine et al. (2014). FEM: Feature Enhanced Map

as_p1_map(map_asu=None)
asu_map_from_fmodel(fmodel, map_type)
box_iterator()
get_p1_map_unscaled(fmodel, map_type)
map_coefficients(filter_noise=True)
noise_filtered_map_coefficients(map_coefficients)
omit_box(start, end)
result_as_sf()
to_asu_map(map_data)

mmtbx.maps.fem module

mmtbx.maps.fem.ccp4_map(cg, file_name, mc=None, map_data=None)
class mmtbx.maps.fem.counter(n1, n2, log)

Bases: object

show()
mmtbx.maps.fem.get_map(mc, cg)
mmtbx.maps.fem.good_atoms_selection(crystal_gridding, map_coeffs, xray_structure)
mmtbx.maps.fem.inner_loop(fmodel, wam, missing, crystal_gridding, n, progress_counter, use_max_map)
mmtbx.maps.fem.low_volume_density_elimination(m, fmodel, selection, end=11)
class mmtbx.maps.fem.run(f_obs, r_free_flags, xray_structure, use_resolve, use_omit, use_max_map, sharp, use_unsharp_masking, resolution_factor, n_inner_loop=10, n_outer_loop=10, log=None)

Bases: object

compute_original_map()
create_fmodel(update_f_part1=True, show=False)
outer_loop()
prepare_f_obs_and_flags()
remove_common_isotropic_adp()
resolve_filter_map()
write_output_files(mtz_file_name, ccp4_map_file_name, fem_label, orig_label)
zero_below_threshold(m)
mmtbx.maps.fem.truncate_with_roots(m, fmodel, c1, c2, cutoff, scale, zero_all_interblob_region=True, as_int=False, average_peak_volume=None, selection=None)

mmtbx.maps.fobs_minus_fobs_map module

class mmtbx.maps.fobs_minus_fobs_map.compute_fo_minus_fo_map(data_arrays, xray_structure, log=None, silent=False, output_file=None, peak_search=False, map_cutoff=None, peak_search_params=None, r_free_arrays=None, write_map=True, multiscale=False, anomalous=False)

Bases: object

write_map_file()
mmtbx.maps.fobs_minus_fobs_map.finish_job(result)
mmtbx.maps.fobs_minus_fobs_map.fo_minus_fo_master_params()
class mmtbx.maps.fobs_minus_fobs_map.launcher(args, file_name, output_dir=None, log_file=None, job_title=None)

Bases: target_with_save_result

run()
mmtbx.maps.fobs_minus_fobs_map.run(args, command_name='phenix.fobs_minus_fobs_map', log=None)
mmtbx.maps.fobs_minus_fobs_map.validate_params(params, callback=None)

mmtbx.maps.kick module

mmtbx.maps.kick.compute_map_and_combine(map_coeffs, crystal_gridding, map_data, thresholds=[0, 0.1, 0.2, 0.3, 0.4, 0.5])
mmtbx.maps.kick.get_map(map_coeffs, crystal_gridding)
mmtbx.maps.kick.kick_fmodel(fmodel, map_type, crystal_gridding, number_of_kicks, macro_cycles, missing=None, kick_completeness=0.95)
mmtbx.maps.kick.kick_map_coeffs(map_coeffs, crystal_gridding, number_of_kicks, macro_cycles, missing=None, kick_completeness=0.95, phases_only=False)
mmtbx.maps.kick.randomize_completeness(map_coeffs, crystal_gridding, n_cycles=10, thresholds=[0, 0.1, 0.2, 0.3, 0.4, 0.5])
mmtbx.maps.kick.randomize_completeness2(map_coeffs, crystal_gridding, n_cycles=10)
mmtbx.maps.kick.randomize_struture_factors(map_coeffs, number_of_kicks, phases_only=False)
class mmtbx.maps.kick.run(fmodel, map_type='2mFo-DFc', mask_data=None, crystal_gridding=None, number_of_kicks=100, macro_cycles=10, kick_completeness=0.95, omit=True)

Bases: object

Note 1: Assume Fmodel has correct scaling already (all f_parts etc). Note 2: More macro_cycles can substantially clean map in extremely bad cases.

convert_to_non_anomalous(fmodel)
map_coeffs_from_map(map_data)
write_mc(file_name='mc.mtz')
class mmtbx.maps.kick.weighted_average(fmodel, map_coefficients)

Bases: object

random_weight_averaged_map_coefficients(random_scale, random_seed, n_cycles, fraction_keep, missing=None)

mmtbx.maps.superpose module

class mmtbx.maps.superpose.common_frame_of_reference(all_sites_cart, lsq_fits, buffer=10.0, move_to_frame_of_reference=True, log=<colorama.ansitowin32.StreamWrapper object>)

Bases: object

apply_new_sites_to_hierarchies(pdb_hierarchies)
apply_new_sites_to_xray_structure(xray_structures)
get_inverse_matrix(i)
inverse_transform_hierarchy(i, pdb_hierarchy)
mmtbx.maps.superpose.exercise()
mmtbx.maps.superpose.generate_p1_box(pdb_hierarchy, buffer=10.0)
mmtbx.maps.superpose.mask_grid(xrs, buffer, map_data, n_real)
mmtbx.maps.superpose.transform_map_by_lsq_fit(fft_map, unit_cell, lsq_fit_obj, pdb_hierarchy, d_min, buffer=10.0, resolution_factor=0.25, file_name=None, log=<colorama.ansitowin32.StreamWrapper object>, mask_map_grid=False, format='xplor')
class mmtbx.maps.superpose.transform_maps(map_coeffs, map_types, pdb_hierarchy, unit_cell, lsq_fit_obj, output_files, resolution_factor=0.25, buffer=10.0, auto_run=True, format='ccp4')

Bases: object

run()

mmtbx.maps.utils module

mmtbx.maps.utils.create_map_from_pdb_and_mtz(pdb_file, mtz_file, output_file, fill=False, out=None, llg_map=False, remove_unknown_scatterering_type=False, assume_pdb_data=False)

Convenience function, used by phenix.fetch_pdb.

Parameters:

remove_unknown_scatterering_type – if True, atoms with unknown scattering type will be removed before extracting the X-ray structure (which would crash otherwise).

mmtbx.maps.utils.decode_resolve_map_coeffs(*args, **kwds)
mmtbx.maps.utils.extract_map_coeffs(*args, **kwds)
mmtbx.maps.utils.extract_phenix_refine_map_coeffs(*args, **kwds)
class mmtbx.maps.utils.fast_maps_from_hkl_file(file_name, xray_structure, scattering_table='wk1995', f_label=None, r_free_label=None, map_out=None, log=<colorama.ansitowin32.StreamWrapper object>, auto_run=True, quiet=False, anomalous_map=False, fill_maps=True, save_fmodel=False, llg_map=False, assume_pdb_data=False)

Bases: object

Automation wrapper for calculating standard 2mFo-DFc, mFo-DFc, and optional anomalous difference map coefficients, starting from an X-ray structure and an MTZ file name. This will attempt to guess the correct data array automatically (giving preference to anomalous amplitudes).

The organization of this class is somewhat messy, as it is designed to be run in parallel for multiple structures.

get_maps_from_fmodel()
run()
mmtbx.maps.utils.format_map_coeffs_for_resolve(*args, **kwds)
class mmtbx.maps.utils.generate_water_omit_map(fmodel, pdb_hierarchy, skip_if_no_waters=False, exclude_free_r_reflections=False, fill_missing_f_obs=False, write_f_model=False, log=None)

Bases: object

Calculate standard map coefficients (with R-free flags omitted, and missing reflections filled) with water occupancies set to zero. Used for ligand fitting after refinement.

write_map_coeffs(file_name)
write_pdb_file(file_name)
mmtbx.maps.utils.get_anomalous_map(fmodel)
mmtbx.maps.utils.get_default_map_coefficients_phil(show_filled=True, disable_filled=False, exclude_free_r_reflections=False, enclosing_scope=None, as_phil_object=True)

Generate the map coefficients PHIL string for phenix.maps, phenix.refine, and related applications, based on several toggles which control the handling of missing reflections and R-free flags. Used by the Phenix GUI.

mmtbx.maps.utils.get_map_coeff_labels(*args, **kwds)
mmtbx.maps.utils.get_map_coeffs_for_build(server)
mmtbx.maps.utils.get_maps_from_fmodel(fmodel, fill_missing_f_obs=True, exclude_free_r_reflections=False)
mmtbx.maps.utils.map_coeffs_from_mtz_file(*args, **kwds)
mmtbx.maps.utils.write_map_coeffs(*args, **kwds)
mmtbx.maps.utils.write_xplor_map(sites_cart, unit_cell, map_data, n_real, file_name, buffer=10)
mmtbx.maps.utils.write_xplor_map_file(coeffs, frac_min, frac_max, file_base)

Module contents

mmtbx.maps.b_factor_sharpening_by_map_kurtosis_maximization(map_coeffs, show=True, b_sharp_best=None, b_only=False, b_min=-100, b_max=100, b_step=5)
mmtbx.maps.cast_map_coeff_params(map_type_obj)
mmtbx.maps.compute_f_calc(fmodel, params)
class mmtbx.maps.compute_map_coefficients(fmodel, params, mtz_dataset=None, post_processing_callback=None, pdb_hierarchy=None, log=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>)

Bases: object

write_mtz_file(file_name, mtz_history_buffer=None, r_free_flags=None)
mmtbx.maps.compute_xplor_maps(fmodel, params, atom_selection_manager=None, file_name_prefix=None, file_name_base=None, post_processing_callback=None, pdb_hierarchy=None)
mmtbx.maps.map_and_map_coeff_master_params()
mmtbx.maps.map_coefficients_from_fmodel(params, fmodel=None, map_calculation_server=None, post_processing_callback=None, pdb_hierarchy=None)
class mmtbx.maps.map_coeffs_mtz_label_manager(map_params)

Bases: object

amplitudes()
phases(root_label, anomalous_sign=None)
mmtbx.maps.maps_including_IO_master_params()
class mmtbx.maps.write_xplor_map_file(params, coeffs, atom_selection_manager=None, xray_structure=None)

Bases: object

atom_iselection()
box_around_selection(iselection, buffer)