Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

VV stress under quench model review #3366

Open
wants to merge 8 commits into
base: main
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
102 changes: 95 additions & 7 deletions process/sctfcoil.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,7 @@
import logging
import copy
import numba
import numpy as np

from process.fortran import rebco_variables
from process.fortran import global_variables
Expand Down Expand Up @@ -2331,13 +2332,14 @@ def vv_stress_on_quench(self):
# Area of the radial plate taken to be the area of steel in the WP
# TODO: value clipped due to #1883
S_rp=numpy.clip(sctfcoil_module.a_tf_steel, 0, None),
# TODO: Does this calculation of Scc exclude the area of the case down the side?
S_cc=sctfcoil_module.a_case_front + sctfcoil_module.a_case_nose,
S_cc=sctfcoil_module.a_case_front
+ sctfcoil_module.a_case_nose
+ 2.0 * sctfcoil_module.t_lat_case_av,
taud=tfcoil_variables.tdmptf,
# TODO: is this the correct current?
I_op=sctfcoil_module.tfc_current / tfcoil_variables.n_tf_turn,
# VV properties
d_vv=build_variables.d_vv_in,
d_vv=build_variables.d_vv_shell_thickness,
)

def tf_field_and_force(self):
Expand Down Expand Up @@ -7185,6 +7187,88 @@ def _inductance_factor(H, Ri, Ro, Rm, theta1):
)


@staticmethod
def lambda_term(tau, omega):
"""Lambda Term

Author: Alexander Pearce, UKAEA

:param tau: tau_{s,k} = (R_{s,k} - R_{c,k}) / R_k
:param omega: omega_k = R_{c,k}/R_k

The lammbda fucntion used inegral in inductance calcuation found
in Y. Itoh et al. The full form of the functions are given in appendix A.
"""
p = 1.0 - omega**2.0

if p < 0:
integral = (1.0 / np.sqrt(np.abs(p))) * np.arcsin(
(1.0 + omega * tau) / (tau + omega)
)
else:
integral = (1.0 / np.sqrt(np.abs(p))) * np.log(
(2.0 * (1.0 + tau * omega - np.sqrt(p * (1 - tau**2.0)))) / (tau + omega)
)

return integral


@staticmethod
def _theta_factor_integral(Ro_vv, Ri_vv, Rm_vv, H_vv, theta1_vv):
"""Theta Factor Integral
Author: Alexander Pearce, UKAEA

:param Ro_vv: the radius of the outboard edge of the VV CCL
:param Ri_vv: the radius of the inboard edge of the VV CCL
:param Rm_vv: the radius where the maximum height of the VV CCL occurs
:param H_vv: the maximum height of the VV CCL
:param theta1_vv: the polar angle of the point at which one circular arc is
joined to another circular arc in the approximation to the VV CCL

The calcuation of the theta factor found in Eq 4 of Y. Itoh et al. The
full form of the integral is given in appendix A.
"""
theta2 = np.pi / 2.0 + theta1_vv
a = (Ro_vv - Ri_vv) / 2.0
Rbar = (Ro_vv + Ri_vv) / 2.0
A = Rbar / a
delta = (Rbar - Rm_vv) / a
kappa = H_vv / a
iota = (1.0 + delta) / kappa

denom = np.cos(theta1_vv) + np.sin(theta1_vv) - 1.0

R1 = H_vv * ((np.cos(theta1_vv) + iota * (np.sin(theta1_vv) - 1.0)) / denom)
R2 = H_vv * ((np.cos(theta1_vv) - 1.0 + iota * np.sin(theta1_vv)) / denom)
R3 = H_vv * (1 - delta) / kappa

Rc1 = (H_vv / kappa) * (A + 1.0) - R1
Rc2 = Rc1 + (R1 - R2) * np.cos(theta1_vv)
Rc3 = Rc2
Zc2 = (R1 - R2) * np.sin(theta1_vv)
Zc3 = Zc2 + R2 - R3

tau = np.array(
[
[np.cos(theta1_vv), np.cos(theta1_vv + theta2), -1.0],
[1.0, np.cos(theta1_vv), np.cos(theta1_vv + theta2)],
]
)

omega = np.array([Rc1 / R1, Rc2 / R2, Rc3 / R3])

# Assume up down symmetry and let Zc6 = - Zc3
chi1 = (Zc3 + np.abs(-Zc3)) / Ri_vv
chi2 = 0.0

for k in range(len(omega)):
chi2 = chi2 + np.abs(
lambda_term(tau[1, k], omega[k]) - lambda_term(tau[0, k], omega[k])
)

return (chi1 + 2.0 * chi2) / (2.0 * np.pi)


def vv_stress_on_quench(
# TF shape
H_coil,
Expand Down Expand Up @@ -7248,7 +7332,7 @@ def vv_stress_on_quench(
:param taud: the discharge time of the TF coil when quench occurs
:param I_op: the 'normal' operating current of the TF coil

:param d_vv: the thickness of the vacuum vessel
:param d_vv: the thickness of the vacuum vessel shell

:returns: the maximum stress experienced by the vacuum vessel

Expand All @@ -7257,18 +7341,22 @@ def vv_stress_on_quench(
The theta1 quantity for the TF coil and VV is not very meaningful. The
impact of it of the inductance is rather small. Generally, the paper seems to
suggest the TF coil is between 40 and 60, as this is the range they calculate
the surrogates over. No range is provided for the VV but the example using
JA DEMO is 1 degree suggesting the quantity will be very small.
the surrogates over. The thickness of the VV considers an ITER like design and
only the outer and inner shells that which act of conductuve structural material.

References
----------
1. ITOH, Yasuyuki & Utoh, Hiroyasu & SAKAMOTO, Yoshiteru & Hiwatari, Ryoji. (2020).
Empirical Formulas for Estimating Self and Mutual Inductances of Toroidal Field Coils and Structures.
Plasma and Fusion Research. 15. 1405078-1405078. 10.1585/pfr.15.1405078.
"""
# Convert angles into radians
theta1_vv_rad = np.pi * (theta1_vv / 180.0)

# Poloidal loop resistance (PLR) in ohms
theta_vv = _theta_factor_integral(Ro_vv, Ri_vv, Rm_vv, H_vv, theta1_vv_rad)
plr_coil = ((0.5 * ccl_length_coil) / (n_tf * (S_cc + S_rp))) * 1e-6
plr_vv = ((0.84 / d_vv) * 0.94) * 1e-6
plr_vv = ((0.84 / d_vv) * theta_vv) * 1e-6

# relevant self-inductances in henry (H)
coil_structure_self_inductance = (
Expand Down
4 changes: 4 additions & 0 deletions source/fortran/build_variables.f90
Original file line number Diff line number Diff line change
Expand Up @@ -78,6 +78,9 @@ module build_variables
real(dp) :: d_vv_bot
!! vacuum vessel underside thickness (TF coil / shield) (m)

real(dp) :: d_vv_shell_thickness
!! vacuum vessel double walled shell thickness (m)

real(dp) :: f_avspace
!! F-value for stellarator radial space check (`constraint equation 83`)

Expand Down Expand Up @@ -321,6 +324,7 @@ subroutine init_build_variables
d_vv_out = 0.07D0
d_vv_top = 0.07D0
d_vv_bot = 0.07D0
d_vv_shell_thickness = 0.12D0
f_avspace = 1.0D0
fcspc = 0.6D0
fseppc = 3.5D8
Expand Down
2 changes: 1 addition & 1 deletion source/fortran/tfcoil_variables.f90
Original file line number Diff line number Diff line change
Expand Up @@ -1054,7 +1054,7 @@ subroutine init_tfcoil_variables
i_cp_joints = -1
cryo_cool_req = 0.0D0
theta1_coil = 45.0D0
theta1_vv = 1.0D0
theta1_vv = 1.0D0 ! 1 Deg
max_vv_stress = 143.0D6
end subroutine init_tfcoil_variables
end module tfcoil_variables
9 changes: 5 additions & 4 deletions tests/unit/test_sctfcoil.py
Original file line number Diff line number Diff line change
Expand Up @@ -14031,7 +14031,7 @@ def test_vv_stress_on_quench():
d_vv=0.12, # for 6.6 restistance -> lambda2 = 2.1
)
)
== 56835032.21809308
== 57045874.69917925
)


Expand Down Expand Up @@ -14059,8 +14059,9 @@ def test_vv_stress_on_quench_integration(sctfcoil, monkeypatch):
monkeypatch.setattr(sctfcoil_module, "a_tf_steel", 0.55) # Section 3

# Sum from Section 3
monkeypatch.setattr(sctfcoil_module, "a_case_front", 0.47)
monkeypatch.setattr(sctfcoil_module, "a_case_nose", 0.47)
monkeypatch.setattr(sctfcoil_module, "a_case_front", 0.42)
monkeypatch.setattr(sctfcoil_module, "a_case_nose", 0.42)
monkeypatch.setattr(sctfcoil_module, "t_lat_case_av", 0.05)

monkeypatch.setattr(build_variables, "vgap_xpoint_divertor", 0.05) # Baseline 2018
monkeypatch.setattr(build_variables, "shldtth", 0.3) # Baseline 2018
Expand Down Expand Up @@ -14094,4 +14095,4 @@ def test_vv_stress_on_quench_integration(sctfcoil, monkeypatch):

sctfcoil.vv_stress_on_quench()

assert pytest.approx(sctfcoil_module.vv_stress_quench) == 56834395.24352395
assert pytest.approx(sctfcoil_module.vv_stress_quench) == 56893800.120420754
Loading