-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathload_Bmns.py
277 lines (204 loc) · 8.38 KB
/
load_Bmns.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
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
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
275
276
277
# -*- coding: utf-8 -*-
"""
load_Bmns
Reads Fourier harmonics of B from M3D-C1 output in netCDF format
Author: Brendan Carrick Lyons
Date created: Wed Jan 20 17:25:40 2016
Date edited:
"""
import numpy as np
from numpy import pi
import xarray as xr
import xarray.ufuncs as xru
from scipy.interpolate import interp1d
from scipy.interpolate import RectBivariateSpline as RBS
from C1py.geofac import geofac
def load_Bmns(folder='./', phasing=0., slice=0,
cur_up=1., cur_low=1., cur_mid=None,
probeg_up=False, probeg_low=False, probeg_mid=False,
machine=None, iplasma=None, code='m3dc1',
ntor=None,phase=False,Jmn=False):
if machine is 'diiid':
conv = 1.0
elif machine is 'aug':
conv = -1.0
elif machine is 'kstar':
conv = 1.0
if phase is not True:
phase = False
phasing = np.asarray(phasing)
if code is 'all':
Bc1 = load_Bmns(folder=folder, phasing=phasing, slice=slice,
cur_up=cur_up, cur_low=cur_low, machine=machine,
iplasma=iplasma, code='m3dc1', ntor=ntor)
Bmars = load_Bmns(folder=folder, phasing=phasing, slice=slice,
cur_up=cur_up, cur_low=cur_low, machine=machine,
iplasma=iplasma, code='mars', ntor=ntor)
Bipec = load_Bmns(folder=folder, phasing=phasing, slice=slice,
cur_up=cur_up, cur_low=cur_low, machine=machine,
iplasma=iplasma, code='ipec', ntor=ntor)
return (Bc1, Bmars, Bipec)
elif code is 'm3dc1':
if Jmn:
file_up = folder+'/jmn3_upper-'+str(slice)+'.cdf'
file_low = folder+'/jmn3_lower-'+str(slice)+'.cdf'
else:
file_up = folder+'/bmn_upper-'+str(slice)+'.cdf'
file_low = folder+'/bmn_lower-'+str(slice)+'.cdf'
ds_up = xr.open_dataset(file_up)
ds_low = xr.open_dataset(file_low)
ntor = ds_up.attrs['ntor']
if probeg_up:
fac = 2.
else:
fac = geofac(machine=machine,ntor=ntor)
cur_up = fac*cur_up
if probeg_low:
fac = 2.
else:
fac = geofac(machine=machine,ntor=ntor)
cur_low = fac*cur_low
m = ds_up.m
try:
print("Using Bmn netCDF version " +str(ds_up.version))
Psi = ds_up.psi_norm
except AttributeError:
print("Using Bmn netCDF version 0")
Psi = ds_up.psi
Bup = ds_up.bmn_real.data + 1j*ds_up.bmn_imag.data
Blow = ds_low.bmn_real.data + 1j*ds_low.bmn_imag.data
if cur_mid is not None:
if Jmn:
file_mid = folder+'/jmn_middle-'+str(slice)+'.cdf'
else:
file_mid = folder+'/bmn_middle-'+str(slice)+'.cdf'
ds_mid = xr.open_dataset(file_mid)
if probeg_mid:
fac = 2.
else:
fac = geofac(machine=machine,ntor=ntor)
cur_mid = fac*cur_mid
Bmid = ds_mid.bmn_real.data + 1j*ds_mid.bmn_imag.data
if (iplasma is not None) and (slice > 0):
if Jmn:
file_up = folder+'/jmn3_upper-0.cdf'
file_low = folder+'/jmn3_lower-0.cdf'
else:
file_up = folder+'/bmn_upper-0.cdf'
file_low = folder+'/bmn_lower-0.cdf'
ds_vu = xr.open_dataset(file_up)
ds_vl = xr.open_dataset(file_low)
Bup -= ds_vu.bmn_real.data
Bup -= 1j*ds_vu.bmn_imag.data
Blow -= ds_vl.bmn_real.data
Blow -= 1j*ds_vl.bmn_imag.data
if cur_mid is not None:
if Jmn:
file_mid= folder+'/jmn3_middle-0.cdf'
else:
file_mid= folder+'/bmn_middle-0.cdf'
ds_vm = xr.open_dataset(file_mid)
Bmid -= ds_vm.bmn_real.data
Bmid -= 1j*ds_vm.bmn_imag.data
q = ds_up.q.data
elif code is 'mars':
file = folder + '/mars_bmn.nc'
dset = xr.open_dataset(file)
m = dset.m
Psi = dset.s**2.0
if slice is 0:
Bup = dset.vacuum_upper_real.data + 1j*dset.vacuum_upper_imag.data
Blow = dset.vacuum_lower_real.data + 1j*dset.vacuum_lower_imag.data
elif slice is 1:
Bup = dset.total_upper_real.data + 1j*dset.total_upper_imag.data
Blow = dset.total_lower_real.data + 1j*dset.total_lower_imag.data
if (iplasma is not None):
Bup -= dset.vacuum_upper_real.data
Bup -= 1j*dset.vacuum_upper_imag.data
Blow -= dset.vacuum_lower_real.data
Blow -= 1j*dset.vacuum_lower_imag.data
else:
return NotImplemented
Psi1 = dset.q_s**2.0
q1 = dset.q
f = interp1d(Psi1,q1,bounds_error=False)
q = f(Psi)
# flip signs because theta is flipped
conv = -conv
m = -m
q = -q
elif code is 'ipec':
if (slice == 0) or (iplasma is not None):
# get the vacuum field
base = folder+'ipec_vbnormal'
sr = 7; # skip this many rows
cr = 3; # column of real Bmn (zero-based)
ci = 4; # column of imaginary Bmn (zero-based)
(Pv,mv,Bvu,Bvl,qv) = B_ipec(base,sr,cr,ci)
if slice == 1:
# get the total field
base = folder+'ipec_xbnormal';
sr = 8; # skip this many rows
cr = 7; # column of real Bmn
ci = 8; # column of imaginary Bmn
(Pt,mt,Btu,Btl,qt) = B_ipec(base,sr,cr,ci)
if slice == 0:
# just the vacuum field
(Psi,m,Bup,Blow,q) = (Pv,mv,Bvu,Bvl,qv)
elif (slice == 1) and (iplasma is None):
# just the total field
(Psi,m,Bup,Blow,q) = (Pt,mt,Btu,Btl,qt)
else:
# we need to subtract the vacuum from the total field
Rv = RBS(Pv,mv,Bvu.real)
Iv = RBS(Pv,mv,Bvu.imag)
Bvu = Rv(Pt,mt) + 1j*Iv(Pt,mt)
Rv = RBS(Pv,mv,Bvl.real)
Iv = RBS(Pv,mv,Bvl.imag)
Bvl = Rv(Pt,mt) + 1j*Iv(Pt,mt)
(Psi,m,Bup,Blow,q) = (Pt,mt,Btu-Bvu,Btl-Bvl,qt)
conv = -conv
m = -m
q = -q
else:
return NotImplemented
# Quantities used for all code results
cur = xr.DataArray(np.cos(pi*phasing/180.) + conv*1j*np.sin(pi*phasing/180.),
coords=[('phasing',phasing)])
Bmn_up = cur_up*xr.DataArray(Bup,coords=[('Psi',Psi),('m',m)])
Bmn_low = cur_low*xr.DataArray(Blow,coords=[('Psi',Psi),('m',m)])
Bmn = cur*Bmn_up + Bmn_low
if cur_mid is not None:
cur2 = xr.DataArray(np.cos(pi*phasing/180.) + conv*1j*np.sin(pi*phasing/180.),
coords=[('phasing2',phasing)])
Bmn_mid = cur_mid*xr.DataArray(Bmid,coords=[('Psi',Psi),('m',m)])
Bmn = Bmn + cur2*Bmn_mid
if phase:
Bmn = xru.angle(Bmn,deg=True) % 360.
else:
Bmn = np.abs(Bmn)
q = xr.DataArray(q,coords=[('Psi',Psi)])
if Jmn:
Bmns = xr.Dataset({'Jmn':Bmn,'q':q})
else:
Bmns = xr.Dataset({'Bmn':Bmn,'q':q})
Bmns.attrs['ntor'] = ntor
Bmns.attrs['phase'] = phase
return Bmns
def B_ipec(base,sr,cr,ci):
file_up = base+'_iu.out'
file_low = base+'_il.out'
Aup = np.loadtxt(file_up,skiprows=sr)
Psi = Aup[:,0]
q = Aup[:,1]
m = Aup[:,2].astype(int)
Bup = Aup[:,cr] + 1j*Aup[:,ci]
Alow = np.loadtxt(file_low,skiprows=sr)
Blow = Alow[:,cr] + 1j*Alow[:,ci]
Mm = max(m) - min(m) + 1
Psi = Psi.reshape(-1,Mm)[:,0]
q = q.reshape(-1,Mm)[:,0]
m = m.reshape(-1,Mm)[0,:]
Bup = 1.e4*Bup.reshape(-1,Mm)
Blow = 1.e4*Blow.reshape(-1,Mm)
return (Psi,m,Bup,Blow,q)