Skip to content

plasma_geometry

PlasmaGeom

Source code in process/models/physics/plasma_geometry.py
 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
class PlasmaGeom:
    def __init__(self):
        self.outfile = constants.NOUT

    def plasma_geometry(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 == 0
        ):  # 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 == 1
        ):  # 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 == 2
        ):  # 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 == 3
        ):  # 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 == 4
        ):  # 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 == 5
        ):  # 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 == 6
        ):  # 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 == 7
        ):  # 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 == 8
        ):  # 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 == 9
        ):  # 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 == 10:
            # physics_variables.kappa95 found from physics_variables.aspect ratio and stabilty 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 == 11:
            # 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

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

        #  Scrape-off layer thicknesses
        if physics_variables.i_plasma_wall_gap == 0:
            build_variables.dr_fw_plasma_gap_outboard = 0.1e0 * physics_variables.rminor
            build_variables.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 == 1
        ):
            (
                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

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

    @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
        tri : 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,
        )

outfile = constants.NOUT instance-attribute

plasma_geometry()

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
 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
def plasma_geometry(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 == 0
    ):  # 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 == 1
    ):  # 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 == 2
    ):  # 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 == 3
    ):  # 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 == 4
    ):  # 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 == 5
    ):  # 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 == 6
    ):  # 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 == 7
    ):  # 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 == 8
    ):  # 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 == 9
    ):  # 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 == 10:
        # physics_variables.kappa95 found from physics_variables.aspect ratio and stabilty 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 == 11:
        # 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

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

    #  Scrape-off layer thicknesses
    if physics_variables.i_plasma_wall_gap == 0:
        build_variables.dr_fw_plasma_gap_outboard = 0.1e0 * physics_variables.rminor
        build_variables.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 == 1
    ):
        (
            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

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
tri 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
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
@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
    tri : 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
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
@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
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
@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
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
@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
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
@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
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
@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,
    )

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
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
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
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
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
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
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
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
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)