Skip to content

profiles

Profile

Bases: ABC

Abstract base class used to create and hold profiles (temperature, density)

Source code in process/models/physics/profiles.py
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
class Profile(ABC):
    """Abstract base class used to create and hold profiles (temperature, density)"""

    def __init__(self, profile_size: int):
        """
        Initialize a Profiles object.

        Parameters:
        - profile_size (int): The size of the profile.

        Attributes:
        - profile_size (int): The size of the profile.
        - profile_x (ndarray): An array of values ranging from 0 to profile_size-1.
        - profile_y (ndarray): An array of zeros with length profile_size.
        - profile_integ (int): The integral of the profile_y array.
        - profile_dx (int): The step size between consecutive values in profile_x.
        """
        self.profile_size = profile_size
        self.profile_x = np.arange(self.profile_size)
        self.profile_y = np.zeros(self.profile_size)
        self.profile_integ = 0
        self.profile_dx = 0

    def normalise_profile_x(self):
        """Normalizes the x-dimension of the profile.

        This method divides the values in the `profile_x` attribute by the maximum value
        in the `profile_x` array, resulting in a normalized version of the x-dimension.

        Example:
            If `profile_x` is [1, 2, 3, 4, 5], after normalization it will become
            [0.2, 0.4, 0.6, 0.8, 1.0].

        Note:
            This method modifies the `profile_x` attribute in-place.


        """
        self.profile_x = self.profile_x / max(self.profile_x)

    def calculate_profile_dx(self):
        """Calculates the differential between points in the profile.

        This method calculates the differential between points in the profile by dividing the maximum x value in the profile
        by the difference in size between the points. The result is stored in the `profile_dx` attribute.
        """
        self.profile_dx = max(self.profile_x) / (self.profile_size - 1)

    @abstractmethod
    def calculate_profile_y(self):
        """Use a profile function to act on self.profile_x to calculate and set the
        values of self.profile_y.
        """

    def integrate_profile_y(self):
        """Integrate profile_y values using scipy.integrate.simpson() function.

        This method calculates the integral of the profile_y values using the Simpson's rule
        provided by the scipy.integrate.simpson() function. The integral is stored in the
        self.profile_integ attribute.
        """
        self.profile_integ = sp.integrate.simpson(
            self.profile_y, x=self.profile_x, dx=self.profile_dx
        )

profile_size = profile_size instance-attribute

profile_x = np.arange(self.profile_size) instance-attribute

profile_y = np.zeros(self.profile_size) instance-attribute

profile_integ = 0 instance-attribute

profile_dx = 0 instance-attribute

normalise_profile_x()

Normalizes the x-dimension of the profile.

This method divides the values in the profile_x attribute by the maximum value in the profile_x array, resulting in a normalized version of the x-dimension.

Example: If profile_x is [1, 2, 3, 4, 5], after normalization it will become [0.2, 0.4, 0.6, 0.8, 1.0].

Note: This method modifies the profile_x attribute in-place.

Source code in process/models/physics/profiles.py
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
def normalise_profile_x(self):
    """Normalizes the x-dimension of the profile.

    This method divides the values in the `profile_x` attribute by the maximum value
    in the `profile_x` array, resulting in a normalized version of the x-dimension.

    Example:
        If `profile_x` is [1, 2, 3, 4, 5], after normalization it will become
        [0.2, 0.4, 0.6, 0.8, 1.0].

    Note:
        This method modifies the `profile_x` attribute in-place.


    """
    self.profile_x = self.profile_x / max(self.profile_x)

calculate_profile_dx()

Calculates the differential between points in the profile.

This method calculates the differential between points in the profile by dividing the maximum x value in the profile by the difference in size between the points. The result is stored in the profile_dx attribute.

Source code in process/models/physics/profiles.py
52
53
54
55
56
57
58
def calculate_profile_dx(self):
    """Calculates the differential between points in the profile.

    This method calculates the differential between points in the profile by dividing the maximum x value in the profile
    by the difference in size between the points. The result is stored in the `profile_dx` attribute.
    """
    self.profile_dx = max(self.profile_x) / (self.profile_size - 1)

calculate_profile_y() abstractmethod

Use a profile function to act on self.profile_x to calculate and set the values of self.profile_y.

Source code in process/models/physics/profiles.py
60
61
62
63
64
@abstractmethod
def calculate_profile_y(self):
    """Use a profile function to act on self.profile_x to calculate and set the
    values of self.profile_y.
    """

integrate_profile_y()

Integrate profile_y values using scipy.integrate.simpson() function.

This method calculates the integral of the profile_y values using the Simpson's rule provided by the scipy.integrate.simpson() function. The integral is stored in the self.profile_integ attribute.

Source code in process/models/physics/profiles.py
66
67
68
69
70
71
72
73
74
75
def integrate_profile_y(self):
    """Integrate profile_y values using scipy.integrate.simpson() function.

    This method calculates the integral of the profile_y values using the Simpson's rule
    provided by the scipy.integrate.simpson() function. The integral is stored in the
    self.profile_integ attribute.
    """
    self.profile_integ = sp.integrate.simpson(
        self.profile_y, x=self.profile_x, dx=self.profile_dx
    )

NeProfile

Bases: Profile

Electron density profile class. Contains a function to calculate the electron density profile and store the data.

Source code in process/models/physics/profiles.py
 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
class NeProfile(Profile):
    """Electron density profile class. Contains a function to calculate the electron density profile and
    store the data.
    """

    def run(self):
        """Subroutine which calls profile functions and stores neprofile data."""
        self.normalise_profile_x()
        self.calculate_profile_dx()
        self.set_physics_variables()
        self.calculate_profile_y(
            self.profile_x,
            physics_variables.radius_plasma_pedestal_density_norm,
            physics_variables.nd_plasma_electron_on_axis,
            physics_variables.nd_plasma_pedestal_electron,
            physics_variables.nd_plasma_separatrix_electron,
            physics_variables.alphan,
        )
        self.integrate_profile_y()

    def calculate_profile_y(
        self,
        rho: np.array,
        radius_plasma_pedestal_density_norm: float,
        n0: float,
        nped: float,
        nsep: float,
        alphan: float,
    ):
        """This routine calculates the density at each normalised minor radius position
        rho for a HELIOS-type density pedestal profile (neprofile).

        Parameters
        ----------
        rho :
            Normalised minor radius vector.
        radius_plasma_pedestal_density_norm :
            Normalised minor radius pedestal position.
        n0 :
            Central density (/m3).
        nped :
            Pedestal density (/m3).
        nsep :
            Separatrix density (/m3)
        alphan :
            Density peaking parameter.
        """

        if physics_variables.i_plasma_pedestal == 0:
            self.profile_y = n0 * (1 - rho**2) ** alphan

        # Input checks

        if n0 < nped:
            logger.info(
                f"NPROFILE: density pedestal is higher than core density. {nped = }, {n0 = }"
            )
        rho_index = rho <= radius_plasma_pedestal_density_norm
        self.profile_y[rho_index] = (
            nped
            + (n0 - nped)
            * (1 - (rho[rho_index] / radius_plasma_pedestal_density_norm) ** 2) ** alphan
        )
        # Invert the rho_index
        self.profile_y[~rho_index] = nsep + (nped - nsep) * (1 - rho[~rho_index]) / (
            1 - radius_plasma_pedestal_density_norm
        )

    @staticmethod
    def ncore(
        radius_plasma_pedestal_density_norm: float,
        nped: float,
        nsep: float,
        nav: float,
        alphan: float,
    ) -> float:
        """This routine calculates the core density of a pedestalised profile.
        The solution comes from integrating and summing the two separate density profiles for the core
        and pedestal region within their bounds. This has to be multiplied by the torus volume element before integration which leads
        to an added rho term in each part of the profile. When dividing by the volume of integration to get the average density
        the simplification leads to a factor of 2 having to be multiplied on to each of the integration results.
        This function for the average density can then be re-arranged to calculate the central plasma density n_0 / ncore.

        Parameters
        ----------
        radius_plasma_pedestal_density_norm :
             The normalised minor radius pedestal position.
        nped :
            The pedestal density (/m3).
        nsep :
            The separatrix density (/m3).
        nav :
            The electron density (/m3).
        alphan :
            The density peaking parameter

        Returns
        -------
        :
            The core density.

        References:
            Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650
        """

        ncore = (
            1
            / (3 * radius_plasma_pedestal_density_norm**2)
            * (
                3 * nav * (1 + alphan)
                + nsep
                * (1 + alphan)
                * (
                    -2
                    + radius_plasma_pedestal_density_norm
                    + radius_plasma_pedestal_density_norm**2
                )
                - nped
                * (
                    (1 + alphan) * (1 + radius_plasma_pedestal_density_norm)
                    + (alphan - 2) * radius_plasma_pedestal_density_norm**2
                )
            )
        )

        if ncore < 0.0:
            # Allows solver to continue and
            # warns the user to raise the lower bound on nd_plasma_electrons_vol_avg if the run did not converge
            logger.error(
                "ncore is going negative when solving. Please raise the value of nd_plasma_electrons_vol_avg and or its lower limit."
            )
            ncore = 1.0e-6
        return ncore

    def set_physics_variables(self):
        """Calculates and sets physics variables required for the profile."""

        if physics_variables.i_plasma_pedestal == 0:
            physics_variables.nd_plasma_electron_on_axis = (
                physics_variables.nd_plasma_electrons_vol_avg
                * (1.0 + physics_variables.alphan)
            )
        elif physics_variables.i_plasma_pedestal == 1:
            physics_variables.nd_plasma_electron_on_axis = self.ncore(
                physics_variables.radius_plasma_pedestal_density_norm,
                physics_variables.nd_plasma_pedestal_electron,
                physics_variables.nd_plasma_separatrix_electron,
                physics_variables.nd_plasma_electrons_vol_avg,
                physics_variables.alphan,
            )
        physics_variables.nd_plasma_ions_on_axis = (
            physics_variables.nd_plasma_ions_total_vol_avg
            / physics_variables.nd_plasma_electrons_vol_avg
            * physics_variables.nd_plasma_electron_on_axis
        )

run()

Subroutine which calls profile functions and stores neprofile data.

Source code in process/models/physics/profiles.py
83
84
85
86
87
88
89
90
91
92
93
94
95
96
def run(self):
    """Subroutine which calls profile functions and stores neprofile data."""
    self.normalise_profile_x()
    self.calculate_profile_dx()
    self.set_physics_variables()
    self.calculate_profile_y(
        self.profile_x,
        physics_variables.radius_plasma_pedestal_density_norm,
        physics_variables.nd_plasma_electron_on_axis,
        physics_variables.nd_plasma_pedestal_electron,
        physics_variables.nd_plasma_separatrix_electron,
        physics_variables.alphan,
    )
    self.integrate_profile_y()

calculate_profile_y(rho, radius_plasma_pedestal_density_norm, n0, nped, nsep, alphan)

This routine calculates the density at each normalised minor radius position rho for a HELIOS-type density pedestal profile (neprofile).

Parameters:

Name Type Description Default
rho array

Normalised minor radius vector.

required
radius_plasma_pedestal_density_norm float

Normalised minor radius pedestal position.

required
n0 float

Central density (/m3).

required
nped float

Pedestal density (/m3).

required
nsep float

Separatrix density (/m3)

required
alphan float

Density peaking parameter.

required
Source code in process/models/physics/profiles.py
 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
def calculate_profile_y(
    self,
    rho: np.array,
    radius_plasma_pedestal_density_norm: float,
    n0: float,
    nped: float,
    nsep: float,
    alphan: float,
):
    """This routine calculates the density at each normalised minor radius position
    rho for a HELIOS-type density pedestal profile (neprofile).

    Parameters
    ----------
    rho :
        Normalised minor radius vector.
    radius_plasma_pedestal_density_norm :
        Normalised minor radius pedestal position.
    n0 :
        Central density (/m3).
    nped :
        Pedestal density (/m3).
    nsep :
        Separatrix density (/m3)
    alphan :
        Density peaking parameter.
    """

    if physics_variables.i_plasma_pedestal == 0:
        self.profile_y = n0 * (1 - rho**2) ** alphan

    # Input checks

    if n0 < nped:
        logger.info(
            f"NPROFILE: density pedestal is higher than core density. {nped = }, {n0 = }"
        )
    rho_index = rho <= radius_plasma_pedestal_density_norm
    self.profile_y[rho_index] = (
        nped
        + (n0 - nped)
        * (1 - (rho[rho_index] / radius_plasma_pedestal_density_norm) ** 2) ** alphan
    )
    # Invert the rho_index
    self.profile_y[~rho_index] = nsep + (nped - nsep) * (1 - rho[~rho_index]) / (
        1 - radius_plasma_pedestal_density_norm
    )

ncore(radius_plasma_pedestal_density_norm, nped, nsep, nav, alphan) staticmethod

This routine calculates the core density of a pedestalised profile. The solution comes from integrating and summing the two separate density profiles for the core and pedestal region within their bounds. This has to be multiplied by the torus volume element before integration which leads to an added rho term in each part of the profile. When dividing by the volume of integration to get the average density the simplification leads to a factor of 2 having to be multiplied on to each of the integration results. This function for the average density can then be re-arranged to calculate the central plasma density n_0 / ncore.

Parameters:

Name Type Description Default
radius_plasma_pedestal_density_norm float

The normalised minor radius pedestal position.

required
nped float

The pedestal density (/m3).

required
nsep float

The separatrix density (/m3).

required
nav float

The electron density (/m3).

required
alphan float

The density peaking parameter

required

Returns:

Name Type Description
float

The core density.

References float

Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650

Source code in process/models/physics/profiles.py
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
@staticmethod
def ncore(
    radius_plasma_pedestal_density_norm: float,
    nped: float,
    nsep: float,
    nav: float,
    alphan: float,
) -> float:
    """This routine calculates the core density of a pedestalised profile.
    The solution comes from integrating and summing the two separate density profiles for the core
    and pedestal region within their bounds. This has to be multiplied by the torus volume element before integration which leads
    to an added rho term in each part of the profile. When dividing by the volume of integration to get the average density
    the simplification leads to a factor of 2 having to be multiplied on to each of the integration results.
    This function for the average density can then be re-arranged to calculate the central plasma density n_0 / ncore.

    Parameters
    ----------
    radius_plasma_pedestal_density_norm :
         The normalised minor radius pedestal position.
    nped :
        The pedestal density (/m3).
    nsep :
        The separatrix density (/m3).
    nav :
        The electron density (/m3).
    alphan :
        The density peaking parameter

    Returns
    -------
    :
        The core density.

    References:
        Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650
    """

    ncore = (
        1
        / (3 * radius_plasma_pedestal_density_norm**2)
        * (
            3 * nav * (1 + alphan)
            + nsep
            * (1 + alphan)
            * (
                -2
                + radius_plasma_pedestal_density_norm
                + radius_plasma_pedestal_density_norm**2
            )
            - nped
            * (
                (1 + alphan) * (1 + radius_plasma_pedestal_density_norm)
                + (alphan - 2) * radius_plasma_pedestal_density_norm**2
            )
        )
    )

    if ncore < 0.0:
        # Allows solver to continue and
        # warns the user to raise the lower bound on nd_plasma_electrons_vol_avg if the run did not converge
        logger.error(
            "ncore is going negative when solving. Please raise the value of nd_plasma_electrons_vol_avg and or its lower limit."
        )
        ncore = 1.0e-6
    return ncore

set_physics_variables()

Calculates and sets physics variables required for the profile.

Source code in process/models/physics/profiles.py
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
def set_physics_variables(self):
    """Calculates and sets physics variables required for the profile."""

    if physics_variables.i_plasma_pedestal == 0:
        physics_variables.nd_plasma_electron_on_axis = (
            physics_variables.nd_plasma_electrons_vol_avg
            * (1.0 + physics_variables.alphan)
        )
    elif physics_variables.i_plasma_pedestal == 1:
        physics_variables.nd_plasma_electron_on_axis = self.ncore(
            physics_variables.radius_plasma_pedestal_density_norm,
            physics_variables.nd_plasma_pedestal_electron,
            physics_variables.nd_plasma_separatrix_electron,
            physics_variables.nd_plasma_electrons_vol_avg,
            physics_variables.alphan,
        )
    physics_variables.nd_plasma_ions_on_axis = (
        physics_variables.nd_plasma_ions_total_vol_avg
        / physics_variables.nd_plasma_electrons_vol_avg
        * physics_variables.nd_plasma_electron_on_axis
    )

TeProfile

Bases: Profile

Electron temperature profile class. Contains a function to calculate the temperature profile and store the data.

Source code in process/models/physics/profiles.py
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
class TeProfile(Profile):
    """Electron temperature profile class. Contains a function to calculate the temperature profile and store the data."""

    def run(self):
        """Subroutine to initialise neprofile and execute calculations."""
        self.normalise_profile_x()
        self.calculate_profile_dx()
        self.set_physics_variables()
        self.calculate_profile_y(
            self.profile_x,
            physics_variables.radius_plasma_pedestal_temp_norm,
            physics_variables.temp_plasma_electron_on_axis_kev,
            physics_variables.temp_plasma_pedestal_kev,
            physics_variables.temp_plasma_separatrix_kev,
            physics_variables.alphat,
            physics_variables.tbeta,
        )
        self.integrate_profile_y()

    def calculate_profile_y(
        self,
        rho: np.array,
        radius_plasma_pedestal_temp_norm: float,
        t0: float,
        temp_plasma_pedestal_kev: float,
        temp_plasma_separatrix_kev: float,
        alphat: float,
        tbeta: float,
    ):
        """Calculates the temperature at a normalised minor radius position rho for a pedestalised profile (teprofile).
        If i_plasma_pedestal = 0 the original parabolic profile form is used instead.

        Parameters
        ----------
        rho : np.array
            Normalised minor radius.
        radius_plasma_pedestal_temp_norm : float
            Normalised minor radius pedestal position.
        t0 : float
            Central temperature (keV).
        temp_plasma_pedestal_kev : float
            Pedestal temperature (keV).
        temp_plasma_separatrix_kev : float
            Separatrix temperature (keV).
        alphat : float
            Temperature peaking parameter.
        tbeta : float
            Second temperature exponent.

        References:
            Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650
        """
        if physics_variables.i_plasma_pedestal == 0:
            # profile values of 0 cause divide by 0 errors so ensure the profile value is at least 1e-8
            # which is small enough that it won't make a difference to any calculations
            self.profile_y = np.maximum(t0 * (1 - rho**2) ** alphat, 1e-8)
            return

        if t0 < temp_plasma_pedestal_kev:
            logger.info(
                f"TPROFILE: temperature pedestal is higher than core temperature. {temp_plasma_pedestal_kev = }, {t0 = }"
            )

        rho_index = rho <= radius_plasma_pedestal_temp_norm
        self.profile_y[rho_index] = (
            temp_plasma_pedestal_kev
            + (t0 - temp_plasma_pedestal_kev)
            * (1 - (rho[rho_index] / radius_plasma_pedestal_temp_norm) ** tbeta)
            ** alphat
        )
        self.profile_y[~rho_index] = temp_plasma_separatrix_kev + (
            temp_plasma_pedestal_kev - temp_plasma_separatrix_kev
        ) * (1 - rho[~rho_index]) / (1 - radius_plasma_pedestal_temp_norm)

    @staticmethod
    def tcore(
        radius_plasma_pedestal_temp_norm: float,
        temp_plasma_pedestal_kev: float,
        temp_plasma_separatrix_kev: float,
        tav: float,
        alphat: float,
        tbeta: float,
    ) -> float:
        """This routine calculates the core temperature (keV)
        of a pedestalised profile. The solution comes from integrating and summing the two seprate temperature profiles for the core
        and pedestal region within their bounds. This has to be multiplied by the torus volume element before integration which leads
        to an added rho term in each part of the profile. When dividing by the volume of integration to get the average temperature
        the simplification leads to a factor of 2 having to be multiplied on to each of the integration results.
        This function for the average temperature can then be re-arranged to calculate the central plasma temeprature T_0 / tcore.
        References:
            Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650

        Parameters
        ----------
        radius_plasma_pedestal_temp_norm : float
            Normalised minor radius pedestal position.
        temp_plasma_pedestal_kev : float
            Pedestal temperature (keV).
        temp_plasma_separatrix_kev : float
            Separatrix temperature (keV).
        tav : float
            Volume average temperature (keV).
        alphat : float
            Temperature peaking parameter.
        tbeta : float
            Second temperature exponent.

        Returns
        -------
        float
            Core temperature.
        """
        #  Calculate core temperature

        return temp_plasma_pedestal_kev + (
            (
                tbeta
                * (
                    3 * tav
                    + temp_plasma_separatrix_kev
                    * (
                        -2.0
                        + radius_plasma_pedestal_temp_norm
                        + radius_plasma_pedestal_temp_norm**2
                    )
                    - temp_plasma_pedestal_kev
                    * (
                        1
                        + radius_plasma_pedestal_temp_norm
                        + radius_plasma_pedestal_temp_norm**2
                    )
                )
            )
            / (
                6
                * radius_plasma_pedestal_temp_norm**2
                * sp.special.beta(1 + alphat, 2 / tbeta)
            )
        )

    def set_physics_variables(self):
        """Calculates and sets physics variables required for the temperature profile."""
        if physics_variables.i_plasma_pedestal == 0:
            physics_variables.temp_plasma_electron_on_axis_kev = (
                physics_variables.temp_plasma_electron_vol_avg_kev
                * (1.0 + physics_variables.alphat)
            )
        elif physics_variables.i_plasma_pedestal == 1:
            physics_variables.temp_plasma_electron_on_axis_kev = self.tcore(
                physics_variables.radius_plasma_pedestal_temp_norm,
                physics_variables.temp_plasma_pedestal_kev,
                physics_variables.temp_plasma_separatrix_kev,
                physics_variables.temp_plasma_electron_vol_avg_kev,
                physics_variables.alphat,
                physics_variables.tbeta,
            )

        physics_variables.temp_plasma_ion_on_axis_kev = (
            physics_variables.temp_plasma_ion_vol_avg_kev
            / physics_variables.temp_plasma_electron_vol_avg_kev
            * physics_variables.temp_plasma_electron_on_axis_kev
        )

run()

Subroutine to initialise neprofile and execute calculations.

Source code in process/models/physics/profiles.py
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
def run(self):
    """Subroutine to initialise neprofile and execute calculations."""
    self.normalise_profile_x()
    self.calculate_profile_dx()
    self.set_physics_variables()
    self.calculate_profile_y(
        self.profile_x,
        physics_variables.radius_plasma_pedestal_temp_norm,
        physics_variables.temp_plasma_electron_on_axis_kev,
        physics_variables.temp_plasma_pedestal_kev,
        physics_variables.temp_plasma_separatrix_kev,
        physics_variables.alphat,
        physics_variables.tbeta,
    )
    self.integrate_profile_y()

calculate_profile_y(rho, radius_plasma_pedestal_temp_norm, t0, temp_plasma_pedestal_kev, temp_plasma_separatrix_kev, alphat, tbeta)

Calculates the temperature at a normalised minor radius position rho for a pedestalised profile (teprofile). If i_plasma_pedestal = 0 the original parabolic profile form is used instead.

Parameters:

Name Type Description Default
rho array

Normalised minor radius.

required
radius_plasma_pedestal_temp_norm float

Normalised minor radius pedestal position.

required
t0 float

Central temperature (keV).

required
temp_plasma_pedestal_kev float

Pedestal temperature (keV).

required
temp_plasma_separatrix_kev float

Separatrix temperature (keV).

required
alphat float

Temperature peaking parameter.

required
tbeta float

Second temperature exponent.

required
References

Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650

required
Source code in process/models/physics/profiles.py
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
def calculate_profile_y(
    self,
    rho: np.array,
    radius_plasma_pedestal_temp_norm: float,
    t0: float,
    temp_plasma_pedestal_kev: float,
    temp_plasma_separatrix_kev: float,
    alphat: float,
    tbeta: float,
):
    """Calculates the temperature at a normalised minor radius position rho for a pedestalised profile (teprofile).
    If i_plasma_pedestal = 0 the original parabolic profile form is used instead.

    Parameters
    ----------
    rho : np.array
        Normalised minor radius.
    radius_plasma_pedestal_temp_norm : float
        Normalised minor radius pedestal position.
    t0 : float
        Central temperature (keV).
    temp_plasma_pedestal_kev : float
        Pedestal temperature (keV).
    temp_plasma_separatrix_kev : float
        Separatrix temperature (keV).
    alphat : float
        Temperature peaking parameter.
    tbeta : float
        Second temperature exponent.

    References:
        Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650
    """
    if physics_variables.i_plasma_pedestal == 0:
        # profile values of 0 cause divide by 0 errors so ensure the profile value is at least 1e-8
        # which is small enough that it won't make a difference to any calculations
        self.profile_y = np.maximum(t0 * (1 - rho**2) ** alphat, 1e-8)
        return

    if t0 < temp_plasma_pedestal_kev:
        logger.info(
            f"TPROFILE: temperature pedestal is higher than core temperature. {temp_plasma_pedestal_kev = }, {t0 = }"
        )

    rho_index = rho <= radius_plasma_pedestal_temp_norm
    self.profile_y[rho_index] = (
        temp_plasma_pedestal_kev
        + (t0 - temp_plasma_pedestal_kev)
        * (1 - (rho[rho_index] / radius_plasma_pedestal_temp_norm) ** tbeta)
        ** alphat
    )
    self.profile_y[~rho_index] = temp_plasma_separatrix_kev + (
        temp_plasma_pedestal_kev - temp_plasma_separatrix_kev
    ) * (1 - rho[~rho_index]) / (1 - radius_plasma_pedestal_temp_norm)

tcore(radius_plasma_pedestal_temp_norm, temp_plasma_pedestal_kev, temp_plasma_separatrix_kev, tav, alphat, tbeta) staticmethod

This routine calculates the core temperature (keV) of a pedestalised profile. The solution comes from integrating and summing the two seprate temperature profiles for the core and pedestal region within their bounds. This has to be multiplied by the torus volume element before integration which leads to an added rho term in each part of the profile. When dividing by the volume of integration to get the average temperature the simplification leads to a factor of 2 having to be multiplied on to each of the integration results. This function for the average temperature can then be re-arranged to calculate the central plasma temeprature T_0 / tcore. References: Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650

Parameters:

Name Type Description Default
radius_plasma_pedestal_temp_norm float

Normalised minor radius pedestal position.

required
temp_plasma_pedestal_kev float

Pedestal temperature (keV).

required
temp_plasma_separatrix_kev float

Separatrix temperature (keV).

required
tav float

Volume average temperature (keV).

required
alphat float

Temperature peaking parameter.

required
tbeta float

Second temperature exponent.

required

Returns:

Type Description
float

Core temperature.

Source code in process/models/physics/profiles.py
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
@staticmethod
def tcore(
    radius_plasma_pedestal_temp_norm: float,
    temp_plasma_pedestal_kev: float,
    temp_plasma_separatrix_kev: float,
    tav: float,
    alphat: float,
    tbeta: float,
) -> float:
    """This routine calculates the core temperature (keV)
    of a pedestalised profile. The solution comes from integrating and summing the two seprate temperature profiles for the core
    and pedestal region within their bounds. This has to be multiplied by the torus volume element before integration which leads
    to an added rho term in each part of the profile. When dividing by the volume of integration to get the average temperature
    the simplification leads to a factor of 2 having to be multiplied on to each of the integration results.
    This function for the average temperature can then be re-arranged to calculate the central plasma temeprature T_0 / tcore.
    References:
        Jean, J. (2011). HELIOS: A Zero-Dimensional Tool for Next Step and Reactor Studies. Fusion Science and Technology, 59(2), 308-349. https://doi.org/10.13182/FST11-A11650

    Parameters
    ----------
    radius_plasma_pedestal_temp_norm : float
        Normalised minor radius pedestal position.
    temp_plasma_pedestal_kev : float
        Pedestal temperature (keV).
    temp_plasma_separatrix_kev : float
        Separatrix temperature (keV).
    tav : float
        Volume average temperature (keV).
    alphat : float
        Temperature peaking parameter.
    tbeta : float
        Second temperature exponent.

    Returns
    -------
    float
        Core temperature.
    """
    #  Calculate core temperature

    return temp_plasma_pedestal_kev + (
        (
            tbeta
            * (
                3 * tav
                + temp_plasma_separatrix_kev
                * (
                    -2.0
                    + radius_plasma_pedestal_temp_norm
                    + radius_plasma_pedestal_temp_norm**2
                )
                - temp_plasma_pedestal_kev
                * (
                    1
                    + radius_plasma_pedestal_temp_norm
                    + radius_plasma_pedestal_temp_norm**2
                )
            )
        )
        / (
            6
            * radius_plasma_pedestal_temp_norm**2
            * sp.special.beta(1 + alphat, 2 / tbeta)
        )
    )

set_physics_variables()

Calculates and sets physics variables required for the temperature profile.

Source code in process/models/physics/profiles.py
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
def set_physics_variables(self):
    """Calculates and sets physics variables required for the temperature profile."""
    if physics_variables.i_plasma_pedestal == 0:
        physics_variables.temp_plasma_electron_on_axis_kev = (
            physics_variables.temp_plasma_electron_vol_avg_kev
            * (1.0 + physics_variables.alphat)
        )
    elif physics_variables.i_plasma_pedestal == 1:
        physics_variables.temp_plasma_electron_on_axis_kev = self.tcore(
            physics_variables.radius_plasma_pedestal_temp_norm,
            physics_variables.temp_plasma_pedestal_kev,
            physics_variables.temp_plasma_separatrix_kev,
            physics_variables.temp_plasma_electron_vol_avg_kev,
            physics_variables.alphat,
            physics_variables.tbeta,
        )

    physics_variables.temp_plasma_ion_on_axis_kev = (
        physics_variables.temp_plasma_ion_vol_avg_kev
        / physics_variables.temp_plasma_electron_vol_avg_kev
        * physics_variables.temp_plasma_electron_on_axis_kev
    )