Skip to content
Merged
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
95 changes: 37 additions & 58 deletions freegs4e/equilibrium.py
Original file line number Diff line number Diff line change
Expand Up @@ -157,9 +157,6 @@ def __init__(
nx, ny, generator, nlevels=1, ncycle=1, niter=2, direct=True
)

# separatrix data not yet calculated
self._separatrix_data_flag = False

def create_psi_plasma_default(
self, adaptive_centre=False, gpars=(0.5, 0.5, 0, 2)
):
Expand Down Expand Up @@ -914,10 +911,8 @@ def separatrix_area(self):
Area of the last closed flux surface (plasma boundary) [m^2].
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return self._sep_area

Expand All @@ -934,10 +929,8 @@ def separatrix_length(self):
Circumference of the last closed flux surface (plasma boundary) [m].
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return self._sep_length

Expand Down Expand Up @@ -1241,10 +1234,8 @@ def geometricAxis(self):
Geometric axis (R,Z) position [m].
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return np.array(
[
Expand Down Expand Up @@ -1296,10 +1287,8 @@ def minorRadius(self):
Minor radius of the plasma [m].
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return (self._sep_Rmax - self._sep_Rmin) / 2

Expand All @@ -1317,9 +1306,8 @@ def aspectRatio(self):
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return (self._sep_Rmax + self._sep_Rmin) / (
self._sep_Rmax - self._sep_Rmin
Expand All @@ -1338,10 +1326,8 @@ def geometricElongation(self):
Geometric elongation of the plasma.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return (self._sep_Zmax - self._sep_Zmin) / (
self._sep_Rmax - self._sep_Rmin
Expand All @@ -1360,10 +1346,8 @@ def geometricElongation_upper(self):
Geometric elongation of the upper part of the plasma.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return (
2
Expand All @@ -1384,10 +1368,8 @@ def geometricElongation_lower(self):
Geometric elongation of the lower part of the plasma.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

return (
2
Expand All @@ -1408,10 +1390,8 @@ def effectiveElongation(self):
Effective elongation of the plasma.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

R_geom = (self._sep_Rmax + self._sep_Rmin) / 2
R_minor = (self._sep_Rmax - self._sep_Rmin) / 2
Expand Down Expand Up @@ -1452,10 +1432,8 @@ def triangularity_upper(self):
Upper triangularity of the plasma.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

R_geom = (self._sep_Rmax + self._sep_Rmin) / 2
R_minor = (self._sep_Rmax - self._sep_Rmin) / 2
Expand All @@ -1475,10 +1453,8 @@ def triangularity_lower(self):
Lower triangularity of the plasma.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

R_geom = (self._sep_Rmax + self._sep_Rmin) / 2
R_minor = (self._sep_Rmax - self._sep_Rmin) / 2
Expand All @@ -1498,10 +1474,8 @@ def triangularity(self):
Triangularity of the plasma.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

R_geom = (self._sep_Rmax + self._sep_Rmin) / 2
R_minor = (self._sep_Rmax - self._sep_Rmin) / 2
Expand Down Expand Up @@ -1531,10 +1505,8 @@ def squareness(self):
Lower inner squareness.
"""

# check if metrics are already calculated
if self._separatrix_data_flag is False:
self._separatrix_metrics() # call function
self._separatrix_data_flag = True # update flag
# calculate separatrix metrics
self._separatrix_metrics()

# create shapely object for plasma core
plasma_boundary = sh.Polygon(self.separatrix())
Expand Down Expand Up @@ -2089,21 +2061,24 @@ def strikepoints(

def dr_sep(
self,
Z_level=None,
sign_flip=False,
):
"""
Computes the radial separation (on the inboard and outboard sides) at the midplane (Z = 0)
Computes the radial separation (on the inboard and outboard sides) at vertical position `Z_level`
between flux surfaces passing through the lower and upper X-points.

The inboard dr_sep is defined as the difference in R (radial position) between
the two flux surfaces intersecting Z = 0 on the inboard side (smaller R values),
the two flux surfaces intersecting `Z_level` on the inboard side (smaller R values),
each corresponding to one of the X-points.

Similarly, the outboard dr_sep is the radial separation of the same flux surfaces
on the outboard side (larger R values).

Parameters
----------
Z_level : float
Vertical height at which to calculate dr_sep [m].
sign_flip : bool
By convention, the function assumes the first X-point corresponds to the lower
X-point, and the second to the upper X-point. So for an upper single null plasma,
Expand All @@ -2117,6 +2092,10 @@ def dr_sep(
[dr_sep_in, dr_sep_out].
"""

# check if Z position specified, else set to Zgeometric
if Z_level is None:
Z_level = self.Zgeometric()

# extract required quantities
psi_bndry = self.psi_bndry
xpts = self.xpt
Expand Down Expand Up @@ -2147,8 +2126,8 @@ def dr_sep(
psi_boundary_lines1.append(item)

# use the shapely package to find where each psi_boundary_line intersects
# the Z = 0 line
R_limits = np.array([[self.Rmin, 0], [self.Rmax, 0]])
# the Z = Z_level line
R_limits = np.array([[self.Rmin, Z_level], [self.Rmax, Z_level]])
R_line = sh.LineString(R_limits)

# set of intersection points for the first flux surface
Expand Down