Skip to content

superconductors

jcrit_rebco(temp_conductor, b_conductor)

Calculate the critical current density for a "REBCO" 2nd generation HTS superconductor.

Parameters:

Name Type Description Default
temp_conductor float

Superconductor temperature in Kelvin (K).

required
b_conductor float

Magnetic field at the superconductor in Tesla (T).

required

Returns:

Type Description
tuple[float, bool]

A tuple containing: - j_critical: Critical current density in the superconductor (A/m²). - validity: A boolean indicating whether the input parameters are within the valid range.

Notes
  • Validity range:
    • Temperature: 4.2 K ≤ temp_conductor ≤ 72.0 K
    • Magnetic field:
      • For temp_conductor < 65 K: 0.0 T ≤ b_conductor ≤ 15.0 T
      • For temp_conductor ≥ 65 K: 0.0 T ≤ b_conductor ≤ 11.5 T
Source code in process/models/superconductors.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
76
77
78
79
80
81
82
83
84
85
86
87
def jcrit_rebco(temp_conductor: float, b_conductor: float) -> tuple[float, bool]:
    """Calculate the critical current density for a "REBCO" 2nd generation HTS superconductor.

    Parameters
    ----------
    temp_conductor : float
        Superconductor temperature in Kelvin (K).
    b_conductor : float
        Magnetic field at the superconductor in Tesla (T).

    Returns
    -------
    tuple[float, bool]
        A tuple containing:
        - j_critical: Critical current density in the superconductor (A/m²).
        - validity: A boolean indicating whether the input parameters are within the valid range.

    Notes
    -----
    - Validity range:
        - Temperature: 4.2 K ≤ temp_conductor ≤ 72.0 K
        - Magnetic field:
            - For temp_conductor < 65 K: 0.0 T ≤ b_conductor ≤ 15.0 T
            - For temp_conductor ≥ 65 K: 0.0 T ≤ b_conductor ≤ 11.5 T

    """
    # Critical temperature (K) at zero field and strain.
    temp_c0max = 90.0
    # Upper critical field (T) for the superconductor at zero temperature and strain.
    b_c20max = 132.5

    C = 1.82962e8  # scaling constant
    p = 0.5875
    q = 1.7
    alpha = 1.54121
    beta = 1.96679
    oneoveralpha = 1 / alpha

    validity = True

    if (temp_conductor < 4.2) or (temp_conductor > 72.0):
        validity = False
    if temp_conductor < 65:
        if (b_conductor < 0.0) or (b_conductor > 15.0):
            validity = False
    else:
        if (b_conductor < 0.0) or (b_conductor > 11.5):
            validity = False

    if not validity:
        logger.error(
            f"""jcrit_rebco: input out of range
            temperature: {temp_conductor}
            Field: {b_conductor}
            """
        )

    if temp_conductor < temp_c0max:
        # Normal case
        birr = b_c20max * (1 - temp_conductor / temp_c0max) ** alpha
    else:
        # If temp is greater than critical temp, ensure result is real but negative.
        birr = b_c20max * (1 - temp_conductor / temp_c0max)

    if b_conductor < birr:
        # Normal case
        factor = (b_conductor / birr) ** p * (1 - b_conductor / birr) ** q
        j_critical = (C / b_conductor) * (birr**beta) * factor
    else:
        # Field is too high
        # Ensure result is real but negative, and varies with temperature.
        # tcb = critical temperature at field b
        tcb = temp_c0max * (1 - (b_conductor / b_c20max) ** oneoveralpha)
        j_critical = -(temp_conductor - tcb)

    return j_critical, validity

current_sharing_rebco(bfield, j)

Current sharing temperature for "REBCO" 2nd generation HTS superconductor

Parameters:

Name Type Description Default
bfield

Magnetic field at superconductor (T)

required
j

Current density in superconductor (A/m2)

required

Returns:

Name Type Description
current_sharing_t

Current sharing temperature (K)

Source code in process/models/superconductors.py
 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
def current_sharing_rebco(bfield, j):
    """Current sharing temperature for "REBCO" 2nd generation HTS superconductor

    Parameters
    ----------
    bfield :
        Magnetic field at superconductor (T)
    j :
        Current density in superconductor (A/m2)

    Returns
    -------
    current_sharing_t :
        Current sharing temperature (K)
    """

    def deltaj_rebco(temperature):
        jcritical, _ = jcrit_rebco(temperature, bfield)
        return jcritical - j

    # No additional arguments are required for deltaj_rebco since it only has one argument.

    estimate = 10.0
    another_estimate = 20.0
    current_sharing_t, _root_result = optimize.newton(
        deltaj_rebco,
        estimate,
        tol=1e-6,
        rtol=1e-6,
        maxiter=50,
        x1=another_estimate,
        full_output=True,
        disp=True,
    )

    return current_sharing_t

itersc(temp_conductor, b_conductor, strain, b_c20max, temp_c0max)

Calculate the critical current density, critical field, and critical temperature for an ITER Nb₃Sn superconductor using the ITER Nb₃Sn critical surface model.

Parameters:

Name Type Description Default
temp_conductor float

Superconductor temperature (K).

required
b_conductor float

Magnetic field at the conductor (T).

required
strain float

Strain in the superconductor.

required
b_c20max float

Upper critical field (T) for the superconductor at zero temperature and strain.

required
temp_c0max float

Critical temperature (K) at zero field and strain.

required

Returns:

Type Description
tuple[float, float, float]

A tuple containing: - j_critical: Critical current density in the superconductor (A/m²). - b_critical: Critical field (T). - temp_critical: Critical temperature (K).

Notes
  • This routine uses the ITER Nb₃Sn critical surface model.
  • The model assumes a strand size of 0.82 mm diameter.
References
Source code in process/models/superconductors.py
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
def itersc(
    temp_conductor: float,
    b_conductor: float,
    strain: float,
    b_c20max: float,
    temp_c0max: float,
) -> tuple[float, float, float]:
    """Calculate the critical current density, critical field, and critical temperature
    for an ITER Nb₃Sn superconductor using the ITER Nb₃Sn critical surface model.

    Parameters
    ----------
    temp_conductor : float
        Superconductor temperature (K).
    b_conductor : float
        Magnetic field at the conductor (T).
    strain : float
        Strain in the superconductor.
    b_c20max : float
        Upper critical field (T) for the superconductor at zero temperature and strain.
    temp_c0max : float
        Critical temperature (K) at zero field and strain.

    Returns
    -------
    tuple[float, float, float]
        A tuple containing:
        - j_critical: Critical current density in the superconductor (A/m²).
        - b_critical: Critical field (T).
        - temp_critical: Critical temperature (K).

    Notes
    -----
    - This routine uses the ITER Nb₃Sn critical surface model.
    - The model assumes a strand size of 0.82 mm diameter.

    References
    ----------
    - ITER DDD 11-7: Magnets - conductors (2NBKXY) (2009),
      https://user.iter.org/?uid=2NBKXY&action=get_document
    """
    # Scaling constant C [AT/mm²]
    csc = 19922.0
    # Low field exponent
    p = 0.63
    # High field exponent
    q = 2.1
    # Strain fitting constant C_{a1}
    ca1 = 44.48
    # Strain fitting constant C_{a2}
    ca2 = 0.0
    # Residual strain component epsilon_{0,a}
    epsilon_0a = 0.00256

    # ITER strand diameter (mm)
    diter = 0.82

    # ITER strand copper fraction
    f_a_strand_copper = 0.5

    j_scaling, b_critical, temp_critical = bottura_scaling(
        csc=csc,
        p=p,
        q=q,
        c_a1=ca1,
        c_a2=ca2,
        epsilon_0a=epsilon_0a,
        temp_conductor=temp_conductor,
        b_conductor=b_conductor,
        epsilon=strain,
        b_c20max=b_c20max,
        temp_c0max=temp_c0max,
    )

    # Critical current density in superconductor (A/m²)
    # ITER parameters are for the current in a single strand,
    # not per unit area, so scalefac converts to current density.
    # Convert from mm² to m².
    scalefac = np.pi * (0.5 * diter) ** 2 * (1.0 - f_a_strand_copper) / 1.0e6

    j_critical = j_scaling / scalefac

    return j_critical, b_critical, temp_critical

jcrit_nbti(temp_conductor, b_conductor, c0, b_c20max, temp_c0max)

Calculate the critical current density and critical temperature for a NbTi superconductor strand using the old empirical Lubell scaling law.

Parameters:

Name Type Description Default
temp_conductor float

Superconductor temperature (K).

required
b_conductor float

Magnetic field at the conductor (T).

required
c0 float

Scaling constant (A/m²).

required
b_c20max float

Upper critical field (T) for the superconductor at zero temperature and strain.

required
temp_c0max float

Critical temperature (K) at zero field and strain.

required

Returns:

Type Description
tuple[float, float]

A tuple containing: - jcrit: Critical current density in the superconductor (A/m²). - tcrit: Critical temperature (K).

Notes
- If the magnetic field exceeds the upper critical field (bmax > bc20max),
  the critical temperature is adjusted to ensure a real (negative) value.
- If the temperature exceeds the critical temperature, the critical surface
  is considered exceeded, and the reduced temperature (tbar) becomes negative.

- This model uses an antiquated scaling law for NbTi, which is not used in the superconductor field for nearly 30 years.
- It is simplistic and linear in J_c (B) (for a fixed temperature), lacking accuracy in both the high and low field regions.
References
- M. Lubell, Empirical scaling formulas for critical current and critical field for commercial NbTi,”
IEEE Transactions on Magnetics, vol. 19, no. 3, pp. 754-757, May 1983,
doi: https://doi.org/10.1109/tmag.1983.1062311.
Source code in process/models/superconductors.py
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
def jcrit_nbti(
    temp_conductor: float,
    b_conductor: float,
    c0: float,
    b_c20max: float,
    temp_c0max: float,
) -> tuple[float, float]:
    """Calculate the critical current density and critical temperature for a NbTi superconductor strand using the old empirical
        Lubell scaling law.

    Parameters
    ----------
    temp_conductor : float
        Superconductor temperature (K).
    b_conductor : float
        Magnetic field at the conductor (T).
    c0 : float
        Scaling constant (A/m²).
    b_c20max : float
        Upper critical field (T) for the superconductor at zero temperature and strain.
    temp_c0max : float
        Critical temperature (K) at zero field and strain.

    Returns
    -------
    tuple[float, float]
        A tuple containing:
        - jcrit: Critical current density in the superconductor (A/m²).
        - tcrit: Critical temperature (K).

    Notes
    -----
        - If the magnetic field exceeds the upper critical field (bmax > bc20max),
          the critical temperature is adjusted to ensure a real (negative) value.
        - If the temperature exceeds the critical temperature, the critical surface
          is considered exceeded, and the reduced temperature (tbar) becomes negative.

        - This model uses an antiquated scaling law for NbTi, which is not used in the superconductor field for nearly 30 years.
        - It is simplistic and linear in J_c (B) (for a fixed temperature), lacking accuracy in both the high and low field regions.

    References
    ----------
        - M. Lubell, “Empirical scaling formulas for critical current and critical field for commercial NbTi,”
        IEEE Transactions on Magnetics, vol. 19, no. 3, pp. 754-757, May 1983,
        doi: https://doi.org/10.1109/tmag.1983.1062311.
    """
    bratio = b_conductor / b_c20max

    # Critical temperature (K) at field
    temp_critical = (
        (temp_c0max * (1.0 - bratio) ** 0.59)
        if bratio < 1
        else (temp_c0max * (1.0 - bratio))
    )

    # Reduced temperature
    tbar = 1.0 - temp_conductor / temp_critical

    # Critical current density (A/m²)
    j_critical = c0 * (1.0 - bratio) * tbar

    return j_critical, temp_critical

bi2212(b_conductor, jstrand, temp_conductor, f_strain)

Fitted parameterization to Bi-2212 superconductor properties.

This function calculates the critical current density and the temperature margin
for Bi-2212 superconductor in the TF coils using a fit by M. Kovari to measurements
described in the reference, specifically from the points shown in Figure 6.

Bi-2212 (Bi₂Sr₂CaCu₂O₈₋ₓ) is a first-generation high-temperature superconductor.
It needs to be operated below about 10K but remains superconducting at much higher
fields at that temperature than Nb₃Sn, etc.

Parameters:

Name Type Description Default
b_conductor float

Magnetic field at conductor (T).

required
jstrand float

Current density in strand (A/m²).

required
temp_conductor float

Superconductor temperature (K).

required
f_strain float

Adjustment factor (≤ 1) to account for strain, radiation damage, fatigue, or AC losses.

required

Returns:

Type Description
tuple[float, float]

A tuple containing: - j_critical: Critical current density in strand (A/m²). - temp_margin: Temperature margin (K).

Notes
-The model's range of validity is:
    T < 20K
    Adjusted field b < 104 T
    B > 6 T
References
- D. C. Larbalestier, J. Jiang, U. P. Trociewitz, F. Kametani, and E. E. Hellstrom,
A transformative superconducting magnet technology for fields well above 30 T using isotropic round wire multifilament Bi2Sr2CaCu2O8-x conductor,
May 06, 2013. https://www.researchgate.net/publication/236627864_A_transformative_superconducting_magnet_technology_for_fields_well_above_30_T_using_isotropic_round_wire_multifilament_Bi2Sr2CaCu2O8-x_conductor

Raises:

Type Description
ProcessValueError

If the input parameters are outside the range of validity.

Source code in process/models/superconductors.py
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
def bi2212(b_conductor, jstrand, temp_conductor, f_strain):
    """Fitted parameterization to Bi-2212 superconductor properties.

        This function calculates the critical current density and the temperature margin
        for Bi-2212 superconductor in the TF coils using a fit by M. Kovari to measurements
        described in the reference, specifically from the points shown in Figure 6.

        Bi-2212 (Bi₂Sr₂CaCu₂O₈₋ₓ) is a first-generation high-temperature superconductor.
        It needs to be operated below about 10K but remains superconducting at much higher
        fields at that temperature than Nb₃Sn, etc.

    Parameters
    ----------
    b_conductor : float
        Magnetic field at conductor (T).
    jstrand : float
        Current density in strand (A/m²).
    temp_conductor : float
        Superconductor temperature (K).
    f_strain : float
        Adjustment factor (≤ 1) to account for strain, radiation damage,
        fatigue, or AC losses.

    Returns
    -------
    tuple[float, float]
        A tuple containing:
        - j_critical: Critical current density in strand (A/m²).
        - temp_margin: Temperature margin (K).

    Notes
    -----
        -The model's range of validity is:
            T < 20K
            Adjusted field b < 104 T
            B > 6 T

    References
    -----------
        - D. C. Larbalestier, J. Jiang, U. P. Trociewitz, F. Kametani, and E. E. Hellstrom,
        “A transformative superconducting magnet technology for fields well above 30 T using isotropic round wire multifilament Bi2Sr2CaCu2O8-x conductor,”
        May 06, 2013. https://www.researchgate.net/publication/236627864_A_transformative_superconducting_magnet_technology_for_fields_well_above_30_T_using_isotropic_round_wire_multifilament_Bi2Sr2CaCu2O8-x_conductor

    Raises
    ------
    ProcessValueError
        If the input parameters are outside the range of validity.
    """

    b = b_conductor / np.exp(-0.168 * (temp_conductor - 4.2))

    #  Engineering (i.e. strand) critical current density (A/m2)

    j_critical = f_strain * (1.175e9 * np.exp(-0.02115 * b) - 1.288e8)

    #  Temperature margin (K)
    #  Simple inversion of above calculation, using actual current density
    #  in strand instead of jcrit

    temp_margin = (
        1.0
        / 0.168
        * np.log(
            np.log(1.175e9 / (jstrand / f_strain + 1.288e8)) / (0.02115 * b_conductor)
        )
        + 4.2
        - temp_conductor
    )

    #  Check if ranges of validity have been violated

    if (temp_conductor > 20.0) or (b_conductor < 6.0) or (b > 104.0):
        raise ProcessValueError(
            "Fit extrapolated outside of range of validity",
            temp_conductor=temp_conductor,
            b_conductor=b_conductor,
            b=b,
        )

    return j_critical, temp_margin

gl_nbti(temp_conductor, b_conductor, strain, b_c20max, t_c0)

Calculate the critical current density, critical field, and critical temperature of an ITER Nb-Ti strand based on the Ginzburg-Landau theory of superconductivity.

Parameters:

Name Type Description Default
temp_conductor float

Superconductor temperature [K].

required
b_conductor float

Magnetic field at the conductor [T].

required
strain float

Intrinsic strain in the superconductor [%].

required
b_c20max float

Strain-dependent upper critical field at zero temperature [T].

required
t_c0 float

Strain-dependent critical temperature at zero strain [K].

required

Returns:

Type Description
tuple[float, float, float]

A tuple containing: - j_critical : Critical current density in the superconductor [A/m²]. - b_critical : Critical magnetic field [T]. - t_critical : Critical temperature [K].

References
  • Model based on: S B L Chislett-Mcdonald, Y. Tsui, E. Surrey, M. Kovari, and D. P. Hampshire, “The magnetic field, temperature, strain and angular dependence of the critical current density for Nb-Ti,” Journal of Physics Conference Series, vol. 1559, no. 1, pp. 012063-012063, Jun. 2020, doi: https://doi.org/10.1088/1742-6596/1559/1/012063.
Source code in process/models/superconductors.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
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
def gl_nbti(
    temp_conductor: float,
    b_conductor: float,
    strain: float,
    b_c20max: float,
    t_c0: float,
) -> tuple[float, float, float]:
    """Calculate the critical current density, critical field, and critical temperature
    of an ITER Nb-Ti strand based on the Ginzburg-Landau theory of superconductivity.

    Parameters
    ----------
    temp_conductor : float
        Superconductor temperature [K].
    b_conductor : float
        Magnetic field at the conductor [T].
    strain : float
        Intrinsic strain in the superconductor [%].
    b_c20max : float
        Strain-dependent upper critical field at zero temperature [T].
    t_c0 : float
        Strain-dependent critical temperature at zero strain [K].

    Returns
    -------
    tuple[float, float, float]
        A tuple containing:
        - j_critical : Critical current density in the superconductor [A/m²].
        - b_critical : Critical magnetic field [T].
        - t_critical : Critical temperature [K].

    References
    ----------
    - Model based on: S B L Chislett-Mcdonald, Y. Tsui, E. Surrey, M. Kovari, and D. P. Hampshire,
    “The magnetic field, temperature, strain and angular dependence of the critical current density for Nb-Ti,”
    Journal of Physics Conference Series, vol. 1559, no. 1, pp. 012063-012063, Jun. 2020, doi:
    https://doi.org/10.1088/1742-6596/1559/1/012063.
    """

    a_0 = 1102e6
    p = 0.49
    q = 0.56
    n = 1.83
    v = 1.42
    u = 0.0
    w = 2.2

    # Strain fitted parameters to a single filament
    c2 = -0.0025
    c3 = -0.0003
    c4 = -0.0001
    epsilon_m = -0.002e-2

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

    epsilon_i = strain - epsilon_m

    strain_func = (
        1 + c2 * (epsilon_i) ** 2 + c3 * (epsilon_i) ** 3 + c4 * (epsilon_i) ** 4
    )

    t_e = t_c0 * strain_func ** (1 / w)

    t_reduced = temp_conductor / t_e

    a_e = a_0 * strain_func ** (u / w)

    # Critical Field
    b_critical = b_c20max * (1 - t_reduced**v) * strain_func

    b_reduced = b_conductor / b_critical

    # Critical temperature (K)
    t_critical = t_e

    # Critical current density (A/m2)
    if b_reduced <= 1.0:
        j_critical = (
            a_e
            * (t_e * (1 - t_reduced**2)) ** 2
            * b_critical ** (n - 3)
            * b_reduced ** (p - 1)
            * (1 - b_reduced) ** q
        )
    else:  # Fudge to yield negative single valued function of Jc for fields above Bc2
        j_critical = (
            a_e
            * (t_e * (1 - t_reduced**2)) ** 2
            * b_critical ** (n - 3)
            * b_reduced ** (p - 1)
            * (1 - b_reduced) ** 1.0
        )

    return j_critical, b_critical, t_critical

gl_rebco(temp_conductor, b_conductor, strain, b_c20max, t_c0)

Calculate the critical current density, critical field, and critical temperature for a SuperPower REBCO tape based on measurements by P. Branch at Durham University and the Ginzburg-Landau theory of superconductivity

Parameters:

Name Type Description Default
temp_conductor float

Coolant/superconductor temperature (K).

required
b_conductor float

Magnetic field at conductor (T).

required
strain float

Intrinsic strain in superconductor (%).

required
b_c20max float

Strain-dependent upper critical field at zero temperature (T).

required
t_c0 float

Strain-dependent critical temperature at zero strain (K).

required

Returns:

Type Description
tuple[float, float, float]

Tuple containing: - j_critical: Critical current density in superconductor (A/m²). - b_critical: Critical magnetic field (T). - temp_critical: Critical temperature (K).

References
  • Model based on: S B L Chislett-Mcdonald, Y. Tsui, E. Surrey, M. Kovari, and D. P. Hampshire, “The magnetic field, temperature, strain and angular dependence of the critical current density for Nb-Ti,” Journal of Physics Conference Series, vol. 1559, no. 1, pp. 012063-012063, Jun. 2020, doi:
https://doi.org/10.1088/1742-6596/1559/1/012063.

-Fit to state-of-the-art measurements at 4.2 K:P. Branch, K. Osamura, and D. Hampshire, “Weak emergence in the angular dependence of the critical current density of the high temperature superconductor coated conductor REBCO,” Superconductor Science and Technology, vol. 33, no. 10, p. 104006, Sep. 2020, doi: 10.1088/1361-6668/abaebe.

Source code in process/models/superconductors.py
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
def gl_rebco(
    temp_conductor: float,
    b_conductor: float,
    strain: float,
    b_c20max: float,
    t_c0: float,
) -> tuple[float, float, float]:
    """Calculate the critical current density, critical field, and critical temperature
    for a SuperPower REBCO tape based on measurements by P. Branch at Durham University and
    the Ginzburg-Landau theory of superconductivity

    Parameters
    ----------
    temp_conductor : float
        Coolant/superconductor temperature (K).
    b_conductor : float
        Magnetic field at conductor (T).
    strain : float
        Intrinsic strain in superconductor (%).
    b_c20max : float
        Strain-dependent upper critical field at zero temperature (T).
    t_c0 : float
        Strain-dependent critical temperature at zero strain (K).

    Returns
    -------
    tuple[float, float, float]
        Tuple containing:
        - j_critical: Critical current density in superconductor (A/m²).
        - b_critical: Critical magnetic field (T).
        - temp_critical: Critical temperature (K).

    References
    ----------
    - Model based on: S B L Chislett-Mcdonald, Y. Tsui, E. Surrey, M. Kovari, and D. P. Hampshire,
    “The magnetic field, temperature, strain and angular dependence of the critical current density for Nb-Ti,”
    Journal of Physics Conference Series, vol. 1559, no. 1, pp. 012063-012063, Jun. 2020, doi:
    https://doi.org/10.1088/1742-6596/1559/1/012063.
    -
    -Fit to state-of-the-art measurements at 4.2 K:P. Branch, K. Osamura, and D. Hampshire,
    “Weak emergence in the angular dependence of the critical current density of the high temperature superconductor coated conductor REBCO,”
    Superconductor Science and Technology, vol. 33, no. 10, p. 104006, Sep. 2020,
    doi: 10.1088/1361-6668/abaebe.
    """
    # Critical current density pre-factor
    a_0 = 2.95e2

    # Flux pinning field scaling parameters
    p = 0.32
    q = 2.50
    n = 3.33
    # temperature scaling parameter
    s = 5.27

    # Strain scaling parameters
    c2 = -0.0191
    c3 = 0.0039
    c4 = 0.00103
    epsilon_m = 0.058

    # Strain conversion parameters
    u = 0.0
    w = 2.2

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

    epsilon_i = strain - epsilon_m

    strain_func = (
        1 + c2 * (epsilon_i) ** 2 + c3 * (epsilon_i) ** 3 + c4 * (epsilon_i) ** 4
    )

    t_e = t_c0 * strain_func ** (1 / w)

    t_reduced = temp_conductor / t_e

    a_e = a_0 * strain_func ** (u / w)

    #  Critical Field
    b_critical = b_c20max * (1 - t_reduced) ** s * strain_func

    b_reduced = b_conductor / b_critical

    #  Critical temperature (K)
    temp_critical = t_e

    #  Critical current density (A/m2)
    j_critical = (
        a_e
        * (t_e * (1 - t_reduced**2)) ** 2
        * b_critical ** (n - 3)
        * b_reduced ** (p - 1)
        * (1 - b_reduced) ** q
    )

    return j_critical, b_critical, temp_critical

hijc_rebco(temp_conductor, b_conductor, b_c20max, t_c0, dr_hts_tape, dx_hts_tape_rebco, dx_hts_tape_total)

Calculates the critical current density, critical field, and critical temperature for a high current density REBCO tape based Wolf et al. parameterization with data from Hazelton and Zhai et al.

Parameters:

Name Type Description Default
temp_conductor float

Superconductor temperature (K).

required
b_conductor float

Magnetic field at conductor (T).

required
b_c20max float

Upper critical field (T) for superconductor at zero temperature and strain.

required
t_c0 float

Critical temperature (K) at zero field and strain.

required
dr_hts_tape float

Width of the tape (m).

required
dx_hts_tape_rebco float

Thickness of the REBCO layer (m).

required
dx_hts_tape_total float

Total thickness of the tape (m).

required

Returns:

Type Description
tuple[float, float, float]

Tuple containing: - j_critical: Critical current density in superconductor (A/m²). - b_critical: Critical field (T). - temp_critical: Critical temperature (K).

Notes
  • The parameter A is transformed into a function A(T) based on a Newton polynomial fit considering A(4.2 K) = 2.2e8, A(20 K) = 2.3e8 and A(65 K) = 3.5e8.

  • A scaling factor of 0.4 was originally applied to j_critical for CORC cables, but is not used here.

References
  • Based in part on the parameterization described in: M. J. Wolf, Nadezda Bagrets, W. H. Fietz, C. Lange, and K.-P. Weiss, “Critical Current Densities of 482 A/mm2 in HTS CrossConductors at 4.2 K and 12 T,” IEEE Transactions on Applied Superconductivity, vol. 28, no. 4, pp. 1-4, Jun. 2018, doi: https://doi.org/10.1109/tasc.2018.2815767.

  • Fit values based on: D. W. Hazelton, “4th Workshop on Accelerator Magnets in HTS (WAMHTS-4) | 2G HTS Wire Development at SuperPower,” Indico, 2017. https://indico.cern.ch/event/588810/contributions/2473740/ (accessed May 20, 2025).

-The high Ic parameterization is a result of modifications based on Ic values observed in: Y. Zhai, D. van der Laan, P. Connolly, and C. Kessel, “Conceptual design of HTS magnets for fusion nuclear science facility,” Fusion Engineering and Design, vol. 168, p. 112611, Jul. 2021, doi: https://doi.org/10.1016/j.fusengdes.2021.112611.

Source code in process/models/superconductors.py
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
def hijc_rebco(
    temp_conductor: float,
    b_conductor: float,
    b_c20max: float,
    t_c0: float,
    dr_hts_tape: float,
    dx_hts_tape_rebco: float,
    dx_hts_tape_total: float,
) -> tuple[float, float, float]:
    """Calculates the critical current density, critical field, and critical temperature
    for a high current density REBCO tape based Wolf et al. parameterization with data from Hazelton
    and Zhai et al.

    Parameters
    ----------
    temp_conductor : float
        Superconductor temperature (K).
    b_conductor : float
        Magnetic field at conductor (T).
    b_c20max : float
        Upper critical field (T) for superconductor at zero temperature and strain.
    t_c0 : float
        Critical temperature (K) at zero field and strain.
    dr_hts_tape : float
        Width of the tape (m).
    dx_hts_tape_rebco : float
        Thickness of the REBCO layer (m).
    dx_hts_tape_total : float
        Total thickness of the tape (m).

    Returns
    -------
    tuple[float, float, float]
        Tuple containing:
        - j_critical: Critical current density in superconductor (A/m²).
        - b_critical: Critical field (T).
        - temp_critical: Critical temperature (K).

    Notes
    -----
    - The parameter A is transformed into a function A(T) based on a Newton polynomial fit
      considering A(4.2 K) = 2.2e8, A(20 K) = 2.3e8 and A(65 K) = 3.5e8.

    - A scaling factor of 0.4 was originally applied to j_critical for CORC cables, but is not used here.

    References
    ----------
    - Based in part on the parameterization described in:
    M. J. Wolf, Nadezda Bagrets, W. H. Fietz, C. Lange, and K.-P. Weiss,
    “Critical Current Densities of 482 A/mm2 in HTS CrossConductors at 4.2 K and 12 T,”
    IEEE Transactions on Applied Superconductivity, vol. 28, no. 4, pp. 1-4, Jun. 2018,
    doi: https://doi.org/10.1109/tasc.2018.2815767.

    - Fit values based on:
    D. W. Hazelton, “4th Workshop on Accelerator Magnets in HTS (WAMHTS-4) | 2G HTS Wire Development at SuperPower,”
    Indico, 2017. https://indico.cern.ch/event/588810/contributions/2473740/ (accessed May 20, 2025).

    -The high Ic parameterization is a result of modifications based on Ic values observed in:
    Y. Zhai, D. van der Laan, P. Connolly, and C. Kessel, “Conceptual design of HTS magnets for fusion nuclear science facility,”
    Fusion Engineering and Design, vol. 168, p. 112611, Jul. 2021,
    doi: https://doi.org/10.1016/j.fusengdes.2021.112611.
    """

    a = 1.4
    b = 2.005
    # critical current density prefactor
    a_0 = 2.2e8
    # flux pinning field scaling parameters
    p = 0.39
    q = 0.9
    # strain conversion parameters
    u = 33450.0
    v = -176577.0

    # Critical Field (T)
    # B_crit(T) calculated using temperature and critical temperature
    b_critical = b_c20max * (1.0 - temp_conductor / t_c0) ** a

    # Critical temperature (K)
    # scaled to match behaviour in GL_REBCO routine,
    # ONLY TO BE USED until a better suggestion is received
    temp_critical = 0.999965 * t_c0

    # finding A(T); constants based on a Newton polynomial fit to pubished data
    a_t = a_0 + (u * temp_conductor**2) + (v * temp_conductor)

    # Critical current density (A/m2)
    # In the original formula bcrit must be > bmax to prevent NaNs.
    # However, negative jcrit is permissible (I think).
    # So when bcrit < bmax, I reverse the sign of the bracket,
    # giving a negative but real value of jcrit.

    if b_critical > b_conductor:
        j_critical = (
            (a_t / b_conductor)
            * b_critical**b
            * (b_conductor / b_critical) ** p
            * (1 - b_conductor / b_critical) ** q
        )
    else:
        j_critical = (
            (a_t / b_conductor)
            * b_critical**b
            * (b_conductor / b_critical) ** p
            * (b_conductor / b_critical - 1) ** q
        )

    # Jc times HTS area: default area is width 4mm times HTS layer thickness 1 um,
    # divided by the tape area to provide engineering Jc per tape,!
    # A scaling factor of 0.4 used to be applied below to assume the difference
    # between tape stacks and CORC cable layouts.

    j_critical = (
        j_critical
        * (dr_hts_tape * dx_hts_tape_rebco)
        / (dr_hts_tape * dx_hts_tape_total)
    )

    return j_critical, b_critical, temp_critical

western_superconducting_nb3sn(temp_conductor, b_conductor, strain, b_c20max, temp_c0max)

Calculate the critical current density, critical field, and critical temperature for a WST Nb₃Sn superconductor using the ITER Nb₃Sn critical surface model.

Parameters:

Name Type Description Default
temp_conductor float

Superconductor temperature (K).

required
b_conductor float

Magnetic field at the superconductor (T).

required
strain float

Strain in the superconductor.

required
b_c20max float

Upper critical field (T) for the superconductor at zero temperature and strain.

required
temp_c0max float

Critical temperature (K) at zero field and strain.

required

Returns:

Type Description
tuple[float, float, float]

A tuple containing: - j_critical: Critical current density in the superconductor (A/m²). - b_crititical: Critical field (T). - temp_critical: Critical temperature (K).

Notes
  • This routine uses the WST Nb3Sn critical surface model.
  • The scaling constants and parameters are based on the reference below.
  • This assumes a strand size of 1.5mm.
  • Compared to the EUTF4 (OST) ( European qualification samples for TF conductor), the performance of the WST at low strain is superior by about 10%.
References
Source code in process/models/superconductors.py
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
def western_superconducting_nb3sn(
    temp_conductor: float,
    b_conductor: float,
    strain: float,
    b_c20max: float,
    temp_c0max: float,
) -> tuple[float, float, float]:
    """Calculate the critical current density, critical field, and critical temperature
    for a WST Nb₃Sn superconductor using the ITER Nb₃Sn critical surface model.

    Parameters
    ----------
    temp_conductor : float
        Superconductor temperature (K).
    b_conductor : float
        Magnetic field at the superconductor (T).
    strain : float
        Strain in the superconductor.
    b_c20max : float
        Upper critical field (T) for the superconductor at zero temperature and strain.
    temp_c0max : float
        Critical temperature (K) at zero field and strain.

    Returns
    -------
    tuple[float, float, float]
        A tuple containing:
        - j_critical: Critical current density in the superconductor (A/m²).
        - b_crititical: Critical field (T).
        - temp_critical: Critical temperature (K).

    Notes
    -----
    - This routine uses the WST Nb3Sn critical surface model.
    - The scaling constants and parameters are based on the reference below.
    - This assumes a strand size of 1.5mm.
    - Compared to the EUTF4 (OST) ( European qualification samples for TF conductor),
    the performance of the WST at low strain is superior by about 10%.

    References
    ----------
    - V. Corato, “EUROFUSION WPMAG-REP(16) 16565 Common operating values for DEMO magnets design for 2016 REPORT.”
    Accessed: May 12, 2025. [Online].
    Available: https://scipub.euro-fusion.org/wp-content/uploads/eurofusion/WPMAGREP16_16565_submitted.pdf

    - “Introduction of WST,” 2015. Accessed: May 12, 2025. [Online].
    Available: https://indico.cern.ch/event/340703/contributions/802232/attachments/668814/919331/WST_INTRO_2015-3_for_FCC_WEEK.pdf
    """
    # Scaling constant C [AT/mm²]
    csc = 83075.0
    # Low field exponent p
    p = 0.593
    # High field exponent q
    q = 2.156
    # Strain fitting constant C_{a1}
    c_a1 = 50.06
    # Strain fitting constant C_{a2}
    c_a2 = 0.0
    # Residual strain component epsilon_{0,a}
    epsilon_0a = 0.00312

    j_scaling, b_critical, t_critical = bottura_scaling(
        csc=csc,
        p=p,
        q=q,
        c_a1=c_a1,
        c_a2=c_a2,
        epsilon_0a=epsilon_0a,
        temp_conductor=temp_conductor,
        b_conductor=b_conductor,
        epsilon=strain,
        b_c20max=b_c20max,
        temp_c0max=temp_c0max,
    )

    # Scale from mm² to m²
    scalefac = 1.0e6
    j_critical = j_scaling * scalefac
    return j_critical, b_critical, t_critical

bottura_scaling(csc, p, q, c_a1, c_a2, epsilon_0a, temp_conductor, b_conductor, epsilon, b_c20max, temp_c0max)

Implements the scaling from: Jc(B,T,epsilon) Parameterization for the ITER Nb₃Sn Production,

Parameters:

Name Type Description Default
csc float

Scaling constant C [AT/mm²].

required
p float

Low field exponent of the pinning force

required
q float

High field exponent of the pinning force

required
c_a1 float

Strain fitting constant C_{a1}.

required
c_a2 float

Strain fitting constant C_{a2}.

required
epsilon_0a float

Residual strain component

required
temp_conductor float

Superconductor temperature (K).

required
b_conductor float

Magnetic field at conductor (T).

required
epsilon float

Strain in superconductor.

required
b_c20max float

Upper critical field (T) for superconductor at zero temperature and strain.

required
temp_c0max float

Critical temperature (K) at zero field and strain.

required

Returns:

Type Description
tuple[float, float, float]

A tuple containing: - jscaling: Critical current density scaling factor (A/m²). - bcrit: Critical field (T). - tcrit: Critical temperature (K).

Notes
  • This is a generic scaling proposed for the characterization and production of ITER Nb₃Sn strands. This is also known as the "ITER-2008 parametrization."

  • Parameter ranges are strain (1.5% to 0.4%), temperature (2.35 to 16 K), and field (0.5 to 19 T). The ITER-2008 parameterization achieves an average accuracy error of 3.8 Amps, with the best at 1.5 Amps and the worst at 7.5 Amps.

  • The strain function is suitable only in the moderate strain region, down to 0.8%.

References
  • L. Bottura and B. Bordini, “J_{C}(B,T, arepsilon) Parameterization for the ITER { m Nb}_{3}{ m Sn} Production,” IEEE Transactions on Applied Superconductivity, vol. 19, no. 3, pp. 1521-1524, Jun. 2009, doi: https://doi.org/10.1109/tasc.2009.2018278.
Source code in process/models/superconductors.py
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
def bottura_scaling(
    csc: float,
    p: float,
    q: float,
    c_a1: float,
    c_a2: float,
    epsilon_0a: float,
    temp_conductor: float,
    b_conductor: float,
    epsilon: float,
    b_c20max: float,
    temp_c0max: float,
) -> tuple[float, float, float]:
    """Implements the scaling from:
    Jc(B,T,epsilon) Parameterization for the ITER Nb₃Sn Production,

    Parameters
    ----------
    csc : float
        Scaling constant C [AT/mm²].
    p : float
        Low field exponent of the pinning force
    q : float
        High field exponent of the pinning force
    c_a1 : float
        Strain fitting constant C_{a1}.
    c_a2 : float
        Strain fitting constant C_{a2}.
    epsilon_0a : float
        Residual strain component
    temp_conductor : float
        Superconductor temperature (K).
    b_conductor : float
        Magnetic field at conductor (T).
    epsilon : float
        Strain in superconductor.
    b_c20max : float
        Upper critical field (T) for superconductor at zero temperature and strain.
    temp_c0max : float
        Critical temperature (K) at zero field and strain.

    Returns
    -------
    tuple[float, float, float]
        A tuple containing:
        - jscaling: Critical current density scaling factor (A/m²).
        - bcrit: Critical field (T).
        - tcrit: Critical temperature (K).

    Notes
    -----
    - This is a generic scaling proposed for the characterization and production of
    ITER Nb₃Sn strands. This is also known as the "ITER-2008 parametrization."

    - Parameter ranges are strain (1.5% to 0.4%), temperature (2.35 to 16 K), and field (0.5 to 19 T).
    The ITER-2008 parameterization achieves an average accuracy error of 3.8 Amps, with the best at 1.5 Amps and the worst at 7.5 Amps.

    - The strain function is suitable only in the moderate strain region, down to 0.8%.

    References
    ----------
    - L. Bottura and B. Bordini, “$J_{C}(B,T,\varepsilon)$ Parameterization for the ITER ${\rm Nb}_{3}{\rm Sn}$ Production,”
    IEEE Transactions on Applied Superconductivity, vol. 19, no. 3, pp. 1521-1524, Jun. 2009,
    doi: https://doi.org/10.1109/tasc.2009.2018278.
    """

    epsilon_sh = (c_a2 * epsilon_0a) / (np.sqrt(c_a1**2 - c_a2**2))

    # Strain function
    # 0.83 < s < 1.0, for -0.005 < strain < 0.005
    strfun = np.sqrt(epsilon_sh**2 + epsilon_0a**2) - np.sqrt(
        (epsilon - epsilon_sh) ** 2 + epsilon_0a**2
    )
    strfun = strfun * c_a1 - (c_a2 * epsilon)
    strfun = 1.0 + (1 / (1.0 - c_a1 * epsilon_0a)) * strfun

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

    # Critical field at zero temperature and current, corrected for strain
    b_c20_eps = b_c20max * strfun

    # Critical temperature at zero field and current, corrected for strain
    temp_c0_eps = temp_c0max * strfun ** (1 / 3)

    # If input temperature is over the strain adjusted critical temperature then report error
    if temp_conductor / temp_c0_eps >= 1.0:
        logger.error(
            f"Reduced temperature t artificially lowered {temp_conductor=} {temp_c0_eps=}"
        )

    # Reduced temperature at zero field, corrected for strain
    # f_temp_conductor_critical > 1 is permitted, indicating the temperature is above the critical value at zero field.
    f_temp_conductor_critical_no_field = temp_conductor / temp_c0_eps

    # If input field is over the strain adjusted critical field then report error
    if b_conductor / b_c20_eps >= 1.0:
        logger.error(
            f"Reduced field bzero artificially lowered {b_conductor=} {b_c20_eps=}"
        )

    # Reduced field at zero temperature, taking account of strain
    f_b_conductor_critical_no_temp = b_conductor / b_c20_eps

    # Critical temperature at given strain and field
    # tcrit is not used in the calculation of jcrit.
    if (
        f_b_conductor_critical_no_temp < 1.0
    ):  # Normal case, field is within critical surface
        temp_critical = temp_c0_eps * (1.0 - f_b_conductor_critical_no_temp) ** (
            1 / 1.52
        )
    else:  # Abnormal case, field is too high.
        temp_critical = -temp_c0_eps * abs(1.0 - f_b_conductor_critical_no_temp) ** (
            1 / 1.52
        )  # Prevents NaNs. Sets tcrit negative

    # Critical field (T) at given strain and temperature
    b_critical = b_c20_eps * (1.0 - f_temp_conductor_critical_no_field**1.52)

    jc1 = (csc / b_conductor) * strfun

    # Check if we are inside the critical surface
    if (
        (f_temp_conductor_critical_no_field > 0)
        and (f_temp_conductor_critical_no_field < 1)
        and (b_conductor > 0)
        and (b_conductor < b_critical)
        and (b_critical > 0)
    ):
        # Reduced field at given strain and temperature
        b_reduced = b_conductor / b_critical

        jc2 = (1.0 - f_temp_conductor_critical_no_field**1.52) * (
            1.0 - f_temp_conductor_critical_no_field**2
        )
        jc3 = b_reduced**p * (1.0 - b_reduced) ** q
        j_scaling = jc1 * jc2 * jc3

    else:
        # Outside the critical surface.
        # We construct a simple function that is always negative and
        # becomes more negative as field and temperature increase.
        jc2 = f_temp_conductor_critical_no_field
        jc3 = b_conductor / max(b_critical, 1.0e-8)
        j_scaling = -abs(jc1 * jc2 * jc3)

    return j_scaling, b_critical, temp_critical

calculate_croco_cable_geometry(dia_croco_strand, dx_croco_strand_copper, dx_hts_tape_rebco, dx_hts_tape_copper, dx_hts_tape_hastelloy)

Calculate geometry and areas for a CroCo cable strand.

Parameters:

Name Type Description Default
dia_croco_strand float

Diameter of CroCo strand (m)

required
dx_croco_strand_copper float

Thickness of copper layer in CroCo strand (m)

required
dx_hts_tape_rebco float

Thickness of REBCO layer in HTS tape (m)

required
dx_hts_tape_copper float

Thickness of copper layer in HTS tape (m)

required
dx_hts_tape_hastelloy float

Thickness of Hastelloy layer in HTS tape (m)

required

Returns:

Type Description
tuple[float, float, float, float, float, float, float, float]

Tuple containing: - dia_croco_strand_tape_region: Inner diameter of CroCo strand tape region (m) - n_croco_strand_hts_tapes: Number of HTS tapes in CroCo strand - a_croco_strand_copper_total: Total copper area in CroCo strand (m²) - a_croco_strand_hastelloy: Total Hastelloy area in CroCo strand (m²) - a_croco_strand_solder: Total solder area in CroCo strand (m²) - a_croco_strand_rebco: Total REBCO area in CroCo strand (m²) - a_croco_strand: Total area of CroCo strand (m²) - dr_hts_tape: Width of the tape (m)

Source code in process/models/superconductors.py
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
def calculate_croco_cable_geometry(
    dia_croco_strand: float,
    dx_croco_strand_copper: float,
    dx_hts_tape_rebco: float,
    dx_hts_tape_copper: float,
    dx_hts_tape_hastelloy: float,
) -> tuple[
    float,  # dia_croco_strand_tape_region
    float,  # n_croco_strand_hts_tapes
    float,  # a_croco_strand_copper_total
    float,  # a_croco_strand_hastelloy
    float,  # a_croco_strand_solder
    float,  # a_croco_strand_rebco
    float,  # a_croco_strand
    float,  # dr_hts_tape
]:
    """Calculate geometry and areas for a CroCo cable strand.

    Parameters
    ----------
    dia_croco_strand : float
        Diameter of CroCo strand (m)
    dx_croco_strand_copper : float
        Thickness of copper layer in CroCo strand (m)
    dx_hts_tape_rebco : float
        Thickness of REBCO layer in HTS tape (m)
    dx_hts_tape_copper : float
        Thickness of copper layer in HTS tape (m)
    dx_hts_tape_hastelloy : float
        Thickness of Hastelloy layer in HTS tape (m)

    Returns
    -------
    tuple[float, float, float, float, float, float, float, float]
        Tuple containing:
        - dia_croco_strand_tape_region: Inner diameter of CroCo strand tape region (m)
        - n_croco_strand_hts_tapes: Number of HTS tapes in CroCo strand
        - a_croco_strand_copper_total: Total copper area in CroCo strand (m²)
        - a_croco_strand_hastelloy: Total Hastelloy area in CroCo strand (m²)
        - a_croco_strand_solder: Total solder area in CroCo strand (m²)
        - a_croco_strand_rebco: Total REBCO area in CroCo strand (m²)
        - a_croco_strand: Total area of CroCo strand (m²)
        - dr_hts_tape: Width of the tape (m)
    """

    # Calculate the inner diameter of the CroCo strand tape region
    dia_croco_strand_tape_region = dia_croco_strand - 2.0 * dx_croco_strand_copper
    if dia_croco_strand_tape_region <= 0.0:
        logger.error("Negative inner CroCo cable diameter")

    # Total thickness of HTS tape
    dx_hts_tape_total = dx_hts_tape_rebco + dx_hts_tape_copper + dx_hts_tape_hastelloy

    scaling = dia_croco_strand_tape_region / 5.4e-3
    dr_hts_tape = scaling * 3.75e-3

    # Calculate the height of HTS tapes in the CroCo strand
    dx_croco_strand_tape_stack = np.sqrt(
        dia_croco_strand_tape_region**2 - dr_hts_tape**2
    )
    # Number of HTS tapes in the CroCo strand
    n_croco_strand_hts_tapes = dx_croco_strand_tape_stack / dx_hts_tape_total

    # Area of copper in the CroCo strand (copper tube + copper in HTS tapes)
    a_croco_strand_copper_total = (
        np.pi * dx_croco_strand_copper * dia_croco_strand
        - np.pi * dx_croco_strand_copper**2
        + dx_hts_tape_copper * dr_hts_tape * n_croco_strand_hts_tapes
    )
    # Area of Hastelloy in the CroCo strand
    a_croco_strand_hastelloy = (
        dx_hts_tape_hastelloy * dr_hts_tape * n_croco_strand_hts_tapes
    )
    # Area of solder in the CroCo strand surrounding the HTS tapes
    a_croco_strand_solder = (
        np.pi / 4.0 * dia_croco_strand_tape_region**2
        - dx_croco_strand_tape_stack * dr_hts_tape
    )

    # Area of REBCO in the CroCo strand
    a_croco_strand_rebco = dx_hts_tape_rebco * dr_hts_tape * n_croco_strand_hts_tapes
    # Total area of the CroCo strand
    a_croco_strand = np.pi / 4.0 * dia_croco_strand**2

    return (
        dia_croco_strand_tape_region,
        n_croco_strand_hts_tapes,
        a_croco_strand_copper_total,
        a_croco_strand_hastelloy,
        a_croco_strand_solder,
        a_croco_strand_rebco,
        a_croco_strand,
        dr_hts_tape,
    )

croco(j_crit_sc, conductor_area, dia_croco_strand, dx_croco_strand_copper)

'CroCo' (cross-conductor) strand and cable design for 'REBCO' 2nd generation HTS superconductor Updated 13/11/18 using data from Lewandowska et al 2018.

Parameters:

Name Type Description Default
j_crit_sc
required
conductor_area
required
dia_croco_strand
required
dx_croco_strand_copper
required
Source code in process/models/superconductors.py
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
def croco(j_crit_sc, conductor_area, dia_croco_strand, dx_croco_strand_copper):
    """'CroCo' (cross-conductor) strand and cable design for
    'REBCO' 2nd generation HTS superconductor
    Updated 13/11/18 using data from Lewandowska et al 2018.

    Parameters
    ----------
    j_crit_sc :

    conductor_area :

    dia_croco_strand :

    dx_croco_strand_copper :

    """

    (
        rebco_variables.dia_croco_strand_tape_region,
        rebco_variables.n_croco_strand_hts_tapes,
        a_croco_strand_copper_total,
        a_croco_strand_hastelloy,
        a_croco_strand_solder,
        a_croco_strand_rebco,
        a_croco_strand,
        rebco_variables.dr_hts_tape,
    ) = calculate_croco_cable_geometry(
        dia_croco_strand,
        dx_croco_strand_copper,
        rebco_variables.dx_hts_tape_rebco,
        rebco_variables.dx_hts_tape_copper,
        rebco_variables.dx_hts_tape_hastelloy,
    )

    rebco_variables.a_croco_strand_copper_total = a_croco_strand_copper_total
    rebco_variables.a_croco_strand_hastelloy = a_croco_strand_hastelloy
    rebco_variables.a_croco_strand_solder = a_croco_strand_solder
    rebco_variables.a_croco_strand_rebco = a_croco_strand_rebco
    rebco_variables.a_croco_strand = a_croco_strand

    croco_strand_critical_current = j_crit_sc * a_croco_strand_rebco

    # Conductor properties
    # conductor%number_croco = conductor%acs*(1.0-cable_helium_fraction-copper_bar)/a_croco_strand
    conductor_critical_current = croco_strand_critical_current * 6.0
    # Area of core = area of strand
    conductor_copper_bar_area = a_croco_strand
    conductor_copper_area = a_croco_strand_copper_total * 6.0 + conductor_copper_bar_area
    conductor_copper_fraction = conductor_copper_area / conductor_area

    # Helium area is set by the user.
    # conductor_helium_area = cable_helium_fraction * conductor_acs
    conductor_helium_area = np.pi / 2.0 * dia_croco_strand**2
    conductor_helium_fraction = conductor_helium_area / conductor_area

    conductor_hastelloy_area = a_croco_strand_hastelloy * 6.0
    conductor_hastelloy_fraction = conductor_hastelloy_area / conductor_area

    conductor_solder_area = a_croco_strand_solder * 6.0
    conductor_solder_fraction = conductor_solder_area / conductor_area

    conductor_rebco_area = a_croco_strand_rebco * 6.0
    conductor_rebco_fraction = conductor_rebco_area / conductor_area

    return (
        a_croco_strand,
        croco_strand_critical_current,
        conductor_copper_area,
        conductor_copper_fraction,
        conductor_copper_bar_area,
        conductor_hastelloy_area,
        conductor_hastelloy_fraction,
        conductor_helium_area,
        conductor_helium_fraction,
        conductor_solder_area,
        conductor_solder_fraction,
        conductor_rebco_area,
        conductor_rebco_fraction,
        conductor_critical_current,
    )

superconductor_current_density_margin(temp_superconductor, i_superconductor_type, j_superconductor, b_superconductor, strain, bc20m, tc0m, c0=0.0)

Calculate the current density margin for a superconductor. The current density margin is the difference between the operating current density and the critical current density of a superconductor at a given temperature and field. It is zero at the current-sharing temperature.

Superconductor material codes: 1: ITER Nb₃Sn 3: NbTi (Lubell scaling) 4: ITER Nb₃Sn (alternate) 5: Western Superconducting Nb₃Sn 7: Ginzburg-Landau NbTi 8: Ginzburg-Landau REBCO 9: High current density REBCO

Parameters:

Name Type Description Default
temp_superconductor float

Superconductor Temperature (K)

required
i_superconductor_type int

Switch for superconductor material

required
j_superconductor float

Superconductor actual current density (A/m²)

required
b_superconductor float

Magnetic field at the superconductor (T)

required
strain float

Superconductor strain

required
bc20m float

Upper critical field (T)

required
tc0m float

Critical temperature (K)

required
c0 float

Scaling constant (A/m²), required for NbTi

0.0

Returns:

Type Description
float

Current density margin (A/m²)

Source code in process/models/superconductors.py
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
def superconductor_current_density_margin(
    temp_superconductor: float,
    i_superconductor_type: int,
    j_superconductor: float,
    b_superconductor: float,
    strain: float,
    bc20m: float,
    tc0m: float,
    c0: float = 0.0,
) -> float:
    """Calculate the current density margin for a superconductor.
    The current density margin is the difference between the operating current density and
    the critical current density of a superconductor at a given temperature and field.
    It is zero at the current-sharing temperature.

    Superconductor material codes:
        1: ITER Nb₃Sn
        3: NbTi (Lubell scaling)
        4: ITER Nb₃Sn (alternate)
        5: Western Superconducting Nb₃Sn
        7: Ginzburg-Landau NbTi
        8: Ginzburg-Landau REBCO
        9: High current density REBCO

    Parameters
    ----------
    temp_superconductor : float
        Superconductor Temperature (K)
    i_superconductor_type : int
        Switch for superconductor material
    j_superconductor : float
        Superconductor actual current density (A/m²)
    b_superconductor : float
        Magnetic field at the superconductor (T)
    strain : float
        Superconductor strain
    bc20m : float
        Upper critical field (T)
    tc0m : float
        Critical temperature (K)
    c0 : float, optional
        Scaling constant (A/m²), required for NbTi

    Returns
    -------
    float
        Current density margin (A/m²)

    """
    material_functions = {
        1: lambda: itersc(temp_superconductor, b_superconductor, strain, bc20m, tc0m)[0],
        3: lambda: jcrit_nbti(temp_superconductor, b_superconductor, c0, bc20m, tc0m)[0],
        4: lambda: itersc(temp_superconductor, b_superconductor, strain, bc20m, tc0m)[0],
        5: lambda: western_superconducting_nb3sn(
            temp_superconductor, b_superconductor, strain, bc20m, tc0m
        )[0],
        7: lambda: gl_nbti(temp_superconductor, b_superconductor, strain, bc20m, tc0m)[
            0
        ],
        8: lambda: gl_rebco(temp_superconductor, b_superconductor, strain, bc20m, tc0m)[
            0
        ],
        9: lambda: hijc_rebco(
            temp_superconductor,
            b_superconductor,
            bc20m,
            tc0m,
            rebco_variables.dr_hts_tape,
            rebco_variables.dx_hts_tape_rebco,
            rebco_variables.dx_hts_tape_total,
        )[0],
    }

    try:
        j_superconductor_critical = material_functions[i_superconductor_type]()
    except KeyError as err:
        raise ValueError(
            f"Unknown superconductor material code: {i_superconductor_type}"
        ) from err

    return j_superconductor_critical - j_superconductor