dxtbx.model¶
dxtbx.model.detector¶
- class dxtbx.model.detector.DetectorFactory¶
Bases:
object
A factory class for detector objects, which will encapsulate standard detector designs to make it a little easier to get started with these. In cases where a CBF image is provided a full description can be used, in other cases assumptions will be made about the experiment configuration. In all cases information is provided in the CBF coordinate frame.
- static complex(sensor, origin, fast, slow, pixel, size, trusted_range=(0.0, 0.0), mask=[], px_mm=None, mu=0.0, gain=None, pedestal=None, identifier='')¶
A complex detector model, where you know exactly where everything is. This is useful for implementation of the Rigaku Saturn header format, as that is exactly what is in there. Origin, fast and slow are vectors in the CBF reference frame, pixel is the dimensions as a tuple as is size.
- static from_dict(d, t=None)¶
Convert the dictionary to a detector model
- Params:
d The dictionary of parameters t The template dictionary to use
- Returns:
The detector model
- static from_phil(params, reference=None, beam=None)¶
Convert phil parameters into detector model
- static generate_from_phil(params, beam=None)¶
Create a new detector model from phil parameters
- static imgCIF(cif_file, sensor)¶
Initialize a detector model from an imgCIF file.
- static imgCIF_H(cbf_handle, sensor)¶
Initialize a detector model from an imgCIF file handle, where it is assumed that the file has already been read.
- static make_detector(stype, fast_axis, slow_axis, origin, pixel_size, image_size, trusted_range=(0.0, 1000000.0), px_mm=None, name='Panel', thickness=0.0, material='', mu=0.0, gain=None, pedestal=None, identifier='')¶
Ensure all types are correct before creating c++ detector class.
- static overwrite_from_phil(params, detector, beam=None)¶
Overwrite from phil parameters
- static sensor(name)¶
Return the correct sensor token for a given name, for example:
ccd, CCD image_plate, IMAGE_PLATE pad, PAD
to the appropriate static token which will be used as a handle everywhere else in this. Also allow existing token to be passed in.
- static simple(sensor, distance, beam_centre, fast_direction, slow_direction, pixel_size, image_size, trusted_range=(0.0, 0.0), mask=[], px_mm=None, mu=0.0, gain=None, pedestal=None, identifier='')¶
Construct a simple detector at a given distance from the sample along the direct beam presumed to be aligned with -z, offset by the beam centre - the directions of which are given by the fast and slow directions, which are themselves given as +x, +y, -x, -y. The pixel size is given in mm in the fast and slow directions and the image size is given in pixels in the same order. Everything else is the same as for the main reference frame.
- static two_theta(sensor, distance, beam_centre, fast_direction, slow_direction, two_theta_direction, two_theta_angle, pixel_size, image_size, trusted_range=(0.0, 0.0), mask=[], px_mm=None, mu=0.0, gain=None, pedestal=None, identifier='')¶
Construct a simple detector at a given distance from the sample along the direct beam presumed to be aligned with -z, offset by the beam centre - the directions of which are given by the fast and slow directions, which are themselves given as +x, +y, -x, -y. The pixel size is given in mm in the fast and slow directions and the image size is given in pixels in the same order. Everything else is the same as for the main reference frame. Also given are the direction of the two-theta axis and the angle in degrees by which the detector is moved.
- dxtbx.model.detector.merge_panel_scope_extracts_by_id(panel_params)¶
dxtbx.model.beam¶
- class dxtbx.model.beam.BeamFactory¶
Bases:
object
A factory class for beam objects, which encapsulate standard beam models. In cases where a full cbf description is available this will be used, otherwise simplified descriptions can be applied.
- static complex(sample_to_source: Tuple[float, float, float], polarization_fraction: float, polarization_plane_normal: Tuple[float, float, float], wavelength: float) Beam ¶
Full access to the constructor for cases where we do know everything that we need…
- static from_dict(dict: dict, template: dict = None) Beam | PolychromaticBeam ¶
Convert the dictionary to a beam model
- static from_phil(params: scope_extract, reference: Beam | PolychromaticBeam | None = None) Beam | PolychromaticBeam ¶
Convert the phil parameters into a beam model
- static imgCIF(cif_file: str) Beam ¶
Initialize a detector model from an imgCIF file. N.B. the definition of the polarization plane is not completely helpful in this - it is the angle between the polarization plane and the +Y laboratory frame vector.
- static imgCIF_H(cbf_handle: cbf_handle_struct) Beam ¶
Initialize a detector model from an imgCIF file. N.B. the definition of the polarization plane is not completely helpful in this - it is the angle between the polarization plane and the +Y laboratory frame vector. This example works from a cbf_handle, which is already configured.
- static make_beam(sample_to_source: Tuple[float, float, float] = None, wavelength: float = None, s0: Tuple[float, float, float] = None, unit_s0: Tuple[float, float, float] = None, divergence: float = None, sigma_divergence: float = None) Beam ¶
- static make_polarized_beam(sample_to_source: Tuple[float, float, float] = None, wavelength: float = None, s0: Tuple[float, float, float] = None, unit_s0: Tuple[float, float, float] = None, polarization: Tuple[float, float, float] = None, polarization_fraction: float = None, divergence: float = None, sigma_divergence: float = None, flux: float = None, transmission: float = None, probe: Probe = dxtbx_model_ext.Probe.xray) Beam ¶
- static make_polychromatic_beam(direction: Tuple[float, float, float], divergence: float = 0.0, sigma_divergence: float = 0.0, polarization_normal: Tuple[float, float, float] = (0.0, 1.0, 0.0), polarization_fraction: float = 0.5, flux: float = 0.0, transmission: float = 1.0, probe: Probe = dxtbx_model_ext.Probe.xray, sample_to_source_distance: float = 0.0, deg: bool = True, wavelength_range: tuple[float, float] = (0.0, 0.0)) PolychromaticBeam ¶
- static simple(wavelength: float) Beam ¶
Construct a beam object on the principle that the beam is aligned with the +z axis, as is quite normal. Also assume the beam has polarization fraction 0.999 and is polarized in the x-z plane, unless it has a wavelength shorter than 0.05 Å in which case we assume electron diffraction and return an unpolarized beam model.
dxtbx.model.scan¶
- class dxtbx.model.scan.ScanFactory¶
Bases:
object
A factory for scan instances, to help with constructing the classes in a set of common circumstances.
- static add(scans)¶
Sum a list of scans wrapping the slightly clumsy idiomatic method: sum(scans[1:], scans[0]).
- static from_dict(d, t=None)¶
Convert the dictionary to a scan model
- Params:
d The dictionary of parameters t The template dictionary to use
- Returns:
The scan model
- static from_phil(params, reference=None)¶
Generate a scan model from phil parameters
- static imgCIF(cif_file)¶
Initialize a scan model from an imgCIF file.
- static imgCIF_H(cif_file, cbf_handle)¶
Initialize a scan model from an imgCIF file handle, where it is assumed that the file has already been read.
- static make_scan(image_range, exposure_times, oscillation, epochs, batch_offset=0, deg=True)¶
- static make_scan_from_properties(image_range, properties, batch_offset=0, deg=True)¶
- static search(filename)¶
Get a list of files which appear to match the template and directory implied by the input filename. This could well be used to get a list of image headers to read and hence construct scans from.
- static single_file(filename, exposure_times, osc_start, osc_width, epoch)¶
Construct an scan instance for a single image.
dxtbx.model.goniometer¶
- class dxtbx.model.goniometer.Goniometer¶
Bases:
GoniometerBase
- static from_dict((dict)arg1) Goniometer : ¶
- C++ signature :
dxtbx::model::Goniometer* from_dict(boost::python::dict)
- get_fixed_rotation((Goniometer)arg1) tuple : ¶
- C++ signature :
scitbx::mat3<double> get_fixed_rotation(dxtbx::model::Goniometer {lvalue})
- get_num_scan_points((Goniometer)arg1) int : ¶
- C++ signature :
unsigned long get_num_scan_points(dxtbx::model::Goniometer {lvalue})
- get_rotation_axis((Goniometer)arg1) tuple : ¶
- C++ signature :
scitbx::vec3<double> get_rotation_axis(dxtbx::model::Goniometer {lvalue})
- get_rotation_axis_datum((Goniometer)arg1) tuple : ¶
- C++ signature :
scitbx::vec3<double> get_rotation_axis_datum(dxtbx::model::Goniometer {lvalue})
- get_setting_rotation((Goniometer)arg1) tuple : ¶
- C++ signature :
scitbx::mat3<double> get_setting_rotation(dxtbx::model::Goniometer {lvalue})
- get_setting_rotation_at_scan_point((Goniometer)arg1, (object)arg2) tuple : ¶
- C++ signature :
scitbx::mat3<double> get_setting_rotation_at_scan_point(dxtbx::model::Goniometer {lvalue},unsigned long)
- get_setting_rotation_at_scan_points((Goniometer)arg1) mat3_double : ¶
- C++ signature :
scitbx::af::shared<scitbx::mat3<double> > get_setting_rotation_at_scan_points(dxtbx::model::Goniometer {lvalue})
- is_similar_to((Goniometer)arg1, (Goniometer)other[, (object)rotation_axis_tolerance=1e-06[, (object)fixed_rotation_tolerance=1e-06[, (object)setting_rotation_tolerance=1e-06]]]) bool : ¶
- C++ signature :
bool is_similar_to(dxtbx::model::Goniometer {lvalue},dxtbx::model::Goniometer [,double=1e-06 [,double=1e-06 [,double=1e-06]]])
- property num_scan_points¶
- reset_scan_points((Goniometer)arg1) None : ¶
- C++ signature :
void reset_scan_points(dxtbx::model::Goniometer {lvalue})
- rotate_around_origin((Goniometer)arg1, (object)axis, (object)angle[, (object)deg=True]) None : ¶
- C++ signature :
void rotate_around_origin(dxtbx::model::Goniometer {lvalue},scitbx::vec3<double>,double [,bool=True])
- set_fixed_rotation((Goniometer)arg1, (object)arg2) None : ¶
- C++ signature :
void set_fixed_rotation(dxtbx::model::Goniometer {lvalue},scitbx::mat3<double>)
- set_rotation_axis((Goniometer)arg1, (object)arg2) None : ¶
- C++ signature :
void set_rotation_axis(dxtbx::model::Goniometer {lvalue},scitbx::vec3<double>)
- set_rotation_axis_datum((Goniometer)arg1, (object)arg2) None : ¶
- C++ signature :
void set_rotation_axis_datum(dxtbx::model::Goniometer {lvalue},scitbx::vec3<double>)
- set_setting_rotation((Goniometer)arg1, (object)arg2) None : ¶
- C++ signature :
void set_setting_rotation(dxtbx::model::Goniometer {lvalue},scitbx::mat3<double>)
- set_setting_rotation_at_scan_points((Goniometer)arg1, (mat3_double)arg2) None : ¶
- C++ signature :
void set_setting_rotation_at_scan_points(dxtbx::model::Goniometer {lvalue},scitbx::af::const_ref<scitbx::mat3<double>, scitbx::af::trivial_accessor>)
set_setting_rotation_at_scan_points( (Goniometer)arg1, (tuple)arg2) -> None :
- C++ signature :
void set_setting_rotation_at_scan_points(dxtbx::model::Goniometer {lvalue},boost::python::tuple)
set_setting_rotation_at_scan_points( (Goniometer)arg1, (list)arg2) -> None :
- C++ signature :
void set_setting_rotation_at_scan_points(dxtbx::model::Goniometer {lvalue},boost::python::list)
- to_dict((Goniometer)arg1) dict : ¶
- C++ signature :
boost::python::dict to_dict(dxtbx::model::Goniometer)
- class dxtbx.model.goniometer.GoniometerFactory¶
Bases:
object
A factory class for goniometer objects, which will encapsulate some standard goniometer designs to make it a little easier to get started with all of this - for cases when we are not using a CBF. When we have a CBF just use that factory method and everything will be peachy.
- static from_dict(d, t=None)¶
Convert the dictionary to a goniometer model
- Params:
d The dictionary of parameters t The template dictionary to use
- Returns:
The goniometer model
- static from_phil(params, reference=None)¶
Convert the phil parameters into a goniometer model
- static imgCIF(cif_file)¶
Initialize a goniometer model from an imgCIF file.
- static imgCIF_H(cbf_handle)¶
Initialize a goniometer model from an imgCIF file handle, where it is assumed that the file has already been read.
- static kappa(alpha, omega, kappa, phi, direction, scan_axis)¶
Return a kappa goniometer where omega is the primary axis (i,e. aligned with X in the CBF coordinate frame) and has the kappa arm with angle alpha attached to it, aligned with -z, +y, +z or -y at omega = 0, that being the direction, which in turn has phi fixed to it which should initially be coincident with omega. We also need to know which axis is being used for the scan i.e. phi or omega. All angles should be given in degrees. This will work by first constructing the rotation axes and then composing them to the scan axis and fixed component of the rotation.
- static known_axis(axis)¶
Return an goniometer instance for a known rotation axis, assuming that nothing is known about the fixed element of the rotation axis.
- static make_goniometer(rotation_axis, fixed_rotation)¶
- static make_kappa_goniometer(alpha, omega, kappa, phi, direction, scan_axis)¶
- static make_multi_axis_goniometer(axes, angles, names, scan_axis)¶
- static multi_axis(axes, angles, names, scan_axis)¶
- static multi_axis_goniometer_from_phil(params, reference=None)¶
- static single_axis()¶
Construct a single axis goniometer which is canonical in the CBF reference frame.
- static single_axis_goniometer_from_phil(params, reference=None)¶
Generate or overwrite a single axis goniometer
- static single_axis_reverse()¶
Construct a single axis goniometer which is canonical in the CBF reference frame, but reversed in rotation.
- class dxtbx.model.goniometer.KappaGoniometer¶
Bases:
Goniometer
- get_alpha_angle((KappaGoniometer)arg1) float : ¶
- C++ signature :
double get_alpha_angle(dxtbx::model::KappaGoniometer {lvalue})
- get_direction((KappaGoniometer)arg1) KappaDirection : ¶
- C++ signature :
dxtbx::model::KappaGoniometer::Direction get_direction(dxtbx::model::KappaGoniometer {lvalue})
- get_kappa_angle((KappaGoniometer)arg1) float : ¶
- C++ signature :
double get_kappa_angle(dxtbx::model::KappaGoniometer {lvalue})
- get_kappa_axis((KappaGoniometer)arg1) tuple : ¶
- C++ signature :
scitbx::vec3<double> get_kappa_axis(dxtbx::model::KappaGoniometer {lvalue})
- get_omega_angle((KappaGoniometer)arg1) float : ¶
- C++ signature :
double get_omega_angle(dxtbx::model::KappaGoniometer {lvalue})
- get_omega_axis((KappaGoniometer)arg1) tuple : ¶
- C++ signature :
scitbx::vec3<double> get_omega_axis(dxtbx::model::KappaGoniometer {lvalue})
- get_phi_angle((KappaGoniometer)arg1) float : ¶
- C++ signature :
double get_phi_angle(dxtbx::model::KappaGoniometer {lvalue})
- get_phi_axis((KappaGoniometer)arg1) tuple : ¶
- C++ signature :
scitbx::vec3<double> get_phi_axis(dxtbx::model::KappaGoniometer {lvalue})
- get_scan_axis((KappaGoniometer)arg1) KappaScanAxis : ¶
- C++ signature :
dxtbx::model::KappaGoniometer::ScanAxis get_scan_axis(dxtbx::model::KappaGoniometer {lvalue})
- class dxtbx.model.goniometer.MultiAxisGoniometer¶
Bases:
Goniometer
- static from_dict((dict)arg1) MultiAxisGoniometer : ¶
- C++ signature :
dxtbx::model::MultiAxisGoniometer* from_dict(boost::python::dict)
- get_angles((MultiAxisGoniometer)arg1) double : ¶
- C++ signature :
scitbx::af::shared<double> get_angles(dxtbx::model::MultiAxisGoniometer {lvalue})
- get_axes((MultiAxisGoniometer)arg1) vec3_double : ¶
- C++ signature :
scitbx::af::shared<scitbx::vec3<double> > get_axes(dxtbx::model::MultiAxisGoniometer {lvalue})
- get_names((MultiAxisGoniometer)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > get_names(dxtbx::model::MultiAxisGoniometer {lvalue})
- get_scan_axis((MultiAxisGoniometer)arg1) int : ¶
- C++ signature :
unsigned long get_scan_axis(dxtbx::model::MultiAxisGoniometer {lvalue})
- set_angles((MultiAxisGoniometer)arg1, (double)arg2) None : ¶
- C++ signature :
void set_angles(dxtbx::model::MultiAxisGoniometer {lvalue},scitbx::af::shared<double>)
- set_axes((MultiAxisGoniometer)arg1, (vec3_double)arg2) None : ¶
- C++ signature :
void set_axes(dxtbx::model::MultiAxisGoniometer {lvalue},scitbx::af::const_ref<scitbx::vec3<double>, scitbx::af::trivial_accessor>)
- to_dict((MultiAxisGoniometer)arg1) dict : ¶
- C++ signature :
boost::python::dict to_dict(dxtbx::model::MultiAxisGoniometer)
dxtbx.model.crystal¶
- class dxtbx.model.crystal.CrystalFactory¶
Bases:
object
- static from_dict(d, t=None)¶
Convert the dictionary to a crystal model
- Params:
d The dictionary of parameters t The template dictionary to use
- Returns:
The crystal model
- static from_mosflm_matrix(mosflm_A_matrix, unit_cell=None, wavelength=None, space_group=None)¶
Create a crystal_model from a Mosflm A matrix (a*, b*, c*); N.B. assumes the mosflm coordinate frame:
/! Y-axis / ! ^ / ! ! / ! ! / ! ! / / Xd ! ! / / * ^ ! ! / ! 3 ! ! !/ X-ray beam ! ! ! /------------------------/--!---->X-axis / ! / *1! <-/- ! / ! / \+ve phi ! Yd / / / ! 2 / / ! * / Z-axis Ys ^ _/ Rotation ! /| Xs axis !/ O
Also assume that the mosaic spread is 0.
- Parameters:
mosflm_A_matrix (tuple of floats) – The A matrix in Mosflm convention.
unit_cell (cctbx.uctbx.unit_cell) – The unit cell parameters which are used to determine the wavelength from the Mosflm A matrix.
wavelength (float) – The wavelength to scale the A matrix
space_group (cctbx.sgtbx.space_group) – If the space group is None then the space_group will be assigned as P1
- Returns:
A crystal model derived from the given Mosflm A matrix
- Return type:
crystal_model