Skip to content

divertor

Divertor

Module containing divertor routines

This module contains routines relevant for calculating the divertor parameters for a fusion power plant.

Source code in process/models/divertor.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
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
class Divertor:
    """Module containing divertor routines

    This module contains routines relevant for calculating the
    divertor parameters for a fusion power plant.
    """

    def __init__(self):
        self.outfile = constants.NOUT  # output file unit

    def run(self, output: bool):
        """Routine to call the divertor model

        This subroutine calls the divertor routine. This routine scales
        dimensions, powers and field levels which are used as input to
        the Harrison divertor model.

        Parameters
        ----------
        output :
            indicate whether output should be written to the output file, or not
        """

        fwbs.p_div_nuclear_heat_total_mw = self.incident_neutron_power(
            p_plasma_neutron_mw=pv.p_plasma_neutron_mw,
            f_ster_div_single=fwbs.f_ster_div_single,
            n_divertors=dv.n_divertors,
        )

        fwbs.p_div_rad_total_mw = self.incident_radiation_power(
            p_plasma_rad_mw=pv.p_plasma_rad_mw,
            f_ster_div_single=fwbs.f_ster_div_single,
            n_divertors=dv.n_divertors,
        )

        if dv.i_div_heat_load == 0 and output:
            po.ovarre(
                self.outfile,
                "Divertor heat load (MW/m2)",
                "(pflux_div_heat_load_mw)",
                dv.pflux_div_heat_load_mw,
            )
            return
        if dv.i_div_heat_load == 1:
            self.divtart(
                pv.rmajor,
                pv.rminor,
                pv.triang,
                bv.dr_fw_plasma_gap_inboard,
                bv.dz_xpoint_divertor,
                pv.p_plasma_separatrix_mw,
                output=output,
                i_single_null=pv.i_single_null,
                dz_divertor=dv.dz_divertor,
            )
            return
        if dv.i_div_heat_load == 2:
            self.divwade(
                pv.rmajor,
                pv.rminor,
                pv.aspect,
                pv.b_plasma_toroidal_on_axis,
                pv.b_plasma_poloidal_average,
                pv.p_plasma_separatrix_mw,
                dv.f_div_flux_expansion,
                pv.nd_plasma_separatrix_electron,
                dv.deg_div_field_plate,
                pv.rad_fraction_sol,
                pv.f_p_div_lower,
                output=output,
            )
            return

    def divtart(
        self,
        rmajor: float,
        rminor: float,
        triang: float,
        dr_fw_plasma_gap_inboard: float,
        dz_xpoint_divertor: float,
        p_plasma_separatrix_mw: float,
        output: bool,
        i_single_null: int,
        dz_divertor: float,
    ) -> float:
        """Tight aspect ratio tokamak divertor model

        This method calculates the divertor heat load for a tight aspect
        ratio machine, assuming that the power is evenly distributed around the
        divertor chamber by the action of a gaseous target. Each divertor is
        modeled as approximately triangular in the R,Z plane.

        Parameters
        ----------
        rmajor : float
            Plasma major radius (m)
        rminor : float
            Plasma minor radius (m)
        triang : float
            Plasma triangularity
        dr_fw_plasma_gap_inboard : float
            Inboard scrape-off width (m)
        dz_xpoint_divertor : float
            Vertical distance from X-point to divertor (m)
        p_plasma_separatrix_mw : float
            Power to the divertor (MW)
        output : bool
            Indicates whether output should be written to the output file
        i_single_null : int
            1 for single null configuration, 0 for double null
        dz_divertor : float
            Vertical height of the divertor (m)

        Returns
        -------
        float
            Divertor heat load for a tight aspect ratio machine (MW/m2)


        Notes
        -----
            - This model assumes a tight aspect ratio tokamak with a gaseous target
              divertor. The divertor chamber is modeled as triangular in the R,Z plane,
              and the heat load is calculated based on the total divertor surface area.
            - The method accounts for both single null and double null configurations.

        References
        ----------
            - Y.-K. M. Peng, J. B. Hicks, AEA Fusion, Culham (UK), "Engineering feasibility of tight aspect ratio Tokamak (spherical torus) reactors".
              1990. https://inis.iaea.org/records/ey2rf-dah04

            - Y.-K. M. Peng, J. B. Hicks, “Engineering feasibility of tight aspect ratio tokamak (spherical torus) reactors,”
              Osti.gov, 1991. https://www.osti.gov/biblio/1022679 (accessed Mar. 24, 2025).
        """

        #  Thickness of centrepost + first wall at divertor height

        r1 = rmajor - rminor * triang - 3.0e0 * dr_fw_plasma_gap_inboard + tfv.drtop

        #  Outer radius of divertor region

        r2 = rmajor + rminor

        #  Angle of diagonal divertor plate from horizontal

        if dz_xpoint_divertor <= 0.0e0:
            raise ProcessValueError(
                "Non-positive dz_xpoint_divertor", dz_xpoint_divertor=dz_xpoint_divertor
            )

        theta = math.atan(dz_divertor / (r2 - r1))

        #  Vertical plate area

        a1 = 2.0e0 * np.pi * r1 * dz_divertor

        #  Horizontal plate area

        a2 = np.pi * (r2 * r2 - r1 * r1)

        #  Diagonal plate area

        a3 = a2 / (math.cos(theta) * math.cos(theta))

        #  Total divertor area

        # Single null case
        if i_single_null == 1:
            areadv = a1 + a2 + a3
        # Double null case
        elif i_single_null == 0:
            areadv = 2.0 * (a1 + a2 + a3)

        if dv.i_div_heat_load == 1:
            dv.pflux_div_heat_load_mw = p_plasma_separatrix_mw / areadv

        if output and dv.i_div_heat_load == 1:
            po.osubhd(self.outfile, "Divertor Heat Load")
            po.ocmmnt(self.outfile, "Assume an expanded divertor with a gaseous target")
            po.oblnkl(self.outfile)
            po.ovarre(
                self.outfile,
                "Power to the divertor (MW)",
                "(p_plasma_separatrix_mw.)",
                p_plasma_separatrix_mw,
            )
            po.ovarre(self.outfile, "Divertor surface area (m2)", "(areadv)", areadv)
            po.ovarre(
                self.outfile,
                "Divertor heat load (MW/m2)",
                "(pflux_div_heat_load_mw)",
                dv.pflux_div_heat_load_mw,
            )

        elif output:
            po.osubhd(self.outfile, "Divertor Heat Load")
            po.ovarre(
                self.outfile,
                "Power to the divertor (MW)",
                "(p_plasma_separatrix_mw.)",
                p_plasma_separatrix_mw,
            )
            po.ovarre(
                self.outfile,
                "Divertor heat load (MW/m2)",
                "(pflux_div_heat_load_mw)",
                dv.pflux_div_heat_load_mw,
            )
        return dv.pflux_div_heat_load_mw

    def divwade(
        self,
        rmajor: float,
        rminor: float,
        aspect: float,
        b_plasma_toroidal_on_axis: float,
        b_plasma_poloidal_average: float,
        p_plasma_separatrix_mw: float,
        f_div_flux_expansion: float,
        nd_plasma_separatrix_electron: float,
        deg_div_field_plate: float,
        rad_fraction_sol: float,
        f_p_div_lower: float,
        output: bool,
    ) -> float:
        """Divertor heat load model (Wade 2020)

        This subroutine calculates the divertor heat flux for any machine,
        with either a single null or double null configuration.
        It uses the Eich scaling (Eich et al. 2013) and spreading factor (Scarabosio et al. 2014)
        to calculate the SOL width. This is then used with a flux expansion factor to calculate
        the wetted area and then the heat flux.

        Parameters
        ----------
        rmajor : float
            plasma major radius (m)
        rminor : float
            plasma minor radius (m)
        aspect : float
            tokamak aspect ratio
        b_plasma_toroidal_on_axis : float
            toroidal field (T)
        b_plasma_poloidal_average : float
            poloidal field (T)
        p_plasma_separatrix_mw : float
            power to divertor (MW)
        f_div_flux_expansion : float
            plasma flux expansion in divertor
        nd_plasma_separatrix_electron : float
            electron density at separatrix (m-3)
        deg_div_field_plate : float
            field line angle wrt divertor target plate (degrees)
        rad_fraction_sol : float
            SOL radiation fraction
        f_p_div_lower : float
            fraction of power to the lower divertor in double null configuration

        Returns
        -------
        float
            divertor heat load for a tight aspect ratio machine
        """

        # Radius on midplane
        r_omp = rmajor + rminor

        # B fields on midplane
        Bp_omp = -b_plasma_poloidal_average * rmajor / r_omp

        Bt_omp = -b_plasma_toroidal_on_axis * rmajor / r_omp

        # Eich scaling for lambda_q
        lambda_eich = (
            1.35
            * p_plasma_separatrix_mw**-0.02
            * rmajor**0.04
            * b_plasma_poloidal_average**-0.92
            * aspect**0.42
        )

        # Spreading factor
        spread_fact = (
            0.12
            * (nd_plasma_separatrix_electron / 1e19) ** -0.02
            * p_plasma_separatrix_mw**-0.21
            * rmajor**0.71
            * b_plasma_poloidal_average**-0.82
        )

        # SOL width
        lambda_int = lambda_eich + 1.64 * spread_fact

        # Flux angle on midplane
        alpha_mid = math.degrees(math.atan(Bp_omp / Bt_omp))

        # Flux angle in the divertor
        alpha_div = f_div_flux_expansion * alpha_mid

        # Tilt of the separatrix relative to the target in the poloidal plane
        theta_div = math.asin(
            (1 + 1 / alpha_div**2) * math.sin(math.radians(deg_div_field_plate))
        )

        # Wetted area
        area_wetted = (
            2 * np.pi * rmajor * lambda_int * f_div_flux_expansion * math.sin(theta_div)
        )

        # Divertor heat load
        hldiv_base = p_plasma_separatrix_mw * (1 - rad_fraction_sol) / area_wetted

        # For double null, calculate heat loads to upper and lower divertors and use the highest
        if dv.n_divertors == 2:
            hldiv_lower = f_p_div_lower * hldiv_base
            hldiv_upper = (1.0 - f_p_div_lower) * hldiv_base
            dv.pflux_div_heat_load_mw = max(hldiv_lower, hldiv_upper)
        else:
            dv.pflux_div_heat_load_mw = hldiv_base

        if output:
            po.osubhd(self.outfile, "Divertor Heat Load")
            po.ocmmnt(self.outfile, "Assume an expanded divertor with a gaseous target")
            po.oblnkl(self.outfile)
            po.ovarre(
                self.outfile,
                "Flux expansion",
                "(f_div_flux_expansion)",
                f_div_flux_expansion,
            )
            po.ovarre(
                self.outfile,
                "Field line angle wrt to target divertor plate (degrees)",
                "(deg_div_field_plate)",
                deg_div_field_plate,
            )
            po.ovarre(
                self.outfile,
                "Divertor heat load (MW/m2)",
                "(pflux_div_heat_load_mw)",
                dv.pflux_div_heat_load_mw,
            )
        return dv.pflux_div_heat_load_mw

    def incident_radiation_power(
        self,
        p_plasma_rad_mw: float,
        f_ster_div_single: float,
        n_divertors: int,
    ) -> float:
        """Calculates the total incident radiation power on the divertor box.

        Parameters
        ----------
        p_plasma_rad_mw : float
            Total plasma radiated power in megawatts (MW).
        f_ster_div_single : float
            Fraction of the solid angle subtended by a single divertor.
        n_divertors : int
            Number of divertors.

        Returns
        -------
        float
            Total incident radiation power on the divertor box in megawatts (MW).
        """

        return p_plasma_rad_mw * f_ster_div_single * n_divertors

    def incident_neutron_power(
        self,
        p_plasma_neutron_mw: float,
        f_ster_div_single: float,
        n_divertors: int,
    ) -> float:
        """Calculates the total incident neutron power on the divertor box.

        Parameters
        ----------
        p_plasma_neutron_mw : float
            Total plasma neutron power in megawatts (MW).
        f_ster_div_single : float
            Fraction of the solid angle subtended by a single divertor.
        n_divertors : int
            Number of divertors.

        Returns
        -------
        float
            Total incident radiation power on the divertor box in megawatts (MW).
        """

        return p_plasma_neutron_mw * f_ster_div_single * n_divertors

outfile = constants.NOUT instance-attribute

run(output)

Routine to call the divertor model

This subroutine calls the divertor routine. This routine scales dimensions, powers and field levels which are used as input to the Harrison divertor model.

Parameters:

Name Type Description Default
output bool

indicate whether output should be written to the output file, or not

required
Source code in process/models/divertor.py
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
def run(self, output: bool):
    """Routine to call the divertor model

    This subroutine calls the divertor routine. This routine scales
    dimensions, powers and field levels which are used as input to
    the Harrison divertor model.

    Parameters
    ----------
    output :
        indicate whether output should be written to the output file, or not
    """

    fwbs.p_div_nuclear_heat_total_mw = self.incident_neutron_power(
        p_plasma_neutron_mw=pv.p_plasma_neutron_mw,
        f_ster_div_single=fwbs.f_ster_div_single,
        n_divertors=dv.n_divertors,
    )

    fwbs.p_div_rad_total_mw = self.incident_radiation_power(
        p_plasma_rad_mw=pv.p_plasma_rad_mw,
        f_ster_div_single=fwbs.f_ster_div_single,
        n_divertors=dv.n_divertors,
    )

    if dv.i_div_heat_load == 0 and output:
        po.ovarre(
            self.outfile,
            "Divertor heat load (MW/m2)",
            "(pflux_div_heat_load_mw)",
            dv.pflux_div_heat_load_mw,
        )
        return
    if dv.i_div_heat_load == 1:
        self.divtart(
            pv.rmajor,
            pv.rminor,
            pv.triang,
            bv.dr_fw_plasma_gap_inboard,
            bv.dz_xpoint_divertor,
            pv.p_plasma_separatrix_mw,
            output=output,
            i_single_null=pv.i_single_null,
            dz_divertor=dv.dz_divertor,
        )
        return
    if dv.i_div_heat_load == 2:
        self.divwade(
            pv.rmajor,
            pv.rminor,
            pv.aspect,
            pv.b_plasma_toroidal_on_axis,
            pv.b_plasma_poloidal_average,
            pv.p_plasma_separatrix_mw,
            dv.f_div_flux_expansion,
            pv.nd_plasma_separatrix_electron,
            dv.deg_div_field_plate,
            pv.rad_fraction_sol,
            pv.f_p_div_lower,
            output=output,
        )
        return

divtart(rmajor, rminor, triang, dr_fw_plasma_gap_inboard, dz_xpoint_divertor, p_plasma_separatrix_mw, output, i_single_null, dz_divertor)

Tight aspect ratio tokamak divertor model

This method calculates the divertor heat load for a tight aspect ratio machine, assuming that the power is evenly distributed around the divertor chamber by the action of a gaseous target. Each divertor is modeled as approximately triangular in the R,Z plane.

Parameters:

Name Type Description Default
rmajor float

Plasma major radius (m)

required
rminor float

Plasma minor radius (m)

required
triang float

Plasma triangularity

required
dr_fw_plasma_gap_inboard float

Inboard scrape-off width (m)

required
dz_xpoint_divertor float

Vertical distance from X-point to divertor (m)

required
p_plasma_separatrix_mw float

Power to the divertor (MW)

required
output bool

Indicates whether output should be written to the output file

required
i_single_null int

1 for single null configuration, 0 for double null

required
dz_divertor float

Vertical height of the divertor (m)

required

Returns:

Type Description
float

Divertor heat load for a tight aspect ratio machine (MW/m2)

Notes
- This model assumes a tight aspect ratio tokamak with a gaseous target
  divertor. The divertor chamber is modeled as triangular in the R,Z plane,
  and the heat load is calculated based on the total divertor surface area.
- The method accounts for both single null and double null configurations.
References
- Y.-K. M. Peng, J. B. Hicks, AEA Fusion, Culham (UK), "Engineering feasibility of tight aspect ratio Tokamak (spherical torus) reactors".
  1990. https://inis.iaea.org/records/ey2rf-dah04

- Y.-K. M. Peng, J. B. Hicks, “Engineering feasibility of tight aspect ratio tokamak (spherical torus) reactors,”
  Osti.gov, 1991. https://www.osti.gov/biblio/1022679 (accessed Mar. 24, 2025).
Source code in process/models/divertor.py
 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
def divtart(
    self,
    rmajor: float,
    rminor: float,
    triang: float,
    dr_fw_plasma_gap_inboard: float,
    dz_xpoint_divertor: float,
    p_plasma_separatrix_mw: float,
    output: bool,
    i_single_null: int,
    dz_divertor: float,
) -> float:
    """Tight aspect ratio tokamak divertor model

    This method calculates the divertor heat load for a tight aspect
    ratio machine, assuming that the power is evenly distributed around the
    divertor chamber by the action of a gaseous target. Each divertor is
    modeled as approximately triangular in the R,Z plane.

    Parameters
    ----------
    rmajor : float
        Plasma major radius (m)
    rminor : float
        Plasma minor radius (m)
    triang : float
        Plasma triangularity
    dr_fw_plasma_gap_inboard : float
        Inboard scrape-off width (m)
    dz_xpoint_divertor : float
        Vertical distance from X-point to divertor (m)
    p_plasma_separatrix_mw : float
        Power to the divertor (MW)
    output : bool
        Indicates whether output should be written to the output file
    i_single_null : int
        1 for single null configuration, 0 for double null
    dz_divertor : float
        Vertical height of the divertor (m)

    Returns
    -------
    float
        Divertor heat load for a tight aspect ratio machine (MW/m2)


    Notes
    -----
        - This model assumes a tight aspect ratio tokamak with a gaseous target
          divertor. The divertor chamber is modeled as triangular in the R,Z plane,
          and the heat load is calculated based on the total divertor surface area.
        - The method accounts for both single null and double null configurations.

    References
    ----------
        - Y.-K. M. Peng, J. B. Hicks, AEA Fusion, Culham (UK), "Engineering feasibility of tight aspect ratio Tokamak (spherical torus) reactors".
          1990. https://inis.iaea.org/records/ey2rf-dah04

        - Y.-K. M. Peng, J. B. Hicks, “Engineering feasibility of tight aspect ratio tokamak (spherical torus) reactors,”
          Osti.gov, 1991. https://www.osti.gov/biblio/1022679 (accessed Mar. 24, 2025).
    """

    #  Thickness of centrepost + first wall at divertor height

    r1 = rmajor - rminor * triang - 3.0e0 * dr_fw_plasma_gap_inboard + tfv.drtop

    #  Outer radius of divertor region

    r2 = rmajor + rminor

    #  Angle of diagonal divertor plate from horizontal

    if dz_xpoint_divertor <= 0.0e0:
        raise ProcessValueError(
            "Non-positive dz_xpoint_divertor", dz_xpoint_divertor=dz_xpoint_divertor
        )

    theta = math.atan(dz_divertor / (r2 - r1))

    #  Vertical plate area

    a1 = 2.0e0 * np.pi * r1 * dz_divertor

    #  Horizontal plate area

    a2 = np.pi * (r2 * r2 - r1 * r1)

    #  Diagonal plate area

    a3 = a2 / (math.cos(theta) * math.cos(theta))

    #  Total divertor area

    # Single null case
    if i_single_null == 1:
        areadv = a1 + a2 + a3
    # Double null case
    elif i_single_null == 0:
        areadv = 2.0 * (a1 + a2 + a3)

    if dv.i_div_heat_load == 1:
        dv.pflux_div_heat_load_mw = p_plasma_separatrix_mw / areadv

    if output and dv.i_div_heat_load == 1:
        po.osubhd(self.outfile, "Divertor Heat Load")
        po.ocmmnt(self.outfile, "Assume an expanded divertor with a gaseous target")
        po.oblnkl(self.outfile)
        po.ovarre(
            self.outfile,
            "Power to the divertor (MW)",
            "(p_plasma_separatrix_mw.)",
            p_plasma_separatrix_mw,
        )
        po.ovarre(self.outfile, "Divertor surface area (m2)", "(areadv)", areadv)
        po.ovarre(
            self.outfile,
            "Divertor heat load (MW/m2)",
            "(pflux_div_heat_load_mw)",
            dv.pflux_div_heat_load_mw,
        )

    elif output:
        po.osubhd(self.outfile, "Divertor Heat Load")
        po.ovarre(
            self.outfile,
            "Power to the divertor (MW)",
            "(p_plasma_separatrix_mw.)",
            p_plasma_separatrix_mw,
        )
        po.ovarre(
            self.outfile,
            "Divertor heat load (MW/m2)",
            "(pflux_div_heat_load_mw)",
            dv.pflux_div_heat_load_mw,
        )
    return dv.pflux_div_heat_load_mw

divwade(rmajor, rminor, aspect, b_plasma_toroidal_on_axis, b_plasma_poloidal_average, p_plasma_separatrix_mw, f_div_flux_expansion, nd_plasma_separatrix_electron, deg_div_field_plate, rad_fraction_sol, f_p_div_lower, output)

Divertor heat load model (Wade 2020)

This subroutine calculates the divertor heat flux for any machine, with either a single null or double null configuration. It uses the Eich scaling (Eich et al. 2013) and spreading factor (Scarabosio et al. 2014) to calculate the SOL width. This is then used with a flux expansion factor to calculate the wetted area and then the heat flux.

Parameters:

Name Type Description Default
rmajor float

plasma major radius (m)

required
rminor float

plasma minor radius (m)

required
aspect float

tokamak aspect ratio

required
b_plasma_toroidal_on_axis float

toroidal field (T)

required
b_plasma_poloidal_average float

poloidal field (T)

required
p_plasma_separatrix_mw float

power to divertor (MW)

required
f_div_flux_expansion float

plasma flux expansion in divertor

required
nd_plasma_separatrix_electron float

electron density at separatrix (m-3)

required
deg_div_field_plate float

field line angle wrt divertor target plate (degrees)

required
rad_fraction_sol float

SOL radiation fraction

required
f_p_div_lower float

fraction of power to the lower divertor in double null configuration

required

Returns:

Type Description
float

divertor heat load for a tight aspect ratio machine

Source code in process/models/divertor.py
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
def divwade(
    self,
    rmajor: float,
    rminor: float,
    aspect: float,
    b_plasma_toroidal_on_axis: float,
    b_plasma_poloidal_average: float,
    p_plasma_separatrix_mw: float,
    f_div_flux_expansion: float,
    nd_plasma_separatrix_electron: float,
    deg_div_field_plate: float,
    rad_fraction_sol: float,
    f_p_div_lower: float,
    output: bool,
) -> float:
    """Divertor heat load model (Wade 2020)

    This subroutine calculates the divertor heat flux for any machine,
    with either a single null or double null configuration.
    It uses the Eich scaling (Eich et al. 2013) and spreading factor (Scarabosio et al. 2014)
    to calculate the SOL width. This is then used with a flux expansion factor to calculate
    the wetted area and then the heat flux.

    Parameters
    ----------
    rmajor : float
        plasma major radius (m)
    rminor : float
        plasma minor radius (m)
    aspect : float
        tokamak aspect ratio
    b_plasma_toroidal_on_axis : float
        toroidal field (T)
    b_plasma_poloidal_average : float
        poloidal field (T)
    p_plasma_separatrix_mw : float
        power to divertor (MW)
    f_div_flux_expansion : float
        plasma flux expansion in divertor
    nd_plasma_separatrix_electron : float
        electron density at separatrix (m-3)
    deg_div_field_plate : float
        field line angle wrt divertor target plate (degrees)
    rad_fraction_sol : float
        SOL radiation fraction
    f_p_div_lower : float
        fraction of power to the lower divertor in double null configuration

    Returns
    -------
    float
        divertor heat load for a tight aspect ratio machine
    """

    # Radius on midplane
    r_omp = rmajor + rminor

    # B fields on midplane
    Bp_omp = -b_plasma_poloidal_average * rmajor / r_omp

    Bt_omp = -b_plasma_toroidal_on_axis * rmajor / r_omp

    # Eich scaling for lambda_q
    lambda_eich = (
        1.35
        * p_plasma_separatrix_mw**-0.02
        * rmajor**0.04
        * b_plasma_poloidal_average**-0.92
        * aspect**0.42
    )

    # Spreading factor
    spread_fact = (
        0.12
        * (nd_plasma_separatrix_electron / 1e19) ** -0.02
        * p_plasma_separatrix_mw**-0.21
        * rmajor**0.71
        * b_plasma_poloidal_average**-0.82
    )

    # SOL width
    lambda_int = lambda_eich + 1.64 * spread_fact

    # Flux angle on midplane
    alpha_mid = math.degrees(math.atan(Bp_omp / Bt_omp))

    # Flux angle in the divertor
    alpha_div = f_div_flux_expansion * alpha_mid

    # Tilt of the separatrix relative to the target in the poloidal plane
    theta_div = math.asin(
        (1 + 1 / alpha_div**2) * math.sin(math.radians(deg_div_field_plate))
    )

    # Wetted area
    area_wetted = (
        2 * np.pi * rmajor * lambda_int * f_div_flux_expansion * math.sin(theta_div)
    )

    # Divertor heat load
    hldiv_base = p_plasma_separatrix_mw * (1 - rad_fraction_sol) / area_wetted

    # For double null, calculate heat loads to upper and lower divertors and use the highest
    if dv.n_divertors == 2:
        hldiv_lower = f_p_div_lower * hldiv_base
        hldiv_upper = (1.0 - f_p_div_lower) * hldiv_base
        dv.pflux_div_heat_load_mw = max(hldiv_lower, hldiv_upper)
    else:
        dv.pflux_div_heat_load_mw = hldiv_base

    if output:
        po.osubhd(self.outfile, "Divertor Heat Load")
        po.ocmmnt(self.outfile, "Assume an expanded divertor with a gaseous target")
        po.oblnkl(self.outfile)
        po.ovarre(
            self.outfile,
            "Flux expansion",
            "(f_div_flux_expansion)",
            f_div_flux_expansion,
        )
        po.ovarre(
            self.outfile,
            "Field line angle wrt to target divertor plate (degrees)",
            "(deg_div_field_plate)",
            deg_div_field_plate,
        )
        po.ovarre(
            self.outfile,
            "Divertor heat load (MW/m2)",
            "(pflux_div_heat_load_mw)",
            dv.pflux_div_heat_load_mw,
        )
    return dv.pflux_div_heat_load_mw

incident_radiation_power(p_plasma_rad_mw, f_ster_div_single, n_divertors)

Calculates the total incident radiation power on the divertor box.

Parameters:

Name Type Description Default
p_plasma_rad_mw float

Total plasma radiated power in megawatts (MW).

required
f_ster_div_single float

Fraction of the solid angle subtended by a single divertor.

required
n_divertors int

Number of divertors.

required

Returns:

Type Description
float

Total incident radiation power on the divertor box in megawatts (MW).

Source code in process/models/divertor.py
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
def incident_radiation_power(
    self,
    p_plasma_rad_mw: float,
    f_ster_div_single: float,
    n_divertors: int,
) -> float:
    """Calculates the total incident radiation power on the divertor box.

    Parameters
    ----------
    p_plasma_rad_mw : float
        Total plasma radiated power in megawatts (MW).
    f_ster_div_single : float
        Fraction of the solid angle subtended by a single divertor.
    n_divertors : int
        Number of divertors.

    Returns
    -------
    float
        Total incident radiation power on the divertor box in megawatts (MW).
    """

    return p_plasma_rad_mw * f_ster_div_single * n_divertors

incident_neutron_power(p_plasma_neutron_mw, f_ster_div_single, n_divertors)

Calculates the total incident neutron power on the divertor box.

Parameters:

Name Type Description Default
p_plasma_neutron_mw float

Total plasma neutron power in megawatts (MW).

required
f_ster_div_single float

Fraction of the solid angle subtended by a single divertor.

required
n_divertors int

Number of divertors.

required

Returns:

Type Description
float

Total incident radiation power on the divertor box in megawatts (MW).

Source code in process/models/divertor.py
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
def incident_neutron_power(
    self,
    p_plasma_neutron_mw: float,
    f_ster_div_single: float,
    n_divertors: int,
) -> float:
    """Calculates the total incident neutron power on the divertor box.

    Parameters
    ----------
    p_plasma_neutron_mw : float
        Total plasma neutron power in megawatts (MW).
    f_ster_div_single : float
        Fraction of the solid angle subtended by a single divertor.
    n_divertors : int
        Number of divertors.

    Returns
    -------
    float
        Total incident radiation power on the divertor box in megawatts (MW).
    """

    return p_plasma_neutron_mw * f_ster_div_single * n_divertors

LowerDivertor

Bases: Divertor

Module containing lower divertor routines

Source code in process/models/divertor.py
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
class LowerDivertor(Divertor):
    """Module containing lower divertor routines"""

    def run(self, output: bool):
        super().run(output=output)

        dv.p_div_lower_nuclear_heat_mw = self.incident_neutron_power(
            p_plasma_neutron_mw=pv.p_plasma_neutron_mw,
            f_ster_div_single=fwbs.f_ster_div_single,
            n_divertors=1,
        )

        dv.p_div_lower_rad_mw = self.incident_radiation_power(
            p_plasma_rad_mw=pv.p_plasma_rad_mw,
            f_ster_div_single=fwbs.f_ster_div_single,
            n_divertors=1,
        )

run(output)

Source code in process/models/divertor.py
413
414
415
416
417
418
419
420
421
422
423
424
425
426
def run(self, output: bool):
    super().run(output=output)

    dv.p_div_lower_nuclear_heat_mw = self.incident_neutron_power(
        p_plasma_neutron_mw=pv.p_plasma_neutron_mw,
        f_ster_div_single=fwbs.f_ster_div_single,
        n_divertors=1,
    )

    dv.p_div_lower_rad_mw = self.incident_radiation_power(
        p_plasma_rad_mw=pv.p_plasma_rad_mw,
        f_ster_div_single=fwbs.f_ster_div_single,
        n_divertors=1,
    )

UpperDivertor

Bases: Divertor

Module containing upper divertor routines

Source code in process/models/divertor.py
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
class UpperDivertor(Divertor):
    """Module containing upper divertor routines"""

    def run(self, output: bool):
        super().run(output=output)

        dv.p_div_upper_nuclear_heat_mw = self.incident_neutron_power(
            p_plasma_neutron_mw=pv.p_plasma_neutron_mw,
            f_ster_div_single=fwbs.f_ster_div_single,
            n_divertors=1,
        )

        dv.p_div_upper_rad_mw = self.incident_radiation_power(
            p_plasma_rad_mw=pv.p_plasma_rad_mw,
            f_ster_div_single=fwbs.f_ster_div_single,
            n_divertors=1,
        )

run(output)

Source code in process/models/divertor.py
432
433
434
435
436
437
438
439
440
441
442
443
444
445
def run(self, output: bool):
    super().run(output=output)

    dv.p_div_upper_nuclear_heat_mw = self.incident_neutron_power(
        p_plasma_neutron_mw=pv.p_plasma_neutron_mw,
        f_ster_div_single=fwbs.f_ster_div_single,
        n_divertors=1,
    )

    dv.p_div_upper_rad_mw = self.incident_radiation_power(
        p_plasma_rad_mw=pv.p_plasma_rad_mw,
        f_ster_div_single=fwbs.f_ster_div_single,
        n_divertors=1,
    )