-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathtok_geqdsk.py
executable file
·193 lines (158 loc) · 5.59 KB
/
tok_geqdsk.py
1
2
3
4
5
6
7
8
9
10
11
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
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
"""
Read and write 'G' formatted equilibria. This is an R-Z free boundary
format.
Format of G-EQDSK file is specified here:
https://fusion.gat.com/THEORY/efit/g_eqdsk.html
"""
import numpy as np
from utils import file_numbers, writef
def read(f):
""" Reads a G-EQDSK file
Parameters
----------
f = Input file. Can either be a file-like object,
or a string. If a string, then treated as a file name
and opened.
Returns
-------
"""
if isinstance(f, str):
# If the input is a string, treat as file name
with open(f) as fh: # Ensure file is closed
return read(fh) # Call again with file object
# Read the first line, which should contain the mesh sizes
desc = f.readline()
if not desc:
raise IOError("Cannot read from input file")
s = desc.split() # Split by whitespace
if len(s) < 3:
raise IOError("First line must contain at least 3 numbers")
idum = int(s[-3])
nxefit = int(s[-2])
nyefit = int(s[-1])
# Use a generator to read numbers
token = file_numbers(f)
xdim = float(token.next())
zdim = float(token.next())
rcentr = float(token.next())
rgrid1 = float(token.next())
zmid = float(token.next())
rmagx = float(token.next())
zmagx = float(token.next())
simagx = float(token.next())
sibdry = float(token.next())
bcentr = float(token.next())
cpasma = float(token.next())
simagx = float(token.next())
xdum = float(token.next())
rmagx = float(token.next())
xdum = float(token.next())
zmagx = float(token.next())
xdum = float(token.next())
sibdry = float(token.next())
xdum = float(token.next())
xdum = float(token.next())
# Read arrays
def read_array(n, name="Unknown"):
data = np.zeros([n])
try:
for i in np.arange(n):
data[i] = float(token.next())
except:
raise IOError("Failed reading array '"+name+"' of size ", n)
return data
def read_2d(nx, ny, name="Unknown"):
data = np.zeros([nx, ny])
for i in np.arange(nx):
data[i,:] = read_array(ny, name+"["+str(i)+"]")
return data
fpol = read_array(nxefit, "fpol")
pres = read_array(nxefit, "pres")
workk1 = read_array(nxefit, "workk1")
workk2 = read_array(nxefit, "workk2")
psi = read_2d(nxefit, nyefit, "psi")
qpsi = read_array(nxefit, "qpsi")
# Read boundary and limiters, if present
nbdry = int(token.next())
nlim = int(token.next())
if nbdry > 0:
rbdry = np.zeros([nbdry])
zbdry = np.zeros([nbdry])
for i in range(nbdry):
rbdry[i] = float(token.next())
zbdry[i] = float(token.next())
else:
rbdry = [0]
zbdry = [0]
if nlim > 0:
xlim = np.zeros([nlim])
ylim = np.zeros([nlim])
for i in range(nlim):
xlim[i] = float(token.next())
ylim[i] = float(token.next())
else:
xlim = [0]
ylim = [0]
# Construct R-Z mesh
r = np.zeros([nxefit, nyefit])
z = r.copy()
for i in range(nxefit):
r[i,:] = rgrid1 + xdim*i/float(nxefit-1)
for j in range(nyefit):
z[:,j] = (zmid-0.5*zdim) + zdim*j/float(nyefit-1)
# Create dictionary of values to return
result = {'nx': nxefit, 'ny':nyefit, # Number of horizontal and vertical points
'r':r, 'z':z, # Location of the grid-poinst
'rdim':xdim, 'zdim':zdim, # Size of the domain in meters
'rcentr':rcentr, 'bcentr':bcentr, # Reference vacuum toroidal field (m, T)
'rgrid1':rgrid1, # R of left side of domain
'zmid':zmid, # Z at the middle of the domain
'rmagx':rmagx, 'zmagx':zmagx, # Location of magnetic axis
'simagx':simagx, # Poloidal flux at the axis (Weber / rad)
'sibdry':sibdry, # Poloidal flux at plasma boundary (Weber / rad)
'cpasma':cpasma,
'psi':psi, # Poloidal flux in Weber/rad on grid points
'fpol':fpol, # Poloidal current function on uniform flux grid
'pressure':pres, # Plasma pressure in nt/m^2 on uniform flux grid
'qpsi':qpsi, # q values on uniform flux grid
'nbdry':nbdry, 'rbdry':rbdry, 'zbdry':zbdry, # Plasma boundary
'nlim':nlim, 'xlim':xlim, 'ylim':ylim} # Wall boundary
return result
def write(f, data):
""" Write a G-EQDSK file
"""
if isinstance(f, str):
# If the input is a string, treat as file name
with open(f, "w") as fh: # Ensure file is closed
return write(fh, data) # Call again with file object
nx = int(data['nx'])
ny = int(data['ny'])
# Write description
f.write("geqdsk 0" + str(nx) + str(ny))
writef(f, data['rdim'])
writef(f, data['zdim'])
writef(f, data['rcentr'])
writef(f, data['rgrid1'])
writef(f, data['zmid'])
writef(f, data['rmagx'])
writef(f, data['zmagx'])
writef(f, data['simagx'])
writef(f, data['sibdry'])
writef(f, data['bcentr'])
writef(f, data['cpasma'])
writef(f, data['simagx'])
writef(f, 0.0)
writef(f, data['rmagx'])
writef(f, 0.0)
writef(f, data['zmagx'])
writef(f, 0.0)
writef(f, data['sibdry'])
writef(f, 0.0)
writef(f, 0.0)
# Write the arrays
writef(f, data['fpol'])
writef(f, data['pressure'])
writef(f, np.zeros(nx)) # Workk1
writef(f, np.zeros(nx)) # Workk2
writef(f, data['psi'])
writef(f, data['qpsi'])