Structure Refinement: LBCO, HRPT¶
This minimalistic example is designed to show how Rietveld refinement can be performed when both the crystal structure and experiment are defined directly in code. Only the experimentally measured data is loaded from an external file. It also shows how to switch calculation engine.
For this example, constant-wavelength neutron powder diffraction data for La0.5Ba0.5CoO3 from HRPT at PSI is used.
It does not contain any advanced features or options, and includes no comments or explanations — these can be found in the other tutorials. Default values are used for all parameters if not specified. Only essential and self-explanatory code is provided.
The example is intended for users who are already familiar with the EasyDiffraction library and want to quickly get started with a simple refinement. It is also useful for those who want to see what a refinement might look like in code. For a more detailed explanation of the code, please refer to the other tutorials.
🛠️ Import Library¶
import easydiffraction as ed
📦 Define Project¶
project = ed.Project(name='lbco_hrpt')
🧩 Define Structure¶
project.structures.create(name='lbco')
structure = project.structures['lbco']
structure.space_group.name_h_m = 'P m -3 m'
structure.space_group.it_coordinate_system_code = '1'
structure.cell.length_a = 3.88
structure.atom_sites.create(
label='La',
type_symbol='La',
fract_x=0,
fract_y=0,
fract_z=0,
adp_iso=0.5,
occupancy=0.5,
)
structure.atom_sites.create(
label='Ba',
type_symbol='Ba',
fract_x=0,
fract_y=0,
fract_z=0,
adp_iso=0.5,
occupancy=0.5,
)
structure.atom_sites.create(
label='Co',
type_symbol='Co',
fract_x=0.5,
fract_y=0.5,
fract_z=0.5,
adp_iso=0.5,
)
structure.atom_sites.create(
label='O',
type_symbol='O',
fract_x=0,
fract_y=0.5,
fract_z=0.5,
adp_iso=0.5,
)
project.display.structure(struct_name='lbco')
Structure 🧩 'lbco' (Atom view type: 'covalent')
wheel = zoom
right-drag = pan
🔬 Define Experiment¶
data_path = ed.download_data(id=3, destination='data')
Getting data...
Data #3: La0.5Ba0.5CoO3, HRPT (PSI), 300 K
✅ Data #3 downloaded to '../../../data/ed-3.xye'
project.experiments.add_from_data_path(
name='hrpt',
data_path=data_path,
sample_form='powder',
beam_mode='constant wavelength',
radiation_probe='neutron',
)
Data loaded successfully
Experiment 🔬 'hrpt'. Number of data points: 3098.
experiment = project.experiments['hrpt']
experiment.instrument.setup_wavelength = 1.494
experiment.instrument.calib_twotheta_offset = 0.6
experiment.peak.broad_gauss_u = 0.1
experiment.peak.broad_gauss_v = -0.1
experiment.peak.broad_gauss_w = 0.1
experiment.peak.broad_lorentz_y = 0.1
experiment.background.create(id='1', x=10, y=170)
experiment.background.create(id='2', x=30, y=170)
experiment.background.create(id='3', x=50, y=170)
experiment.background.create(id='4', x=110, y=170)
experiment.background.create(id='5', x=165, y=170)
experiment.excluded_regions.create(id='1', start=0, end=5)
experiment.excluded_regions.create(id='2', start=165, end=180)
experiment.linked_phases.create(id='lbco', scale=10.0)
🚀 Perform Analysis¶
Without Constraints¶
structure.cell.length_a.free = True
structure.atom_sites['La'].adp_iso.free = True
structure.atom_sites['Ba'].adp_iso.free = True
structure.atom_sites['Co'].adp_iso.free = True
structure.atom_sites['O'].adp_iso.free = True
experiment.instrument.calib_twotheta_offset.free = True
experiment.peak.broad_gauss_u.free = True
experiment.peak.broad_gauss_v.free = True
experiment.peak.broad_gauss_w.free = True
experiment.peak.broad_lorentz_y.free = True
experiment.background['1'].y.free = True
experiment.background['2'].y.free = True
experiment.background['3'].y.free = True
experiment.background['4'].y.free = True
experiment.background['5'].y.free = True
experiment.linked_phases['lbco'].scale.free = True
project.analysis.fit()
Standard fitting
📋 Using experiment 🔬 'hrpt' for 'single' fitting
🚀 Starting fit process with 'lmfit (leastsq)'...
📈 Goodness-of-fit progress:
| iteration | time (s) | χ² | change / status | |
|---|---|---|---|---|
| 1 | 1 | 0.29 | 165.02 | |
| 2 | 28 | 1.52 | 33.58 | 79.7% ↓ |
| 3 | 45 | 2.28 | 10.82 | 67.8% ↓ |
| 4 | 63 | 3.10 | 6.49 | 40.0% ↓ |
| 5 | 81 | 3.95 | 3.35 | 48.4% ↓ |
| 6 | 98 | 4.72 | 2.24 | 33.1% ↓ |
| 7 | 116 | 5.77 | 1.91 | 14.7% ↓ |
| 8 | 133 | 6.55 | 1.50 | 21.3% ↓ |
| 9 | 150 | 7.32 | 1.45 | 3.6% ↓ |
| 10 | 167 | 8.09 | 1.34 | 7.8% ↓ |
| 11 | 185 | 8.92 | 1.29 | 3.4% ↓ |
| 12 | 276 | 13.17 | 1.29 |
🏆 Best goodness-of-fit (reduced χ²) is 1.29 at iteration 275
✅ Fitting complete.
project.display.fit.results()
⚙️ Settings used:
| Name | Value | Description | |
|---|---|---|---|
| 1 | max_iterations | 1000 | Maximum solver iterations. |
📋 Least-squares fit results:
| Metric | Value | |
|---|---|---|
| 1 | 🧪 Minimizer | lmfit (leastsq) |
| 2 | ✅ Overall status | success |
| 3 | ⏱️ Fitting time (seconds) | 13.17 |
| 4 | 🔁 Iterations | 273 |
| 5 | 📏 Goodness-of-fit (reduced χ²) | 1.29 |
| 6 | 📏 R-factor (Rf, %) | 5.63 |
| 7 | 📏 R-factor squared (Rf², %) | 5.27 |
| 8 | 📏 Weighted R-factor (wR, %) | 4.41 |
📈 Refined parameters:
| datablock | category | entry | parameter | units | start | value | s.u. | change | |
|---|---|---|---|---|---|---|---|---|---|
| 1 | lbco | cell | length_a | Å | 3.8800 | 3.8909 | 0.0000 | 0.28 % ↑ | |
| 2 | lbco | atom_site | La | adp_iso | Ų | 0.5000 | 0.5052 | 3231.7297 | 1.04 % ↑ |
| 3 | lbco | atom_site | Ba | adp_iso | Ų | 0.5000 | 0.5049 | 5252.5594 | 0.97 % ↑ |
| 4 | lbco | atom_site | Co | adp_iso | Ų | 0.5000 | 0.2370 | 0.0611 | 52.60 % ↓ |
| 5 | lbco | atom_site | O | adp_iso | Ų | 0.5000 | 1.3935 | 0.0167 | 178.71 % ↑ |
| 6 | hrpt | linked_phases | lbco | scale | 10.0000 | 9.1351 | 0.0642 | 8.65 % ↓ | |
| 7 | hrpt | peak | broad_gauss_u | deg² | 0.1000 | 0.0816 | 0.0031 | 18.43 % ↓ | |
| 8 | hrpt | peak | broad_gauss_v | deg² | -0.1000 | -0.1159 | 0.0067 | 15.93 % ↑ | |
| 9 | hrpt | peak | broad_gauss_w | deg² | 0.1000 | 0.1204 | 0.0033 | 20.45 % ↑ | |
| 10 | hrpt | peak | broad_lorentz_y | deg | 0.1000 | 0.0844 | 0.0021 | 15.55 % ↓ | |
| 11 | hrpt | instrument | twotheta_offset | deg | 0.6000 | 0.6226 | 0.0010 | 3.76 % ↑ | |
| 12 | hrpt | background | 1 | y | 170.0000 | 168.5585 | 1.3671 | 0.85 % ↓ | |
| 13 | hrpt | background | 2 | y | 170.0000 | 164.3357 | 0.9992 | 3.33 % ↓ | |
| 14 | hrpt | background | 3 | y | 170.0000 | 166.8881 | 0.7388 | 1.83 % ↓ | |
| 15 | hrpt | background | 4 | y | 170.0000 | 175.4006 | 0.6578 | 3.18 % ↑ | |
| 16 | hrpt | background | 5 | y | 170.0000 | 174.2813 | 0.9105 | 2.52 % ↑ |
• value = refined value from least-squares minimization
• s.u. = standard uncertainty (one sigma), from the covariance matrix
• change = relative change from start, in %; ↑ = increase, ↓ = decrease
project.display.fit.correlations()
project.display.pattern(expt_name='hrpt')
With Constraints¶
# As can be seen from the parameter-correlation plot, the isotropic
# displacement parameters of La and Ba are highly correlated. Because
# La and Ba share the same mixed-occupancy site, their contributions to
# the neutron diffraction pattern are difficult to separate, especially
# since their coherent scattering lengths are not very different.
# Therefore, it is necessary to constrain them to be equal. First we
# define aliases and then use them to create a constraint.
project.analysis.aliases.create(
label='biso_La',
param=project.structures['lbco'].atom_sites['La'].adp_iso,
)
project.analysis.aliases.create(
label='biso_Ba',
param=project.structures['lbco'].atom_sites['Ba'].adp_iso,
)
project.analysis.constraints.create(expression='biso_Ba = biso_La')
project.analysis.minimizer.show_supported()
project.analysis.minimizer.type = 'lmfit'
Minimizer types
| Type | Description | ||
|---|---|---|---|
| 1 | bumps | BUMPS library using the default Levenberg-Marquardt method | |
| 2 | bumps (amoeba) | BUMPS library with Nelder-Mead simplex method | |
| 3 | bumps (de) | BUMPS library with differential evolution method | |
| 4 | bumps (dream) | BUMPS library with DREAM Bayesian sampling | |
| 5 | bumps (lm) | BUMPS library with Levenberg-Marquardt method | |
| 6 | dfols | DFO-LS library for derivative-free least-squares optimization | |
| 7 | emcee | emcee affine-invariant ensemble Bayesian sampling | |
| 8 | lmfit | LMFIT library using the default Levenberg-Marquardt method | |
| 9 | lmfit (least_squares) | LMFIT library with SciPy's trust region reflective algorithm | |
| 10 | * | lmfit (leastsq) | LMFIT library with Levenberg-Marquardt least squares method |
Current minimizer changed to
lmfit
project.analysis.fit()
Standard fitting
📋 Using experiment 🔬 'hrpt' for 'single' fitting
🚀 Starting fit process with 'lmfit (leastsq)'...
📈 Goodness-of-fit progress:
| iteration | time (s) | χ² | change / status | |
|---|---|---|---|---|
| 1 | 1 | 0.05 | 1.29 | |
| 2 | 36 | 2.00 | 1.29 |
🏆 Best goodness-of-fit (reduced χ²) is 1.29 at iteration 35
✅ Fitting complete.
project.display.fit.results()
⚙️ Settings used:
| Name | Value | Description | |
|---|---|---|---|
| 1 | max_iterations | 1000 | Maximum solver iterations. |
📋 Least-squares fit results:
| Metric | Value | |
|---|---|---|
| 1 | 🧪 Minimizer | lmfit |
| 2 | ✅ Overall status | success |
| 3 | ⏱️ Fitting time (seconds) | 2.00 |
| 4 | 🔁 Iterations | 33 |
| 5 | 📏 Goodness-of-fit (reduced χ²) | 1.29 |
| 6 | 📏 R-factor (Rf, %) | 5.63 |
| 7 | 📏 R-factor squared (Rf², %) | 5.27 |
| 8 | 📏 Weighted R-factor (wR, %) | 4.41 |
📈 Refined parameters:
| datablock | category | entry | parameter | units | start | value | s.u. | change | |
|---|---|---|---|---|---|---|---|---|---|
| 1 | lbco | cell | length_a | Å | 3.8909 | 3.8909 | 0.0000 | 0.00 % ↑ | |
| 2 | lbco | atom_site | La | adp_iso | Ų | 0.5052 | 0.5051 | 0.0278 | 0.03 % ↓ |
| 3 | lbco | atom_site | Co | adp_iso | Ų | 0.2370 | 0.2370 | 0.0564 | 0.00 % ↑ |
| 4 | lbco | atom_site | O | adp_iso | Ų | 1.3935 | 1.3935 | 0.0160 | 0.00 % ↓ |
| 5 | hrpt | linked_phases | lbco | scale | 9.1351 | 9.1351 | 0.0538 | 0.00 % ↓ | |
| 6 | hrpt | peak | broad_gauss_u | deg² | 0.0816 | 0.0816 | 0.0031 | 0.00 % ↑ | |
| 7 | hrpt | peak | broad_gauss_v | deg² | -0.1159 | -0.1159 | 0.0066 | 0.00 % ↑ | |
| 8 | hrpt | peak | broad_gauss_w | deg² | 0.1204 | 0.1204 | 0.0032 | 0.00 % ↑ | |
| 9 | hrpt | peak | broad_lorentz_y | deg | 0.0844 | 0.0844 | 0.0021 | 0.00 % ↓ | |
| 10 | hrpt | instrument | twotheta_offset | deg | 0.6226 | 0.6226 | 0.0010 | 0.00 % ↑ | |
| 11 | hrpt | background | 1 | y | 168.5585 | 168.5585 | 1.3669 | 0.00 % ↓ | |
| 12 | hrpt | background | 2 | y | 164.3357 | 164.3357 | 0.9990 | 0.00 % ↑ | |
| 13 | hrpt | background | 3 | y | 166.8881 | 166.8881 | 0.7386 | 0.00 % ↑ | |
| 14 | hrpt | background | 4 | y | 175.4006 | 175.4006 | 0.6488 | 0.00 % ↑ | |
| 15 | hrpt | background | 5 | y | 174.2813 | 174.2813 | 0.8944 | 0.00 % ↓ |
• value = refined value from least-squares minimization
• s.u. = standard uncertainty (one sigma), from the covariance matrix
• change = relative change from start, in %; ↑ = increase, ↓ = decrease
project.display.fit.correlations()
project.display.pattern(expt_name='hrpt')
Switch Calculator¶
experiment.calculator.show_supported()
Calculator types
| Type | Description | ||
|---|---|---|---|
| 1 | crysfml | CrysFML library for crystallographic calculations | |
| 2 | * | cryspy | CrysPy library for crystallographic calculations |
experiment.calculator.type = 'crysfml'
Calculator for experiment 'hrpt' changed to
crysfml
project.analysis.fit()
Standard fitting
📋 Using experiment 🔬 'hrpt' for 'single' fitting
🚀 Starting fit process with 'lmfit (leastsq)'...
📈 Goodness-of-fit progress:
| iteration | time (s) | χ² | change / status | |
|---|---|---|---|---|
| 1 | 1 | 0.05 | 1.29 | |
| 2 | 117 | 5.07 | 1.29 | |
| 3 | 122 | 5.36 | 1.29 |
🏆 Best goodness-of-fit (reduced χ²) is 1.29 at iteration 112
✅ Fitting complete.
project.display.fit.results()
⚙️ Settings used:
| Name | Value | Description | |
|---|---|---|---|
| 1 | max_iterations | 1000 | Maximum solver iterations. |
📋 Least-squares fit results:
| Metric | Value | |
|---|---|---|
| 1 | 🧪 Minimizer | lmfit |
| 2 | ✅ Overall status | success |
| 3 | ⏱️ Fitting time (seconds) | 5.36 |
| 4 | 🔁 Iterations | 119 |
| 5 | 📏 Goodness-of-fit (reduced χ²) | 1.29 |
| 6 | 📏 R-factor (Rf, %) | 5.62 |
| 7 | 📏 R-factor squared (Rf², %) | 5.27 |
| 8 | 📏 Weighted R-factor (wR, %) | 4.44 |
📈 Refined parameters:
| datablock | category | entry | parameter | units | start | value | s.u. | change | |
|---|---|---|---|---|---|---|---|---|---|
| 1 | lbco | cell | length_a | Å | 3.8909 | 3.8908 | 0.0000 | 0.00 % ↓ | |
| 2 | lbco | atom_site | La | adp_iso | Ų | 0.5051 | 0.5192 | 0.0255 | 2.79 % ↑ |
| 3 | lbco | atom_site | Co | adp_iso | Ų | 0.2370 | 0.2163 | 0.0550 | 8.74 % ↓ |
| 4 | lbco | atom_site | O | adp_iso | Ų | 1.3935 | 1.3955 | 0.0155 | 0.14 % ↑ |
| 5 | hrpt | linked_phases | lbco | scale | 9.1351 | 9.1561 | 0.0516 | 0.23 % ↑ | |
| 6 | hrpt | peak | broad_gauss_u | deg² | 0.0816 | 0.0815 | 0.0029 | 0.10 % ↓ | |
| 7 | hrpt | peak | broad_gauss_v | deg² | -0.1159 | -0.1159 | 0.0063 | 0.06 % ↓ | |
| 8 | hrpt | peak | broad_gauss_w | deg² | 0.1204 | 0.1204 | 0.0030 | 0.05 % ↓ | |
| 9 | hrpt | peak | broad_lorentz_y | deg | 0.0844 | 0.0846 | 0.0021 | 0.13 % ↑ | |
| 10 | hrpt | instrument | twotheta_offset | deg | 0.6226 | 0.6200 | 0.0007 | 0.41 % ↓ | |
| 11 | hrpt | background | 1 | y | 168.5585 | 169.6910 | 1.3655 | 0.67 % ↑ | |
| 12 | hrpt | background | 2 | y | 164.3357 | 164.7858 | 0.9957 | 0.27 % ↑ | |
| 13 | hrpt | background | 3 | y | 166.8881 | 167.0204 | 0.7365 | 0.08 % ↑ | |
| 14 | hrpt | background | 4 | y | 175.4006 | 175.2816 | 0.6506 | 0.07 % ↓ | |
| 15 | hrpt | background | 5 | y | 174.2813 | 175.1232 | 0.8851 | 0.48 % ↑ |
• value = refined value from least-squares minimization
• s.u. = standard uncertainty (one sigma), from the covariance matrix
• change = relative change from start, in %; ↑ = increase, ↓ = decrease
project.display.fit.correlations()
project.display.pattern(expt_name='hrpt')
💾 Save Project¶
project.save_as(dir_path='projects/ed_2_lbco_hrpt')
Saving project 📦 'lbco_hrpt' to '../../../projects/ed_2_lbco_hrpt'
├── 📄 project.cif
├── 📁 structures/
│ └── 📄 lbco.cif
├── 📁 experiments/
│ └── 📄 hrpt.cif
├── 📁 analysis/
│ └── 📄 analysis.cif
└── 📁 reports/
└── 📄 lbco_hrpt.html