Coverage for phenomena/tmm.py: 78%

125 statements  

« prev     ^ index     » next       coverage.py v7.11.0, created at 2025-11-16 22:49 +1300

1""" 

2*Transfer Matrix Method* (TMM) for Multilayer Optical Calculations 

3================================================================== 

4 

5Implements the *Transfer Matrix Method* for computing optical properties of 

6multilayer planar stacks. 

7 

8The *TMM* calculations handle two orthogonal polarisation states: 

9 

10- **s-polarisation (Transverse Electric)**: The electric field vector 

11 :math:`\\vec{E}_s` oscillates perpendicular to the plane of incidence. 

12 This is sometimes called TE-mode because the electric field is 

13 transverse to the plane of incidence. 

14 

15- **p-polarisation (TM - Transverse Magnetic)**: The electric field vector 

16 :math:`\\vec{E}_p` oscillates parallel to the plane of incidence. 

17 This is sometimes called TM polarisation because the magnetic field is 

18 transverse to the plane of incidence. 

19 

20The plane of incidence is defined by the incident light ray and the surface 

21normal vector. For unpolarised light, the reflectance and transmittance are 

22the average of the s and p components: :math:`R = (R_s + R_p)/2`. 

23 

24- :func:`colour.phenomena.snell_law` 

25- :func:`colour.phenomena.polarised_light_magnitude_elements` 

26- :func:`colour.phenomena.polarised_light_reflection_amplitude` 

27- :func:`colour.phenomena.polarised_light_reflection_coefficient` 

28- :func:`colour.phenomena.polarised_light_transmission_amplitude` 

29- :func:`colour.phenomena.polarised_light_transmission_coefficient` 

30- :class:`colour.phenomena.TransferMatrixResult_TMM` 

31- :func:`colour.phenomena.matrix_transfer_tmm` 

32 

33References 

34---------- 

35- :cite:`Byrnes2016` : Byrnes, S. J. (2016). Multilayer optical 

36 calculations. arXiv:1603.02720 [Physics]. 

37 http://arxiv.org/abs/1603.02720 

38""" 

39 

40from __future__ import annotations 

41 

42from dataclasses import dataclass 

43from typing import TYPE_CHECKING 

44 

45import numpy as np 

46 

47from colour.constants import DTYPE_COMPLEX_DEFAULT 

48from colour.utilities import ( 

49 MixinDataclassArithmetic, 

50 as_complex_array, 

51 as_float_array, 

52 tsplit, 

53 tstack, 

54 zeros, 

55) 

56 

57if TYPE_CHECKING: 

58 from colour.hints import ArrayLike, NDArrayComplex, NDArrayFloat, Tuple 

59 

60__author__ = "Colour Developers" 

61__copyright__ = "Copyright 2013 Colour Developers" 

62__license__ = "BSD-3-Clause - https://opensource.org/licenses/BSD-3-Clause" 

63__maintainer__ = "Colour Developers" 

64__email__ = "colour-developers@colour-science.org" 

65__status__ = "Production" 

66 

67__all__ = [ 

68 "snell_law", 

69 "polarised_light_magnitude_elements", 

70 "polarised_light_reflection_amplitude", 

71 "polarised_light_reflection_coefficient", 

72 "polarised_light_transmission_amplitude", 

73 "polarised_light_transmission_coefficient", 

74 "TransferMatrixResult", 

75 "matrix_transfer_tmm", 

76] 

77 

78 

79def _tsplit_complex(a: ArrayLike) -> NDArrayComplex: 

80 """ 

81 Split the specified stacked array along the last axis (tail) 

82 to produce an array of complex arrays. 

83 

84 Convenience wrapper around :func:`colour.utilities.tsplit` that 

85 automatically uses ``DTYPE_COMPLEX_DEFAULT`` for complex number 

86 operations in *Transfer Matrix Method* calculations. 

87 

88 Parameters 

89 ---------- 

90 a 

91 Stacked array to split. 

92 

93 Returns 

94 ------- 

95 :class:`numpy.ndarray` 

96 Array of complex arrays. 

97 """ 

98 

99 return tsplit(a, dtype=DTYPE_COMPLEX_DEFAULT) # type: ignore[arg-type] 

100 

101 

102def _tstack_complex(a: ArrayLike) -> NDArrayComplex: 

103 """ 

104 Stack the specified array of arrays along the last axis (tail) 

105 to produce a stacked complex array. 

106 

107 Convenience wrapper around :func:`colour.utilities.tstack` that 

108 automatically uses ``DTYPE_COMPLEX_DEFAULT`` for complex number 

109 operations in *Transfer Matrix Method* calculations. 

110 

111 Parameters 

112 ---------- 

113 a 

114 Array of arrays to stack along the last axis. 

115 

116 Returns 

117 ------- 

118 :class:`numpy.ndarray` 

119 Stacked complex array. 

120 

121 References 

122 ---------- 

123 :cite:`Byrnes2016` 

124 """ 

125 

126 return tstack(a, dtype=DTYPE_COMPLEX_DEFAULT) # type: ignore[arg-type] 

127 

128 

129def snell_law( 

130 n_1: ArrayLike, 

131 n_2: ArrayLike, 

132 theta_i: ArrayLike, 

133) -> NDArrayFloat: 

134 """ 

135 Compute the refraction angle using *Snell's Law*. 

136 

137 Parameters 

138 ---------- 

139 n_1 

140 Refractive index of the incident medium :math:`n_1`. 

141 n_2 

142 Refractive index of the refracting medium :math:`n_2`. 

143 theta_i 

144 Incident angle :math:`\\theta_i` in degrees. 

145 

146 Returns 

147 ------- 

148 :class:`numpy.ndarray` 

149 Refracted angle in degrees. 

150 

151 Notes 

152 ----- 

153 - *Snell's Law* relates the angles of incidence and refraction when light 

154 passes through a boundary between two different media (*Equation 3* 

155 from :cite:`Byrnes2016`): 

156 

157 .. math:: 

158 

159 n_i \\sin \\theta_i = n_j \\sin \\theta_j 

160 

161 Where: 

162 

163 - :math:`n_i, n_j`: Refractive indices of the incident and refracting media 

164 - :math:`\\theta_i, \\theta_j`: Angles of incidence and refraction 

165 

166 References 

167 ---------- 

168 :cite:`Byrnes2016` 

169 

170 Examples 

171 -------- 

172 >>> snell_law(1.0, 1.5, 30.0) # doctest: +ELLIPSIS 

173 19.4712206... 

174 """ 

175 

176 n_1 = np.real(as_complex_array(n_1)) 

177 n_2 = np.real(as_complex_array(n_2)) 

178 theta_i = np.radians(as_float_array(theta_i)) 

179 

180 # Apply Snell's law: n_i * sin(theta_i) = n_j * sin(theta_j) (Byrnes Eq. 3) 

181 return np.degrees(np.arcsin(n_1 * np.sin(theta_i) / n_2)) 

182 

183 

184def polarised_light_magnitude_elements( 

185 n_1: ArrayLike, 

186 n_2: ArrayLike, 

187 theta_i: ArrayLike, 

188 theta_t: ArrayLike, 

189) -> Tuple[NDArrayComplex, NDArrayComplex, NDArrayComplex, NDArrayComplex]: 

190 """ 

191 Compute common magnitude elements for *Fresnel* equations. 

192 

193 This function computes the common terms used in the *Fresnel* equations 

194 for both s-polarisation (perpendicular) and p-polarisation (parallel) 

195 components of light at a dielectric interface. 

196 

197 Parameters 

198 ---------- 

199 n_1 

200 Refractive index of the incident medium :math:`n_1`. 

201 n_2 

202 Refractive index of the transmitted medium :math:`n_2`. 

203 theta_i 

204 Incident angle :math:`\\theta_i` in degrees. 

205 theta_t 

206 Transmitted angle :math:`\\theta_t` in degrees. 

207 

208 Returns 

209 ------- 

210 :class:`tuple` 

211 Tuple of precomputed magnitude elements: 

212 :math:`(n_1 \\cos \\theta_i, n_1 \\cos \\theta_t, n_2 \\cos \\theta_i, 

213 n_2 \\cos \\theta_t)` 

214 

215 Notes 

216 ----- 

217 These magnitude elements are fundamental components in the *Fresnel* equations: 

218 

219 - :math:`n_1 \\cos \\theta_i`: Incident medium magnitude 

220 - :math:`n_1 \\cos \\theta_t`: Incident medium magnitude (transmitted angle) 

221 - :math:`n_2 \\cos \\theta_i`: Transmitted medium magnitude (incident angle) 

222 - :math:`n_2 \\cos \\theta_t`: Transmitted medium magnitude 

223 

224 These terms appear in all *Fresnel* amplitude and power coefficients for 

225 both reflection and transmission at dielectric interfaces. 

226 

227 Examples 

228 -------- 

229 >>> polarised_light_magnitude_elements(1.0, 1.5, 0.0, 0.0) 

230 ((1+0j), (1+0j), (1.5+0j), (1.5+0j)) 

231 """ 

232 

233 n_1 = as_complex_array(n_1) 

234 n_2 = as_complex_array(n_2) 

235 

236 cos_theta_i = np.cos(np.radians(as_float_array(theta_i))) 

237 cos_theta_t = np.cos(np.radians(as_float_array(theta_t))) 

238 

239 n_1_cos_theta_i = n_1 * cos_theta_i 

240 n_1_cos_theta_t = n_1 * cos_theta_t 

241 n_2_cos_theta_i = n_2 * cos_theta_i 

242 n_2_cos_theta_t = n_2 * cos_theta_t 

243 

244 return n_1_cos_theta_i, n_1_cos_theta_t, n_2_cos_theta_i, n_2_cos_theta_t 

245 

246 

247def polarised_light_reflection_amplitude( 

248 n_1: ArrayLike, 

249 n_2: ArrayLike, 

250 theta_i: ArrayLike, 

251 theta_t: ArrayLike, 

252) -> NDArrayComplex: 

253 """ 

254 Compute *Fresnel* reflection amplitude coefficients. 

255 

256 This function computes the complex reflection amplitude coefficients for 

257 both s-polarisation (perpendicular) and p-polarisation (parallel) components 

258 of electromagnetic waves at a dielectric interface. 

259 

260 Parameters 

261 ---------- 

262 n_1 

263 Refractive index of the incident medium :math:`n_1`. 

264 n_2 

265 Refractive index of the transmitted medium :math:`n_2`. 

266 theta_i 

267 Incident angle :math:`\\theta_i` in degrees. 

268 theta_t 

269 Transmitted angle :math:`\\theta_t` in degrees. 

270 

271 Returns 

272 ------- 

273 :class:`numpy.ndarray` 

274 *Fresnel* reflection amplitude coefficients for s and p polarisations 

275 stacked along the last axis. The array contains :math:`[r_s, r_p]` where 

276 :math:`r_s` and :math:`r_p` are the complex reflection coefficients. 

277 

278 Notes 

279 ----- 

280 The *Fresnel* reflection amplitude coefficients are given by (*Equation 6* 

281 from :cite:`Byrnes2016`): 

282 

283 .. math:: 

284 

285 r_s &= \\frac{n_1 \\cos \\theta_1 - n_2 \\cos \\theta_2}{n_1 \\cos \ 

286\\theta_1 + n_2 \\cos \\theta_2} \\\\ 

287 r_p &= \\frac{n_2 \\cos \\theta_1 - n_1 \\cos \\theta_2}{n_2 \\cos \ 

288\\theta_1 + n_1 \\cos \\theta_2} 

289 

290 Where: 

291 

292 - :math:`r_s`: s-polarisation reflection amplitude (electric field perpendicular 

293 to the plane of incidence) 

294 - :math:`r_p`: p-polarisation reflection amplitude (electric field parallel 

295 to the plane of incidence) 

296 - :math:`n_1, n_2`: Refractive indices of incident and transmitted media 

297 - :math:`\\theta_1, \\theta_2`: Incident and transmitted angles 

298 

299 Examples 

300 -------- 

301 >>> polarised_light_reflection_amplitude(1.0, 1.5, 0.0, 0.0) 

302 array([-0.2+0.j, -0.2+0.j]) 

303 """ 

304 

305 n_1_cos_theta_i, n_1_cos_theta_t, n_2_cos_theta_i, n_2_cos_theta_t = ( 

306 polarised_light_magnitude_elements(n_1, n_2, theta_i, theta_t) 

307 ) 

308 

309 # Fresnel reflection amplitudes (Byrnes Eq. 6) 

310 r_s = (n_1_cos_theta_i - n_2_cos_theta_t) / (n_1_cos_theta_i + n_2_cos_theta_t) 

311 r_p = (n_1_cos_theta_t - n_2_cos_theta_i) / (n_1_cos_theta_t + n_2_cos_theta_i) 

312 

313 return _tstack_complex([r_s, r_p]) 

314 

315 

316def polarised_light_reflection_coefficient( 

317 n_1: ArrayLike, 

318 n_2: ArrayLike, 

319 theta_i: ArrayLike, 

320 theta_t: ArrayLike, 

321) -> NDArrayComplex: 

322 """ 

323 Compute *Fresnel* reflection power coefficients (reflectance). 

324 

325 This function computes the reflection power coefficients, which represent 

326 the fraction of incident power that is reflected at a dielectric interface 

327 for both s-polarisation (perpendicular) and p-polarisation (parallel) components. 

328 

329 Parameters 

330 ---------- 

331 n_1 

332 Refractive index of the incident medium :math:`n_1`. 

333 n_2 

334 Refractive index of the transmitted medium :math:`n_2`. 

335 theta_i 

336 Incident angle :math:`\\theta_i` in degrees. 

337 theta_t 

338 Transmitted angle :math:`\\theta_t` in degrees. 

339 

340 Returns 

341 ------- 

342 :class:`numpy.ndarray` 

343 *Fresnel* reflection power coefficients (reflectance) for s and p 

344 polarisations stacked along the last axis. The array contains 

345 :math:`[R_s, R_p]`. 

346 

347 Notes 

348 ----- 

349 The *Fresnel* reflection power coefficients (reflectance) are given by: 

350 

351 .. math:: 

352 

353 R_s &= |r_s|^2 = \\left|\\frac{n_1 \\cos \\theta_i - n_2 \\cos \ 

354\\theta_t}{n_1 \\cos \\theta_i + n_2 \\cos \\theta_t}\\right|^2 \\\\ 

355 R_p &= |r_p|^2 = \\left|\\frac{n_1 \\cos \\theta_t - n_2 \\cos \ 

356\\theta_i}{n_1 \\cos \\theta_t + n_2 \\cos \\theta_i}\\right|^2 

357 

358 Where: 

359 

360 - :math:`R_s`: s-polarisation reflectance (fraction of incident power reflected) 

361 - :math:`R_p`: p-polarisation reflectance (fraction of incident power reflected) 

362 - :math:`r_s, r_p`: complex reflection amplitude coefficients 

363 - The s-polarisation electric field is perpendicular to the plane of incidence 

364 - The p-polarisation electric field is parallel to the plane of incidence 

365 

366 The reflectance values satisfy: :math:`0 \\leq R_s, R_p \\leq 1`. 

367 

368 References 

369 ---------- 

370 :cite:`Byrnes2016` 

371 

372 Examples 

373 -------- 

374 >>> result = polarised_light_reflection_coefficient(1.0, 1.5, 0.0, 0.0) 

375 >>> result.real 

376 array([ 0.04, 0.04]) 

377 """ 

378 

379 # Reflectance: R = |r|^2 (Byrnes Eq. 23) 

380 R = np.abs(polarised_light_reflection_amplitude(n_1, n_2, theta_i, theta_t)) ** 2 

381 

382 return as_complex_array(R) 

383 

384 

385def polarised_light_transmission_amplitude( 

386 n_1: ArrayLike, 

387 n_2: ArrayLike, 

388 theta_i: ArrayLike, 

389 theta_t: ArrayLike, 

390) -> NDArrayComplex: 

391 """ 

392 Compute *Fresnel* transmission amplitude coefficients. 

393 

394 This function computes the complex transmission amplitude coefficients for 

395 both s-polarisation (perpendicular) and p-polarisation (parallel) components 

396 of electromagnetic waves at a dielectric interface. 

397 

398 Parameters 

399 ---------- 

400 n_1 

401 Refractive index of the incident medium :math:`n_1`. 

402 n_2 

403 Refractive index of the transmitted medium :math:`n_2`. 

404 theta_i 

405 Incident angle :math:`\\theta_i` in degrees. 

406 theta_t 

407 Transmitted angle :math:`\\theta_t` in degrees. 

408 

409 Returns 

410 ------- 

411 :class:`numpy.ndarray` 

412 *Fresnel* transmission amplitude coefficients for s and p polarisations 

413 stacked along the last axis. The array contains :math:`[t_s, t_p]` where 

414 :math:`t_s` and :math:`t_p` are the complex transmission coefficients. 

415 

416 Notes 

417 ----- 

418 The *Fresnel* transmission amplitude coefficients are given by (*Equation 6* 

419 from :cite:`Byrnes2016`): 

420 

421 .. math:: 

422 

423 t_s &= \\frac{2n_1 \\cos \\theta_1}{n_1 \\cos \\theta_1 + n_2 \\cos \ 

424\\theta_2} \\\\ 

425 t_p &= \\frac{2n_1 \\cos \\theta_1}{n_2 \\cos \\theta_1 + n_1 \\cos \ 

426\\theta_2} 

427 

428 Where: 

429 

430 - :math:`t_s`: s-polarisation transmission amplitude (electric field perpendicular 

431 to the plane of incidence) 

432 - :math:`t_p`: p-polarisation transmission amplitude (electric field parallel 

433 to the plane of incidence) 

434 - :math:`n_1, n_2`: Refractive indices of incident and transmitted media 

435 - :math:`\\theta_1, \\theta_2`: Incident and transmitted angles 

436 

437 Examples 

438 -------- 

439 >>> polarised_light_transmission_amplitude(1.0, 1.5, 0.0, 0.0) 

440 array([ 0.8+0.j, 0.8+0.j]) 

441 """ 

442 

443 n_1_cos_theta_i, n_1_cos_theta_t, n_2_cos_theta_i, n_2_cos_theta_t = ( 

444 polarised_light_magnitude_elements(n_1, n_2, theta_i, theta_t) 

445 ) 

446 

447 two_n_1_cos_theta_i = 2 * n_1_cos_theta_i 

448 

449 # Fresnel transmission amplitudes (Byrnes Eq. 6) 

450 t_s = two_n_1_cos_theta_i / (n_1_cos_theta_i + n_2_cos_theta_t) 

451 t_p = two_n_1_cos_theta_i / (n_2_cos_theta_i + n_1_cos_theta_t) 

452 

453 return _tstack_complex([t_s, t_p]) 

454 

455 

456def polarised_light_transmission_coefficient( 

457 n_1: ArrayLike, 

458 n_2: ArrayLike, 

459 theta_i: ArrayLike, 

460 theta_t: ArrayLike, 

461) -> NDArrayComplex: 

462 """ 

463 Compute *Fresnel* transmission power coefficients (transmittance). 

464 

465 This function computes the transmission power coefficients, which represent 

466 the fraction of incident power that is transmitted through a dielectric interface 

467 for both s-polarisation (perpendicular) and p-polarisation (parallel) components. 

468 

469 Parameters 

470 ---------- 

471 n_1 

472 Refractive index of the incident medium :math:`n_1`. 

473 n_2 

474 Refractive index of the transmitted medium :math:`n_2`. 

475 theta_i 

476 Incident angle :math:`\\theta_i` in degrees. 

477 theta_t 

478 Transmitted angle :math:`\\theta_t` in degrees. 

479 

480 Returns 

481 ------- 

482 :class:`numpy.ndarray` 

483 *Fresnel* transmission power coefficients (transmittance) for s and p 

484 polarisations stacked along the last axis. The array contains 

485 :math:`[T_s, T_p]`. 

486 

487 Notes 

488 ----- 

489 The *Fresnel* transmission power coefficients (transmittance) are given by: 

490 

491 .. math:: 

492 

493 T_s &= \\frac{n_2 \\cos \\theta_t}{n_1 \\cos \\theta_i} |t_s|^2 = \ 

494\\frac{n_2 \\cos \\theta_t}{n_1 \\cos \\theta_i} \\left|\\frac{2n_1 \\cos \ 

495\\theta_i}{n_1 \\cos \\theta_i + n_2 \\cos \\theta_t}\\right|^2 \\\\ 

496 T_p &= \\frac{n_2 \\cos \\theta_t}{n_1 \\cos \\theta_i} |t_p|^2 = \ 

497\\frac{n_2 \\cos \\theta_t}{n_1 \\cos \\theta_i} \\left|\\frac{2n_1 \\cos \ 

498\\theta_i}{n_2 \\cos \\theta_i + n_1 \\cos \\theta_t}\\right|^2 

499 

500 Where: 

501 

502 - :math:`T_s`: s-polarisation transmittance (fraction of incident power 

503 transmitted) 

504 - :math:`T_p`: p-polarisation transmittance (fraction of incident power 

505 transmitted) 

506 - :math:`t_s, t_p`: complex transmission amplitude coefficients 

507 - The s-polarisation electric field is perpendicular to the plane of incidence 

508 - The p-polarisation electric field is parallel to the plane of incidence 

509 

510 The refractive index factor 

511 :math:`\\frac{n_2 \\cos \\theta_t}{n_1 \\cos \\theta_i}` accounts for the 

512 change in beam cross-section and energy density in the transmission medium. 

513 

514 **Energy Conservation**: For non-absorbing media: 

515 :math:`R_s + T_s = 1` and :math:`R_p + T_p = 1`, where :math:`R_s, R_p` are the 

516 corresponding reflectance coefficients. 

517 

518 The transmittance values satisfy: :math:`0 \\leq T_s, T_p \\leq 1`. 

519 

520 References 

521 ---------- 

522 :cite:`Byrnes2016` 

523 

524 Examples 

525 -------- 

526 >>> polarised_light_transmission_coefficient(1.0, 1.5, 0.0, 0.0) 

527 array([ 0.96+0.j, 0.96+0.j]) 

528 """ 

529 

530 n_1 = as_complex_array(n_1) 

531 n_2 = as_complex_array(n_2) 

532 

533 n_1_cos_theta_i, _n_1_cos_theta_t, _n_2_cos_theta_i, n_2_cos_theta_t = ( 

534 polarised_light_magnitude_elements(n_1, n_2, theta_i, theta_t) 

535 ) 

536 

537 # Transmittance with beam cross-section correction (Byrnes Eq. 21-22) 

538 T = (n_2_cos_theta_t / n_1_cos_theta_i)[..., None] * np.abs( 

539 polarised_light_transmission_amplitude(n_1, n_2, theta_i, theta_t) 

540 ) ** 2 

541 return as_complex_array(T) 

542 

543 

544@dataclass 

545class TransferMatrixResult(MixinDataclassArithmetic): 

546 """ 

547 Define the *Transfer Matrix Method* calculation results. 

548 

549 Parameters 

550 ---------- 

551 M_s 

552 Transfer matrix for s-polarisation :math:`M_s`, shape 

553 (..., wavelengths_count, 2, 2). 

554 M_p 

555 Transfer matrix for p-polarisation :math:`M_p`, shape 

556 (..., wavelengths_count, 2, 2). 

557 theta 

558 Propagation angles in each layer :math:`\\theta_j` (degrees), shape 

559 (..., n_layers+2). Includes [incident, layer_1, ..., layer_n, substrate]. 

560 n 

561 Complete multilayer stack :math:`n_j`, shape 

562 (..., n_layers+2, wavelengths_count). Includes 

563 [n_incident, n_layer_1, ..., n_layer_n, n_substrate]. 

564 

565 References 

566 ---------- 

567 :cite:`Byrnes2016` 

568 """ 

569 

570 M_s: NDArrayComplex 

571 M_p: NDArrayComplex 

572 theta: NDArrayFloat 

573 n: NDArrayComplex 

574 

575 

576def matrix_transfer_tmm( 

577 n: ArrayLike, 

578 t: ArrayLike, 

579 theta: ArrayLike, 

580 wavelength: ArrayLike, 

581) -> TransferMatrixResult: 

582 """ 

583 Calculate transfer matrices for multilayer thin film structures using the 

584 *Transfer Matrix Method*. 

585 

586 This function constructs the transfer matrices for s-polarised and 

587 p-polarised light propagating through a multilayer structure. The transfer 

588 matrices encode the optical properties of the structure and are used to 

589 calculate reflectance and transmittance. 

590 

591 Parameters 

592 ---------- 

593 n 

594 Complete refractive index stack :math:`n_j` including incident medium, 

595 layers, and substrate. Shape: (media_count,) or 

596 (media_count, wavelengths_count). Can be complex for absorbing 

597 materials. The array should contain [n_incident, n_layer_1, ..., 

598 n_layer_n, n_substrate]. 

599 t 

600 Thicknesses of each layer :math:`t_j` in nanometers (excluding incident 

601 and substrate). Shape: (layers_count,) or (thickness_count, layers_count). 

602 

603 - **1D array** ``[t1, t2, ...]``: One thickness per layer for a single 

604 multilayer configuration. Shape: ``(layers_count,)`` 

605 - **2D array** ``[[t1, t2, ...], [t1', t2', ...]]``: Multiple thickness 

606 configurations for outer product broadcasting. Shape: 

607 ``(thickness_count, layers_count)`` 

608 

609 Most users should use :func:`thin_film_tmm` or :func:`multilayer_tmm` 

610 instead, which provide simpler interfaces. 

611 theta 

612 Incident angle :math:`\\theta` in degrees. Scalar or array of shape 

613 (angles_count,) for angle broadcasting. 

614 wavelength 

615 Vacuum wavelength values :math:`\\lambda` in nanometers. 

616 

617 Returns 

618 ------- 

619 :class:`colour.TransferMatrixResult` 

620 Transfer matrix calculation results containing M_s, M_p, theta, 

621 and n arrays. 

622 

623 Examples 

624 -------- 

625 Single layer at one wavelength: 

626 

627 >>> result = matrix_transfer_tmm( 

628 ... n=[1.0, 1.5, 1.0], 

629 ... t=[250], 

630 ... theta=0, 

631 ... wavelength=550, 

632 ... ) 

633 >>> result.M_s.shape 

634 (1, 1, 1, 2, 2) 

635 >>> result.theta.shape 

636 (1, 3) 

637 

638 Multiple wavelengths: 

639 

640 >>> result = matrix_transfer_tmm( 

641 ... n=[1.0, 1.5, 1.0], 

642 ... t=[250], 

643 ... theta=0, 

644 ... wavelength=[400, 500, 600], 

645 ... ) 

646 >>> result.M_s.shape 

647 (3, 1, 1, 2, 2) 

648 

649 Multiple angles (angle broadcasting): 

650 

651 >>> result = matrix_transfer_tmm( 

652 ... n=[1.0, 1.5, 1.0], 

653 ... t=[250], 

654 ... theta=[0, 30, 45, 60], 

655 ... wavelength=[400, 500, 600], 

656 ... ) 

657 >>> result.M_s.shape 

658 (3, 4, 1, 2, 2) 

659 >>> result.theta.shape 

660 (4, 3) 

661 

662 Notes 

663 ----- 

664 - The *Transfer Matrix Method* relates the field amplitudes across the entire 

665 multilayer structure (*Equations 10-15* from :cite:`Byrnes2016`): 

666 

667 .. math:: 

668 

669 \\begin{pmatrix} v_n \\\\ w_n \\end{pmatrix} = M_n \\begin{pmatrix} \ 

670v_{n+1} \\\\ w_{n+1} \\end{pmatrix} 

671 

672 Where :math:`M_n` combines the layer propagation and interface matrices: 

673 

674 .. math:: 

675 

676 M_n = L_n \\cdot I_{n,n+1} = 

677 \\begin{pmatrix} 

678 e^{-i\\delta_n} & 0 \\\\ 0 & e^{i\\delta_n} 

679 \\end{pmatrix} 

680 \\frac{1}{t_{n,n+1}} 

681 \\begin{pmatrix} 

682 1 & r_{n,n+1} \\\\ r_{n,n+1} & 1 

683 \\end{pmatrix} 

684 

685 The overall transfer matrix :math:`\\tilde{M}` for the complete structure is: 

686 

687 .. math:: 

688 

689 \\tilde{M} = \\frac{1}{t_{0,1}} 

690 \\begin{pmatrix} 

691 1 & r_{0,1} \\\\ r_{0,1} & 1 

692 \\end{pmatrix} 

693 M_1 M_2 \\cdots M_{N-2} 

694 

695 From which the overall reflection and transmission coefficients are extracted: 

696 

697 .. math:: 

698 

699 \\begin{pmatrix} 1 \\\\ r \\end{pmatrix} = \\tilde{M} \ 

700\\begin{pmatrix} t \\\\ 0 \\end{pmatrix} 

701 

702 .. math:: 

703 

704 t = \\frac{1}{\\tilde{M}_{00}}, \\quad r = \ 

705\\frac{\\tilde{M}_{10}}{\\tilde{M}_{00}} 

706 

707 Where: 

708 

709 - :math:`v_n, w_n`: Forward and backward field amplitudes in layer :math:`n` 

710 - :math:`M_n`: Transfer matrix for layer :math:`n` 

711 - :math:`L_n`: Layer propagation matrix 

712 - :math:`I_{n,n+1}`: Interface matrix between layers :math:`n` and :math:`n+1` 

713 - :math:`\\tilde{M}`: Overall transfer matrix 

714 - :math:`r, t`: Overall reflection and transmission amplitude coefficients 

715 

716 - Supports complex refractive indices for absorbing materials. 

717 - **Angle broadcasting**: All computations are vectorized across angles. 

718 The output always includes the angle dimension. 

719 - The transfer matrices always have shape (angles_count, wavelengths_count, 2, 2), 

720 even for scalar theta (angles_count=1). 

721 

722 References 

723 ---------- 

724 :cite:`Byrnes2016` 

725 """ 

726 

727 n = as_complex_array(n) 

728 t = as_float_array(t) 

729 theta = np.atleast_1d(as_float_array(theta)) 

730 wavelength = np.atleast_1d(as_float_array(wavelength)) 

731 

732 angles_count = theta.shape[0] 

733 wavelengths_count = wavelength.shape[0] 

734 

735 # Convert 1D n to column vector and tile across wavelengths 

736 # (M,) -> (M, 1) -> (M, W) 

737 if n.ndim == 1: 

738 n = np.transpose(np.atleast_2d(n)) 

739 n = np.tile(n, (1, wavelengths_count)) 

740 

741 # (1, layers_count) 

742 if t.ndim == 1: 

743 t = t[np.newaxis, :] 

744 

745 media_count = n.shape[0] 

746 layers_count = media_count - 2 

747 

748 thickness_count = t.shape[0] 

749 

750 n_0 = n[0, 0] if n.ndim == 2 else n[0] 

751 

752 # Snell's law: n_i * sin(theta_i) = n_j * sin(theta_j) (Byrnes Eq. 3) 

753 # Broadcasting: theta (A,) → theta_media (A, M) 

754 theta_media = snell_law( 

755 n_0, (n[:, 0] if n.ndim == 2 else n)[:, None], theta[None, :] 

756 ).T 

757 

758 # Fresnel coefficients (Byrnes Eq. 6) 

759 # Broadcasting: n (M, W), theta_media (A, M) → coefficients (A, M-1, W) 

760 n_1 = n[:-1, :] # (M-1, W) 

761 n_2 = n[1:, :] # (M-1, W) 

762 theta_1 = theta_media[:, :-1] # (A, M-1) 

763 theta_2 = theta_media[:, 1:] # (A, M-1) 

764 

765 r_media_s, r_media_p = _tsplit_complex( 

766 polarised_light_reflection_amplitude( 

767 n_1[None, :, :], # (1, M-1, W) 

768 n_2[None, :, :], # (1, M-1, W) 

769 theta_1[:, :, None], # (A, M-1, 1) 

770 theta_2[:, :, None], # (A, M-1, 1) 

771 ) 

772 ) # Output: (A, M-1, W) 

773 

774 t_media_s, t_media_p = _tsplit_complex( 

775 polarised_light_transmission_amplitude( 

776 n_1[None, :, :], # (1, M-1, W) 

777 n_2[None, :, :], # (1, M-1, W) 

778 theta_1[:, :, None], # (A, M-1, 1) 

779 theta_2[:, :, None], # (A, M-1, 1) 

780 ) 

781 ) # Output: (A, M-1, W) 

782 

783 # Phase accumulation: delta = d * k_z (Byrnes Eq. 8) 

784 # Broadcasting directly in (W, A, T, L) order 

785 n_previous = n[0:layers_count, :] # (L, W) - Media before each layer 

786 n_layer = n[1 : layers_count + 1, :] # (L, W) - Each layer's refractive index 

787 theta_layer = theta_media[:, 0:layers_count] # (A, L) 

788 

789 theta_radians = np.radians(theta_layer)[:, :, None] # (A, L, 1) 

790 k_z_layers = np.sqrt( 

791 n_layer[None, :, :] ** 2 

792 - n_previous[None, :, :] ** 2 * np.sin(theta_radians) ** 2 

793 ) # (A, L, W) 

794 

795 # Compute phase: delta = (2π/λ) * d * k_z 

796 phase_factor = 2 * np.pi / wavelength[:, None, None, None] # (W, 1, 1, 1) 

797 # Reshape k_z from (A, L, W) to (W, A, 1, L) for broadcasting with thickness 

798 k_z = np.transpose(k_z_layers, (2, 0, 1))[:, :, None, :] # (W, A, 1, L) 

799 delta = phase_factor * t[None, None, :, :] * k_z # (W, A, T, L) 

800 

801 A = np.exp(1j * delta) # (W, A, T, L) 

802 

803 # Layer matrices: M_n = L_n * I_{n,n+1} (Byrnes Eq. 10-11) 

804 # (W, A, T, L, 2, 2, 2) for [wavelengths, angles, thickness, layers, 2x2, pol] 

805 M = zeros( 

806 (wavelengths_count, angles_count, thickness_count, layers_count, 2, 2, 2), 

807 dtype=DTYPE_COMPLEX_DEFAULT, # pyright: ignore 

808 ) 

809 

810 r_s = r_media_s[:, 1 : layers_count + 1, :] # (A, L, W) 

811 r_p = r_media_p[:, 1 : layers_count + 1, :] # (A, L, W) 

812 t_s = t_media_s[:, 1 : layers_count + 1, :] # (A, L, W) 

813 t_p = t_media_p[:, 1 : layers_count + 1, :] # (A, L, W) 

814 

815 # Broadcast Fresnel coefficients from (A, L, W) to (W, A, 1, L) 

816 # (A,L,W) -> (W,A,L) -> (W,A,1,L) 

817 r_s_b = np.transpose(r_s, (2, 0, 1))[:, :, None, :] 

818 r_p_b = np.transpose(r_p, (2, 0, 1))[:, :, None, :] 

819 t_s_b = np.transpose(t_s, (2, 0, 1))[:, :, None, :] 

820 t_p_b = np.transpose(t_p, (2, 0, 1))[:, :, None, :] 

821 

822 M[:, :, :, :, 0, 0, 0] = 1 / (A * t_s_b) 

823 M[:, :, :, :, 0, 1, 0] = r_s_b / (A * t_s_b) 

824 M[:, :, :, :, 1, 0, 0] = A * r_s_b / t_s_b 

825 M[:, :, :, :, 1, 1, 0] = A / t_s_b 

826 

827 M[:, :, :, :, 0, 0, 1] = 1 / (A * t_p_b) 

828 M[:, :, :, :, 0, 1, 1] = r_p_b / (A * t_p_b) 

829 M[:, :, :, :, 1, 0, 1] = A * r_p_b / t_p_b 

830 M[:, :, :, :, 1, 1, 1] = A / t_p_b 

831 

832 # Initial interface matrix (Byrnes Eq. 11) 

833 # Shape: (W, A, T, 2, 2) 

834 M_s = zeros( 

835 (wavelengths_count, angles_count, thickness_count, 2, 2), 

836 dtype=DTYPE_COMPLEX_DEFAULT, # pyright: ignore 

837 ) 

838 # Fresnel coefficients at incident → first layer interface 

839 t_s_01 = t_media_s[:, 0, :] # (A, W) 

840 r_s_01 = r_media_s[:, 0, :] # (A, W) 

841 M_s[:, :, :, 0, 0] = (1 / t_s_01).T[:, :, None] # (W, A, 1) 

842 M_s[:, :, :, 0, 1] = (r_s_01 / t_s_01).T[:, :, None] 

843 M_s[:, :, :, 1, 0] = (r_s_01 / t_s_01).T[:, :, None] 

844 M_s[:, :, :, 1, 1] = (1 / t_s_01).T[:, :, None] 

845 

846 M_p = zeros( 

847 (wavelengths_count, angles_count, thickness_count, 2, 2), 

848 dtype=DTYPE_COMPLEX_DEFAULT, # pyright: ignore 

849 ) 

850 t_p_01 = t_media_p[:, 0, :] # (A, W) 

851 r_p_01 = r_media_p[:, 0, :] # (A, W) 

852 M_p[:, :, :, 0, 0] = (1 / t_p_01).T[:, :, None] 

853 M_p[:, :, :, 0, 1] = (r_p_01 / t_p_01).T[:, :, None] 

854 M_p[:, :, :, 1, 0] = (r_p_01 / t_p_01).T[:, :, None] 

855 M_p[:, :, :, 1, 1] = (1 / t_p_01).T[:, :, None] 

856 

857 # Overall transfer matrix: M_tilde = I_01 @ M_1 @ M_2 @ ... (Byrnes Eq. 12) 

858 for i in range(layers_count): 

859 M_s = np.matmul(M_s, M[:, :, :, i, :, :, 0]) 

860 M_p = np.matmul(M_p, M[:, :, :, i, :, :, 1]) 

861 

862 return TransferMatrixResult( 

863 M_s=M_s, 

864 M_p=M_p, 

865 theta=theta_media, 

866 n=n, 

867 )