Public API reference¶
All symbols are exported from tldecpy and available as tl.<name> after
import tldecpy as tl.
Fitting¶
tldecpy.fit.multi.fit_multi ¶
Fit a TL glow curve with multiple kinetic components and optional background.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
x
|
ndarray
|
Temperature grid :math: |
required |
y
|
ndarray
|
Observed thermoluminescence intensity :math: |
required |
peaks
|
list[PeakSpec]
|
One :class: |
required |
bg
|
BackgroundSpec | None
|
Background model specification. Pass |
None
|
beta
|
float
|
Heating rate :math: |
1.0
|
robust
|
RobustOptions | None
|
Robust-loss configuration for residual minimisation. Defaults to linear (ordinary least squares) with no Poisson weighting. |
None
|
options
|
FitOptions | None
|
Solver options (local optimizer, tolerances) and optional
:class: |
None
|
strategy
|
('local', 'global_hybrid', 'global_hybrid_pso')
|
Optimisation strategy. |
"local"
|
Returns:
| Type | Description |
|---|---|
MultiFitResult
|
Structured fit output. Key fields:
|
Raises:
| Type | Description |
|---|---|
ValueError
|
If |
ValueError
|
If |
Notes
The Levenberg-Marquardt optimizer ("lm") does not support box bounds
natively. Bounds are enforced via soft penalty terms appended to the
residual vector.
Examples:
>>> import tldecpy as tl
>>> T, I = tl.load_refglow("x001")
>>> peaks, bg = tl.autoinit_multi(T, I, max_peaks=2)
>>> result = tl.fit_multi(T, I, peaks=peaks, bg=bg, beta=1.0)
>>> print(result.converged, f"{result.metrics.FOM:.2f}%")
References
.. [1] Kitis, G., et al. (1998). Thermoluminescence glow-curve deconvolution functions for first, second and general orders of kinetics. J. Phys. D 31, 2636. .. [2] Virtanen, P., et al. (2020). SciPy 1.0: Fundamental Algorithms. Nature Methods 17, 261.
Source code in tldecpy/fit/multi.py
1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 | |
tldecpy.fit.solvers.fit_single_peak ¶
Fit a single TL component using the multi-component engine.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
x
|
ndarray
|
Temperature grid :math: |
required |
y
|
ndarray
|
Measured intensity :math: |
required |
model
|
str
|
Model key/family accepted by :mod: |
"fo"
|
init
|
dict[str, float] | None
|
Initial parameter guesses (for example |
None
|
bounds
|
dict[str, tuple[float, float]] | None
|
Optional parameter bounds as |
None
|
beta
|
float
|
Heating rate :math: |
1.0
|
robust
|
RobustOptions | None
|
Robust fitting configuration. |
None
|
Returns:
| Type | Description |
|---|---|
FitResult
|
Backward-compatible single-peak result schema. |
Source code in tldecpy/fit/solvers.py
tldecpy.fit.automator.iterative_deconvolution ¶
iterative_deconvolution(x, y, *, max_peaks=6, allow_models=('fo', 'go', 'otor_lw'), bg_mode='auto', sensitivity=1.0, min_snr=1.0, widths=None, residual_sigma_threshold=3.0, beta=1.0, robust=None, options=None, strategy='local')
Run automatic TL deconvolution with heuristic seeding and one-shot fitting.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
x
|
ndarray
|
Temperature grid :math: |
required |
y
|
ndarray
|
Measured glow-curve intensity :math: |
required |
max_peaks
|
int
|
Maximum number of components used in the automatic model. |
6
|
allow_models
|
tuple[str, ...]
|
Allowed model families/aliases for automatic seed-to-model mapping. |
("fo", "go", "otor_lw")
|
bg_mode
|
('linear', 'exponential', 'none', 'auto')
|
Background inference mode. |
"linear"
|
sensitivity
|
float
|
Sensitivity used by fallback local-peak detector. |
1.0
|
min_snr
|
float
|
Minimum SNR used by CWT-based detector. |
1.0
|
widths
|
ndarray | None
|
Optional CWT width grid in sample units. |
None
|
residual_sigma_threshold
|
float
|
Kept for backward compatibility with historical iterative API. |
3.0
|
beta
|
float
|
Heating rate :math: |
1.0
|
robust
|
RobustOptions | None
|
Robust loss/weighting settings. |
None
|
options
|
FitOptions | None
|
Local optimizer and uncertainty options. |
None
|
strategy
|
('local', 'global_hybrid', 'global_hybrid_pso')
|
Optimization strategy used by the underlying fitter. |
"local"
|
Returns:
| Type | Description |
|---|---|
MultiFitResult
|
Full multi-peak fit result including diagnostics and uncertainty payloads. |
Source code in tldecpy/fit/automator.py
104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 | |
Initialization¶
tldecpy.fit.init.autoinit_multi ¶
autoinit_multi(x, y, max_peaks=6, allow_models=('fo', 'go', 'otor_lw'), bg_mode='auto', sensitivity=1.0)
Build heuristic peak/background initialisation for multi-peak deconvolution.
The function executes three stages:
- Preprocess — Savitzky-Golay smoothing + robust outlier removal.
- Detect — CWT-based peak finder returns candidate positions, FWHM
and asymmetry coefficient :math:
\mu_g. - Seed — assigns a kinetic model per peak using :math:
\mu_g, estimates activation energy :math:Ewith Chen-style FWHM heuristics, and constructs :class:~tldecpy.schemas.PeakSpecobjects with automatic bounds.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
x
|
ndarray
|
Temperature grid :math: |
required |
y
|
ndarray
|
Measured TL intensity :math: |
required |
max_peaks
|
int
|
Maximum number of peaks to include in the returned list. If more
candidates are detected, the |
6
|
allow_models
|
tuple[str, ...]
|
Allowed model families or canonical keys. The seeding algorithm
selects the best match from this set based on peak asymmetry.
Continuous models ( |
("fo", "go", "otor_lw")
|
bg_mode
|
('linear', 'exponential', 'none', 'auto')
|
Background initialisation mode. |
"linear"
|
sensitivity
|
float
|
Peak detection sensitivity multiplier passed to :func: |
1.0
|
Returns:
| Type | Description |
|---|---|
tuple[list[PeakSpec], BackgroundSpec | None]
|
|
Raises:
| Type | Description |
|---|---|
ValueError
|
If |
Notes
Activation-energy initialisation uses Chen-style peak-shape heuristics:
.. math::
E \approx c_w \frac{k T_m^2}{\omega} - b_w (2 k T_m)
where :math:\omega is the FWHM and the constants :math:c_w, :math:b_w
depend on the kinetic order.
Examples:
>>> import tldecpy as tl
>>> T, I = tl.load_refglow("x002")
>>> peaks, bg = tl.autoinit_multi(T, I, max_peaks=4, allow_models=("fo_rq",))
>>> print(len(peaks), "peaks seeded")
References
.. [1] Chen, R., and McKeever, S. W. S. (1997). Theory of thermoluminescence and related phenomena. World Scientific.
Source code in tldecpy/fit/init.py
227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 | |
tldecpy.fit.init.pick_peaks ¶
Detect candidate TL peaks and return initialization seeds.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
x
|
ndarray
|
Temperature grid in kelvin. |
required |
y
|
ndarray
|
Preprocessed intensity signal. |
required |
sensitivity
|
float
|
Detection sensitivity multiplier. Higher values lower effective prominence/height thresholds and can recover weaker peaks. |
1.0
|
Returns:
| Type | Description |
|---|---|
list[PeakSeed]
|
Ordered list of peak seeds including :math: |
Source code in tldecpy/fit/init.py
tldecpy.fit.init.preprocess ¶
Preprocess a TL glow curve before peak detection.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
x
|
ndarray
|
Temperature grid in kelvin. |
required |
y
|
ndarray
|
Measured intensity values for each temperature sample. |
required |
sg_window
|
int
|
Savitzky-Golay window length (odd integer in samples). |
7
|
sg_poly
|
int
|
Savitzky-Golay polynomial order. |
3
|
remove_outliers
|
bool
|
If |
True
|
Returns:
| Type | Description |
|---|---|
ndarray
|
Smoothed, non-negative intensity signal used by automatic initialization. |
Source code in tldecpy/fit/init.py
tldecpy.fit.detection.detect_peaks_cwt ¶
Detect TL peak candidates using CWT ridge responses (Ricker/Mexican-hat).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
x
|
ndarray
|
Temperature grid :math: |
required |
y
|
ndarray
|
Intensity signal :math: |
required |
min_snr
|
float
|
Minimum signal-to-noise threshold used to retain ridge-consistent candidates. |
1.0
|
widths
|
ndarray | Iterable[float] | None
|
Wavelet widths in sample units. When omitted, an adaptive range based on data length is used. |
None
|
Returns:
| Type | Description |
|---|---|
list[PeakSeed]
|
Peak candidates with position, intensity and shape descriptors compatible with the auto-initialization pipeline. |
Notes
The implementation uses find_peaks_cwt when available and falls back to a
local ridge-threshold detector otherwise.
Source code in tldecpy/fit/detection.py
76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 | |
Simulation¶
tldecpy.simulate.core.simulate ¶
simulate(rhs_func, params, *, T0, T_end, beta, y0, state_keys, method='LSODA', points=1000, noise_config=None)
Integrate TL kinetic ODEs under linear heating and return a typed result.
The time variable is related to temperature by the linear heating law
:math:T(t) = T_0 + \beta t. The ODE is integrated in time and the
output is re-expressed on a temperature grid.
TL intensity is computed as minus the time-derivative of the first state
variable: :math:I(T) = -\dot{y}_0(t).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
rhs_func
|
Callable
|
ODE right-hand side callable. Must accept positional arguments
|
required |
params
|
tuple[float, ...]
|
Model-specific kinetic parameters passed as positional args to
|
required |
T0
|
float
|
Initial temperature in kelvin. Must be < |
required |
T_end
|
float
|
Final temperature in kelvin. Must be > |
required |
beta
|
float
|
Linear heating rate :math: |
required |
y0
|
list[float]
|
Initial state vector passed to |
required |
state_keys
|
list[str]
|
Human-readable names for each state variable. Used as keys in
|
required |
method
|
str
|
Integration algorithm accepted by |
"LSODA"
|
points
|
int
|
Number of evenly-spaced output temperature samples. |
1000
|
noise_config
|
dict | None
|
If provided, additive noise is applied to the intensity array. Recognised keys:
|
None
|
Returns:
| Type | Description |
|---|---|
SimulationResult
|
Typed result with fields:
|
Notes
A RuntimeWarning is emitted if the ODE integrator does not converge.
The result is still returned with whatever data was produced up to the
failure point.
Examples:
>>> import numpy as np
>>> import tldecpy as tl
>>> from tldecpy.simulate.fo import ode_fo, infer_n0_s
>>> T0, T_end, beta = 300.0, 600.0, 1.0
>>> E, Tm = 1.2, 450.0
>>> n0, s = infer_n0_s(Im=1.0, Tm=Tm, E=E, beta=beta)
>>> result = tl.simulate(
... ode_fo, params=(s, E),
... T0=T0, T_end=T_end, beta=beta,
... y0=[n0], state_keys=["n"],
... )
>>> print(result.T.shape, result.I.max())
Source code in tldecpy/simulate/core.py
16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 | |
Data¶
tldecpy.data.refglow.load_refglow ¶
Load one Refglow/GLOCANIN benchmark curve as temperature and intensity arrays.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
key
|
str
|
Dataset identifier. Accepted forms:
Available datasets ( .. list-table:: :header-rows: 1
|
required |
Returns:
| Type | Description |
|---|---|
tuple[ndarray, ndarray]
|
|
Raises:
| Type | Description |
|---|---|
ValueError
|
If |
FileNotFoundError
|
If the dataset file cannot be located inside the package data directory. |
Examples:
References
.. [1] Bos, A. J. J., et al. (1993). Intercomparison of glow curve analysis computer programs: the GLOCANIN project. Radiat. Prot. Dosim. 47, 473. .. [2] Bos, A. J. J., et al. (1994). GLOCANIN: A data set for the intercomparison of glow-curve analysis programs. Radiat. Meas. 23, 393.
Source code in tldecpy/data/refglow.py
tldecpy.data.refglow.list_refglow ¶
List available Refglow dataset identifiers.
Returns:
| Type | Description |
|---|---|
list[str]
|
Available IDs from |
Model registry¶
tldecpy.models.registry.list_models ¶
List available model keys for TL kinetic families.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
order
|
str | None
|
Optional family filter ( |
None
|
include_aliases
|
bool
|
When |
False
|
Returns:
| Type | Description |
|---|---|
list[str]
|
Canonical model keys (and optional aliases). |
Source code in tldecpy/models/registry.py
tldecpy.models.registry.get_model ¶
Resolve and return a callable model implementation.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
key
|
str
|
Canonical key, alias, or family key. |
required |
Returns:
| Type | Description |
|---|---|
Callable
|
Numerical model function that maps TL parameters to intensity. |
Source code in tldecpy/models/registry.py
tldecpy.version.get_version_info ¶
Return package and API version metadata.
Returns:
| Type | Description |
|---|---|
VersionInfo
|
Structured version payload containing package and API versions. |