1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 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 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 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 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597
|
<!DOCTYPE html>
<html class="writer-html5" lang="en" data-content_root="./">
<head>
<meta charset="utf-8" /><meta name="viewport" content="width=device-width, initial-scale=1" />
<meta name="viewport" content="width=device-width, initial-scale=1.0" />
<title>pymatgen.analysis.interfaces package — pymatgen 2025.6.14 documentation</title>
<link rel="stylesheet" type="text/css" href="assets/pygments.css?v=03e43079" />
<link rel="stylesheet" type="text/css" href="assets/css/theme.css?v=e59714d7" />
<link rel="stylesheet" type="text/css" href="assets/css/custom.css" />
<link rel="canonical" href="https://pymatgen.orgpymatgen.analysis.interfaces.html"/>
<script src="assets/jquery.js?v=5d32c60e"></script>
<script src="assets/_sphinx_javascript_frameworks_compat.js?v=2cd50e6c"></script>
<script src="assets/documentation_options.js?v=03bec786"></script>
<script src="assets/doctools.js?v=9bcbadda"></script>
<script src="assets/sphinx_highlight.js?v=dc90522c"></script>
<script src="assets/js/theme.js"></script>
<link rel="index" title="Index" href="genindex.html" />
<link rel="search" title="Search" href="search.html" />
</head>
<body class="wy-body-for-nav">
<div class="wy-grid-for-nav">
<nav data-toggle="wy-nav-shift" class="wy-nav-side">
<div class="wy-side-scroll">
<div class="wy-side-nav-search" style="background: linear-gradient(0deg, rgba(23,63,162,1) 0%, rgba(0,70,192,1) 100%)" >
<a href="index.html">
</a>
<div role="search">
<form id="rtd-search-form" class="wy-form" action="search.html" method="get">
<input type="text" name="q" placeholder="Search docs" aria-label="Search docs" />
<input type="hidden" name="check_keywords" value="yes" />
<input type="hidden" name="area" value="default" />
</form>
</div>
</div><div class="wy-menu wy-menu-vertical" data-spy="affix" role="navigation" aria-label="Navigation menu">
<!-- Local TOC -->
<div class="local-toc"><ul>
<li><a class="reference internal" href="#">pymatgen.analysis.interfaces package</a><ul>
<li><a class="reference internal" href="#submodules">Submodules</a></li>
<li><a class="reference internal" href="#module-pymatgen.analysis.interfaces.coherent_interfaces">pymatgen.analysis.interfaces.coherent_interfaces module</a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.coherent_interfaces.CoherentInterfaceBuilder"><code class="docutils literal notranslate"><span class="pre">CoherentInterfaceBuilder</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.coherent_interfaces.CoherentInterfaceBuilder.get_interfaces"><code class="docutils literal notranslate"><span class="pre">CoherentInterfaceBuilder.get_interfaces()</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.coherent_interfaces.from_2d_to_3d"><code class="docutils literal notranslate"><span class="pre">from_2d_to_3d()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.coherent_interfaces.get_2d_transform"><code class="docutils literal notranslate"><span class="pre">get_2d_transform()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.coherent_interfaces.get_rot_3d_for_2d"><code class="docutils literal notranslate"><span class="pre">get_rot_3d_for_2d()</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#module-pymatgen.analysis.interfaces.substrate_analyzer">pymatgen.analysis.interfaces.substrate_analyzer module</a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer"><code class="docutils literal notranslate"><span class="pre">SubstrateAnalyzer</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer.calculate"><code class="docutils literal notranslate"><span class="pre">SubstrateAnalyzer.calculate()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer.generate_surface_vectors"><code class="docutils literal notranslate"><span class="pre">SubstrateAnalyzer.generate_surface_vectors()</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.elastic_energy"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.elastic_energy</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.film_miller"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.film_miller</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.from_zsl"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.from_zsl()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.ground_state_energy"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.ground_state_energy</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.strain"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.strain</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.substrate_miller"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.substrate_miller</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.total_energy"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.total_energy</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.von_mises_strain"><code class="docutils literal notranslate"><span class="pre">SubstrateMatch.von_mises_strain</span></code></a></li>
</ul>
</li>
</ul>
</li>
<li><a class="reference internal" href="#module-pymatgen.analysis.interfaces.zsl">pymatgen.analysis.interfaces.zsl module</a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator"><code class="docutils literal notranslate"><span class="pre">ZSLGenerator</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator.generate_sl_transformation_sets"><code class="docutils literal notranslate"><span class="pre">ZSLGenerator.generate_sl_transformation_sets()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator.get_equiv_transformations"><code class="docutils literal notranslate"><span class="pre">ZSLGenerator.get_equiv_transformations()</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch"><code class="docutils literal notranslate"><span class="pre">ZSLMatch</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.film_sl_vectors"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.film_sl_vectors</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.film_transformation"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.film_transformation</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.film_vectors"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.film_vectors</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.match_area"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.match_area</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.match_transformation"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.match_transformation</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_sl_vectors"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.substrate_sl_vectors</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_transformation"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.substrate_transformation</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_vectors"><code class="docutils literal notranslate"><span class="pre">ZSLMatch.substrate_vectors</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.fast_norm"><code class="docutils literal notranslate"><span class="pre">fast_norm()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.gen_sl_transform_matrices"><code class="docutils literal notranslate"><span class="pre">gen_sl_transform_matrices()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.get_factors"><code class="docutils literal notranslate"><span class="pre">get_factors()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.is_same_vectors"><code class="docutils literal notranslate"><span class="pre">is_same_vectors()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.reduce_vectors"><code class="docutils literal notranslate"><span class="pre">reduce_vectors()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.rel_angle"><code class="docutils literal notranslate"><span class="pre">rel_angle()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.rel_strain"><code class="docutils literal notranslate"><span class="pre">rel_strain()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.vec_angle"><code class="docutils literal notranslate"><span class="pre">vec_angle()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.vec_area"><code class="docutils literal notranslate"><span class="pre">vec_area()</span></code></a></li>
</ul>
</li>
</ul>
</li>
</ul>
</div>
</div>
</div>
</nav>
<section data-toggle="wy-nav-shift" class="wy-nav-content-wrap"><nav class="wy-nav-top" aria-label="Mobile navigation menu" style="background: linear-gradient(0deg, rgba(23,63,162,1) 0%, rgba(0,70,192,1) 100%)" >
<i data-toggle="wy-nav-top" class="fa fa-bars"></i>
<a href="index.html">pymatgen</a>
</nav>
<div class="wy-nav-content">
<div class="rst-content style-external-links">
<div role="navigation" aria-label="Page navigation">
<ul class="wy-breadcrumbs">
<li><a href="index.html" class="icon icon-home" aria-label="Home"></a></li>
<li class="breadcrumb-item active">pymatgen.analysis.interfaces package</li>
<li class="wy-breadcrumbs-aside">
<a href="https://github.com/materialsproject/pymatgen/blob/master/docs_rst/pymatgen.analysis.interfaces.rst" class="fa fa-github"> Edit on GitHub</a>
</li>
</ul>
<hr/>
</div>
<div role="main" class="document" itemscope="itemscope" itemtype="http://schema.org/Article">
<div itemprop="articleBody">
<section id="module-pymatgen.analysis.interfaces">
<span id="pymatgen-analysis-interfaces-package"></span><h1>pymatgen.analysis.interfaces package<a class="headerlink" href="#module-pymatgen.analysis.interfaces" title="Link to this heading"></a></h1>
<p>Module that implements various algorithms related to interface construction and analysis.</p>
<section id="submodules">
<h2>Submodules<a class="headerlink" href="#submodules" title="Link to this heading"></a></h2>
</section>
<section id="module-pymatgen.analysis.interfaces.coherent_interfaces">
<span id="pymatgen-analysis-interfaces-coherent-interfaces-module"></span><h2>pymatgen.analysis.interfaces.coherent_interfaces module<a class="headerlink" href="#module-pymatgen.analysis.interfaces.coherent_interfaces" title="Link to this heading"></a></h2>
<p>This module provides classes to store, generate, and manipulate material interfaces.</p>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.coherent_interfaces.CoherentInterfaceBuilder">
<em class="property"><span class="k"><span class="pre">class</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">CoherentInterfaceBuilder</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">substrate_structure</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><span class="pre">Structure</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_structure</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><span class="pre">Structure</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_miller</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">tuple</span><span class="p"><span class="pre">[</span></span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">]</span></span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_miller</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">tuple</span><span class="p"><span class="pre">[</span></span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">]</span></span></span></em>, <em class="sig-param"><span class="n"><span class="pre">zslgen</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator" title="pymatgen.analysis.interfaces.zsl.ZSLGenerator"><span class="pre">ZSLGenerator</span></a><span class="w"> </span><span class="p"><span class="pre">|</span></span><span class="w"> </span><span class="pre">None</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">None</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">termination_ftol</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">0.25</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">label_index</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">bool</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">False</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">filter_out_sym_slabs</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">bool</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">True</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/coherent_interfaces.py#L22-L258"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.coherent_interfaces.CoherentInterfaceBuilder" title="Link to this definition"></a></dt>
<dd><p>Bases: <code class="xref py py-class docutils literal notranslate"><span class="pre">object</span></code></p>
<p>This class constructs the coherent interfaces between two crystalline slabs
Coherency is defined by matching lattices not sub-planes.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>substrate_structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – substrate structure</p></li>
<li><p><strong>film_structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – film structure</p></li>
<li><p><strong>film_miller</strong> (<em>tuple</em><em>[</em><em>int</em><em>, </em><em>int</em><em>, </em><em>int</em><em>]</em>) – miller index for the film layer</p></li>
<li><p><strong>substrate_miller</strong> (<em>tuple</em><em>[</em><em>int</em><em>, </em><em>int</em><em>, </em><em>int</em><em>]</em>) – miller index for the substrate layer</p></li>
<li><p><strong>zslgen</strong> (<a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator" title="pymatgen.analysis.interfaces.zsl.ZSLGenerator"><em>ZSLGenerator</em></a><em> | </em><em>None</em>) – BiDirectionalZSL if you want custom lattice matching tolerances for coherency.</p></li>
<li><p><strong>termination_ftol</strong> (<em>float</em>) – tolerance to distinguish different terminating atomic planes.</p></li>
<li><p><strong>label_index</strong> (<em>bool</em>) – If True add an extra index at the beginning of the termination label.</p></li>
<li><p><strong>filter_out_sym_slabs</strong> (<em>bool</em>) – If True filter out identical slabs with different terminations.
This might need to be set as False to find more non-identical terminations because slab
identity separately does not mean combinational identity.</p></li>
</ul>
</dd>
</dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.coherent_interfaces.CoherentInterfaceBuilder.get_interfaces">
<span class="sig-name descname"><span class="pre">get_interfaces</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">termination</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">tuple</span><span class="p"><span class="pre">[</span></span><span class="pre">str</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">str</span><span class="p"><span class="pre">]</span></span></span></em>, <em class="sig-param"><span class="n"><span class="pre">gap</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">2.0</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">vacuum_over_film</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">20.0</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_thickness</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">1</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_thickness</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">1</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">in_layers</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">bool</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">True</span></span></em><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">Iterator</span><span class="p"><span class="pre">[</span></span><a class="reference internal" href="pymatgen.core.html#pymatgen.core.interface.Interface" title="pymatgen.core.interface.Interface"><span class="pre">Interface</span></a><span class="p"><span class="pre">]</span></span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/coherent_interfaces.py#L164-L258"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.coherent_interfaces.CoherentInterfaceBuilder.get_interfaces" title="Link to this definition"></a></dt>
<dd><p>Generate interface structures given the film and substrate structure
as well as the desired terminations.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>termination</strong> (<em>tuple</em><em>[</em><em>str</em><em>, </em><em>str</em><em>]</em>) – termination from self.termination list</p></li>
<li><p><strong>gap</strong> (<em>float</em><em>, </em><em>optional</em>) – gap between film and substrate. Defaults to 2.0.</p></li>
<li><p><strong>vacuum_over_film</strong> (<em>float</em><em>, </em><em>optional</em>) – vacuum over the top of the film. Defaults to 20.0.</p></li>
<li><p><strong>film_thickness</strong> (<em>float</em><em>, </em><em>optional</em>) – the film thickness. Defaults to 1.</p></li>
<li><p><strong>substrate_thickness</strong> (<em>float</em><em>, </em><em>optional</em>) – substrate thickness. Defaults to 1.</p></li>
<li><p><strong>in_layers</strong> (<em>bool</em><em>, </em><em>optional</em>) – set the thickness in layer units. Defaults to True.</p></li>
</ul>
</dd>
<dt class="field-even">Yields<span class="colon">:</span></dt>
<dd class="field-even"><p><em>Iterator[Interface]</em> – interfaces from slabs</p>
</dd>
</dl>
</dd></dl>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.coherent_interfaces.from_2d_to_3d">
<span class="sig-name descname"><span class="pre">from_2d_to_3d</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">mat</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">ndarray</span></span></em><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">ndarray</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/coherent_interfaces.py#L292-L296"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.coherent_interfaces.from_2d_to_3d" title="Link to this definition"></a></dt>
<dd><p>Convert a 2D matrix to a 3D matrix.</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.coherent_interfaces.get_2d_transform">
<span class="sig-name descname"><span class="pre">get_2d_transform</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">start</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">Sequence</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">end</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">Sequence</span></span></em><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">np.ndarray</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/coherent_interfaces.py#L284-L289"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.coherent_interfaces.get_2d_transform" title="Link to this definition"></a></dt>
<dd><p>Gets a 2d transformation matrix
that converts start to end.</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.coherent_interfaces.get_rot_3d_for_2d">
<span class="sig-name descname"><span class="pre">get_rot_3d_for_2d</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">film_matrix</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">sub_matrix</span></span></em><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">ndarray</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/coherent_interfaces.py#L261-L281"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.coherent_interfaces.get_rot_3d_for_2d" title="Link to this definition"></a></dt>
<dd><p>Find transformation matrix that will rotate and strain the film to the substrate while preserving the c-axis.</p>
</dd></dl>
</section>
<section id="module-pymatgen.analysis.interfaces.substrate_analyzer">
<span id="pymatgen-analysis-interfaces-substrate-analyzer-module"></span><h2>pymatgen.analysis.interfaces.substrate_analyzer module<a class="headerlink" href="#module-pymatgen.analysis.interfaces.substrate_analyzer" title="Link to this heading"></a></h2>
<p>This module provides classes to identify optimal substrates for film growth.</p>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer">
<em class="property"><span class="k"><span class="pre">class</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">SubstrateAnalyzer</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">film_max_miller</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">1</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_max_miller</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">1</span></span></em>, <em class="sig-param"><span class="o"><span class="pre">**</span></span><span class="n"><span class="pre">kwargs</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/substrate_analyzer.py#L86-L203"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer" title="Link to this definition"></a></dt>
<dd><p>Bases: <a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator" title="pymatgen.analysis.interfaces.zsl.ZSLGenerator"><code class="xref py py-class docutils literal notranslate"><span class="pre">ZSLGenerator</span></code></a></p>
<p>This class applies a set of search criteria to identify suitable
substrates for film growth. It first uses a topological search by Zur
and McGill to identify matching super-lattices on various faces of the
two materials. Additional criteria can then be used to identify the most
suitable substrate. Currently, the only additional criteria is the
elastic strain energy of the super-lattices.</p>
<p>Initialize the substrate analyzer.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>zslgen</strong> (<a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator" title="pymatgen.analysis.interfaces.zsl.ZSLGenerator"><em>ZSLGenerator</em></a>) – Defaults to a ZSLGenerator with standard
tolerances, but can be fed one with custom tolerances</p></li>
<li><p><strong>film_max_miller</strong> (<em>int</em>) – maximum miller index to generate for film
surfaces</p></li>
<li><p><strong>substrate_max_miller</strong> (<em>int</em>) – maximum miller index to generate for
substrate surfaces.</p></li>
</ul>
</dd>
</dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer.calculate">
<span class="sig-name descname"><span class="pre">calculate</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">film</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><span class="pre">Structure</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><span class="pre">Structure</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">elasticity_tensor</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">None</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_millers</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">ArrayLike</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">None</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_millers</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">ArrayLike</span></span><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">None</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">ground_state_energy</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">0</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">lowest</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">False</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/substrate_analyzer.py#L150-L203"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer.calculate" title="Link to this definition"></a></dt>
<dd><p>Find all topological matches for the substrate and calculates elastic
strain energy and total energy for the film if elasticity tensor and
ground state energy are provided:</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>film</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – conventional standard structure for the film</p></li>
<li><p><strong>substrate</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – conventional standard structure for the
substrate</p></li>
<li><p><strong>elasticity_tensor</strong> (<a class="reference internal" href="pymatgen.analysis.elasticity.html#pymatgen.analysis.elasticity.elastic.ElasticTensor" title="pymatgen.analysis.elasticity.elastic.ElasticTensor"><em>ElasticTensor</em></a>) – elasticity tensor for the film
in the IEEE orientation</p></li>
<li><p><strong>film_millers</strong> (<em>array</em>) – film facets to consider in search as defined by
miller indices</p></li>
<li><p><strong>substrate_millers</strong> (<em>array</em>) – substrate facets to consider in search as
defined by miller indices</p></li>
<li><p><strong>ground_state_energy</strong> (<em>float</em>) – ground state energy for the film</p></li>
<li><p><strong>lowest</strong> (<em>bool</em>) – only consider lowest matching area for each surface</p></li>
</ul>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer.generate_surface_vectors">
<span class="sig-name descname"><span class="pre">generate_surface_vectors</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">film</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><span class="pre">Structure</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><span class="pre">Structure</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_millers</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">ArrayLike</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_millers</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">ArrayLike</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/substrate_analyzer.py#L112-L148"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateAnalyzer.generate_surface_vectors" title="Link to this definition"></a></dt>
<dd><p>Generate the film/substrate slab combinations for a set of given
miller indices.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>film</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – film structure</p></li>
<li><p><strong>substrate</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – substrate structure</p></li>
<li><p><strong>film_millers</strong> (<em>array</em>) – all miller indices to generate slabs for
film</p></li>
<li><p><strong>substrate_millers</strong> (<em>array</em>) – all miller indices to generate slabs
for substrate</p></li>
</ul>
</dd>
</dl>
</dd></dl>
</dd></dl>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch">
<em class="property"><span class="k"><span class="pre">class</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">SubstrateMatch</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">film_sl_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_sl_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_transformation</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_transformation</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_miller</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">tuple</span><span class="p"><span class="pre">[</span></span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">]</span></span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_miller</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">tuple</span><span class="p"><span class="pre">[</span></span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">]</span></span></span></em>, <em class="sig-param"><span class="n"><span class="pre">strain</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.analysis.elasticity.html#pymatgen.analysis.elasticity.strain.Strain" title="pymatgen.analysis.elasticity.strain.Strain"><span class="pre">Strain</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">von_mises_strain</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">ground_state_energy</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">elastic_energy</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/substrate_analyzer.py#L19-L83"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch" title="Link to this definition"></a></dt>
<dd><p>Bases: <a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch" title="pymatgen.analysis.interfaces.zsl.ZSLMatch"><code class="xref py py-class docutils literal notranslate"><span class="pre">ZSLMatch</span></code></a></p>
<p>A substrate match building on the Zur and McGill algorithm. This match class includes the miller
planes of the film and substrate the full strain tensor, the Von Mises strain, the ground state
energy if provided, and the elastic energy.</p>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.elastic_energy">
<span class="sig-name descname"><span class="pre">elastic_energy</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">float</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/substrate_analyzer.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.elastic_energy" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.film_miller">
<span class="sig-name descname"><span class="pre">film_miller</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">tuple</span><span class="p"><span class="pre">[</span></span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">]</span></span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/substrate_analyzer.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.film_miller" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.from_zsl">
<em class="property"><span class="k"><span class="pre">classmethod</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">from_zsl</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">match</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch" title="pymatgen.analysis.interfaces.zsl.ZSLMatch"><span class="pre">ZSLMatch</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">film</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><span class="pre">Structure</span></a></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_miller</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_miller</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">elasticity_tensor</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">None</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">ground_state_energy</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">0</span></span></em><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">Self</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/substrate_analyzer.py#L34-L78"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.from_zsl" title="Link to this definition"></a></dt>
<dd><p>Generate a substrate match from a ZSL match plus metadata.</p>
</dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.ground_state_energy">
<span class="sig-name descname"><span class="pre">ground_state_energy</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">float</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/substrate_analyzer.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.ground_state_energy" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.strain">
<span class="sig-name descname"><span class="pre">strain</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><a class="reference internal" href="pymatgen.analysis.elasticity.html#pymatgen.analysis.elasticity.strain.Strain" title="pymatgen.analysis.elasticity.strain.Strain"><span class="pre">Strain</span></a></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/substrate_analyzer.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.strain" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.substrate_miller">
<span class="sig-name descname"><span class="pre">substrate_miller</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">tuple</span><span class="p"><span class="pre">[</span></span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">int</span><span class="p"><span class="pre">]</span></span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/substrate_analyzer.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.substrate_miller" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py property">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.total_energy">
<em class="property"><span class="k"><span class="pre">property</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">total_energy</span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/substrate_analyzer.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.total_energy" title="Link to this definition"></a></dt>
<dd><p>Total energy of this match.</p>
</dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.von_mises_strain">
<span class="sig-name descname"><span class="pre">von_mises_strain</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">float</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/substrate_analyzer.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.substrate_analyzer.SubstrateMatch.von_mises_strain" title="Link to this definition"></a></dt>
<dd></dd></dl>
</dd></dl>
</section>
<section id="module-pymatgen.analysis.interfaces.zsl">
<span id="pymatgen-analysis-interfaces-zsl-module"></span><h2>pymatgen.analysis.interfaces.zsl module<a class="headerlink" href="#module-pymatgen.analysis.interfaces.zsl" title="Link to this heading"></a></h2>
<p>This module implements the Zur and McGill lattice matching algorithm.</p>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLGenerator">
<em class="property"><span class="k"><span class="pre">class</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">ZSLGenerator</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">max_area_ratio_tol</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">0.09</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">max_area</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">400</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">max_length_tol</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">0.03</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">max_angle_tol</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">0.01</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">bidirectional</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">False</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L62-L217"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator" title="Link to this definition"></a></dt>
<dd><p>Bases: <code class="xref py py-class docutils literal notranslate"><span class="pre">MSONable</span></code></p>
<p>This class generate matching interface super lattices based on the methodology
of lattice vector matching for heterostructural interfaces proposed by
Zur and McGill:
Journal of Applied Physics 55 (1984), 378 ; doi: 10.1063/1.333084
The process of generating all possible matching super lattices is:
1.) Reduce the surface lattice vectors and calculate area for the surfaces
2.) Generate all super lattice transformations within a maximum allowed area</p>
<blockquote>
<div><p>limit that give nearly equal area super-lattices for the two
surfaces - generate_sl_transformation_sets</p>
</div></blockquote>
<dl>
<dt>3.) For each superlattice set:</dt><dd><p>1.) Reduce super lattice vectors
2.) Check length and angle between film and substrate super lattice</p>
<blockquote>
<div><p>vectors to determine if the super lattices are the nearly same
and therefore coincident - get_equiv_transformations.</p>
</div></blockquote>
</dd>
<dt>Initialize a Zur Super Lattice Generator for a specific film and</dt><dd><p>substrate.</p>
</dd>
</dl>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>max_area_ratio_tol</strong> (<em>float</em>) – Max tolerance on ratio of
super-lattices to consider equal</p></li>
<li><p><strong>max_area</strong> (<em>float</em>) – max super lattice area to generate in search</p></li>
<li><p><strong>max_length_tol</strong> – maximum length tolerance in checking if two
vectors are of nearly the same length</p></li>
<li><p><strong>max_angle_tol</strong> – maximum angle tolerance in checking of two sets
of vectors have nearly the same angle between them.</p></li>
</ul>
</dd>
</dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLGenerator.generate_sl_transformation_sets">
<span class="sig-name descname"><span class="pre">generate_sl_transformation_sets</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">film_area</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">int</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_area</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">int</span></span></em><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">Iterator</span><span class="p"><span class="pre">[</span></span><span class="pre">tuple</span><span class="p"><span class="pre">]</span></span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L111-L144"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator.generate_sl_transformation_sets" title="Link to this definition"></a></dt>
<dd><p>Generate transformation sets for film/substrate pair given the
area of the unit cell area for the film and substrate. The
transformation sets map the film and substrate unit cells to super
lattices with a maximum area.</p>
<dl class="field-list">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>film_area</strong> (<em>int</em>) – the unit cell area for the film</p></li>
<li><p><strong>substrate_area</strong> (<em>int</em>) – the unit cell area for the substrate</p></li>
</ul>
</dd>
<dt class="field-even">Yields<span class="colon">:</span></dt>
<dd class="field-even"><p><em>transformation_sets</em> –</p>
<dl class="simple">
<dt>a set of transformation_sets defined as:</dt><dd><p>1.) the transformation matrices for the film to create a
super lattice of area i*film area
2.) the transformation matrices for the substrate to create
a super lattice of area j*film area.</p>
</dd>
</dl>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLGenerator.get_equiv_transformations">
<span class="sig-name descname"><span class="pre">get_equiv_transformations</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">transformation_sets</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_vectors</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_vectors</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L146-L188"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLGenerator.get_equiv_transformations" title="Link to this definition"></a></dt>
<dd><p>Applies the transformation_sets to the film and substrate vectors
to generate super-lattices and checks if they matches.
Returns all matching vectors sets.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>transformation_sets</strong> (<em>array</em>) – an array of transformation sets:
each transformation set is an array with the (i,j)
indicating the area multiples of the film and substrate it
corresponds to, an array with all possible transformations
for the film area multiple i and another array for the
substrate area multiple j.</p></li>
<li><p><strong>film_vectors</strong> (<em>array</em>) – film vectors to generate super lattices</p></li>
<li><p><strong>substrate_vectors</strong> (<em>array</em>) – substrate vectors to generate super
lattices</p></li>
</ul>
</dd>
</dl>
</dd></dl>
</dd></dl>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch">
<em class="property"><span class="k"><span class="pre">class</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">ZSLMatch</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">film_sl_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_sl_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_vectors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">film_transformation</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">substrate_transformation</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L19-L59"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch" title="Link to this definition"></a></dt>
<dd><p>Bases: <code class="xref py py-class docutils literal notranslate"><span class="pre">MSONable</span></code></p>
<p>A match from the Zur and McGill Algorithm. The super_lattice vectors are listed
as _sl_vectors. These are reduced according to the algorithm in the paper which
effectively a rotation in 3D space. Use the match_transformation property to get
the appropriate transformation matrix.</p>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.film_sl_vectors">
<span class="sig-name descname"><span class="pre">film_sl_vectors</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">list</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.film_sl_vectors" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.film_transformation">
<span class="sig-name descname"><span class="pre">film_transformation</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">list</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.film_transformation" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.film_vectors">
<span class="sig-name descname"><span class="pre">film_vectors</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">list</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.film_vectors" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py property">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.match_area">
<em class="property"><span class="k"><span class="pre">property</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">match_area</span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.match_area" title="Link to this definition"></a></dt>
<dd><p>The area of the match between the substrate and film super lattice vectors.</p>
</dd></dl>
<dl class="py property">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.match_transformation">
<em class="property"><span class="k"><span class="pre">property</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">match_transformation</span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.match_transformation" title="Link to this definition"></a></dt>
<dd><p>The transformation matrix to convert the film super lattice vectors to the substrate.</p>
</dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_sl_vectors">
<span class="sig-name descname"><span class="pre">substrate_sl_vectors</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">list</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_sl_vectors" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_transformation">
<span class="sig-name descname"><span class="pre">substrate_transformation</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">list</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_transformation" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_vectors">
<span class="sig-name descname"><span class="pre">substrate_vectors</span></span><em class="property"><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="pre">list</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/interfaces/zsl.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.ZSLMatch.substrate_vectors" title="Link to this definition"></a></dt>
<dd></dd></dl>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.fast_norm">
<span class="sig-name descname"><span class="pre">fast_norm</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">a</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L263-L270"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.fast_norm" title="Link to this definition"></a></dt>
<dd><p>Much faster variant of numpy linalg norm.</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.gen_sl_transform_matrices">
<span class="sig-name descname"><span class="pre">gen_sl_transform_matrices</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">area_multiple</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L220-L242"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.gen_sl_transform_matrices" title="Link to this definition"></a></dt>
<dd><p>Generates the transformation matrices that convert a set of 2D
vectors into a super lattice of integer area multiple as proven
in Cassels:</p>
<p>Cassels, John William Scott. An introduction to the geometry of
numbers. Springer Science & Business Media, 2012.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>area_multiple</strong> (<em>int</em>) – integer multiple of unit cell area for super</p></li>
<li><p><strong>area</strong> (<em>lattice</em>)</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>transformation matrices to convert unit vectors to
super lattice vectors</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p>matrix_list</p>
</dd>
</dl>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.get_factors">
<span class="sig-name descname"><span class="pre">get_factors</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">n</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L313-L318"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.get_factors" title="Link to this definition"></a></dt>
<dd><p>Generate all factors of n.</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.is_same_vectors">
<span class="sig-name descname"><span class="pre">is_same_vectors</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">vec_set1</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">vec_set2</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">bidirectional</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">False</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">max_length_tol</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">0.03</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">max_angle_tol</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">0.01</span></span></em><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">bool</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L345-L364"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.is_same_vectors" title="Link to this definition"></a></dt>
<dd><p>Determine if two sets of vectors are the same within length and angle
tolerances
:param vec_set1: an array of two vectors
:type vec_set1: array[array]
:param vec_set2: second array of two vectors.
:type vec_set2: array[array]</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.reduce_vectors">
<span class="sig-name descname"><span class="pre">reduce_vectors</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">a</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">b</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L287-L310"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.reduce_vectors" title="Link to this definition"></a></dt>
<dd><p>Generate independent and unique basis vectors based on the
methodology of Zur and McGill.</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.rel_angle">
<span class="sig-name descname"><span class="pre">rel_angle</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">vec_set1</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">vec_set2</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L251-L260"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.rel_angle" title="Link to this definition"></a></dt>
<dd><p>Calculate the relative angle between two vector sets.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>vec_set1</strong> (<em>array</em><em>[</em><em>array</em><em>]</em>) – an array of two vectors</p></li>
<li><p><strong>vec_set2</strong> (<em>array</em><em>[</em><em>array</em><em>]</em>) – second array of two vectors</p></li>
</ul>
</dd>
</dl>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.rel_strain">
<span class="sig-name descname"><span class="pre">rel_strain</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">vec1</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">vec2</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L245-L248"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.rel_strain" title="Link to this definition"></a></dt>
<dd><p>Calculate relative strain between two vectors.</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.vec_angle">
<span class="sig-name descname"><span class="pre">vec_angle</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">a</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">b</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L273-L278"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.vec_angle" title="Link to this definition"></a></dt>
<dd><p>Calculate angle between two vectors.</p>
</dd></dl>
<dl class="py function">
<dt class="sig sig-object py" id="pymatgen.analysis.interfaces.zsl.vec_area">
<span class="sig-name descname"><span class="pre">vec_area</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">a</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">b</span></span></em><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/interfaces/zsl.py#L281-L284"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.interfaces.zsl.vec_area" title="Link to this definition"></a></dt>
<dd><p>Area of lattice plane defined by two vectors.</p>
</dd></dl>
</section>
</section>
</div>
</div>
<footer>
<hr/>
<div role="contentinfo">
<p>© Copyright 2011, Pymatgen Development Team.</p>
</div>
Built with <a href="https://www.sphinx-doc.org/">Sphinx</a> using a
<a href="https://github.com/readthedocs/sphinx_rtd_theme">theme</a>
provided by <a href="https://readthedocs.org">Read the Docs</a>.
</footer>
</div>
</div>
</section>
</div>
<script>
jQuery(function () {
SphinxRtdTheme.Navigation.enable(true);
});
</script>
</body>
</html>
|