Skip to content

plasma_geometry

Plasma geometry models and shape definitions for PROCESS.

PlasmaShapeModelType

Bases: IntEnum

Enum for plasma shape model types.

Source code in process/models/physics/plasma_geometry.py
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
class PlasmaShapeModelType(IntEnum):
    """Enum for plasma shape model types."""

    PROCESS_ORIGINAL = (0, "PROCESS Original Double Arc")
    SAUTER = (1, "Sauter")

    def __new__(cls, value: int, full_name: str):
        """Create a new PlasmaShapeModelType instance."""
        obj = int.__new__(cls, value)
        obj._value_ = value
        obj._full_name_ = full_name
        return obj

    @DynamicClassAttribute
    def full_name(self):
        """Get the full name of the plasma shape model."""
        return self._full_name_

PROCESS_ORIGINAL = (0, 'PROCESS Original Double Arc') class-attribute instance-attribute

SAUTER = (1, 'Sauter') class-attribute instance-attribute

full_name()

Get the full name of the plasma shape model.

Source code in process/models/physics/plasma_geometry.py
33
34
35
36
@DynamicClassAttribute
def full_name(self):
    """Get the full name of the plasma shape model."""
    return self._full_name_

PlasmaGeometryModels

Bases: IntEnum

Enum for plasma geometry model types.

Source code in process/models/physics/plasma_geometry.py
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
class PlasmaGeometryModels(IntEnum):
    """Enum for plasma geometry model types."""

    USER_INPUT = (0, "User Input")
    IPDG89 = (1, "IPDG89")
    STAR_CODE = (2, "STAR Code")
    ZOHM_ITER = (3, "Zohm ITER Scaling")
    MAST_DATA = (4, "Fit to MAST data")
    FIESTA_RUNS = (5, "Fiesta Runs")
    CREATE_DATA_EU_DEMO = (6, "CREATE Data EU DEMO")
    MENARD_2016 = (7, "Menard 2016 ST Scaling")
    UNKNOWN = (8, "Unknown")
    MENARD_1997 = (9, "Menard 1997 ST Scaling")

    def __new__(cls, value: int, description: str):
        """Create a new PlasmaGeometryModels instance."""
        obj = int.__new__(cls, value)
        obj._value_ = value
        obj._description_ = description
        return obj

    @DynamicClassAttribute
    def description(self):
        """Get the description of the plasma geometry model."""
        return self._description_

USER_INPUT = (0, 'User Input') class-attribute instance-attribute

IPDG89 = (1, 'IPDG89') class-attribute instance-attribute

STAR_CODE = (2, 'STAR Code') class-attribute instance-attribute

ZOHM_ITER = (3, 'Zohm ITER Scaling') class-attribute instance-attribute

MAST_DATA = (4, 'Fit to MAST data') class-attribute instance-attribute

FIESTA_RUNS = (5, 'Fiesta Runs') class-attribute instance-attribute

CREATE_DATA_EU_DEMO = (6, 'CREATE Data EU DEMO') class-attribute instance-attribute

MENARD_2016 = (7, 'Menard 2016 ST Scaling') class-attribute instance-attribute

UNKNOWN = (8, 'Unknown') class-attribute instance-attribute

MENARD_1997 = (9, 'Menard 1997 ST Scaling') class-attribute instance-attribute

description()

Get the description of the plasma geometry model.

Source code in process/models/physics/plasma_geometry.py
60
61
62
63
@DynamicClassAttribute
def description(self):
    """Get the description of the plasma geometry model."""
    return self._description_

PlasmaGeometryModelType

Bases: IntEnum

Enum for i_plasma_geometry plasma geometry model types.

Source code in process/models/physics/plasma_geometry.py
 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
class PlasmaGeometryModelType(IntEnum):
    """Enum for i_plasma_geometry plasma geometry model types."""

    IPDG89_X_POINT = (
        0,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.IPDG89,
    )
    STAR_FIESTA = (
        1,
        PlasmaGeometryModels.STAR_CODE,
        PlasmaGeometryModels.STAR_CODE,
        PlasmaGeometryModels.FIESTA_RUNS,
        PlasmaGeometryModels.FIESTA_RUNS,
    )
    ZOHM_ITER_X_POINT = (
        2,
        PlasmaGeometryModels.ZOHM_ITER,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.IPDG89,
    )
    ZOHM_ITER_95 = (
        3,
        PlasmaGeometryModels.ZOHM_ITER,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.USER_INPUT,
    )
    IPDG89_95 = (
        4,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.USER_INPUT,
    )
    MAST_DATA_95 = (
        5,
        PlasmaGeometryModels.MAST_DATA,
        PlasmaGeometryModels.MAST_DATA,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.USER_INPUT,
    )
    MAST_DATA_X_POINT = (
        6,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.MAST_DATA,
        PlasmaGeometryModels.MAST_DATA,
    )
    FIESTA_RUNS_95 = (
        7,
        PlasmaGeometryModels.FIESTA_RUNS,
        PlasmaGeometryModels.FIESTA_RUNS,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.USER_INPUT,
    )
    FIESTA_RUNS_X_POINT = (
        8,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.FIESTA_RUNS,
        PlasmaGeometryModels.FIESTA_RUNS,
    )
    INDUCTANCE_SCALING_X_POINT = (
        9,
        PlasmaGeometryModels.UNKNOWN,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.IPDG89,
    )
    CREATE_DATA_EU_DEMO_X_POINT = (
        10,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.CREATE_DATA_EU_DEMO,
        PlasmaGeometryModels.IPDG89,
    )
    MENARD_2016_X_POINT = (
        11,
        PlasmaGeometryModels.MENARD_2016,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.IPDG89,
    )
    MENARD_1997_X_POINT = (
        12,
        PlasmaGeometryModels.MENARD_1997,
        PlasmaGeometryModels.USER_INPUT,
        PlasmaGeometryModels.IPDG89,
        PlasmaGeometryModels.IPDG89,
    )

    def __new__(
        cls,
        value: int,
        kappa_model: PlasmaGeometryModels,
        triang_model: PlasmaGeometryModels,
        kappa95_model: PlasmaGeometryModels,
        triang95_model: PlasmaGeometryModels,
    ):
        """Create a new PlasmaGeometryModelType instance."""
        obj = int.__new__(cls, value)
        obj._value_ = value
        obj._kappa_model_ = kappa_model
        obj._triang_model_ = triang_model
        obj._kappa95_model_ = kappa95_model
        obj._triang95_model_ = triang95_model
        return obj

    @DynamicClassAttribute
    def kappa_model(self):
        """Get the kappa model."""
        return self._kappa_model_

    @DynamicClassAttribute
    def triang_model(self):
        """Get the triangularity model."""
        return self._triang_model_

    @DynamicClassAttribute
    def kappa95_model(self):
        """Get the kappa95 model."""
        return self._kappa95_model_

    @DynamicClassAttribute
    def triang95_model(self):
        """Get the triangularity95 model."""
        return self._triang95_model_

IPDG89_X_POINT = (0, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.IPDG89) class-attribute instance-attribute

STAR_FIESTA = (1, PlasmaGeometryModels.STAR_CODE, PlasmaGeometryModels.STAR_CODE, PlasmaGeometryModels.FIESTA_RUNS, PlasmaGeometryModels.FIESTA_RUNS) class-attribute instance-attribute

ZOHM_ITER_X_POINT = (2, PlasmaGeometryModels.ZOHM_ITER, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.IPDG89) class-attribute instance-attribute

ZOHM_ITER_95 = (3, PlasmaGeometryModels.ZOHM_ITER, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.USER_INPUT) class-attribute instance-attribute

IPDG89_95 = (4, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.USER_INPUT) class-attribute instance-attribute

MAST_DATA_95 = (5, PlasmaGeometryModels.MAST_DATA, PlasmaGeometryModels.MAST_DATA, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.USER_INPUT) class-attribute instance-attribute

MAST_DATA_X_POINT = (6, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.MAST_DATA, PlasmaGeometryModels.MAST_DATA) class-attribute instance-attribute

FIESTA_RUNS_95 = (7, PlasmaGeometryModels.FIESTA_RUNS, PlasmaGeometryModels.FIESTA_RUNS, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.USER_INPUT) class-attribute instance-attribute

FIESTA_RUNS_X_POINT = (8, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.FIESTA_RUNS, PlasmaGeometryModels.FIESTA_RUNS) class-attribute instance-attribute

INDUCTANCE_SCALING_X_POINT = (9, PlasmaGeometryModels.UNKNOWN, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.IPDG89) class-attribute instance-attribute

CREATE_DATA_EU_DEMO_X_POINT = (10, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.CREATE_DATA_EU_DEMO, PlasmaGeometryModels.IPDG89) class-attribute instance-attribute

MENARD_2016_X_POINT = (11, PlasmaGeometryModels.MENARD_2016, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.IPDG89) class-attribute instance-attribute

MENARD_1997_X_POINT = (12, PlasmaGeometryModels.MENARD_1997, PlasmaGeometryModels.USER_INPUT, PlasmaGeometryModels.IPDG89, PlasmaGeometryModels.IPDG89) class-attribute instance-attribute

kappa_model()

Get the kappa model.

Source code in process/models/physics/plasma_geometry.py
178
179
180
181
@DynamicClassAttribute
def kappa_model(self):
    """Get the kappa model."""
    return self._kappa_model_

triang_model()

Get the triangularity model.

Source code in process/models/physics/plasma_geometry.py
183
184
185
186
@DynamicClassAttribute
def triang_model(self):
    """Get the triangularity model."""
    return self._triang_model_

kappa95_model()

Get the kappa95 model.

Source code in process/models/physics/plasma_geometry.py
188
189
190
191
@DynamicClassAttribute
def kappa95_model(self):
    """Get the kappa95 model."""
    return self._kappa95_model_

triang95_model()

Get the triangularity95 model.

Source code in process/models/physics/plasma_geometry.py
193
194
195
196
@DynamicClassAttribute
def triang95_model(self):
    """Get the triangularity95 model."""
    return self._triang95_model_

PlasmaGeom

Bases: Model

Class for calculating plasma geometry parameters.

Source code in process/models/physics/plasma_geometry.py
 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
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
class PlasmaGeom(Model):
    """Class for calculating plasma geometry parameters."""

    def __init__(self):
        self.outfile = constants.NOUT

    def run(self):
        """Plasma geometry parameters

        This method calculates the plasma geometry parameters based on various shaping
        terms and input values. It updates the `physics_variables` with calculated
        values for kappa, triangularity, surface area, volume, etc.

        References
        ----------
            - J D Galambos, STAR Code : Spherical Tokamak Analysis and Reactor Code,
              unpublished internal Oak Ridge document
            - H. Zohm et al, On the Physics Guidelines for a Tokamak DEMO,
              FTP/3-3, Proc. IAEA Fusion Energy Conference, October 2012, San Diego
        """
        xsi = 0.0e0
        xso = 0.0e0
        thetai = 0.0e0
        thetao = 0.0e0
        xi = 0.0e0
        xo = 0.0e0

        # Define plasma minor radius from major radius and aspect ratio
        physics_variables.rminor = physics_variables.rmajor / physics_variables.aspect

        # Define the inverse aspect ratio
        physics_variables.eps = 1.0e0 / physics_variables.aspect

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry == PlasmaGeometryModelType.IPDG89_X_POINT
        ):  # Use input kappa, physics_variables.triang values
            #  Rough estimate of 95% values
            #  ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
            #  (close to previous estimate of (physics_variables.kappa - 0.04) / 1.1
            #  over a large physics_variables.kappa range)

            physics_variables.kappa95 = physics_variables.kappa / 1.12e0
            physics_variables.triang95 = physics_variables.triang / 1.50e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry == PlasmaGeometryModelType.STAR_FIESTA
        ):  # ST scaling with physics_variables.aspect ratio [STAR Code]
            physics_variables.q95_min = 3.0e0 * (
                1.0e0 + 2.6e0 * physics_variables.eps**2.8e0
            )

            physics_variables.kappa = 2.05e0 * (
                1.0e0 + 0.44e0 * physics_variables.eps**2.1e0
            )
            physics_variables.triang = 0.53e0 * (
                1.0e0 + 0.77e0 * physics_variables.eps**3
            )

            # SIM 10/09/2020: Switched to FIESTA ST scaling  from IPDG89
            physics_variables.kappa95 = (
                physics_variables.kappa - 0.39467e0
            ) / 0.90698e0  # Fit to FIESTA (Issue #1086)
            physics_variables.triang95 = (
                physics_variables.triang - 0.048306e0
            ) / 1.3799e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry
            == PlasmaGeometryModelType.ZOHM_ITER_X_POINT
        ):  # Zohm et al. ITER scaling for elongation, input physics_variables.triang
            physics_variables.kappa = physics_variables.fkzohm * min(
                2.0e0, 1.5e0 + 0.5e0 / (physics_variables.aspect - 1.0e0)
            )

            # ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
            physics_variables.kappa95 = physics_variables.kappa / 1.12e0
            physics_variables.triang95 = physics_variables.triang / 1.50e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry == PlasmaGeometryModelType.ZOHM_ITER_95
        ):  # Zohm et al. ITER scaling for elongation, input physics_variables.triang95
            physics_variables.kappa = physics_variables.fkzohm * min(
                2.0e0, 1.5e0 + 0.5e0 / (physics_variables.aspect - 1.0e0)
            )

            # ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
            physics_variables.triang = 1.5e0 * physics_variables.triang95

            physics_variables.kappa95 = physics_variables.kappa / 1.12e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry == PlasmaGeometryModelType.IPDG89_95
        ):  # Use input kappa95, physics_variables.triang95 values
            # ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
            physics_variables.kappa = 1.12e0 * physics_variables.kappa95
            physics_variables.triang = 1.5e0 * physics_variables.triang95

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry == PlasmaGeometryModelType.MAST_DATA_95
        ):  # Use input kappa95, physics_variables.triang95 values
            # Fit to MAST data (Issue #1086)
            physics_variables.kappa = 0.91300e0 * physics_variables.kappa95 + 0.38654e0
            physics_variables.triang = 0.77394e0 * physics_variables.triang95 + 0.18515e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry
            == PlasmaGeometryModelType.MAST_DATA_X_POINT
        ):  # Use input kappa, physics_variables.triang values
            # Fit to MAST data (Issue #1086)
            physics_variables.kappa95 = (physics_variables.kappa - 0.38654e0) / 0.91300e0
            physics_variables.triang95 = (
                physics_variables.triang - 0.18515e0
            ) / 0.77394e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry == PlasmaGeometryModelType.FIESTA_RUNS_95
        ):  # Use input kappa95, physics_variables.triang95 values
            # Fit to FIESTA (Issue #1086)
            physics_variables.kappa = 0.90698e0 * physics_variables.kappa95 + 0.39467e0
            physics_variables.triang = 1.3799e0 * physics_variables.triang95 + 0.048306e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry
            == PlasmaGeometryModelType.FIESTA_RUNS_X_POINT
        ):  # Use input kappa, physics_variables.triang values
            # Fit to FIESTA (Issue #1086)
            physics_variables.kappa95 = (physics_variables.kappa - 0.39467e0) / 0.90698e0
            physics_variables.triang95 = (
                physics_variables.triang - 0.048306e0
            ) / 1.3799e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry
            == PlasmaGeometryModelType.INDUCTANCE_SCALING_X_POINT
        ):  # Use input triang, physics_variables.ind_plasma_internal_norm values
            # physics_variables.kappa found from physics_variables.aspect ratio and
            # plasma internal inductance li(3)
            physics_variables.kappa = (
                1.09e0 + 0.26e0 / physics_variables.ind_plasma_internal_norm
            ) * (1.5e0 / physics_variables.aspect) ** 0.4e0

            physics_variables.kappa95 = physics_variables.kappa / 1.12e0
            physics_variables.triang95 = physics_variables.triang / 1.50e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry
            == PlasmaGeometryModelType.CREATE_DATA_EU_DEMO_X_POINT
        ):
            # physics_variables.kappa95 found from physics_variables.aspect ratio and
            # stability margin Based on fit to CREATE data. ref Issue #1399
            # valid for EU-DEMO like machine - physics_variables.aspect ratio 2.6 - 3.6
            # Model updated see Issue #1648
            a = 3.68436807e0
            b = -0.27706527e0
            c = 0.87040251e0
            d = -18.83740952e0
            e = -0.27267618e0
            f = 20.5141261e0

            physics_variables.kappa95 = (
                -d
                - c * physics_variables.aspect
                - np.sqrt(
                    (c**2.0e0 - 4.0e0 * a * b) * physics_variables.aspect**2.0e0
                    + (2.0e0 * d * c - 4.0e0 * a * e) * physics_variables.aspect
                    + d**2.0e0
                    - 4.0e0 * a * f
                    + 4.0e0 * a * physics_variables.m_s_limit
                )
            ) / (2.0e0 * a)

            if physics_variables.kappa95 > 1.77:
                ratio = 1.77 / physics_variables.kappa95
                corner_fudge = 0.3 * (physics_variables.kappa95 - 1.77) / ratio
                physics_variables.kappa95 = (
                    physics_variables.kappa95 ** (ratio) + corner_fudge
                )

            physics_variables.kappa = 1.12e0 * physics_variables.kappa95
            physics_variables.triang95 = physics_variables.triang / 1.50e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry
            == PlasmaGeometryModelType.MENARD_2016_X_POINT
        ):
            # See Issue #1439
            # physics_variables.triang is an input
            # physics_variables.kappa found from physics_variables.aspect ratio scaling
            # on p32 of Menard: Menard, et al. "Fusion Nuclear Science Facilities
            # and Pilot Plants Based on the Spherical Tokamak." Nucl. Fusion, 2016, 44.

            physics_variables.kappa = 0.95e0 * (
                1.9e0 + 1.9e0 / physics_variables.aspect**1.4e0
            )
            physics_variables.kappa95 = physics_variables.kappa / 1.12e0
            physics_variables.triang95 = physics_variables.triang / 1.50e0

        # ======================================================================

        if (
            physics_variables.i_plasma_geometry
            == PlasmaGeometryModelType.MENARD_1997_X_POINT
        ):
            # physics_variables.triang is an input
            # physics_variables.kappa found from physics_variables.aspect ratio scaling from
            # J.E. Menard et al 1997 Nucl. Fusion 37 595 and assume max controllable kappa
            # and assume lᵢ(3) is held constant

            physics_variables.kappa = (
                2.93e0 * (1.8e0 / physics_variables.aspect) ** 0.4e0
            )

            physics_variables.kappa95 = physics_variables.kappa / 1.12e0
            physics_variables.triang95 = physics_variables.triang / 1.50e0

        # ======================================================================

        #  Scrape-off layer thicknesses
        if physics_variables.i_plasma_wall_gap == 0:
            self.data.build.dr_fw_plasma_gap_outboard = 0.1e0 * physics_variables.rminor
            self.data.build.dr_fw_plasma_gap_inboard = 0.1e0 * physics_variables.rminor

        # ======================================================================

        # Find parameters of arcs describing plasma surfaces
        xi, thetai, xo, thetao = self.plasma_angles_arcs(
            physics_variables.rminor,
            physics_variables.kappa,
            physics_variables.triang,
        )

        #  Surface area - inboard and outboard.  These are not given by Sauter but
        #  the outboard area is required by DCLL and divertor
        xsi, xso = self.plasma_surface_area(
            physics_variables.rmajor,
            physics_variables.rminor,
            xi,
            thetai,
            xo,
            thetao,
        )
        physics_variables.a_plasma_surface_outboard = xso

        # ======================================================================

        # i_plasma_current = 8 specifies use of the Sauter geometry as well as plasma
        # current.
        if (
            physics_variables.i_plasma_current == 8
            or physics_variables.i_plasma_shape == PlasmaShapeModelType.SAUTER
        ):
            (
                physics_variables.len_plasma_poloidal,
                physics_variables.a_plasma_surface,
                physics_variables.a_plasma_poloidal,
                physics_variables.vol_plasma,
            ) = self.sauter_geometry(
                physics_variables.rminor,
                physics_variables.rmajor,
                physics_variables.kappa,
                physics_variables.triang,
                physics_variables.plasma_square,
            )

        else:
            #  Poloidal perimeter
            physics_variables.len_plasma_poloidal = self.plasma_poloidal_perimeter(
                xi, thetai, xo, thetao
            )

            #  Volume
            physics_variables.vol_plasma = (
                physics_variables.f_vol_plasma
                * self.plasma_volume(
                    physics_variables.rmajor,
                    physics_variables.rminor,
                    xi,
                    thetai,
                    xo,
                    thetao,
                )
            )

            #  Cross-sectional area
            physics_variables.a_plasma_poloidal = self.plasma_cross_section(
                xi, thetai, xo, thetao
            )

            #  Surface area - sum of inboard and outboard.
            physics_variables.a_plasma_surface = xsi + xso

        # ======================================================================

    def output(self):
        """Output plasma geometry parameters to file.

        Raises
        ------
        ProcessValueError
            If `n_divertors` has an illegal value.
        """
        po.oheadr(self.outfile, "Plasma Geometry")

        if self.data.stellarator.istell == 0:
            if self.data.divertor.n_divertors == 0:
                po.ocmmnt(self.outfile, "Plasma configuration = limiter")
            elif self.data.divertor.n_divertors == 1:
                po.ocmmnt(self.outfile, "Plasma configuration = single null divertor")
            elif self.data.divertor.n_divertors == 2:
                po.ocmmnt(self.outfile, "Plasma configuration = double null divertor")
            else:
                raise ProcessValueError(
                    "Illegal value of n_divertors",
                    n_divertors=self.data.divertor.n_divertors,
                )
        else:
            po.ocmmnt(self.outfile, "Plasma configuration = stellarator")

        if self.data.stellarator.istell == 0:
            if physics_variables.itart == 0:
                physics_variables.itart_r = physics_variables.itart
                po.ovarin(
                    self.outfile,
                    "Tokamak aspect ratio = Conventional, itart = 0",
                    "(itart)",
                    physics_variables.itart_r,
                )
            elif physics_variables.itart == 1:
                physics_variables.itart_r = physics_variables.itart
                po.ovarin(
                    self.outfile,
                    "Tokamak aspect ratio = Spherical, itart = 1",
                    "(itart)",
                    physics_variables.itart_r,
                )

        po.ovarin(
            self.outfile,
            "Plasma shaping model",
            "(i_plasma_shape)",
            physics_variables.i_plasma_shape,
        )

        po.osubhd(
            self.outfile,
            f"{PlasmaShapeModelType(physics_variables.i_plasma_shape).full_name} "
            "plasma shape model is used :",
        )

        po.ovarrf(
            self.outfile, "Major radius (R₀) (m)", "(rmajor)", physics_variables.rmajor
        )
        po.ovarrf(
            self.outfile,
            "Minor radius (a) (m)",
            "(rminor)",
            physics_variables.rminor,
            "OP ",
        )
        po.ovarrf(self.outfile, "Aspect ratio (A)", "(aspect)", physics_variables.aspect)
        po.ovarrf(
            self.outfile,
            "Plasma squareness (ζ)",
            "(plasma_square)",
            physics_variables.plasma_square,
            "IP",
        )
        po.oblnkl(self.outfile)

        po.ovarin(
            self.outfile,
            "Plasma geometry model",
            "(i_plasma_geometry)",
            physics_variables.i_plasma_geometry,
        )
        po.oblnkl(self.outfile)
        po.ovarre(
            self.outfile,
            "ITER Physics Basis definition of elongation (κₐ)",
            "(kappa_ipb)",
            physics_variables.kappa_ipb,
            "OP ",
        )
        po.oblnkl(self.outfile)
        geom_type = PlasmaGeometryModelType(physics_variables.i_plasma_geometry)

        if self.data.stellarator.istell == 0:
            po.ocmmnt(
                self.outfile,
                f"X-Point Elongation set from: {geom_type.kappa_model.description}",
            )

            po.ovarrf(
                self.outfile,
                "Elongation, X-point (κₐ)",
                "(kappa)",
                physics_variables.kappa,
                "IP"
                if geom_type.kappa_model == PlasmaGeometryModels.USER_INPUT
                else "OP",
            )
            if geom_type.kappa_model == PlasmaGeometryModels.ZOHM_ITER:
                po.ovarrf(
                    self.outfile,
                    "Zohm scaling adjustment factor",
                    "(fkzohm)",
                    physics_variables.fkzohm,
                )

            po.oblnkl(self.outfile)

            po.ocmmnt(
                self.outfile,
                f"95% Elongation set from: {geom_type.kappa95_model.description}",
            )
            po.ovarrf(
                self.outfile,
                "Elongation, 95% surface (κ₉₅)",
                "(kappa95)",
                physics_variables.kappa95,
                "IP"
                if geom_type.kappa95_model == PlasmaGeometryModels.USER_INPUT
                else "OP",
            )

            po.oblnkl(self.outfile)

            po.ocmmnt(
                self.outfile,
                f"X-Point Triangularity set from: {geom_type.triang_model.description}",
            )

            po.ovarrf(
                self.outfile,
                "Triangularity, X-point (δ)",
                "(triang)",
                physics_variables.triang,
                "IP "
                if geom_type.triang_model == PlasmaGeometryModels.USER_INPUT
                else "OP",
            )
            po.oblnkl(self.outfile)
            po.ocmmnt(
                self.outfile,
                f"95% Triangularity set from: {geom_type.triang95_model.description}",
            )

            po.ovarrf(
                self.outfile,
                "Triangularity, 95% surface (δ₉₅)",
                "(triang95)",
                physics_variables.triang95,
                "IP "
                if geom_type.triang95_model == PlasmaGeometryModels.USER_INPUT
                else "OP",
            )

            po.oblnkl(self.outfile)

            po.ovarrf(
                self.outfile,
                "Plasma poloidal perimeter (m)",
                "(len_plasma_poloidal)",
                physics_variables.len_plasma_poloidal,
                "OP ",
            )

            po.ovarre(
                self.outfile,
                "Plasma cross-sectional area (m²)",
                "(a_plasma_poloidal)",
                physics_variables.a_plasma_poloidal,
                "OP ",
            )
            po.ovarre(
                self.outfile,
                "Plasma surface area (m²)",
                "(a_plasma_surface)",
                physics_variables.a_plasma_surface,
                "OP ",
            )
            po.ovarre(
                self.outfile,
                "Plasma volume (m³)",
                "(vol_plasma)",
                physics_variables.vol_plasma,
                "OP ",
            )

        po.oblnkl(self.outfile)

    @staticmethod
    def plasma_angles_arcs(
        a: float, kappa: float, triang: float
    ) -> tuple[float, float, float, float]:
        """Routine to find parameters used for calculating geometrical
        properties for double-null plasmas.

        This function finds plasma geometrical parameters, using the
        revolution of two intersecting arcs around the device centreline.
        This calculation is appropriate for plasmas with a separatrix.


        Parameters
        ----------
        a : float
            Plasma minor radius (m)
        kappa : float
            Plasma separatrix elongation
        triang : float
            Plasma separatrix triangularity

        Returns
        -------
        tuple
            A tuple containing:
            - xi (float): Radius of arc describing inboard surface (m)
            - thetai (float): Half-angle of arc describing inboard surface
            - xo (float): Radius of arc describing outboard surface (m)
            - thetao (float): Half-angle of arc describing outboard surface


        References
        ----------
        - F/MI/PJK/LOGBOOK14, p.42
        - F/PL/PJK/PROCESS/CODE/047
        """
        t = 1.0e0 - triang
        denomi = (kappa**2 - t**2) / (2.0e0 * t)
        thetai = np.arctan(kappa / denomi)
        xi = a * (denomi + 1.0e0 - triang)

        #  Find radius and half-angle of outboard arc

        n = 1.0e0 + triang
        denomo = (kappa**2 - n**2) / (2.0e0 * n)
        thetao = np.arctan(kappa / denomo)
        xo = a * (denomo + 1.0e0 + triang)

        return xi, thetai, xo, thetao

    @staticmethod
    def plasma_poloidal_perimeter(
        xi: float, thetai: float, xo: float, thetao: float
    ) -> float:
        """Calculate the poloidal perimeter of the plasma.

        Parameters
        ----------
        xi : float
            Radius of arc describing inboard surface (m)
        thetai : float
            Half-angle of arc describing inboard surface
        xo : float
            Radius of arc describing outboard surface (m)
        thetao : float
            Half-angle of arc describing outboard surface

        Returns
        -------
        float
            Poloidal perimeter (m)
        """
        return 2.0e0 * (xo * thetao + xi * thetai)

    @staticmethod
    def plasma_surface_area(
        rmajor: float, rminor: float, xi: float, thetai: float, xo: float, thetao: float
    ) -> tuple[float, float]:
        """Plasma surface area calculation

        This function finds the plasma surface area, using the
        revolution of two intersecting arcs around the device centreline.
        This calculation is appropriate for plasmas with a separatrix.

        Parameters
        ----------
        rmajor : float
            Plasma major radius (m)
        rminor : float
            Plasma minor radius (m)
        xi : float
            Radius of arc describing inboard surface (m)
        thetai : float
            Half-angle of arc describing inboard surface
        xo : float
            Radius of arc describing outboard surface (m)
        thetao : float
            Half-angle of arc describing outboard surface

        Returns
        -------
        tuple
            A tuple containing:
            - xsi (float): Inboard surface area (m^2)
            - xso (float): Outboard surface area (m^2)


        References
        ----------
        - F/MI/PJK/LOGBOOK14, p.43
        """
        fourpi = 4.0e0 * np.pi

        rc = rmajor - rminor + xi
        xsi = fourpi * xi * (rc * thetai - xi * np.sin(thetai))

        rc = rmajor + rminor - xo
        xso = fourpi * xo * (rc * thetao + xo * np.sin(thetao))

        return xsi, xso

    @staticmethod
    def plasma_volume(
        rmajor: float, rminor: float, xi: float, thetai: float, xo: float, thetao: float
    ) -> float:
        """Plasma volume calculation
        This function finds the plasma volume, using the
        revolution of two intersecting arcs around the device centreline.
        This calculation is appropriate for plasmas with a separatrix.

        Parameters
        ----------
        rmajor : float
            Plasma major radius (m)
        rminor : float
            Plasma minor radius (m)
        xi : float
            Radius of arc describing inboard surface (m)
        thetai : float
            Half-angle of arc describing inboard surface
        xo : float
            Radius of arc describing outboard surface (m)
        thetao : float
            Half-angle of arc describing outboard surface

        Returns
        -------
        float
            Plasma volume (m^3)


        References
        ----------
        - F/MI/PJK/LOGBOOK14, p.43
        """
        third = 1.0e0 / 3.0e0

        rc = rmajor - rminor + xi
        vin = (
            2.0
            * np.pi
            * xi
            * (
                rc**2 * np.sin(thetai)
                - rc * xi * thetai
                - 0.5e0 * rc * xi * np.sin(2.0e0 * thetai)
                + xi * xi * np.sin(thetai)
                - third * xi * xi * (np.sin(thetai)) ** 3
            )
        )

        rc = rmajor + rminor - xo
        vout = (
            2.0
            * np.pi
            * xo
            * (
                rc**2 * np.sin(thetao)
                + rc * xo * thetao
                + 0.5e0 * rc * xo * np.sin(2.0e0 * thetao)
                + xo * xo * np.sin(thetao)
                - third * xo * xo * (np.sin(thetao)) ** 3
            )
        )

        return vout - vin

    @staticmethod
    def plasma_cross_section(
        xi: float, thetai: float, xo: float, thetao: float
    ) -> float:
        """Plasma cross-sectional area calculation

        This function finds the plasma cross-sectional area, using the
        revolution of two intersecting arcs around the device centreline.
        This calculation is appropriate for plasmas with a separatrix.

        Parameters
        ----------
        xi : float
            Radius of arc describing inboard surface (m)
        thetai : float
            Half-angle of arc describing inboard surface
        xo : float
            Radius of arc describing outboard surface (m)
        thetao : float
            Half-angle of arc describing outboard surface

        Returns
        -------
        float
            Plasma cross-sectional area (m^2)


        References
        ----------
        - F/MI/PJK/LOGBOOK14, p.41
        """
        return xo**2 * (thetao - np.cos(thetao) * np.sin(thetao)) + xi**2 * (
            thetai - np.cos(thetai) * np.sin(thetai)
        )

    @staticmethod
    def sauter_geometry(
        a: float, r0: float, kappa: float, triang: float, square: float
    ) -> tuple[float, float, float, float, float]:
        """Calculate the plasma geometry parameters using the Sauter geometry model.

        Parameters
        ----------
        a : float
            Plasma minor radius (m)
        r0 : float
            Plasma major radius (m)
        kappa : float
            Plasma separatrix elongation
        triang : float
            Plasma separatrix triangularity
        square : float
            Plasma squareness

        Returns
        -------
        tuple
            A tuple containing:
            - len_plasma_poloidal (float): Poloidal perimeter
            - a_plasma_surface (float): Surface area
            - a_plasma_poloidal (float): Cross-section area
            - vol_plasma (float): Plasma volume


        References
        ----------
            - O. Sauter, “Geometric formulas for system codes including the effect of
              negative triangularity,” Fusion Engineering and Design, vol. 112,
              pp. 633-645, Nov. 2016,
              doi: https://doi.org/10.1016/j.fusengdes.2016.04.033.
        """
        # Calculate w07 parameter from paper from squareness assuming top-down symmetry
        w07 = square + 1

        # Inverse aspect ratio
        eps = a / r0

        # Poloidal perimeter (named Lp in Sauter)
        len_plasma_poloidal = (
            2.0e0
            * np.pi
            * a
            * (1 + 0.55 * (kappa - 1))
            * (1 + 0.08 * triang**2)
            * (1 + 0.2 * (w07 - 1))
        )

        # Surface area (named Ap in Sauter)
        a_plasma_surface = (
            2.0e0 * np.pi * r0 * (1 - 0.32 * triang * eps) * len_plasma_poloidal
        )

        # Cross-section area (named S_phi in Sauter)
        a_plasma_poloidal = np.pi * a**2 * kappa * (1 + 0.52 * (w07 - 1))

        # Volume
        vol_plasma = 2.0e0 * np.pi * r0 * (1 - 0.25 * triang * eps) * a_plasma_poloidal

        return (
            len_plasma_poloidal,
            a_plasma_surface,
            a_plasma_poloidal,
            vol_plasma,
        )

    @staticmethod
    def calculate_iter_physics_basis_elongation(
        vol_plasma: float, rmajor: float, rminor: float
    ) -> float:
        """Calculate the ITER physics basis elongation.

        This method calculates the ITER physics basis elongation (κₐ) based on the
        aspect ratio and a scaling factor.

        Parameters
        ----------
        vol_plasma :
            Plasma volume (m³)
        rmajor :
            Plasma major radius (m)
        rminor :
            Plasma minor radius (m)

        Returns
        -------
        kappa_ipb:
            The calculated ITER physics basis elongation.

        References
        ----------
            - Otto Kardaun, N. K. Thomsen, and Alexander Chudnovskiy,
              “Corrections to a sequence of papers in Nuclear Fusion,” Nuclear Fusion,
              vol. 48, no. 9, pp. 099801099801, Aug. 2008, doi: https://doi.org/10.1088/0029-5515/48/9/099801.
        """
        return (vol_plasma / (2.0 * np.pi * rmajor)) / (np.pi * rminor**2)

outfile = constants.NOUT instance-attribute

run()

Plasma geometry parameters

This method calculates the plasma geometry parameters based on various shaping terms and input values. It updates the physics_variables with calculated values for kappa, triangularity, surface area, volume, etc.

References
- J D Galambos, STAR Code : Spherical Tokamak Analysis and Reactor Code,
  unpublished internal Oak Ridge document
- H. Zohm et al, On the Physics Guidelines for a Tokamak DEMO,
  FTP/3-3, Proc. IAEA Fusion Energy Conference, October 2012, San Diego
Source code in process/models/physics/plasma_geometry.py
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
def run(self):
    """Plasma geometry parameters

    This method calculates the plasma geometry parameters based on various shaping
    terms and input values. It updates the `physics_variables` with calculated
    values for kappa, triangularity, surface area, volume, etc.

    References
    ----------
        - J D Galambos, STAR Code : Spherical Tokamak Analysis and Reactor Code,
          unpublished internal Oak Ridge document
        - H. Zohm et al, On the Physics Guidelines for a Tokamak DEMO,
          FTP/3-3, Proc. IAEA Fusion Energy Conference, October 2012, San Diego
    """
    xsi = 0.0e0
    xso = 0.0e0
    thetai = 0.0e0
    thetao = 0.0e0
    xi = 0.0e0
    xo = 0.0e0

    # Define plasma minor radius from major radius and aspect ratio
    physics_variables.rminor = physics_variables.rmajor / physics_variables.aspect

    # Define the inverse aspect ratio
    physics_variables.eps = 1.0e0 / physics_variables.aspect

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry == PlasmaGeometryModelType.IPDG89_X_POINT
    ):  # Use input kappa, physics_variables.triang values
        #  Rough estimate of 95% values
        #  ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
        #  (close to previous estimate of (physics_variables.kappa - 0.04) / 1.1
        #  over a large physics_variables.kappa range)

        physics_variables.kappa95 = physics_variables.kappa / 1.12e0
        physics_variables.triang95 = physics_variables.triang / 1.50e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry == PlasmaGeometryModelType.STAR_FIESTA
    ):  # ST scaling with physics_variables.aspect ratio [STAR Code]
        physics_variables.q95_min = 3.0e0 * (
            1.0e0 + 2.6e0 * physics_variables.eps**2.8e0
        )

        physics_variables.kappa = 2.05e0 * (
            1.0e0 + 0.44e0 * physics_variables.eps**2.1e0
        )
        physics_variables.triang = 0.53e0 * (
            1.0e0 + 0.77e0 * physics_variables.eps**3
        )

        # SIM 10/09/2020: Switched to FIESTA ST scaling  from IPDG89
        physics_variables.kappa95 = (
            physics_variables.kappa - 0.39467e0
        ) / 0.90698e0  # Fit to FIESTA (Issue #1086)
        physics_variables.triang95 = (
            physics_variables.triang - 0.048306e0
        ) / 1.3799e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry
        == PlasmaGeometryModelType.ZOHM_ITER_X_POINT
    ):  # Zohm et al. ITER scaling for elongation, input physics_variables.triang
        physics_variables.kappa = physics_variables.fkzohm * min(
            2.0e0, 1.5e0 + 0.5e0 / (physics_variables.aspect - 1.0e0)
        )

        # ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
        physics_variables.kappa95 = physics_variables.kappa / 1.12e0
        physics_variables.triang95 = physics_variables.triang / 1.50e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry == PlasmaGeometryModelType.ZOHM_ITER_95
    ):  # Zohm et al. ITER scaling for elongation, input physics_variables.triang95
        physics_variables.kappa = physics_variables.fkzohm * min(
            2.0e0, 1.5e0 + 0.5e0 / (physics_variables.aspect - 1.0e0)
        )

        # ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
        physics_variables.triang = 1.5e0 * physics_variables.triang95

        physics_variables.kappa95 = physics_variables.kappa / 1.12e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry == PlasmaGeometryModelType.IPDG89_95
    ):  # Use input kappa95, physics_variables.triang95 values
        # ITER Physics Design Guidlines: 1989 (Uckan et al. 1990)
        physics_variables.kappa = 1.12e0 * physics_variables.kappa95
        physics_variables.triang = 1.5e0 * physics_variables.triang95

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry == PlasmaGeometryModelType.MAST_DATA_95
    ):  # Use input kappa95, physics_variables.triang95 values
        # Fit to MAST data (Issue #1086)
        physics_variables.kappa = 0.91300e0 * physics_variables.kappa95 + 0.38654e0
        physics_variables.triang = 0.77394e0 * physics_variables.triang95 + 0.18515e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry
        == PlasmaGeometryModelType.MAST_DATA_X_POINT
    ):  # Use input kappa, physics_variables.triang values
        # Fit to MAST data (Issue #1086)
        physics_variables.kappa95 = (physics_variables.kappa - 0.38654e0) / 0.91300e0
        physics_variables.triang95 = (
            physics_variables.triang - 0.18515e0
        ) / 0.77394e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry == PlasmaGeometryModelType.FIESTA_RUNS_95
    ):  # Use input kappa95, physics_variables.triang95 values
        # Fit to FIESTA (Issue #1086)
        physics_variables.kappa = 0.90698e0 * physics_variables.kappa95 + 0.39467e0
        physics_variables.triang = 1.3799e0 * physics_variables.triang95 + 0.048306e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry
        == PlasmaGeometryModelType.FIESTA_RUNS_X_POINT
    ):  # Use input kappa, physics_variables.triang values
        # Fit to FIESTA (Issue #1086)
        physics_variables.kappa95 = (physics_variables.kappa - 0.39467e0) / 0.90698e0
        physics_variables.triang95 = (
            physics_variables.triang - 0.048306e0
        ) / 1.3799e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry
        == PlasmaGeometryModelType.INDUCTANCE_SCALING_X_POINT
    ):  # Use input triang, physics_variables.ind_plasma_internal_norm values
        # physics_variables.kappa found from physics_variables.aspect ratio and
        # plasma internal inductance li(3)
        physics_variables.kappa = (
            1.09e0 + 0.26e0 / physics_variables.ind_plasma_internal_norm
        ) * (1.5e0 / physics_variables.aspect) ** 0.4e0

        physics_variables.kappa95 = physics_variables.kappa / 1.12e0
        physics_variables.triang95 = physics_variables.triang / 1.50e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry
        == PlasmaGeometryModelType.CREATE_DATA_EU_DEMO_X_POINT
    ):
        # physics_variables.kappa95 found from physics_variables.aspect ratio and
        # stability margin Based on fit to CREATE data. ref Issue #1399
        # valid for EU-DEMO like machine - physics_variables.aspect ratio 2.6 - 3.6
        # Model updated see Issue #1648
        a = 3.68436807e0
        b = -0.27706527e0
        c = 0.87040251e0
        d = -18.83740952e0
        e = -0.27267618e0
        f = 20.5141261e0

        physics_variables.kappa95 = (
            -d
            - c * physics_variables.aspect
            - np.sqrt(
                (c**2.0e0 - 4.0e0 * a * b) * physics_variables.aspect**2.0e0
                + (2.0e0 * d * c - 4.0e0 * a * e) * physics_variables.aspect
                + d**2.0e0
                - 4.0e0 * a * f
                + 4.0e0 * a * physics_variables.m_s_limit
            )
        ) / (2.0e0 * a)

        if physics_variables.kappa95 > 1.77:
            ratio = 1.77 / physics_variables.kappa95
            corner_fudge = 0.3 * (physics_variables.kappa95 - 1.77) / ratio
            physics_variables.kappa95 = (
                physics_variables.kappa95 ** (ratio) + corner_fudge
            )

        physics_variables.kappa = 1.12e0 * physics_variables.kappa95
        physics_variables.triang95 = physics_variables.triang / 1.50e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry
        == PlasmaGeometryModelType.MENARD_2016_X_POINT
    ):
        # See Issue #1439
        # physics_variables.triang is an input
        # physics_variables.kappa found from physics_variables.aspect ratio scaling
        # on p32 of Menard: Menard, et al. "Fusion Nuclear Science Facilities
        # and Pilot Plants Based on the Spherical Tokamak." Nucl. Fusion, 2016, 44.

        physics_variables.kappa = 0.95e0 * (
            1.9e0 + 1.9e0 / physics_variables.aspect**1.4e0
        )
        physics_variables.kappa95 = physics_variables.kappa / 1.12e0
        physics_variables.triang95 = physics_variables.triang / 1.50e0

    # ======================================================================

    if (
        physics_variables.i_plasma_geometry
        == PlasmaGeometryModelType.MENARD_1997_X_POINT
    ):
        # physics_variables.triang is an input
        # physics_variables.kappa found from physics_variables.aspect ratio scaling from
        # J.E. Menard et al 1997 Nucl. Fusion 37 595 and assume max controllable kappa
        # and assume lᵢ(3) is held constant

        physics_variables.kappa = (
            2.93e0 * (1.8e0 / physics_variables.aspect) ** 0.4e0
        )

        physics_variables.kappa95 = physics_variables.kappa / 1.12e0
        physics_variables.triang95 = physics_variables.triang / 1.50e0

    # ======================================================================

    #  Scrape-off layer thicknesses
    if physics_variables.i_plasma_wall_gap == 0:
        self.data.build.dr_fw_plasma_gap_outboard = 0.1e0 * physics_variables.rminor
        self.data.build.dr_fw_plasma_gap_inboard = 0.1e0 * physics_variables.rminor

    # ======================================================================

    # Find parameters of arcs describing plasma surfaces
    xi, thetai, xo, thetao = self.plasma_angles_arcs(
        physics_variables.rminor,
        physics_variables.kappa,
        physics_variables.triang,
    )

    #  Surface area - inboard and outboard.  These are not given by Sauter but
    #  the outboard area is required by DCLL and divertor
    xsi, xso = self.plasma_surface_area(
        physics_variables.rmajor,
        physics_variables.rminor,
        xi,
        thetai,
        xo,
        thetao,
    )
    physics_variables.a_plasma_surface_outboard = xso

    # ======================================================================

    # i_plasma_current = 8 specifies use of the Sauter geometry as well as plasma
    # current.
    if (
        physics_variables.i_plasma_current == 8
        or physics_variables.i_plasma_shape == PlasmaShapeModelType.SAUTER
    ):
        (
            physics_variables.len_plasma_poloidal,
            physics_variables.a_plasma_surface,
            physics_variables.a_plasma_poloidal,
            physics_variables.vol_plasma,
        ) = self.sauter_geometry(
            physics_variables.rminor,
            physics_variables.rmajor,
            physics_variables.kappa,
            physics_variables.triang,
            physics_variables.plasma_square,
        )

    else:
        #  Poloidal perimeter
        physics_variables.len_plasma_poloidal = self.plasma_poloidal_perimeter(
            xi, thetai, xo, thetao
        )

        #  Volume
        physics_variables.vol_plasma = (
            physics_variables.f_vol_plasma
            * self.plasma_volume(
                physics_variables.rmajor,
                physics_variables.rminor,
                xi,
                thetai,
                xo,
                thetao,
            )
        )

        #  Cross-sectional area
        physics_variables.a_plasma_poloidal = self.plasma_cross_section(
            xi, thetai, xo, thetao
        )

        #  Surface area - sum of inboard and outboard.
        physics_variables.a_plasma_surface = xsi + xso

output()

Output plasma geometry parameters to file.

Raises:

Type Description
ProcessValueError

If n_divertors has an illegal value.

Source code in process/models/physics/plasma_geometry.py
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
def output(self):
    """Output plasma geometry parameters to file.

    Raises
    ------
    ProcessValueError
        If `n_divertors` has an illegal value.
    """
    po.oheadr(self.outfile, "Plasma Geometry")

    if self.data.stellarator.istell == 0:
        if self.data.divertor.n_divertors == 0:
            po.ocmmnt(self.outfile, "Plasma configuration = limiter")
        elif self.data.divertor.n_divertors == 1:
            po.ocmmnt(self.outfile, "Plasma configuration = single null divertor")
        elif self.data.divertor.n_divertors == 2:
            po.ocmmnt(self.outfile, "Plasma configuration = double null divertor")
        else:
            raise ProcessValueError(
                "Illegal value of n_divertors",
                n_divertors=self.data.divertor.n_divertors,
            )
    else:
        po.ocmmnt(self.outfile, "Plasma configuration = stellarator")

    if self.data.stellarator.istell == 0:
        if physics_variables.itart == 0:
            physics_variables.itart_r = physics_variables.itart
            po.ovarin(
                self.outfile,
                "Tokamak aspect ratio = Conventional, itart = 0",
                "(itart)",
                physics_variables.itart_r,
            )
        elif physics_variables.itart == 1:
            physics_variables.itart_r = physics_variables.itart
            po.ovarin(
                self.outfile,
                "Tokamak aspect ratio = Spherical, itart = 1",
                "(itart)",
                physics_variables.itart_r,
            )

    po.ovarin(
        self.outfile,
        "Plasma shaping model",
        "(i_plasma_shape)",
        physics_variables.i_plasma_shape,
    )

    po.osubhd(
        self.outfile,
        f"{PlasmaShapeModelType(physics_variables.i_plasma_shape).full_name} "
        "plasma shape model is used :",
    )

    po.ovarrf(
        self.outfile, "Major radius (R₀) (m)", "(rmajor)", physics_variables.rmajor
    )
    po.ovarrf(
        self.outfile,
        "Minor radius (a) (m)",
        "(rminor)",
        physics_variables.rminor,
        "OP ",
    )
    po.ovarrf(self.outfile, "Aspect ratio (A)", "(aspect)", physics_variables.aspect)
    po.ovarrf(
        self.outfile,
        "Plasma squareness (ζ)",
        "(plasma_square)",
        physics_variables.plasma_square,
        "IP",
    )
    po.oblnkl(self.outfile)

    po.ovarin(
        self.outfile,
        "Plasma geometry model",
        "(i_plasma_geometry)",
        physics_variables.i_plasma_geometry,
    )
    po.oblnkl(self.outfile)
    po.ovarre(
        self.outfile,
        "ITER Physics Basis definition of elongation (κₐ)",
        "(kappa_ipb)",
        physics_variables.kappa_ipb,
        "OP ",
    )
    po.oblnkl(self.outfile)
    geom_type = PlasmaGeometryModelType(physics_variables.i_plasma_geometry)

    if self.data.stellarator.istell == 0:
        po.ocmmnt(
            self.outfile,
            f"X-Point Elongation set from: {geom_type.kappa_model.description}",
        )

        po.ovarrf(
            self.outfile,
            "Elongation, X-point (κₐ)",
            "(kappa)",
            physics_variables.kappa,
            "IP"
            if geom_type.kappa_model == PlasmaGeometryModels.USER_INPUT
            else "OP",
        )
        if geom_type.kappa_model == PlasmaGeometryModels.ZOHM_ITER:
            po.ovarrf(
                self.outfile,
                "Zohm scaling adjustment factor",
                "(fkzohm)",
                physics_variables.fkzohm,
            )

        po.oblnkl(self.outfile)

        po.ocmmnt(
            self.outfile,
            f"95% Elongation set from: {geom_type.kappa95_model.description}",
        )
        po.ovarrf(
            self.outfile,
            "Elongation, 95% surface (κ₉₅)",
            "(kappa95)",
            physics_variables.kappa95,
            "IP"
            if geom_type.kappa95_model == PlasmaGeometryModels.USER_INPUT
            else "OP",
        )

        po.oblnkl(self.outfile)

        po.ocmmnt(
            self.outfile,
            f"X-Point Triangularity set from: {geom_type.triang_model.description}",
        )

        po.ovarrf(
            self.outfile,
            "Triangularity, X-point (δ)",
            "(triang)",
            physics_variables.triang,
            "IP "
            if geom_type.triang_model == PlasmaGeometryModels.USER_INPUT
            else "OP",
        )
        po.oblnkl(self.outfile)
        po.ocmmnt(
            self.outfile,
            f"95% Triangularity set from: {geom_type.triang95_model.description}",
        )

        po.ovarrf(
            self.outfile,
            "Triangularity, 95% surface (δ₉₅)",
            "(triang95)",
            physics_variables.triang95,
            "IP "
            if geom_type.triang95_model == PlasmaGeometryModels.USER_INPUT
            else "OP",
        )

        po.oblnkl(self.outfile)

        po.ovarrf(
            self.outfile,
            "Plasma poloidal perimeter (m)",
            "(len_plasma_poloidal)",
            physics_variables.len_plasma_poloidal,
            "OP ",
        )

        po.ovarre(
            self.outfile,
            "Plasma cross-sectional area (m²)",
            "(a_plasma_poloidal)",
            physics_variables.a_plasma_poloidal,
            "OP ",
        )
        po.ovarre(
            self.outfile,
            "Plasma surface area (m²)",
            "(a_plasma_surface)",
            physics_variables.a_plasma_surface,
            "OP ",
        )
        po.ovarre(
            self.outfile,
            "Plasma volume (m³)",
            "(vol_plasma)",
            physics_variables.vol_plasma,
            "OP ",
        )

    po.oblnkl(self.outfile)

plasma_angles_arcs(a, kappa, triang) staticmethod

Routine to find parameters used for calculating geometrical properties for double-null plasmas.

This function finds plasma geometrical parameters, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix.

Parameters:

Name Type Description Default
a float

Plasma minor radius (m)

required
kappa float

Plasma separatrix elongation

required
triang float

Plasma separatrix triangularity

required

Returns:

Type Description
tuple

A tuple containing: - xi (float): Radius of arc describing inboard surface (m) - thetai (float): Half-angle of arc describing inboard surface - xo (float): Radius of arc describing outboard surface (m) - thetao (float): Half-angle of arc describing outboard surface

References
  • F/MI/PJK/LOGBOOK14, p.42
  • F/PL/PJK/PROCESS/CODE/047
Source code in process/models/physics/plasma_geometry.py
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
@staticmethod
def plasma_angles_arcs(
    a: float, kappa: float, triang: float
) -> tuple[float, float, float, float]:
    """Routine to find parameters used for calculating geometrical
    properties for double-null plasmas.

    This function finds plasma geometrical parameters, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.


    Parameters
    ----------
    a : float
        Plasma minor radius (m)
    kappa : float
        Plasma separatrix elongation
    triang : float
        Plasma separatrix triangularity

    Returns
    -------
    tuple
        A tuple containing:
        - xi (float): Radius of arc describing inboard surface (m)
        - thetai (float): Half-angle of arc describing inboard surface
        - xo (float): Radius of arc describing outboard surface (m)
        - thetao (float): Half-angle of arc describing outboard surface


    References
    ----------
    - F/MI/PJK/LOGBOOK14, p.42
    - F/PL/PJK/PROCESS/CODE/047
    """
    t = 1.0e0 - triang
    denomi = (kappa**2 - t**2) / (2.0e0 * t)
    thetai = np.arctan(kappa / denomi)
    xi = a * (denomi + 1.0e0 - triang)

    #  Find radius and half-angle of outboard arc

    n = 1.0e0 + triang
    denomo = (kappa**2 - n**2) / (2.0e0 * n)
    thetao = np.arctan(kappa / denomo)
    xo = a * (denomo + 1.0e0 + triang)

    return xi, thetai, xo, thetao

plasma_poloidal_perimeter(xi, thetai, xo, thetao) staticmethod

Calculate the poloidal perimeter of the plasma.

Parameters:

Name Type Description Default
xi float

Radius of arc describing inboard surface (m)

required
thetai float

Half-angle of arc describing inboard surface

required
xo float

Radius of arc describing outboard surface (m)

required
thetao float

Half-angle of arc describing outboard surface

required

Returns:

Type Description
float

Poloidal perimeter (m)

Source code in process/models/physics/plasma_geometry.py
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
@staticmethod
def plasma_poloidal_perimeter(
    xi: float, thetai: float, xo: float, thetao: float
) -> float:
    """Calculate the poloidal perimeter of the plasma.

    Parameters
    ----------
    xi : float
        Radius of arc describing inboard surface (m)
    thetai : float
        Half-angle of arc describing inboard surface
    xo : float
        Radius of arc describing outboard surface (m)
    thetao : float
        Half-angle of arc describing outboard surface

    Returns
    -------
    float
        Poloidal perimeter (m)
    """
    return 2.0e0 * (xo * thetao + xi * thetai)

plasma_surface_area(rmajor, rminor, xi, thetai, xo, thetao) staticmethod

Plasma surface area calculation

This function finds the plasma surface area, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix.

Parameters:

Name Type Description Default
rmajor float

Plasma major radius (m)

required
rminor float

Plasma minor radius (m)

required
xi float

Radius of arc describing inboard surface (m)

required
thetai float

Half-angle of arc describing inboard surface

required
xo float

Radius of arc describing outboard surface (m)

required
thetao float

Half-angle of arc describing outboard surface

required

Returns:

Type Description
tuple

A tuple containing: - xsi (float): Inboard surface area (m^2) - xso (float): Outboard surface area (m^2)

References
  • F/MI/PJK/LOGBOOK14, p.43
Source code in process/models/physics/plasma_geometry.py
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
@staticmethod
def plasma_surface_area(
    rmajor: float, rminor: float, xi: float, thetai: float, xo: float, thetao: float
) -> tuple[float, float]:
    """Plasma surface area calculation

    This function finds the plasma surface area, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.

    Parameters
    ----------
    rmajor : float
        Plasma major radius (m)
    rminor : float
        Plasma minor radius (m)
    xi : float
        Radius of arc describing inboard surface (m)
    thetai : float
        Half-angle of arc describing inboard surface
    xo : float
        Radius of arc describing outboard surface (m)
    thetao : float
        Half-angle of arc describing outboard surface

    Returns
    -------
    tuple
        A tuple containing:
        - xsi (float): Inboard surface area (m^2)
        - xso (float): Outboard surface area (m^2)


    References
    ----------
    - F/MI/PJK/LOGBOOK14, p.43
    """
    fourpi = 4.0e0 * np.pi

    rc = rmajor - rminor + xi
    xsi = fourpi * xi * (rc * thetai - xi * np.sin(thetai))

    rc = rmajor + rminor - xo
    xso = fourpi * xo * (rc * thetao + xo * np.sin(thetao))

    return xsi, xso

plasma_volume(rmajor, rminor, xi, thetai, xo, thetao) staticmethod

Plasma volume calculation This function finds the plasma volume, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix.

Parameters:

Name Type Description Default
rmajor float

Plasma major radius (m)

required
rminor float

Plasma minor radius (m)

required
xi float

Radius of arc describing inboard surface (m)

required
thetai float

Half-angle of arc describing inboard surface

required
xo float

Radius of arc describing outboard surface (m)

required
thetao float

Half-angle of arc describing outboard surface

required

Returns:

Type Description
float

Plasma volume (m^3)

References
  • F/MI/PJK/LOGBOOK14, p.43
Source code in process/models/physics/plasma_geometry.py
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
@staticmethod
def plasma_volume(
    rmajor: float, rminor: float, xi: float, thetai: float, xo: float, thetao: float
) -> float:
    """Plasma volume calculation
    This function finds the plasma volume, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.

    Parameters
    ----------
    rmajor : float
        Plasma major radius (m)
    rminor : float
        Plasma minor radius (m)
    xi : float
        Radius of arc describing inboard surface (m)
    thetai : float
        Half-angle of arc describing inboard surface
    xo : float
        Radius of arc describing outboard surface (m)
    thetao : float
        Half-angle of arc describing outboard surface

    Returns
    -------
    float
        Plasma volume (m^3)


    References
    ----------
    - F/MI/PJK/LOGBOOK14, p.43
    """
    third = 1.0e0 / 3.0e0

    rc = rmajor - rminor + xi
    vin = (
        2.0
        * np.pi
        * xi
        * (
            rc**2 * np.sin(thetai)
            - rc * xi * thetai
            - 0.5e0 * rc * xi * np.sin(2.0e0 * thetai)
            + xi * xi * np.sin(thetai)
            - third * xi * xi * (np.sin(thetai)) ** 3
        )
    )

    rc = rmajor + rminor - xo
    vout = (
        2.0
        * np.pi
        * xo
        * (
            rc**2 * np.sin(thetao)
            + rc * xo * thetao
            + 0.5e0 * rc * xo * np.sin(2.0e0 * thetao)
            + xo * xo * np.sin(thetao)
            - third * xo * xo * (np.sin(thetao)) ** 3
        )
    )

    return vout - vin

plasma_cross_section(xi, thetai, xo, thetao) staticmethod

Plasma cross-sectional area calculation

This function finds the plasma cross-sectional area, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix.

Parameters:

Name Type Description Default
xi float

Radius of arc describing inboard surface (m)

required
thetai float

Half-angle of arc describing inboard surface

required
xo float

Radius of arc describing outboard surface (m)

required
thetao float

Half-angle of arc describing outboard surface

required

Returns:

Type Description
float

Plasma cross-sectional area (m^2)

References
  • F/MI/PJK/LOGBOOK14, p.41
Source code in process/models/physics/plasma_geometry.py
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
934
@staticmethod
def plasma_cross_section(
    xi: float, thetai: float, xo: float, thetao: float
) -> float:
    """Plasma cross-sectional area calculation

    This function finds the plasma cross-sectional area, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.

    Parameters
    ----------
    xi : float
        Radius of arc describing inboard surface (m)
    thetai : float
        Half-angle of arc describing inboard surface
    xo : float
        Radius of arc describing outboard surface (m)
    thetao : float
        Half-angle of arc describing outboard surface

    Returns
    -------
    float
        Plasma cross-sectional area (m^2)


    References
    ----------
    - F/MI/PJK/LOGBOOK14, p.41
    """
    return xo**2 * (thetao - np.cos(thetao) * np.sin(thetao)) + xi**2 * (
        thetai - np.cos(thetai) * np.sin(thetai)
    )

sauter_geometry(a, r0, kappa, triang, square) staticmethod

Calculate the plasma geometry parameters using the Sauter geometry model.

Parameters:

Name Type Description Default
a float

Plasma minor radius (m)

required
r0 float

Plasma major radius (m)

required
kappa float

Plasma separatrix elongation

required
triang float

Plasma separatrix triangularity

required
square float

Plasma squareness

required

Returns:

Type Description
tuple

A tuple containing: - len_plasma_poloidal (float): Poloidal perimeter - a_plasma_surface (float): Surface area - a_plasma_poloidal (float): Cross-section area - vol_plasma (float): Plasma volume

References
- O. Sauter, Geometric formulas for system codes including the effect of
  negative triangularity,” Fusion Engineering and Design, vol. 112,
  pp. 633-645, Nov. 2016,
  doi: https://doi.org/10.1016/j.fusengdes.2016.04.033.
Source code in process/models/physics/plasma_geometry.py
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
@staticmethod
def sauter_geometry(
    a: float, r0: float, kappa: float, triang: float, square: float
) -> tuple[float, float, float, float, float]:
    """Calculate the plasma geometry parameters using the Sauter geometry model.

    Parameters
    ----------
    a : float
        Plasma minor radius (m)
    r0 : float
        Plasma major radius (m)
    kappa : float
        Plasma separatrix elongation
    triang : float
        Plasma separatrix triangularity
    square : float
        Plasma squareness

    Returns
    -------
    tuple
        A tuple containing:
        - len_plasma_poloidal (float): Poloidal perimeter
        - a_plasma_surface (float): Surface area
        - a_plasma_poloidal (float): Cross-section area
        - vol_plasma (float): Plasma volume


    References
    ----------
        - O. Sauter, “Geometric formulas for system codes including the effect of
          negative triangularity,” Fusion Engineering and Design, vol. 112,
          pp. 633-645, Nov. 2016,
          doi: https://doi.org/10.1016/j.fusengdes.2016.04.033.
    """
    # Calculate w07 parameter from paper from squareness assuming top-down symmetry
    w07 = square + 1

    # Inverse aspect ratio
    eps = a / r0

    # Poloidal perimeter (named Lp in Sauter)
    len_plasma_poloidal = (
        2.0e0
        * np.pi
        * a
        * (1 + 0.55 * (kappa - 1))
        * (1 + 0.08 * triang**2)
        * (1 + 0.2 * (w07 - 1))
    )

    # Surface area (named Ap in Sauter)
    a_plasma_surface = (
        2.0e0 * np.pi * r0 * (1 - 0.32 * triang * eps) * len_plasma_poloidal
    )

    # Cross-section area (named S_phi in Sauter)
    a_plasma_poloidal = np.pi * a**2 * kappa * (1 + 0.52 * (w07 - 1))

    # Volume
    vol_plasma = 2.0e0 * np.pi * r0 * (1 - 0.25 * triang * eps) * a_plasma_poloidal

    return (
        len_plasma_poloidal,
        a_plasma_surface,
        a_plasma_poloidal,
        vol_plasma,
    )

calculate_iter_physics_basis_elongation(vol_plasma, rmajor, rminor) staticmethod

Calculate the ITER physics basis elongation.

This method calculates the ITER physics basis elongation (κₐ) based on the aspect ratio and a scaling factor.

Parameters:

Name Type Description Default
vol_plasma float

Plasma volume (m³)

required
rmajor float

Plasma major radius (m)

required
rminor float

Plasma minor radius (m)

required

Returns:

Name Type Description
kappa_ipb float

The calculated ITER physics basis elongation.

References
- Otto Kardaun, N. K. Thomsen, and Alexander Chudnovskiy,
  “Corrections to a sequence of papers in Nuclear Fusion,” Nuclear Fusion,
  vol. 48, no. 9, pp. 099801099801, Aug. 2008, doi: https://doi.org/10.1088/0029-5515/48/9/099801.
Source code in process/models/physics/plasma_geometry.py
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
@staticmethod
def calculate_iter_physics_basis_elongation(
    vol_plasma: float, rmajor: float, rminor: float
) -> float:
    """Calculate the ITER physics basis elongation.

    This method calculates the ITER physics basis elongation (κₐ) based on the
    aspect ratio and a scaling factor.

    Parameters
    ----------
    vol_plasma :
        Plasma volume (m³)
    rmajor :
        Plasma major radius (m)
    rminor :
        Plasma minor radius (m)

    Returns
    -------
    kappa_ipb:
        The calculated ITER physics basis elongation.

    References
    ----------
        - Otto Kardaun, N. K. Thomsen, and Alexander Chudnovskiy,
          “Corrections to a sequence of papers in Nuclear Fusion,” Nuclear Fusion,
          vol. 48, no. 9, pp. 099801099801, Aug. 2008, doi: https://doi.org/10.1088/0029-5515/48/9/099801.
    """
    return (vol_plasma / (2.0 * np.pi * rmajor)) / (np.pi * rminor**2)

surfa(a, r, k, d)

Plasma surface area calculation

This function finds the plasma surface area, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix. It was the original method in PROCESS.

Parameters:

Name Type Description Default
a float

Plasma minor radius (m)

required
r float

Plasma major radius (m)

required
k float

Plasma separatrix elongation

required
d float

Plasma separatrix triangularity

required

Returns:

Type Description
tuple

A tuple containing: - sa (float): Plasma total surface area (m^2) - so (float): Plasma outboard surface area (m^2)

References
  • J D Galambos, STAR Code : Spherical Tokamak Analysis and Reactor Code, unpublished internal Oak Ridge document
Source code in process/models/physics/plasma_geometry.py
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
def surfa(a: float, r: float, k: float, d: float) -> tuple[float, float]:
    """Plasma surface area calculation

    This function finds the plasma surface area, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.
    It was the original method in PROCESS.

    Parameters
    ----------
    a : float
        Plasma minor radius (m)
    r : float
        Plasma major radius (m)
    k : float
        Plasma separatrix elongation
    d : float
        Plasma separatrix triangularity

    Returns
    -------
    tuple
        A tuple containing:
        - sa (float): Plasma total surface area (m^2)
        - so (float): Plasma outboard surface area (m^2)


    References
    ----------
    - J D Galambos, STAR Code : Spherical Tokamak Analysis and Reactor Code,
    unpublished internal Oak Ridge document
    """
    radco = a * (1.0e0 + (k**2 + d**2 - 1.0e0) / (2.0e0 * (1.0e0 + d)))
    b = k * a
    thto = np.arcsin(b / radco)
    so = 4.0e0 * np.pi * radco * ((r + a - radco) * thto + b)

    #  Inboard side

    radci = a * (1.0e0 + (k**2 + d**2 - 1.0e0) / (2.0e0 * (1.0e0 - d)))
    thti = np.arcsin(b / radci)
    si = 4.0e0 * np.pi * radci * ((r - a + radci) * thti - b)

    sa = so + si

    return sa, so

perim(a, kap, tri)

Plasma poloidal perimeter calculation

This function finds the plasma poloidal perimeter, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix.

Parameters:

Name Type Description Default
a float

Plasma minor radius (m)

required
kap float

Plasma separatrix elongation

required
tri float

Plasma separatrix triangularity

required

Returns:

Type Description
float

Plasma poloidal perimeter (m)

References
  • F/PL/PJK/PROCESS/CODE/047
Source code in process/models/physics/plasma_geometry.py
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
def perim(a: float, kap: float, tri: float) -> float:
    """Plasma poloidal perimeter calculation

    This function finds the plasma poloidal perimeter, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.

    Parameters
    ----------
    a : float
        Plasma minor radius (m)
    kap : float
        Plasma separatrix elongation
    tri : float
        Plasma separatrix triangularity

    Returns
    -------
    float
        Plasma poloidal perimeter (m)


    References
    ----------
    - F/PL/PJK/PROCESS/CODE/047
    """
    #  Inboard arc

    denomi = (tri**2 + kap**2 - 1.0e0) / (2.0e0 * (1.0e0 - tri)) + tri
    thetai = np.arctan(kap / denomi)
    xli = a * (denomi + 1.0e0 - tri)

    #  Outboard arc

    denomo = (tri**2 + kap**2 - 1.0e0) / (2.0e0 * (1.0e0 + tri)) - tri
    thetao = np.arctan(kap / denomo)
    xlo = a * (denomo + 1.0e0 + tri)

    return 2.0e0 * (xlo * thetao + xli * thetai)

fvol(r, a, kap, tri)

Plasma volume calculation

This function finds the plasma volume, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix.

Parameters:

Name Type Description Default
r float

Plasma major radius (m)

required
a float

Plasma minor radius (m)

required
kap float

Plasma separatrix elongation

required
tri float

Plasma separatrix triangularity

required

Returns:

Type Description
float

Plasma volume (m^3)

References
  • F/MI/PJK/LOGBOOK14, p.41
  • F/PL/PJK/PROCESS/CODE/047
Source code in process/models/physics/plasma_geometry.py
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
def fvol(r: float, a: float, kap: float, tri: float) -> float:
    """Plasma volume calculation

    This function finds the plasma volume, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.

    Parameters
    ----------
    r : float
        Plasma major radius (m)
    a : float
        Plasma minor radius (m)
    kap : float
        Plasma separatrix elongation
    tri : float
        Plasma separatrix triangularity

    Returns
    -------
    float
        Plasma volume (m^3)


    References
    ----------
    - F/MI/PJK/LOGBOOK14, p.41
    - F/PL/PJK/PROCESS/CODE/047
    """
    zn = kap * a

    c1 = ((r + a) ** 2 - (r - tri * a) ** 2 - zn**2) / (2.0e0 * (1.0e0 + tri) * a)
    rc1 = r + a - c1
    vout = (
        -0.66666666e0 * np.pi * zn**3
        + 2.0e0 * np.pi * zn * (c1**2 + rc1**2)
        + 2.0e0
        * np.pi
        * c1
        * (zn * np.sqrt(rc1**2 - zn**2) + rc1**2 * np.arcsin(zn / rc1))
    )

    c2 = (-((r - a) ** 2) + (r - tri * a) ** 2 + zn**2) / (2.0e0 * (1.0e0 - tri) * a)
    rc2 = c2 - r + a
    vin = (
        -0.66666e0 * np.pi * zn**3
        + 2.0e0 * np.pi * zn * (rc2**2 + c2**2)
        - 2.0e0
        * np.pi
        * c2
        * (zn * np.sqrt(rc2**2 - zn**2) + rc2**2 * np.arcsin(zn / rc2))
    )

    return vout - vin

xsect0(a, kap, tri)

Plasma cross-sectional area calculation

This function finds the plasma cross-sectional area, using the revolution of two intersecting arcs around the device centreline. This calculation is appropriate for plasmas with a separatrix.

Parameters:

Name Type Description Default
a float

Plasma minor radius (m)

required
kap float

Plasma separatrix elongation

required
tri float

Plasma separatrix triangularity

required

Returns:

Type Description
float

Plasma cross-sectional area (m^2)

References
  • F/MI/PJK/LOGBOOK14, p.41
  • F/PL/PJK/PROCESS/CODE/047
Source code in process/models/physics/plasma_geometry.py
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
def xsect0(a: float, kap: float, tri: float) -> float:
    """Plasma cross-sectional area calculation

    This function finds the plasma cross-sectional area, using the
    revolution of two intersecting arcs around the device centreline.
    This calculation is appropriate for plasmas with a separatrix.

    Parameters
    ----------
    a : float
        Plasma minor radius (m)
    kap : float
        Plasma separatrix elongation
    tri : float
        Plasma separatrix triangularity

    Returns
    -------
    float
        Plasma cross-sectional area (m^2)

    References
    ----------
    - F/MI/PJK/LOGBOOK14, p.41
    - F/PL/PJK/PROCESS/CODE/047
    """
    denomi = (tri**2 + kap**2 - 1.0e0) / (2.0e0 * (1.0e0 - tri)) + tri
    thetai = np.arctan(kap / denomi)
    xli = a * (denomi + 1.0e0 - tri)

    cti = np.cos(thetai)
    sti = np.sin(thetai)

    #  Find radius and half-angle of outboard arc

    denomo = (tri**2 + kap**2 - 1.0e0) / (2.0e0 * (1.0e0 + tri)) - tri
    thetao = np.arctan(kap / denomo)
    xlo = a * (denomo + 1.0e0 + tri)

    cto = np.cos(thetao)
    sto = np.sin(thetao)

    #  Find cross-sectional area

    return xlo**2 * (thetao - cto * sto) + xli**2 * (thetai - cti * sti)