MX Calibrate#

Calibrate a translation table from a set of powder diffraction images taken at various sample-detector distances.

This is a notebook replacement of the MX-calibrate tool from pyFAI with advanced features.

Start with some constant definition:#

[1]:
calibrant_name = "CeO2"
detector_name = "Pilatus 2M"
file_pattern = "massif1/test-powder*.cbf"
result_file = "MX-calibrate.json"
wavelength = None # set a value to override the one in the headers
[2]:
!wget http://www.silx.org/pub/pyFAI/massif1.tar.bz2
!tar -xvjf massif1.tar.bz2
--2025-11-16 11:22:24--  http://www.silx.org/pub/pyFAI/massif1.tar.bz2
Resolving www.silx.org (www.silx.org)... 195.154.237.27
Connecting to www.silx.org (www.silx.org)|195.154.237.27|:80... connected.
HTTP request sent, awaiting response... 200 OK
Length: 6784503 (6.5M) [application/x-bzip2]
Saving to: ‘massif1.tar.bz2.10’

massif1.tar.bz2.10  100%[===================>]   6.47M  22.9MB/s    in 0.3s

2025-11-16 11:22:25 (22.9 MB/s) - ‘massif1.tar.bz2.10’ saved [6784503/6784503]

massif1/
massif1/test-powder_5_0001.poni
massif1/test-powder_4_0001.cbf
massif1/test-powder_6_0001.poni
massif1/test-powder_8_0001.cbf
massif1/test-powder_7_0001.cbf
massif1/test-powder_7_0001.poni
massif1/test-powder_3_0001.cbf
massif1/test-powder_1_0001.poni
massif1/test-powder_4_0001.poni
massif1/test-powder_8_0001.poni
massif1/test-powder_2_0001.poni
massif1/test-powder_6_0001.cbf
massif1/test-powder_2_0001.cbf
massif1/test-powder_5_0001.cbf
massif1/test-powder_3_0001.poni
massif1/test-powder_1_0001.cbf
[3]:
%matplotlib widget
#inline2
[4]:
import os
import glob
import logging
import numpy
from matplotlib.pyplot import subplots
from scipy.stats import linregress
import fabio
import pyFAI
from pyFAI.gui import jupyter
import pyFAI.calibrant
from pyFAI.gui.jupyter.calib import Calibration
from pyFAI.goniometer import GeometryTransformation, GoniometerRefinement
from pyFAI.gui.cli_calibration import AbstractCalibration
import pyFAI.gui.mpl_calib
pyFAI.gui.mpl_calib.logger.setLevel(logging.ERROR)
print(f"Running pyFAI version {pyFAI.version}")
WARNING:pyFAI.gui.matplotlib:Matplotlib already loaded with backend `widget`, setting its backend to `QtAgg` may not work!
Running pyFAI version 2025.11.0-dev0
[5]:
detector = pyFAI.detector_factory(detector_name)
calibrant = pyFAI.calibrant.get_calibrant(calibrant_name)
files = sorted(glob.glob(file_pattern))
print("input files: "+" ".join(files))
input files: massif1/test-powder_1_0001.cbf massif1/test-powder_2_0001.cbf massif1/test-powder_3_0001.cbf massif1/test-powder_4_0001.cbf massif1/test-powder_5_0001.cbf massif1/test-powder_6_0001.cbf massif1/test-powder_7_0001.cbf massif1/test-powder_8_0001.cbf
[6]:
first = fabio.open(files[0])

def get_dectris_headers(fimg):
    """return the dectris headers from a Pilatus detector"""
    res = {}
    for line in fimg.header.get("_array_data.header_contents", "").split("\n"):
        words = line.split()
        if len(words)>=3:
            key = words[1]
            for v in words[2:]:
                try:
                    vf = float(v)
                except:
                    continue
                if not("." in v or "e" in v):
                    vf = int(v)
                res[key] = vf
    return res

get_dectris_headers(first)
[6]:
{'Silicon': 0.00045,
 'Pixel_size': 0.000172,
 'N_oscillations': 1,
 'Chi': 0.0,
 'Phi': 0.0,
 'Kappa': 0.0,
 'Alpha': 0.0,
 'Polarization': 0.99,
 'Detector_2theta': 0.0,
 'Angle_increment': 1.0,
 'Transmission': 100.0,
 'Flux': 436215830143.2828,
 'Detector_Voffset': 0.0,
 'Detector_distance': 0.126474,
 'Wavelength': 0.965459,
 'N_excluded_pixels:': 321,
 'Threshold_setting': 6421,
 'Count_cutoff': 1048500,
 'Tau': 0,
 'Exposure_period': 0.02115,
 'Exposure_time': 0.02,
 'Start_angle': 0.0}
[7]:
if wavelength is None:
    wavelength = get_dectris_headers(first)["Wavelength"] * 1e-10
calibrant.wavelength = wavelength
[8]:
#apply mask to the detector
mask = numpy.logical_or(detector.mask, first.data<0)
detector.mask = mask

Manual calibration of the first image#

[9]:
# Important: select the ring number before right-click on the ring. Finally click on the refine button
calib = Calibration(img=first.data,
                    mask=mask,
                    detector=detector,
                    wavelength=wavelength,
                    calibrant=calibrant)
[10]:
input("Break automatic execution ... press enter after picking at least 2 rings!")
---------------------------------------------------------------------------
KeyboardInterrupt                         Traceback (most recent call last)
Cell In[10], line 1
----> 1 input("Break automatic execution ... press enter after picking at least 2 rings!")

File ~/.venv/py313/lib/python3.13/site-packages/ipykernel/kernelbase.py:1282, in Kernel.raw_input(self, prompt)
   1280     msg = "raw_input was called, but this frontend does not support input requests."
   1281     raise StdinNotImplementedError(msg)
-> 1282 return self._input_request(
   1283     str(prompt),
   1284     self._parent_ident["shell"],
   1285     self.get_parent("shell"),
   1286     password=False,
   1287 )

File ~/.venv/py313/lib/python3.13/site-packages/ipykernel/kernelbase.py:1325, in Kernel._input_request(self, prompt, ident, parent, password)
   1322 except KeyboardInterrupt:
   1323     # re-raise KeyboardInterrupt, to truncate traceback
   1324     msg = "Interrupted by user"
-> 1325     raise KeyboardInterrupt(msg) from None
   1326 except Exception:
   1327     self.log.warning("Invalid Message:", exc_info=True)

KeyboardInterrupt: Interrupted by user
Break automatic execution ... press enter after picking at least 2 rings!
[11]:
calib.extract_cpt()
# calib.geoRef.rot1 = calib.geoRef.rot2 = calib.geoRef.rot3 = 0
calib.refine(fixed=["wavelength", "rot3"])
Before refinement, the geometry is:
Detector Pilatus 2M      PixelSize= 172µm, 172µm         BottomRight (3)
Wavelength= 0.965459 Å
SampleDetDist= 1.262744e-01 m   PONI= 1.471304e-01, 1.208236e-01 m      rot1=0.003736  rot2=-0.003450  rot3=0.000000 rad
DirectBeamDist= 126.276 mm      Center: x=699.720, y=852.877 pix        Tilt= 0.291° tiltPlanRotation= -137.277° λ= 0.965Å
Detector Pilatus 2M      PixelSize= 172µm, 172µm         BottomRight (3)
Wavelength= 0.965459 Å
SampleDetDist= 1.263307e-01 m   PONI= 1.468959e-01, 1.209421e-01 m      rot1=0.004413  rot2=-0.002006  rot3=0.000000 rad
DirectBeamDist= 126.332 mm      Center: x=699.910, y=852.573 pix        Tilt= 0.278° tiltPlanRotation= -155.555° λ= 0.965Å
Detector Pilatus 2M      PixelSize= 172µm, 172µm         BottomRight (3)
Wavelength= 0.965459 Å
SampleDetDist= 1.263307e-01 m   PONI= 1.468959e-01, 1.209421e-01 m      rot1=0.004413  rot2=-0.002006  rot3=0.000000 rad
DirectBeamDist= 126.332 mm      Center: x=699.910, y=852.573 pix        Tilt= 0.278° tiltPlanRotation= -155.555° λ= 0.965Å
Detector Pilatus 2M      PixelSize= 172µm, 172µm         BottomRight (3)
Wavelength= 0.965459 Å
SampleDetDist= 1.263307e-01 m   PONI= 1.468959e-01, 1.209421e-01 m      rot1=0.004413  rot2=-0.002006  rot3=0.000000 rad
DirectBeamDist= 126.332 mm      Center: x=699.910, y=852.573 pix        Tilt= 0.278° tiltPlanRotation= -155.555° λ= 0.965Å
[12]:
#Return to `inline` mode
%matplotlib inline
calib.peakPicker.widget.fig.show()

Check that the beam-center and the distance is correct and how much they are off.

Calibration of all frames in automatic mode:#

[13]:
# Definition of the geometry translation function:
get_distance = lambda fimg: get_dectris_headers(fimg)["Detector_distance"]

geotrans = GeometryTransformation(param_names = ["dist_offset",
                                                 "poni1", "poni2", "rot1","rot2",
                                                "dist_scale", "poni1_scale", "poni2_scale"],
                                  dist_expr="pos * dist_scale + dist_offset",
                                  poni1_expr="pos * poni1_scale + poni1",
                                  poni2_expr="pos * poni2_scale + poni2",
                                  rot1_expr="rot1",
                                  rot2_expr="rot2",
                                  rot3_expr="0.0")

param = {
         "dist_offset": calib.geoRef.dist-get_distance(first),
         "poni1": calib.geoRef.poni1,
         "poni2": calib.geoRef.poni2,
         "rot1": calib.geoRef.rot1,
         "rot2": calib.geoRef.rot2,
         "dist_scale": 1.0,
         "poni1_scale": 0.0,
         "poni2_scale": 0.0,
}

print(param)
{'dist_offset': np.float64(-0.00014332792633586777), 'poni1': np.float64(0.14689593764120265), 'poni2': np.float64(0.12094209262835975), 'rot1': np.float64(0.0044133715540227905), 'rot2': np.float64(-0.0020062119420559553), 'dist_scale': 1.0, 'poni1_scale': 0.0, 'poni2_scale': 0.0}
[14]:
# Definition of the geometry refinement: the parameter order is the same as the param_names



gonioref = GoniometerRefinement(param, #initial guess
                                pos_function=get_distance,
                                trans_function=geotrans,
                                detector=detector,
                                wavelength=wavelength)
print("Empty refinement object:")
print(gonioref)
Empty refinement object:
GoniometerRefinement with 0 geometries labeled: .
[15]:
# Let's populate the goniometer refinement object with all geometries:
for fn in files:
    base = os.path.splitext(fn)[0]
    fimg = fabio.open(fn)
    local_calib = AbstractCalibration(img=fimg.data, mask=mask,
                                      detector=detector,
                                      wavelength=wavelength,
                                      calibrant=calibrant)
    local_calib.preprocess()
    local_calib.fixed = ["wavelength", "rot3"]
    local_calib.ai = gonioref.get_ai(get_distance(fimg))
    local_calib.extract_cpt()
    sg = gonioref.new_geometry(os.path.basename(base), image=fimg.data, metadata=fimg,
                              control_points=local_calib.peakPicker.points,
                              geometry=local_calib.ai,
                              calibrant=calibrant)

print("Filled refinement object:")
print(gonioref)
print(os.linesep+"\tLabel \t Distance")
for k, v in gonioref.single_geometries.items():
    print(k,v.get_position())
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
ERROR:root:No diffraction image available => not showing the contour
Filled refinement object:
GoniometerRefinement with 8 geometries labeled: test-powder_1_0001, test-powder_2_0001, test-powder_3_0001, test-powder_4_0001, test-powder_5_0001, test-powder_6_0001, test-powder_7_0001, test-powder_8_0001.

        Label    Distance
test-powder_1_0001 0.126474
test-powder_2_0001 0.141749
test-powder_3_0001 0.199249
test-powder_4_0001 0.171074
test-powder_5_0001 0.226674
test-powder_6_0001 0.293162
test-powder_7_0001 0.357899
test-powder_8_0001 0.484611
[16]:
fig, ax = subplots(len(files), figsize=(10, 10*len(files)))
for sp, sg in zip(ax, gonioref.single_geometries.values()):
    jupyter.display(sg=sg, ax=sp, label=sg.label)
../../../../_images/usage_tutorial_Goniometer_MX-calibrate_MX-calibrate_18_0.png
../../../../_images/usage_tutorial_Goniometer_MX-calibrate_MX-calibrate_18_1.png
[17]:
gonioref.refine3(fix=["dist_scale", "poni1_scale", "poni2_scale"])
Free parameters: ['dist_offset', 'poni1', 'poni2', 'rot1', 'rot2']
Fixed: {'dist_scale': 1.0, 'poni1_scale': 0.0, 'poni2_scale': 0.0}
 message: Optimization terminated successfully
 success: True
  status: 0
     fun: 2.4377579977009444e-07
       x: [-8.682e-05  1.468e-01  1.208e-01  3.353e-03 -1.419e-03]
     nit: 10
     jac: [ 2.276e-06 -1.443e-06 -3.045e-07  1.175e-07 -3.283e-07]
    nfev: 64
    njev: 10
Constrained Least square 3.8224925736402864e-07 --> 2.4377579977009444e-07
maxdelta on rot1: 0.0044133715540227905 --> 0.0033527499190415286
[17]:
np.float64(2.4377579977009444e-07)
[18]:
gonioref.refine3(fix=[])
Free parameters: ['dist_offset', 'poni1', 'poni2', 'rot1', 'rot2', 'dist_scale', 'poni1_scale', 'poni2_scale']
Fixed: {}
 message: Optimization terminated successfully
 success: True
  status: 0
     fun: 1.8046633143817243e-07
       x: [-3.340e-04  1.468e-01  1.207e-01  4.827e-03 -1.658e-03
            1.002e+00  4.483e-04  2.257e-03]
     nit: 16
     jac: [-2.004e-09 -6.613e-09 -2.448e-10 -3.260e-08 -1.599e-07
           -2.526e-09  2.587e-07 -4.498e-08]
    nfev: 147
    njev: 16
Constrained Least square 2.4377579977009444e-07 --> 1.8046633143817243e-07
maxdelta on poni2_scale: 0.0 --> 0.002256950050514632
[18]:
np.float64(1.8046633143817243e-07)

Interpretation of this fit:

[19]:
gonioref.get_ai(0.2)
[19]:
Detector Pilatus 2M      PixelSize= 172µm, 172µm         BottomRight (3)
Wavelength= 0.965459 Å
SampleDetDist= 1.999847e-01 m   PONI= 1.468775e-01, 1.211879e-01 m      rot1=0.004827  rot2=-0.001658  rot3=0.000000 rad
DirectBeamDist= 199.987 mm      Center: x=698.969, y=852.010 pix        Tilt= 0.292° tiltPlanRotation= -161.037° λ= 0.965Å
[20]:
gonioref.get_ai(0.3)
[20]:
Detector Pilatus 2M      PixelSize= 172µm, 172µm         BottomRight (3)
Wavelength= 0.965459 Å
SampleDetDist= 3.001440e-01 m   PONI= 1.469223e-01, 1.214136e-01 m      rot1=0.004827  rot2=-0.001658  rot3=0.000000 rad
DirectBeamDist= 300.148 mm      Center: x=697.470, y=851.305 pix        Tilt= 0.292° tiltPlanRotation= -161.037° λ= 0.965Å

Persistence of this fit

[21]:
gonioref.save(result_file)
[22]:
with open(result_file) as r:
    print(r.read())
{
  "content": "Goniometer calibration v2",
  "detector": "Pilatus 2M",
  "detector_config": {
    "pixel1": 0.000172,
    "pixel2": 0.000172,
    "orientation": 3
  },
  "wavelength": 9.65459e-11,
  "param": [
    -0.0003339888308120894,
    0.14678782013623856,
    0.12073652734814777,
    0.004826714200712452,
    -0.0016584992025804749,
    1.001593356164191,
    0.00044831771332554327,
    0.002256950050514632
  ],
  "param_names": [
    "dist_offset",
    "poni1",
    "poni2",
    "rot1",
    "rot2",
    "dist_scale",
    "poni1_scale",
    "poni2_scale"
  ],
  "pos_names": [
    "pos"
  ],
  "trans_function": {
    "content": "GeometryTransformation",
    "param_names": [
      "dist_offset",
      "poni1",
      "poni2",
      "rot1",
      "rot2",
      "dist_scale",
      "poni1_scale",
      "poni2_scale"
    ],
    "pos_names": [
      "pos"
    ],
    "dist_expr": "pos * dist_scale + dist_offset",
    "poni1_expr": "pos * poni1_scale + poni1",
    "poni2_expr": "pos * poni2_scale + poni2",
    "rot1_expr": "rot1",
    "rot2_expr": "rot2",
    "rot3_expr": "0.0",
    "constants": {
      "pi": 3.141592653589793
    }
  }
}

Interpretation of the fit:

[23]:
distances = []
f_distances = []
f_poni1 = []
f_poni2 = []
g_distances = []
g_poni1 = []
g_poni2 = []
for sg in gonioref.single_geometries.values():
    distance = sg.get_position()
    distances.append(distance)
    sg.geometry_refinement.refine3(fix=["wavelength", "rot3"])
    f_distances.append(sg.geometry_refinement.dist)
    f_poni1.append(sg.geometry_refinement.poni1)
    f_poni2.append(sg.geometry_refinement.poni2)
    ai = gonioref.get_ai(distance)
    g_distances.append(ai.dist)
    g_poni1.append(ai.poni1)
    g_poni2.append(ai.poni2)

order = numpy.argsort(distances)
distances = numpy.take(distances, order)
f_distances = numpy.take(f_distances, order)
f_poni1 = numpy.take(f_poni1, order)
f_poni2 = numpy.take(f_poni2, order)
g_distances = numpy.take(g_distances, order)
g_poni1 = numpy.take(g_poni1, order)
g_poni1 = numpy.take(g_poni1, order)

[24]:
fig,ax = subplots(2, figsize=(8,10))
ax[0].plot(distances, f_distances, "1", label="fitted individually")
ax[0].plot(distances, g_distances, label="fitted table")
ax[0].set_title("Observed deviations:")
ax[1].set_xlabel("Encoder distance (m)")
ax[1].plot(distances, f_poni1, "1", label="poni1 fitted individually")
ax[1].plot(distances, g_poni1, label="poni1 fitted table")
ax[1].plot(distances, f_poni2, "1", label="poni2 fitted individually")
ax[1].plot(distances, g_poni2, label="poni2 fitted table")
ax[0].set_ylabel("Fitted distance (m)")
ax[1].set_ylabel("Fitted PONIs (m)")
ax[0].legend()
ax[1].legend();
../../../../_images/usage_tutorial_Goniometer_MX-calibrate_MX-calibrate_29_0.png
[25]:
obtained_x = []
obtained_y = []

for dst in distances:
    geo = gonioref.get_ai(dst)
    fit2d = geo.getFit2D()
    obtained_x.append(fit2d["centerX"])
    obtained_y.append(fit2d["centerY"])
fig,ax = subplots()
ax.plot(distances, obtained_x, label="X")
ax.plot(distances, obtained_y, label="Y")
ax.set_title("Observed variation:")
ax.set_ylabel("Beam center (pixels)")
ax.set_xlabel("Encoder distance (m)")
ax.legend();
../../../../_images/usage_tutorial_Goniometer_MX-calibrate_MX-calibrate_30_0.png
[26]:
#Simply print out the result:
lrx = linregress(distances, obtained_x)
lry = linregress(distances, obtained_y)
print(f"Beam-center X: {lrx}")
print(f"Beam-center y: {lry}")
print()
print(f"beamcenter_x = {lrx.intercept:.3f} {lrx.slope:+.3f} * encoder_value")
print(f"beamcenter_y = {lry.intercept:.3f} {lry.slope:+.3f} * encoder_value")
Beam-center X: LinregressResult(slope=np.float64(-14.98542074431889), intercept=np.float64(701.9659269145309), rvalue=np.float64(-0.9999999999999999), pvalue=np.float64(3.422063358560521e-48), stderr=np.float64(9.116191269084659e-08), intercept_stderr=np.float64(2.5076144821866313e-08))
Beam-center y: LinregressResult(slope=np.float64(-7.0514240837066335), intercept=np.float64(853.4207794384212), rvalue=np.float64(-0.9999999999999999), pvalue=np.float64(3.422063358560521e-48), stderr=np.float64(4.289644699557044e-08), intercept_stderr=np.float64(1.1799637430297612e-08))

beamcenter_x = 701.966 -14.985 * encoder_value
beamcenter_y = 853.421 -7.051 * encoder_value

Nota:

The degradation between 0.3 and 0.5m correspond to the image 6->7 and is related to the disparition of the third ring!

Conclusion:#

This notebook demonstrates: * The usage of the geometry calibration in Jupyter-lab to calibrate the first image * The creation of a goniometer-refinement * The population of this goniometer-refinement with automatic control-point extraction * The fit of the table, first with the constrains of a perfecty aligned table, then with a mis-aligned table

In our case the table is miss-aligned in the horizontal direction by 2.3mm/meter (i.e. 2.3 mradian). This should be taken into account when calculating the beam-center at different distances.