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 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933
|
<!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.diffraction 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.diffraction.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.diffraction package</a><ul>
<li><a class="reference internal" href="#submodules">Submodules</a></li>
<li><a class="reference internal" href="#module-pymatgen.analysis.diffraction.core">pymatgen.analysis.diffraction.core module</a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator"><code class="docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.SCALED_INTENSITY_TOL"><code class="docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator.SCALED_INTENSITY_TOL</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.TWO_THETA_TOL"><code class="docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator.TWO_THETA_TOL</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.get_pattern"><code class="docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator.get_pattern()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.get_plot"><code class="docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator.get_plot()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.plot_structures"><code class="docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator.plot_structures()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.show_plot"><code class="docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator.show_plot()</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.DiffractionPattern"><code class="docutils literal notranslate"><span class="pre">DiffractionPattern</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.DiffractionPattern.XLABEL"><code class="docutils literal notranslate"><span class="pre">DiffractionPattern.XLABEL</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.DiffractionPattern.YLABEL"><code class="docutils literal notranslate"><span class="pre">DiffractionPattern.YLABEL</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.core.get_unique_families"><code class="docutils literal notranslate"><span class="pre">get_unique_families()</span></code></a></li>
</ul>
</li>
<li><a class="reference internal" href="#module-pymatgen.analysis.diffraction.neutron">pymatgen.analysis.diffraction.neutron module</a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.neutron.NDCalculator"><code class="docutils literal notranslate"><span class="pre">NDCalculator</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.neutron.NDCalculator.get_pattern"><code class="docutils literal notranslate"><span class="pre">NDCalculator.get_pattern()</span></code></a></li>
</ul>
</li>
</ul>
</li>
<li><a class="reference internal" href="#module-pymatgen.analysis.diffraction.tem">pymatgen.analysis.diffraction.tem module</a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator"><code class="docutils literal notranslate"><span class="pre">TEMCalculator</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.bragg_angles"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.bragg_angles()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.cell_intensity"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.cell_intensity()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.cell_scattering_factors"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.cell_scattering_factors()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.electron_scattering_factors"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.electron_scattering_factors()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.generate_points"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.generate_points()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_first_point"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_first_point()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_interplanar_angle"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_interplanar_angle()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_interplanar_spacings"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_interplanar_spacings()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_pattern"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_pattern()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_2d"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_plot_2d()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_2d_concise"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_plot_2d_concise()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_coeffs"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_plot_coeffs()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_positions"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_positions()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_s2"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.get_s2()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.is_parallel"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.is_parallel()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.normalized_cell_intensity"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.normalized_cell_intensity()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.tem_dots"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.tem_dots()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.wavelength_rel"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.wavelength_rel()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.x_ray_factors"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.x_ray_factors()</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.zone_axis_filter"><code class="docutils literal notranslate"><span class="pre">TEMCalculator.zone_axis_filter()</span></code></a></li>
</ul>
</li>
</ul>
</li>
<li><a class="reference internal" href="#module-pymatgen.analysis.diffraction.xrd">pymatgen.analysis.diffraction.xrd module</a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.xrd.XRDCalculator"><code class="docutils literal notranslate"><span class="pre">XRDCalculator</span></code></a><ul>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.xrd.XRDCalculator.AVAILABLE_RADIATION"><code class="docutils literal notranslate"><span class="pre">XRDCalculator.AVAILABLE_RADIATION</span></code></a></li>
<li><a class="reference internal" href="#pymatgen.analysis.diffraction.xrd.XRDCalculator.get_pattern"><code class="docutils literal notranslate"><span class="pre">XRDCalculator.get_pattern()</span></code></a></li>
</ul>
</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.diffraction package</li>
<li class="wy-breadcrumbs-aside">
<a href="https://github.com/materialsproject/pymatgen/blob/master/docs_rst/pymatgen.analysis.diffraction.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.diffraction">
<span id="pymatgen-analysis-diffraction-package"></span><h1>pymatgen.analysis.diffraction package<a class="headerlink" href="#module-pymatgen.analysis.diffraction" title="Link to this heading"></a></h1>
<p>This package implements various diffraction analyses.</p>
<section id="submodules">
<h2>Submodules<a class="headerlink" href="#submodules" title="Link to this heading"></a></h2>
</section>
<section id="module-pymatgen.analysis.diffraction.core">
<span id="pymatgen-analysis-diffraction-core-module"></span><h2>pymatgen.analysis.diffraction.core module<a class="headerlink" href="#module-pymatgen.analysis.diffraction.core" title="Link to this heading"></a></h2>
<p>This module implements core classes for calculation of diffraction patterns.</p>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator">
<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">AbstractDiffractionPatternCalculator</span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/core.py#L42-L202"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator" title="Link to this definition"></a></dt>
<dd><p>Bases: <code class="xref py py-class docutils literal notranslate"><span class="pre">ABC</span></code></p>
<p>Abstract base class for computing the diffraction pattern of a crystal.</p>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.SCALED_INTENSITY_TOL">
<span class="sig-name descname"><span class="pre">SCALED_INTENSITY_TOL</span></span><em class="property"><span class="w"> </span><span class="p"><span class="pre">=</span></span><span class="w"> </span><span class="pre">0.001</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/diffraction/core.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.SCALED_INTENSITY_TOL" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.TWO_THETA_TOL">
<span class="sig-name descname"><span class="pre">TWO_THETA_TOL</span></span><em class="property"><span class="w"> </span><span class="p"><span class="pre">=</span></span><span class="w"> </span><span class="pre">1e-05</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/diffraction/core.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.TWO_THETA_TOL" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.get_pattern">
<em class="property"><span class="k"><span class="pre">abstractmethod</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">get_pattern</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">scaled</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">True</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">two_theta_range</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">(0,</span> <span class="pre">90)</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/diffraction/core.py#L55-L73"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.get_pattern" title="Link to this definition"></a></dt>
<dd><p>Calculates the diffraction pattern for a structure.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – Input structure</p></li>
<li><p><strong>scaled</strong> (<em>bool</em>) – Whether to return scaled intensities. The maximum
peak is set to a value of 100. Defaults to True. Use False if
you need the absolute values to combine XRD plots.</p></li>
<li><p><strong>two_theta_range</strong> (<em>[</em><em>float</em><em> of </em><em>length 2</em><em>]</em>) – Tuple for range of
two_thetas to calculate in degrees. Defaults to (0, 90). Set to
None if you want all diffracted beams within the limiting
sphere of radius 2 / wavelength.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>DiffractionPattern</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.get_plot">
<span class="sig-name descname"><span class="pre">get_plot</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">two_theta_range</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">float</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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,</span> <span class="pre">90)</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">annotate_peaks</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">'compact'</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">ax</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">plt.Axes</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">with_labels</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">True</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">fontsize</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">16</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">plt.Axes</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/core.py#L75-L153"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.get_plot" title="Link to this definition"></a></dt>
<dd><p>Get the diffraction plot as a matplotlib Axes.</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>structure</strong> – Input structure</p></li>
<li><p><strong>two_theta_range</strong> (<em>tuple</em><em>[</em><em>float</em><em>, </em><em>float</em><em>]</em>) – Range of two_thetas to calculate in degrees.
Defaults to (0, 90). Set to None if you want all diffracted beams within the limiting
sphere of radius 2 / wavelength.</p></li>
<li><p><strong>annotate_peaks</strong> (<em>str</em><em> | </em><em>None</em>) – Whether and how to annotate the peaks
with hkl indices. Default is ‘compact’, i.e. show short
version (oriented vertically), e.g. 100. If ‘full’, show
long version, e.g. (1, 0, 0). If None, do not show anything.</p></li>
<li><p><strong>ax</strong> – matplotlib Axes or None if a new figure should be
created.</p></li>
<li><p><strong>with_labels</strong> – True to add xlabels and ylabels to the plot.</p></li>
<li><p><strong>fontsize</strong> – (int) fontsize for peak labels.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>matplotlib Axes object</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p>plt.Axes</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.plot_structures">
<span class="sig-name descname"><span class="pre">plot_structures</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">structures</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">fontsize</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">6</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/../util/plotting.py#L171-L202"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.plot_structures" title="Link to this definition"></a></dt>
<dd><p>Plot diffraction patterns for multiple structures on the same figure.</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>structures</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – List of structures</p></li>
<li><p><strong>two_theta_range</strong> (<em>[</em><em>float</em><em> of </em><em>length 2</em><em>]</em>) – Tuple for range of
two_thetas to calculate in degrees. Defaults to (0, 90). Set to
None if you want all diffracted beams within the limiting
sphere of radius 2 / wavelength.</p></li>
<li><p><strong>annotate_peaks</strong> (<em>str</em><em> | </em><em>None</em>) – Whether and how to annotate the peaks
with hkl indices. Default is ‘compact’, i.e. show short
version (oriented vertically), e.g. 100. If ‘full’, show
long version, e.g. (1, 0, 0). If None, do not show anything.</p></li>
<li><p><strong>fontsize</strong> – (int) fontsize for peak labels.</p></li>
</ul>
</dd>
</dl>
<p>Keyword arguments controlling the display of the figure:</p>
<table class="docutils align-default">
<thead>
<tr class="row-odd"><th class="head"><p>kwargs</p></th>
<th class="head"><p>Meaning</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>title</p></td>
<td><p>Title of the plot (Default: None).</p></td>
</tr>
<tr class="row-odd"><td><p>show</p></td>
<td><p>True to show the figure (default: True).</p></td>
</tr>
<tr class="row-even"><td><p>savefig</p></td>
<td><p>“abc.png” or “abc.eps” to save the figure to a file.</p></td>
</tr>
<tr class="row-odd"><td><p>size_kwargs</p></td>
<td><p>Dictionary with options passed to fig.set_size_inches
e.g. size_kwargs=dict(w=3, h=4)</p></td>
</tr>
<tr class="row-even"><td><p>tight_layout</p></td>
<td><p>True to call fig.tight_layout (default: False)</p></td>
</tr>
<tr class="row-odd"><td><p>ax_grid</p></td>
<td><p>True (False) to add (remove) grid from all axes in fig.
Default: None i.e. fig is left unchanged.</p></td>
</tr>
<tr class="row-even"><td><p>ax_annotate</p></td>
<td><p>Add labels to subplots e.g. (a), (b).
Default: False</p></td>
</tr>
<tr class="row-odd"><td><p>fig_close</p></td>
<td><p>Close figure. Default: False.</p></td>
</tr>
</tbody>
</table>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.show_plot">
<span class="sig-name descname"><span class="pre">show_plot</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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="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/diffraction/core.py#L155-L169"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator.show_plot" title="Link to this definition"></a></dt>
<dd><p>Show the diffraction plot.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – Input structure</p></li>
<li><p><strong>two_theta_range</strong> (<em>[</em><em>float</em><em> of </em><em>length 2</em><em>]</em>) – Tuple for range of
two_thetas to calculate in degrees. Defaults to (0, 90). Set to
None if you want all diffracted beams within the limiting
sphere of radius 2 / wavelength.</p></li>
<li><p><strong>annotate_peaks</strong> (<em>str</em><em> | </em><em>None</em>) – Whether and how to annotate the peaks
with hkl indices. Default is ‘compact’, i.e. show short
version (oriented vertically), e.g. 100. If ‘full’, show
long version, e.g. (1, 0, 0). If None, do not show anything.</p></li>
</ul>
</dd>
</dl>
</dd></dl>
</dd></dl>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.DiffractionPattern">
<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">DiffractionPattern</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">x</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">y</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">hkls</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">d_hkls</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/diffraction/core.py#L19-L39"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.DiffractionPattern" title="Link to this definition"></a></dt>
<dd><p>Bases: <a class="reference internal" href="pymatgen.core.html#pymatgen.core.spectrum.Spectrum" title="pymatgen.core.spectrum.Spectrum"><code class="xref py py-class docutils literal notranslate"><span class="pre">Spectrum</span></code></a></p>
<p>A representation of a diffraction pattern.</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>x</strong> – Two theta angles.</p></li>
<li><p><strong>y</strong> – Intensities</p></li>
<li><p><strong>hkls</strong> – [{“hkl”: (h, k, l), “multiplicity”: mult}],
where {“hkl”: (h, k, l), “multiplicity”: mult}
is a dict of Miller
indices for all diffracted lattice facets contributing to each
intensity.</p></li>
<li><p><strong>d_hkls</strong> – List of interplanar spacings.</p></li>
</ul>
</dd>
</dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.DiffractionPattern.XLABEL">
<span class="sig-name descname"><span class="pre">XLABEL</span></span><em class="property"><span class="w"> </span><span class="p"><span class="pre">=</span></span><span class="w"> </span><span class="pre">'$2\\Theta$'</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/diffraction/core.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.DiffractionPattern.XLABEL" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.core.DiffractionPattern.YLABEL">
<span class="sig-name descname"><span class="pre">YLABEL</span></span><em class="property"><span class="w"> </span><span class="p"><span class="pre">=</span></span><span class="w"> </span><span class="pre">'Intensity'</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/diffraction/core.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.DiffractionPattern.YLABEL" 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.diffraction.core.get_unique_families">
<span class="sig-name descname"><span class="pre">get_unique_families</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">hkls</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/diffraction/core.py#L205-L237"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.core.get_unique_families" title="Link to this definition"></a></dt>
<dd><p>Get unique families of Miller indices. Families must be permutations
of each other.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><p><strong>hkls</strong> (<em>[</em><em>h</em><em>, </em><em>k</em><em>, </em><em>l</em><em>]</em>) – List of Miller indices.</p>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>multiplicity}: A dict with unique hkl and multiplicity.</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p>{hkl</p>
</dd>
</dl>
</dd></dl>
</section>
<section id="module-pymatgen.analysis.diffraction.neutron">
<span id="pymatgen-analysis-diffraction-neutron-module"></span><h2>pymatgen.analysis.diffraction.neutron module<a class="headerlink" href="#module-pymatgen.analysis.diffraction.neutron" title="Link to this heading"></a></h2>
<p>This module implements a neutron diffraction (ND) pattern calculator.</p>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.neutron.NDCalculator">
<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">NDCalculator</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">wavelength</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">1.54184</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">symprec</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</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">debye_waller_factors</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">None</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/diffraction/neutron.py#L37-L200"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.neutron.NDCalculator" title="Link to this definition"></a></dt>
<dd><p>Bases: <a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator" title="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator"><code class="xref py py-class docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator</span></code></a></p>
<p>Computes the powder neutron diffraction pattern of a crystal structure.
This code is a slight modification of XRDCalculator in
pymatgen.analysis.diffraction.xrd. See it for details of the algorithm.
Main changes by using neutron instead of X-ray are as follows:</p>
<ol class="arabic simple">
<li><p>Atomic scattering length is a constant.</p></li>
<li><p>Polarization correction term of Lorentz factor is unnecessary.</p></li>
</ol>
<p>Reference:
Marc De Graef and Michael E. McHenry, Structure of Materials 2nd ed,
Chapter13, Cambridge University Press 2003.</p>
<p>Initialize the ND calculator with a given radiation.</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>wavelength</strong> (<em>float</em>) – The wavelength of neutron in angstroms.
Defaults to 1.54, corresponds to Cu K_alpha x-ray radiation.</p></li>
<li><p><strong>symprec</strong> (<em>float</em>) – Symmetry precision for structure refinement. If
set to 0, no refinement is done. Otherwise, refinement is
performed using spglib with provided precision.</p></li>
<li><p><strong>symbol</strong> (<em>debye_waller_factors</em><em> (</em><em>{element</em>) – float}): Allows the
specification of Debye-Waller factors. Note that these
factors are temperature dependent.</p></li>
</ul>
</dd>
</dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.neutron.NDCalculator.get_pattern">
<span class="sig-name descname"><span class="pre">get_pattern</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">scaled</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">True</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">two_theta_range</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">(0,</span> <span class="pre">90)</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/diffraction/neutron.py#L69-L200"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.neutron.NDCalculator.get_pattern" title="Link to this definition"></a></dt>
<dd><p>Calculates the powder neutron diffraction pattern for a structure.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – Input structure</p></li>
<li><p><strong>scaled</strong> (<em>bool</em>) – Whether to return scaled intensities. The maximum
peak is set to a value of 100. Defaults to True. Use False if
you need the absolute values to combine ND plots.</p></li>
<li><p><strong>two_theta_range</strong> (<em>[</em><em>float</em><em> of </em><em>length 2</em><em>]</em>) – Tuple for range of
two_thetas to calculate in degrees. Defaults to (0, 90). Set to
None if you want all diffracted beams within the limiting
sphere of radius 2 / wavelength.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>ND pattern</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p><a class="reference internal" href="#pymatgen.analysis.diffraction.core.DiffractionPattern" title="pymatgen.analysis.diffraction.core.DiffractionPattern">DiffractionPattern</a></p>
</dd>
</dl>
</dd></dl>
</dd></dl>
</section>
<section id="module-pymatgen.analysis.diffraction.tem">
<span id="pymatgen-analysis-diffraction-tem-module"></span><h2>pymatgen.analysis.diffraction.tem module<a class="headerlink" href="#module-pymatgen.analysis.diffraction.tem" title="Link to this heading"></a></h2>
<p>TEM pattern calculator.</p>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator">
<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">TEMCalculator</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">symprec</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">float</span><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">voltage</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">200</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">beam_direction</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><span class="w"> </span><span class="o"><span class="pre">=</span></span><span class="w"> </span><span class="default_value"><span class="pre">(0,</span> <span class="pre">0,</span> <span class="pre">1)</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">camera_length</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">int</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">160</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">debye_waller_factors</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</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">float</span><span class="p"><span class="pre">]</span></span><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">cs</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><span class="sig-paren">)</span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L37-L687"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator" title="Link to this definition"></a></dt>
<dd><p>Bases: <a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator" title="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator"><code class="xref py py-class docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator</span></code></a></p>
<p>Compute the TEM pattern of a crystal structure for multiple Laue zones.
Code partially inspired from XRD calculation implementation. X-ray factor to electron factor</p>
<blockquote>
<div><p>conversion based on the International Table of Crystallography.</p>
</div></blockquote>
<dl class="simple">
<dt>#TODO: Could add “number of iterations”, “magnification”, “critical value of beam”,</dt><dd><p>“twin direction” for certain materials, “sample thickness”, and “excitation error s”.</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>symprec</strong> (<em>float</em>) – Symmetry precision for structure refinement. If
set to 0, no refinement is done. Otherwise, refinement is
performed using spglib with provided precision.</p></li>
<li><p><strong>voltage</strong> (<em>float</em>) – The wavelength is a function of the TEM microscope’s
voltage (in kV). Defaults to 200.</p></li>
<li><p><strong>beam_direction</strong> (<em>tuple</em>) – The direction of the electron beam fired onto the sample.
By default, set to [0,0,1], which corresponds to the normal direction
of the sample plane.</p></li>
<li><p><strong>camera_length</strong> (<em>int</em>) – The distance from the sample to the projected diffraction pattern.
By default, set to 160 cm. Units in cm.</p></li>
<li><p><strong>symbol</strong> (<em>debye_waller_factors</em><em> (</em><em>{element</em>) – float}): Allows the
specification of Debye-Waller factors. Note that these
factors are temperature dependent.</p></li>
<li><p><strong>cs</strong> (<em>float</em>) – The chromatic aberration coefficient (in mm). Defaults to 1.</p></li>
</ul>
</dd>
</dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.bragg_angles">
<span class="sig-name descname"><span class="pre">bragg_angles</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">interplanar_spacings</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</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/diffraction/tem.py#L145-L159"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.bragg_angles" title="Link to this definition"></a></dt>
<dd><p>Get the Bragg angles for every hkl point passed in (where n = 1).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><p><strong>interplanar_spacings</strong> (<em>dict</em>) – dictionary of hkl to interplanar spacing</p>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>hkl planes mapped to Bragg angles [radians]</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p>dict[tuple[int, int, int], float]</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.cell_intensity">
<span class="sig-name descname"><span class="pre">cell_intensity</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">bragg_angles</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</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/diffraction/tem.py#L261-L277"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.cell_intensity" title="Link to this definition"></a></dt>
<dd><p>Calculates cell intensity for each hkl plane. For simplicity’s sake, take I = <a href="#id1"><span class="problematic" id="id2">|</span></a>F|**2.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>bragg_angles</strong> (<em>dict</em><em> of </em><em>3-tuple to float</em>) – The Bragg angles for each hkl plane.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>dict of hkl plane to cell intensity</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.cell_scattering_factors">
<span class="sig-name descname"><span class="pre">cell_scattering_factors</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">bragg_angles</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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">dict</span><span class="p"><span class="pre">[</span></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><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></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L234-L259"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.cell_scattering_factors" title="Link to this definition"></a></dt>
<dd><p>Calculates the scattering factor for the whole cell.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>bragg_angles</strong> (<em>dict</em><em> of </em><em>3-tuple to float</em>) – The Bragg angles for each hkl plane.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>dict of hkl plane (3-tuple) to scattering factor (in angstroms).</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.electron_scattering_factors">
<span class="sig-name descname"><span class="pre">electron_scattering_factors</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">bragg_angles</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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">dict</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">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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/diffraction/tem.py#L207-L232"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.electron_scattering_factors" title="Link to this definition"></a></dt>
<dd><p>Calculates atomic scattering factors for electrons using the Mott-Bethe formula (1st order Born approximation).</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>bragg_angles</strong> (<em>dict</em><em> of </em><em>3-tuple to float</em>) – The Bragg angles for each hkl plane.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>dict from atomic symbol to another dict of hkl plane to factor (in angstroms)</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.generate_points">
<em class="property"><span class="k"><span class="pre">static</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">generate_points</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">coord_left</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">int</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">-10</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">coord_right</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">int</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">10</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/diffraction/tem.py#L90-L105"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.generate_points" title="Link to this definition"></a></dt>
<dd><p>Generate a bunch of 3D points that span a cube.</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>coord_left</strong> (<em>int</em>) – The minimum coordinate value.</p></li>
<li><p><strong>coord_right</strong> (<em>int</em>) – The maximum coordinate value.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>2d array</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p>np.array</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_first_point">
<span class="sig-name descname"><span class="pre">get_first_point</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">points</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> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</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/diffraction/tem.py#L359-L377"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_first_point" title="Link to this definition"></a></dt>
<dd><p>Get the first point to be plotted in the 2D DP, corresponding to maximum d/minimum R.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>points</strong> (<em>list</em>) – All points to be checked.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>dict of a hkl plane to max interplanar distance.</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_interplanar_angle">
<em class="property"><span class="k"><span class="pre">static</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">get_interplanar_angle</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">p1</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">p2</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><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">float</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L379-L430"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_interplanar_angle" title="Link to this definition"></a></dt>
<dd><p>Get the interplanar angle (in degrees) between the normal of two crystal planes.
Formulas from International Tables for Crystallography Volume C pp. 2-9.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>p1</strong> (<em>3-tuple</em>) – plane 1</p></li>
<li><p><strong>p2</strong> (<em>3-tuple</em>) – plane 2</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>float</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_interplanar_spacings">
<span class="sig-name descname"><span class="pre">get_interplanar_spacings</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">points</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">]</span></span><span class="w"> </span><span class="p"><span class="pre">|</span></span><span class="w"> </span><span class="pre">np.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">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</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/diffraction/tem.py#L127-L143"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_interplanar_spacings" title="Link to this definition"></a></dt>
<dd><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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – the input structure.</p></li>
<li><p><strong>points</strong> (<em>tuple</em>) – the desired hkl indices.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p><dl class="simple">
<dt>hkl planes mapped to</dt><dd><p>interplanar spacings, in angstroms (float).</p>
</dd>
</dl>
</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p>dict[tuple[int, int, int], float]</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_pattern">
<span class="sig-name descname"><span class="pre">get_pattern</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">scaled</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">bool</span><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">two_theta_range</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">float</span><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></span><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><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">pd.DataFrame</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L279-L317"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_pattern" title="Link to this definition"></a></dt>
<dd><p>Get all relevant TEM DP info in a pandas dataframe.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>scaled</strong> (<em>bool</em>) – Required value for inheritance, does nothing in TEM pattern</p></li>
<li><p><strong>two_theta_range</strong> (<em>tuple</em><em>[</em><em>float</em><em>, </em><em>float</em><em>]</em>) – Required value for inheritance, does nothing in TEM pattern</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>pd.DataFrame</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_2d">
<span class="sig-name descname"><span class="pre">get_plot_2d</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">go.Figure</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L544-L619"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_2d" title="Link to this definition"></a></dt>
<dd><p>Generate the 2D diffraction pattern of the input structure.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><p><strong>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>Figure</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_2d_concise">
<span class="sig-name descname"><span class="pre">get_plot_2d_concise</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">go.Figure</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L621-L687"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_2d_concise" title="Link to this definition"></a></dt>
<dd><p>Generate the concise 2D diffraction pattern of the input structure of a smaller size and without layout.
Does not display.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><p><strong>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>Figure</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_coeffs">
<em class="property"><span class="k"><span class="pre">static</span></span><span class="w"> </span></em><span class="sig-name descname"><span class="pre">get_plot_coeffs</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">p1</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">p2</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">p3</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><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/diffraction/tem.py#L432-L454"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_plot_coeffs" title="Link to this definition"></a></dt>
<dd><p>Calculates coefficients of the vector addition required to generate positions for each DP point
by the Moore-Penrose inverse method.</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>p1</strong> (<em>3-tuple</em>) – The first point. Fixed.</p></li>
<li><p><strong>p2</strong> (<em>3-tuple</em>) – The second point. Fixed.</p></li>
<li><p><strong>p3</strong> (<em>3-tuple</em>) – The point whose coefficients are to be calculated.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>Numpy array</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_positions">
<span class="sig-name descname"><span class="pre">get_positions</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">points</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> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">np.ndarray</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/diffraction/tem.py#L456-L511"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_positions" title="Link to this definition"></a></dt>
<dd><p>Calculates all the positions of each hkl point in the 2D diffraction pattern by vector addition.
Distance in centimeters.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>points</strong> (<em>list</em>) – All points to be checked.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>dict of hkl plane to xy-coordinates.</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.get_s2">
<span class="sig-name descname"><span class="pre">get_s2</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">bragg_angles</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</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/diffraction/tem.py#L161-L175"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.get_s2" title="Link to this definition"></a></dt>
<dd><p>Calculates the s squared parameter (= square of sin theta over lambda) for each hkl plane.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters<span class="colon">:</span></dt>
<dd class="field-odd"><p><strong>bragg_angles</strong> (<em>dict</em>) – The bragg angles for each hkl plane.</p>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p><dl class="simple">
<dt>Dict of hkl plane to s2 parameter, calculates the s squared parameter</dt><dd><p>(= square of sin theta over lambda).</p>
</dd>
</dl>
</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.is_parallel">
<span class="sig-name descname"><span class="pre">is_parallel</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">plane</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">other_plane</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><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/diffraction/tem.py#L340-L357"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.is_parallel" title="Link to this definition"></a></dt>
<dd><p>Checks if two hkl planes are parallel in reciprocal space.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>plane</strong> (<em>3-tuple</em>) – The first plane to be compared.</p></li>
<li><p><strong>other_plane</strong> (<em>3-tuple</em>) – The other plane to be compared.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>True if the planes are parallel, False otherwise.</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p>bool</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.normalized_cell_intensity">
<span class="sig-name descname"><span class="pre">normalized_cell_intensity</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">bragg_angles</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</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/diffraction/tem.py#L319-L338"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.normalized_cell_intensity" title="Link to this definition"></a></dt>
<dd><p>Normalizes the cell_intensity dict to 1, for use in plotting.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>bragg_angles</strong> (<em>dict</em><em> of </em><em>3-tuple to float</em>) – The Bragg angles for each hkl plane.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>dict of hkl plane to normalized cell intensity</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.tem_dots">
<span class="sig-name descname"><span class="pre">tem_dots</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">points</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">list</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L513-L542"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.tem_dots" title="Link to this definition"></a></dt>
<dd><p>Generate all TEM_dot as named tuples that will appear on the 2D diffraction pattern.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>points</strong> (<em>list</em>) – All points to be checked.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>list of TEM_dots</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.wavelength_rel">
<span class="sig-name descname"><span class="pre">wavelength_rel</span></span><span class="sig-paren">(</span><span class="sig-paren">)</span> <span class="sig-return"><span class="sig-return-icon">→</span> <span class="sig-return-typehint"><span class="pre">float</span></span></span><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/core/../analysis/diffraction/tem.py#L79-L88"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.wavelength_rel" title="Link to this definition"></a></dt>
<dd><dl class="simple">
<dt>Calculates the wavelength of the electron beam with relativistic kinematic effects taken</dt><dd><p>into account.</p>
</dd>
</dl>
<dl class="field-list simple">
<dt class="field-odd">Returns<span class="colon">:</span></dt>
<dd class="field-odd"><p>Relativistic Wavelength (in angstroms)</p>
</dd>
<dt class="field-even">Return type<span class="colon">:</span></dt>
<dd class="field-even"><p>float</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.x_ray_factors">
<span class="sig-name descname"><span class="pre">x_ray_factors</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">bragg_angles</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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">dict</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">dict</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">,</span></span><span class="w"> </span><span class="pre">float</span><span class="p"><span class="pre">]</span></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/diffraction/tem.py#L177-L205"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.x_ray_factors" title="Link to this definition"></a></dt>
<dd><p>Calculates x-ray factors, which are required to calculate atomic scattering factors. Method partially inspired
by the equivalent process in the xrd module.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – The input structure.</p></li>
<li><p><strong>bragg_angles</strong> (<em>dict</em>) – Dictionary of hkl plane to Bragg angle.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>dict of atomic symbol to another dict of hkl plane to x-ray factor (in angstroms).</p>
</dd>
</dl>
</dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.tem.TEMCalculator.zone_axis_filter">
<span class="sig-name descname"><span class="pre">zone_axis_filter</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">points</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">list</span><span class="p"><span class="pre">[</span></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><span class="p"><span class="pre">]</span></span><span class="w"> </span><span class="p"><span class="pre">|</span></span><span class="w"> </span><span class="pre">ndarray</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">laue_zone</span></span><span class="p"><span class="pre">:</span></span><span class="w"> </span><span class="n"><span class="pre">int</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</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">list</span><span class="p"><span class="pre">[</span></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><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/diffraction/tem.py#L107-L125"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.tem.TEMCalculator.zone_axis_filter" title="Link to this definition"></a></dt>
<dd><p>Filter out all points that exist within the specified Laue zone according to the zone axis rule.</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>points</strong> (<em>np.ndarray</em>) – The list of points to be filtered.</p></li>
<li><p><strong>laue_zone</strong> (<em>int</em>) – The desired Laue zone.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>list of 3-tuples</p>
</dd>
</dl>
</dd></dl>
</dd></dl>
</section>
<section id="module-pymatgen.analysis.diffraction.xrd">
<span id="pymatgen-analysis-diffraction-xrd-module"></span><h2>pymatgen.analysis.diffraction.xrd module<a class="headerlink" href="#module-pymatgen.analysis.diffraction.xrd" title="Link to this heading"></a></h2>
<p>This module implements an XRD pattern calculator.</p>
<dl class="py class">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.xrd.XRDCalculator">
<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">XRDCalculator</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">wavelength</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">'CuKa'</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">symprec</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</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">debye_waller_factors</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">None</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/diffraction/xrd.py#L57-L280"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.xrd.XRDCalculator" title="Link to this definition"></a></dt>
<dd><p>Bases: <a class="reference internal" href="#pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator" title="pymatgen.analysis.diffraction.core.AbstractDiffractionPatternCalculator"><code class="xref py py-class docutils literal notranslate"><span class="pre">AbstractDiffractionPatternCalculator</span></code></a></p>
<p>Computes the XRD pattern of a crystal structure.</p>
<p>This code is implemented by Shyue Ping Ong as part of UCSD’s NANO106 -
Crystallography of Materials. The formalism for this code is based on
that given in Chapters 11 and 12 of Structure of Materials by Marc De
Graef and Michael E. McHenry. This takes into account the atomic
scattering factors and the Lorentz polarization factor, but not
the Debye-Waller (temperature) factor (for which data is typically not
available). Note that the multiplicity correction is not needed since
this code simply goes through all reciprocal points within the limiting
sphere, which includes all symmetrically equivalent facets. The algorithm
is as follows</p>
<ol class="arabic">
<li><p>Calculate reciprocal lattice of structure. Find all reciprocal points
within the limiting sphere given by frac{2}{lambda}.</p></li>
<li><p>For each reciprocal point mathbf{g_{hkl}} corresponding to
lattice plane (hkl), compute the Bragg condition
sin(theta) = frac{ lambda}{2d_{hkl}}</p></li>
<li><p>Compute the structure factor as the sum of the atomic scattering
factors. The atomic scattering factors are given by</p>
<blockquote>
<div><p>f(s) = Z - 41.78214 times s^2 times sum limits_{i=1}^n a_i exp(-b_is^2)</p>
</div></blockquote>
<p>where s = frac{sin(theta)}{lambda} and a_i
and b_i are the fitted parameters for each element. The
structure factor is then given by</p>
<blockquote>
<div><p>F_{hkl} = sum limits_{j=1}^N f_j exp(2 pi i mathbf{g_{hkl}} cdot mathbf{r})</p>
</div></blockquote>
</li>
<li><p>The intensity is then given by the modulus square of the structure factor.</p>
<blockquote>
<div><p>I_{hkl} = F_{hkl}F_{hkl}^*</p>
</div></blockquote>
</li>
<li><p>Finally, the Lorentz polarization correction factor is applied. This
factor is given by:</p>
<blockquote>
<div><p>P(theta) = frac{1 + cos^2(2 theta)}{sin^2(theta) cos(theta)}</p>
</div></blockquote>
</li>
</ol>
<p>Initialize the XRD calculator with a given radiation.</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>wavelength</strong> (<em>str</em><em> | </em><em>float</em>) – The wavelength can be specified as either a
float or a string. If it is a string, it must be one of the
supported definitions in the AVAILABLE_RADIATION class
variable, which provides useful commonly used wavelengths.
If it is a float, it is interpreted as a wavelength in
angstroms. Defaults to “CuKa”, i.e, Cu K_alpha radiation.</p></li>
<li><p><strong>symprec</strong> (<em>float</em>) – Symmetry precision for structure refinement. If
set to 0, no refinement is done. Otherwise, refinement is
performed using spglib with provided precision.</p></li>
<li><p><strong>symbol</strong> (<em>debye_waller_factors</em><em> (</em><em>{element</em>) – float}): Allows the
specification of Debye-Waller factors. Note that these
factors are temperature dependent.</p></li>
</ul>
</dd>
</dl>
<dl class="py attribute">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.xrd.XRDCalculator.AVAILABLE_RADIATION">
<span class="sig-name descname"><span class="pre">AVAILABLE_RADIATION</span></span><em class="property"><span class="w"> </span><span class="p"><span class="pre">=</span></span><span class="w"> </span><span class="pre">('CuKa',</span> <span class="pre">'CuKa2',</span> <span class="pre">'CuKa1',</span> <span class="pre">'CuKb1',</span> <span class="pre">'MoKa',</span> <span class="pre">'MoKa2',</span> <span class="pre">'MoKa1',</span> <span class="pre">'MoKb1',</span> <span class="pre">'CrKa',</span> <span class="pre">'CrKa2',</span> <span class="pre">'CrKa1',</span> <span class="pre">'CrKb1',</span> <span class="pre">'FeKa',</span> <span class="pre">'FeKa2',</span> <span class="pre">'FeKa1',</span> <span class="pre">'FeKb1',</span> <span class="pre">'CoKa',</span> <span class="pre">'CoKa2',</span> <span class="pre">'CoKa1',</span> <span class="pre">'CoKb1',</span> <span class="pre">'AgKa',</span> <span class="pre">'AgKa2',</span> <span class="pre">'AgKa1',</span> <span class="pre">'AgKb1')</span></em><a class="reference external" href="https://github.com/materialsproject/pymatgen/blob/v2025.6.14/src/pymatgen/analysis/diffraction/xrd.py"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.xrd.XRDCalculator.AVAILABLE_RADIATION" title="Link to this definition"></a></dt>
<dd></dd></dl>
<dl class="py method">
<dt class="sig sig-object py" id="pymatgen.analysis.diffraction.xrd.XRDCalculator.get_pattern">
<span class="sig-name descname"><span class="pre">get_pattern</span></span><span class="sig-paren">(</span><em class="sig-param"><span class="n"><span class="pre">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">scaled</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">True</span></span></em>, <em class="sig-param"><span class="n"><span class="pre">two_theta_range</span></span><span class="o"><span class="pre">=</span></span><span class="default_value"><span class="pre">(0,</span> <span class="pre">90)</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/diffraction/xrd.py#L131-L280"><span class="viewcode-link"><span class="pre">[source]</span></span></a><a class="headerlink" href="#pymatgen.analysis.diffraction.xrd.XRDCalculator.get_pattern" title="Link to this definition"></a></dt>
<dd><p>Calculates the diffraction pattern for a structure.</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>structure</strong> (<a class="reference internal" href="pymatgen.core.html#pymatgen.core.structure.Structure" title="pymatgen.core.structure.Structure"><em>Structure</em></a>) – Input structure</p></li>
<li><p><strong>scaled</strong> (<em>bool</em>) – Whether to return scaled intensities. The maximum
peak is set to a value of 100. Defaults to True. Use False if
you need the absolute values to combine XRD plots.</p></li>
<li><p><strong>two_theta_range</strong> (<em>[</em><em>float</em><em> of </em><em>length 2</em><em>]</em>) – Tuple for range of
two_thetas to calculate in degrees. Defaults to (0, 90). Set to
None if you want all diffracted beams within the limiting
sphere of radius 2 / wavelength.</p></li>
</ul>
</dd>
<dt class="field-even">Returns<span class="colon">:</span></dt>
<dd class="field-even"><p>XRD pattern</p>
</dd>
<dt class="field-odd">Return type<span class="colon">:</span></dt>
<dd class="field-odd"><p><a class="reference internal" href="#pymatgen.analysis.diffraction.core.DiffractionPattern" title="pymatgen.analysis.diffraction.core.DiffractionPattern">DiffractionPattern</a></p>
</dd>
</dl>
</dd></dl>
</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>
|