-
Notifications
You must be signed in to change notification settings - Fork 14
/
vmec.py
1876 lines (1663 loc) · 69.8 KB
/
vmec.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
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
357
358
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
453
454
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
551
552
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
672
673
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
753
754
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
902
903
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
998
999
1000
"""Functions and classes for interfacing with VMEC equilibria."""
import io
import os
import warnings
import matplotlib.pyplot as plt
import numpy as np
from netCDF4 import Dataset, stringtochar
from scipy import integrate, interpolate, optimize
from scipy.constants import mu_0
from desc.basis import DoubleFourierSeries
from desc.compat import ensure_positive_jacobian
from desc.compute.utils import surface_averages
from desc.equilibrium import Equilibrium
from desc.geometry import FourierRZToroidalSurface
from desc.grid import Grid, LinearGrid
from desc.objectives import (
ObjectiveFunction,
get_fixed_axis_constraints,
get_fixed_boundary_constraints,
maybe_add_self_consistency,
)
from desc.objectives.utils import factorize_linear_constraints
from desc.profiles import PowerSeriesProfile, SplineProfile
from desc.transform import Transform
from desc.utils import Timer
from desc.vmec_utils import (
fourier_to_zernike,
ptolemy_identity_fwd,
ptolemy_identity_rev,
zernike_to_fourier,
)
class VMECIO:
"""Performs input from VMEC netCDF files to DESC Equilibrium and vice-versa."""
@classmethod
def load(
cls, path, L=None, M=None, N=None, spectral_indexing="ansi", profile="iota"
):
"""Load a VMEC netCDF file as an Equilibrium.
Parameters
----------
path : str
File path of input data.
L : int, optional
Radial resolution. Default determined by index.
M : int, optional
Poloidal resolution. Default = MPOL-1 from VMEC solution.
N : int, optional
Toroidal resolution. Default = NTOR from VMEC solution.
spectral_indexing : str, optional
Type of Zernike indexing scheme to use. (Default = ``'ansi'``)
profile : {"iota", "current"}
Which profile to use as the equilibrium constraint. (Default = ``'iota'``)
Returns
-------
eq: Equilibrium
Equilibrium that resembles the VMEC data.
Notes
-----
To ensure compatibility with different profile representations in VMEC,
the resulting equilibrium will always have SplineProfile types for all profiles
"""
assert profile in ["iota", "current"]
file = Dataset(path, mode="r")
inputs = {}
version = file.variables["version_"][0]
if version < 9:
warnings.warn(
"VMEC output appears to be from version {}".format(str(version))
+ " while DESC is only designed for compatibility with VMEC version"
+ " 9. Some data may not be loaded correctly."
)
# parameters
inputs["Psi"] = float(file.variables["phi"][-1])
inputs["NFP"] = int(file.variables["nfp"][0])
inputs["M"] = M if M is not None else int(file.variables["mpol"][0] - 1)
inputs["N"] = N if N is not None else int(file.variables["ntor"][0])
inputs["spectral_indexing"] = spectral_indexing
default_L = {
"ansi": inputs["M"],
"fringe": 2 * inputs["M"],
}
inputs["L"] = L if L is not None else default_L[inputs["spectral_indexing"]]
# data
xm = file.variables["xm"][:].filled()
xn = file.variables["xn"][:].filled() / inputs["NFP"]
rmnc = file.variables["rmnc"][:].filled()
zmns = file.variables["zmns"][:].filled()
lmns = file.variables["lmns"][:].filled()
try:
rmns = file.variables["rmns"][:].filled()
zmnc = file.variables["zmnc"][:].filled()
lmnc = file.variables["lmnc"][:].filled()
inputs["sym"] = False
except KeyError:
rmns = np.zeros_like(rmnc)
zmnc = np.zeros_like(zmns)
lmnc = np.zeros_like(lmns)
inputs["sym"] = True
# profiles
r = np.sqrt(np.linspace(0, 1, file.variables["ns"][:].filled()))
pres = file.variables["presf"][:].filled()
inputs["pressure"] = SplineProfile(pres, r, name="pressure")
if profile == "iota":
iota = file.variables["iotaf"][:].filled()
inputs["iota"] = SplineProfile(iota, r, name="iota")
inputs["current"] = None
if profile == "current":
curr = 2 * np.pi / mu_0 * file.variables["buco"][:].filled()
inputs["current"] = SplineProfile(curr, r, name="current")
inputs["iota"] = None
# boundary
m, n, Rb_lmn = ptolemy_identity_fwd(xm, xn, s=rmns[-1, :], c=rmnc[-1, :])
m, n, Zb_lmn = ptolemy_identity_fwd(xm, xn, s=zmns[-1, :], c=zmnc[-1, :])
surface = np.vstack((np.zeros_like(m), m, n, Rb_lmn, Zb_lmn)).T
# need to create surface object here so we can tell it not to flip the
# orientation yet. If we did it here, it would mess up the self-consistency
# stuff later
inputs["surface"] = FourierRZToroidalSurface(
surface[:, 3],
surface[:, 4],
surface[:, 1:3].astype(int),
surface[:, 1:3].astype(int),
inputs["NFP"],
inputs["sym"],
check_orientation=False,
)
# axis
rax_cc = file.variables["raxis_cc"][:].filled()
zax_cs = file.variables["zaxis_cs"][:].filled()
try:
rax_cs = file.variables["raxis_cs"][:].filled()
rax_cc = file.variables["zaxis_cc"][:].filled()
except KeyError:
rax_cs = np.zeros_like(rax_cc)
zax_cc = np.zeros_like(zax_cs)
rax = np.concatenate([-rax_cs[1:][::-1], rax_cc])
zax = np.concatenate([-zax_cs[1:][::-1], zax_cc])
nax = len(rax_cc) - 1
nax = np.arange(-nax, nax + 1)
inputs["axis"] = np.vstack([nax, rax, zax]).T
file.close()
# initialize Equilibrium
eq = Equilibrium(**inputs, check_orientation=False)
# R
m, n, R_mn = ptolemy_identity_fwd(xm, xn, s=rmns, c=rmnc)
eq.R_lmn = fourier_to_zernike(m, n, R_mn, eq.R_basis)
# Z
m, n, Z_mn = ptolemy_identity_fwd(xm, xn, s=zmns, c=zmnc)
eq.Z_lmn = fourier_to_zernike(m, n, Z_mn, eq.Z_basis)
# lambda
m, n, L_mn = ptolemy_identity_fwd(xm, xn, s=lmns, c=lmnc)
eq.L_lmn = fourier_to_zernike(m, n, L_mn, eq.L_basis)
# apply boundary conditions
constraints = get_fixed_axis_constraints(
profiles=False, eq=eq
) + get_fixed_boundary_constraints(eq=eq)
constraints = maybe_add_self_consistency(eq, constraints)
objective = ObjectiveFunction(constraints)
objective.build(verbose=0)
_, _, _, _, _, project, recover = factorize_linear_constraints(
objective, objective
)
args = objective.unpack_state(recover(project(objective.x(eq))), False)[0]
eq.params_dict = args
# now we flip the orientation at the very end
with warnings.catch_warnings():
warnings.filterwarnings(
"ignore", message="Left handed coordinates detected"
)
eq = ensure_positive_jacobian(eq)
return eq
@classmethod
def save(cls, eq, path, surfs=128, verbose=1): # noqa: C901 - FIXME - simplify
"""Save an Equilibrium as a netCDF file in the VMEC format.
Parameters
----------
eq : Equilibrium
Equilibrium to save.
path : str
File path of output data.
surfs: int
Number of flux surfaces to interpolate at (Default = 128).
verbose: int
Level of output (Default = 1).
* 0: no output
* 1: status of quantities computed
* 2: as above plus timing information
Returns
-------
None
"""
timer = Timer()
timer.start("Total time")
""" VMEC netCDF file is generated in VMEC2000/Sources/Input_Output/wrout.f
see lines 300+ for full list of included variables
"""
file = Dataset(path, mode="w", format="NETCDF3_64BIT_OFFSET")
Psi = eq.Psi
NFP = eq.NFP
M = eq.M
N = eq.N
M_nyq = M + 4
N_nyq = N + 2 if N > 0 else 0
# VMEC radial coordinate: s = rho^2 = Psi / Psi(LCFS)
s_full = np.linspace(0, 1, surfs)
s_half = s_full[0:-1] + 0.5 / (surfs - 1)
r_full = np.sqrt(s_full)
r_half = np.sqrt(s_half)
# dimensions
file.createDimension("radius", surfs) # number of flux surfaces
file.createDimension(
"mn_mode", (2 * N + 1) * M + N + 1
) # number of Fourier modes
file.createDimension(
"mn_mode_nyq", (2 * N_nyq + 1) * M_nyq + N_nyq + 1
) # number of Nyquist Fourier modes
file.createDimension("n_tor", N + 1) # number of toroidal Fourier modes
file.createDimension("preset", 21) # dimension of profile inputs
file.createDimension("ndfmax", 101) # used for am_aux & ai_aux
file.createDimension("time", 100) # used for fsq* & wdot
file.createDimension("dim_00001", 1)
file.createDimension("dim_00020", 20)
file.createDimension("dim_00100", 100)
file.createDimension("dim_00200", 200)
# compute data
timer.start("compute")
if verbose > 0:
print("Computing data")
grid_axis = LinearGrid(M=M_nyq, N=N_nyq, rho=np.array([0.0]), NFP=NFP)
grid_lcfs = LinearGrid(M=M_nyq, N=N_nyq, rho=np.array([1.0]), NFP=NFP)
grid_half = LinearGrid(M=M_nyq, N=N_nyq, NFP=NFP, rho=r_half)
grid_full = LinearGrid(M=M_nyq, N=N_nyq, NFP=NFP, rho=r_full)
data_quad = eq.compute(
["R0/a", "V", "<|B|>_rms", "<beta>_vol", "<beta_pol>_vol", "<beta_tor>_vol"]
)
data_axis = eq.compute(["G", "p", "R", "<|B|^2>"], grid=grid_axis)
data_lcfs = eq.compute(["G", "I", "R", "Z"], grid=grid_lcfs)
data_half = eq.compute(
[
"B_rho",
"B_theta",
"B_zeta",
"G",
"I",
"J",
"iota",
"p",
"sqrt(g)",
"V_r(r)",
"|B|",
"<|B|^2>",
],
grid=grid_half,
)
data_full = eq.compute(
[
"B_rho",
"B_theta",
"B_zeta",
"current",
"D_Mercier",
"iota",
"J",
"J^theta*sqrt(g)",
"J^zeta",
"p",
"sqrt(g)",
"<|B|^2>",
"<J*B>",
],
grid=grid_full,
)
timer.stop("compute")
if verbose > 1:
timer.disp("compute")
# parameters
timer.start("parameters")
if verbose > 0:
print("Saving parameters")
version_ = file.createVariable("version_", np.float64)
version_[:] = 9 # VMEC version at the time of this writing
input_extension = file.createVariable("input_extension", "S1", ("dim_00100",))
input_extension[:] = stringtochar(
np.array([" " * 100], "S" + str(file.dimensions["dim_00100"].size))
) # VMEC input filename: input.[input_extension]
# TODO: instead of hard-coding for fixed-boundary, also allow for free-boundary?
mgrid_mode = file.createVariable("mgrid_mode", "S1", ("dim_00001",))
mgrid_mode[:] = stringtochar(
np.array([""], "S" + str(file.dimensions["dim_00001"].size))
)
mgrid_file = file.createVariable("mgrid_file", "S1", ("dim_00200",))
mgrid_file[:] = stringtochar(
np.array(["none" + " " * 196], "S" + str(file.dimensions["dim_00200"].size))
)
ier_flag = file.createVariable("ier_flag", np.int32)
ier_flag.long_name = (
"error flag (DESC always outputs 0; "
+ "manually check for a good equilibrium solution)"
)
ier_flag[:] = 0
lfreeb = file.createVariable("lfreeb__logical__", np.int32)
lfreeb.long_name = "free boundary logical (0 = fixed boundary)"
lfreeb[:] = 0
lrecon = file.createVariable("lrecon__logical__", np.int32)
lrecon.long_name = "reconstruction logical (0 = no reconstruction)"
lrecon[:] = 0
lrfp = file.createVariable("lrfp__logical__", np.int32)
lrfp.long_name = "reverse-field pinch logical (0 = not an RFP)"
lrfp[:] = 0
lasym = file.createVariable("lasym__logical__", np.int32)
lasym.long_name = "asymmetry logical (0 = stellarator symmetry)"
lasym[:] = int(not eq.sym)
nfp = file.createVariable("nfp", np.int32)
nfp.long_name = "number of field periods"
nfp[:] = NFP
ns = file.createVariable("ns", np.int32)
ns.long_name = "number of flux surfaces"
ns[:] = surfs
mpol = file.createVariable("mpol", np.int32)
mpol.long_name = "number of poloidal Fourier modes"
mpol[:] = M + 1
ntor = file.createVariable("ntor", np.int32)
ntor.long_name = "number of positive toroidal Fourier modes"
ntor[:] = N
mnmax = file.createVariable("mnmax", np.int32)
mnmax.long_name = "total number of Fourier modes"
mnmax[:] = file.dimensions["mn_mode"].size
xm = file.createVariable("xm", np.float64, ("mn_mode",))
xm.long_name = "poloidal mode numbers"
xm[:] = np.tile(np.linspace(0, M, M + 1), (2 * N + 1, 1)).T.flatten()[
-file.dimensions["mn_mode"].size :
]
xn = file.createVariable("xn", np.float64, ("mn_mode",))
xn.long_name = "toroidal mode numbers"
xn[:] = np.tile(np.linspace(-N, N, 2 * N + 1) * NFP, M + 1)[
-file.dimensions["mn_mode"].size :
]
mnmax_nyq = file.createVariable("mnmax_nyq", np.int32)
mnmax_nyq.long_name = "total number of Nyquist Fourier modes"
mnmax_nyq[:] = file.dimensions["mn_mode_nyq"].size
xm_nyq = file.createVariable("xm_nyq", np.float64, ("mn_mode_nyq",))
xm_nyq.long_name = "poloidal Nyquist mode numbers"
xm_nyq[:] = np.tile(
np.linspace(0, M_nyq, M_nyq + 1), (2 * N_nyq + 1, 1)
).T.flatten()[-file.dimensions["mn_mode_nyq"].size :]
xn_nyq = file.createVariable("xn_nyq", np.float64, ("mn_mode_nyq",))
xn_nyq.long_name = "toroidal Nyquist mode numbers"
xn_nyq[:] = np.tile(np.linspace(-N_nyq, N_nyq, 2 * N_nyq + 1) * NFP, M_nyq + 1)[
-file.dimensions["mn_mode_nyq"].size :
]
signgs = file.createVariable("signgs", np.int32)
signgs.long_name = "sign of coordinate system Jacobian"
signgs[:] = -1 # VMEC always uses a negative Jacobian
gamma = file.createVariable("gamma", np.float64)
gamma.long_name = "compressibility index (0 = pressure prescribed)"
gamma[:] = 0
nextcur = file.createVariable("nextcur", np.int32)
nextcur.long_name = "number of coils (external currents)"
nextcur[:] = 0 # hard-coded assuming fixed-boundary solve
# TODO: add option for saving spline profiles
power_series = stringtochar(
np.array(
["power_series" + " " * 8], "S" + str(file.dimensions["dim_00020"].size)
)
)
pmass_type = file.createVariable("pmass_type", "S1", ("dim_00020",))
pmass_type.long_name = "parameterization of pressure function"
pmass_type[:] = power_series
piota_type = file.createVariable("piota_type", "S1", ("dim_00020",))
piota_type.long_name = "parameterization of rotational transform function"
piota_type[:] = power_series
pcurr_type = file.createVariable("pcurr_type", "S1", ("dim_00020",))
pcurr_type.long_name = "parameterization of current density function"
pcurr_type[:] = power_series
# scalars computed on a Quadrature grid
Rmajor_p = file.createVariable("Rmajor_p", np.float64)
Rmajor_p.long_name = "major radius"
Rmajor_p.units = "m"
Rmajor_p[:] = data_quad["R0"]
Aminor_p = file.createVariable("Aminor_p", np.float64)
Aminor_p.long_name = "minor radius"
Aminor_p.units = "m"
Aminor_p[:] = data_quad["a"]
aspect = file.createVariable("aspect", np.float64)
aspect.long_name = "aspect ratio = R_major / A_minor"
aspect.units = "None"
aspect[:] = data_quad["R0/a"]
volume_p = file.createVariable("volume_p", np.float64)
volume_p.long_name = "plasma volume"
volume_p.units = "m^3"
volume_p[:] = data_quad["V"]
volavgB = file.createVariable("volavgB", np.float64)
volavgB.long_name = "volume average magnetic field, root mean square"
volavgB.units = "T"
volavgB[:] = data_quad["<|B|>_rms"]
betatotal = file.createVariable("betatotal", np.float64)
betatotal.long_name = "normalized plasma pressure"
betatotal.units = "None"
betatotal[:] = data_quad["<beta>_vol"]
betapol = file.createVariable("betapol", np.float64)
betapol.long_name = "normalized poloidal plasma pressure"
betapol.units = "None"
betapol[:] = data_quad["<beta_pol>_vol"]
betator = file.createVariable("betator", np.float64)
betator.long_name = "normalized toroidal plasma pressure"
betator.units = "None"
betator[:] = data_quad["<beta_tor>_vol"]
# scalars computed at the magnetic axis
rbtor0 = file.createVariable("rbtor0", np.float64)
rbtor0.long_name = "<R*B_tor> on axis"
rbtor0.units = "T*m"
rbtor0[:] = data_axis["G"][0]
b0 = file.createVariable("b0", np.float64)
b0.long_name = "average B_tor on axis"
b0.units = "T"
b0[:] = data_axis["G"][0] / data_axis["R"][0]
betaxis = file.createVariable("betaxis", np.float64)
betaxis.long_name = "2 * mu_0 * pressure / <|B|^2> on the magnetic axis"
betaxis.units = "None"
betaxis[:] = 2 * mu_0 * data_axis["p"][0] / data_axis["<|B|^2>"][0]
# scalars computed at the last closed flux surface
ctor = file.createVariable("ctor", np.float64)
ctor.long_name = "total toroidal plasma current"
ctor.units = "A"
ctor[:] = data_lcfs["I"][0] * 2 * np.pi / mu_0
rbtor = file.createVariable("rbtor", np.float64)
rbtor.long_name = "<R*B_tor> on last closed flux surface"
rbtor.units = "T*m"
rbtor[:] = data_lcfs["G"][0]
# half mesh quantities
pres = file.createVariable("pres", np.float64, ("radius",))
pres.long_name = "pressure on half mesh"
pres.units = "Pa"
pres[0] = 0
pres[1:] = grid_half.compress(data_half["p"])
mass = file.createVariable("mass", np.float64, ("radius",))
mass.long_name = "mass on half mesh"
mass.units = "Pa"
mass[:] = pres[:]
iotas = file.createVariable("iotas", np.float64, ("radius",))
iotas.long_name = "rotational transform on half mesh"
iotas.units = "None"
iotas[0] = 0
iotas[1:] = -grid_half.compress(data_half["iota"]) # - for negative Jacobian
phips = file.createVariable("phips", np.float64, ("radius",))
phips.long_name = "d(phi)/ds * -1/2pi: toroidal flux derivative, on half mesh"
phips[0] = 0
phips[1:] = -Psi * np.ones((surfs - 1,)) / (2 * np.pi)
buco = file.createVariable("buco", np.float64, ("radius",))
buco.long_name = "Boozer toroidal current I, on half mesh"
buco.units = "T*m"
buco[1:] = -grid_half.compress(data_half["I"]) # - for negative Jacobian
buco[0] = 0
bvco = file.createVariable("bvco", np.float64, ("radius",))
bvco.long_name = "Boozer poloidal current G, on half mesh"
bvco.units = "T*m"
bvco[1:] = grid_half.compress(data_half["G"])
bvco[0] = 0
beta_vol = file.createVariable("beta_vol", np.float64, ("radius",))
beta_vol.long_name = "2 * mu_0 * pressure / <|B|^2>, on half mesh"
beta_vol.units = "None"
beta_vol[1:] = (
2
* mu_0
* grid_half.compress(data_half["p"])
/ grid_half.compress(data_half["<|B|^2>"])
)
beta_vol[0] = 0.0
vp = file.createVariable("vp", np.float64, ("radius",))
vp.long_name = "dV/ds normalized by 4*pi^2, on half mesh"
vp.units = "m^3"
vp[1:] = grid_half.compress(data_half["V_r(r)"]) / (
8 * np.pi**2 * grid_half.compress(data_half["rho"])
)
vp[0] = 0
# full mesh quantities
presf = file.createVariable("presf", np.float64, ("radius",))
presf.long_name = "pressure on full mesh"
presf.units = "Pa"
presf[:] = grid_full.compress(data_full["p"])
iotaf = file.createVariable("iotaf", np.float64, ("radius",))
iotaf.long_name = "rotational transform on full mesh"
iotaf.units = "None"
iotaf[:] = -grid_full.compress(data_full["iota"]) # - for negative Jacobian
q_factor = file.createVariable("q_factor", np.float64, ("radius",))
q_factor.long_name = "inverse rotational transform on full mesh"
q_factor.units = "None"
q_factor[:] = 1 / iotaf[:]
phi = file.createVariable("phi", np.float64, ("radius",))
phi.long_name = "toroidal flux, on full mesh"
phi.units = "Wb"
phi[:] = np.linspace(0, Psi, surfs)
phipf = file.createVariable("phipf", np.float64, ("radius",))
phipf.long_name = "d(phi)/ds: toroidal flux derivative, on full mesh"
phipf[:] = Psi * np.ones((surfs,))
chi = file.createVariable("chi", np.float64, ("radius",))
chi.long_name = "poloidal flux, on full mesh"
chi.units = "Wb"
chi[:] = (
-2 # - for negative Jacobian
* Psi
* integrate.cumulative_trapezoid(r_full * iotaf[:], r_full, initial=0)
)
chipf = file.createVariable("chipf", np.float64, ("radius",))
chipf.long_name = "d(chi)/ds: poloidal flux derivative, on full mesh"
chipf[:] = phipf[:] * iotaf[:]
am = file.createVariable("am", np.float64, ("preset",))
am.long_name = "pressure coefficients"
am.units = "Pa"
am[:] = np.zeros((file.dimensions["preset"].size,))
# only using up to 10th order to avoid poor conditioning
am[:11] = PowerSeriesProfile.from_values(
s_full, grid_full.compress(data_full["p"]), order=10, sym=False
).params
ai = file.createVariable("ai", np.float64, ("preset",))
ai.long_name = "rotational transform coefficients"
ai[:] = np.zeros((file.dimensions["preset"].size,))
if eq.iota is not None:
# only using up to 10th order to avoid poor conditioning
ai[:11] = -PowerSeriesProfile.from_values( # - for negative Jacobian
s_full, grid_full.compress(data_full["iota"]), order=10, sym=False
).params
ac = file.createVariable("ac", np.float64, ("preset",))
ac.long_name = "normalized toroidal current density coefficients"
ac[:] = np.zeros((file.dimensions["preset"].size,))
if eq.current is not None:
# only using up to 10th order to avoid poor conditioning
ac[:11] = PowerSeriesProfile.from_values(
s_full, grid_full.compress(data_full["current"]), order=10, sym=False
).params
bdotb = file.createVariable("bdotb", np.float64, ("radius",))
bdotb.long_name = "flux surface average of magnetic field squared, on full mesh"
bdotb.units = "T^2"
bdotb[:] = grid_full.compress(data_full["<|B|^2>"])
bdotb[0] = 0
jdotb = file.createVariable("jdotb", np.float64, ("radius",))
jdotb.long_name = "flux surface average of J*B, on full mesh"
jdotb.units = "N/m^3"
jdotb[:] = grid_full.compress(data_full["<J*B>"])
jdotb[0] = 0
jcuru = file.createVariable("jcuru", np.float64, ("radius",))
jcuru.long_name = "flux surface average of sqrt(g)*J^theta, on full mesh"
jcuru.units = "A/m^3"
jcuru[:] = -surface_averages( # - for negative Jacobian
grid_full,
data_full["J^theta*sqrt(g)"] / (2 * data_full["rho"]),
sqrt_g=data_full["sqrt(g)"],
expand_out=False,
)
jcuru[0] = 0
jcurv = file.createVariable("jcurv", np.float64, ("radius",))
jcuru.long_name = "flux surface average of sqrt(g)*J^zeta, on full mesh"
jcurv.units = "A/m^3"
jcurv[:] = surface_averages(
grid_full,
data_full["sqrt(g)"] * data_full["J^zeta"] / (2 * data_full["rho"]),
sqrt_g=data_full["sqrt(g)"],
expand_out=False,
)
jcurv[0] = 0
DShear = file.createVariable("DShear", np.float64, ("radius",))
DShear.long_name = (
"Mercier stability criterion magnetic shear term, on full mesh"
)
DShear.units = "1/Wb^2"
DShear[:] = grid_full.compress(data_full["D_shear"])
DShear[0] = 0
DCurr = file.createVariable("DCurr", np.float64, ("radius",))
DCurr.long_name = (
"Mercier stability criterion toroidal current term, on full mesh"
)
DCurr.units = "1/Wb^2"
DCurr[:] = grid_full.compress(data_full["D_current"])
DCurr[0] = 0
DWell = file.createVariable("DWell", np.float64, ("radius",))
DWell.long_name = "Mercier stability criterion magnetic well term, on full mesh"
DWell.units = "1/Wb^2"
DWell[:] = grid_full.compress(data_full["D_well"])
DWell[0] = 0
DGeod = file.createVariable("DGeod", np.float64, ("radius",))
DGeod.long_name = (
"Mercier stability criterion geodesic curvature term, on full mesh"
)
DGeod.units = "1/Wb^2"
DGeod[:] = grid_full.compress(data_full["D_geodesic"])
DGeod[0] = 0
DMerc = file.createVariable("DMerc", np.float64, ("radius",))
DMerc.long_name = "Mercier stability criterion, on full mesh"
DMerc.units = "1/Wb^2"
DMerc[:] = grid_full.compress(data_full["D_Mercier"])
DMerc[0] = 0
timer.stop("parameters")
if verbose > 1:
timer.disp("parameters")
# independent variables (exact conversion)
# R axis
idx = np.where(eq.R_basis.modes[:, 1] == 0)[0]
R0_n = np.zeros((2 * N + 1,))
for k in idx:
(l, m, n) = eq.R_basis.modes[k, :]
R0_n[n + N] += (-2 * (l // 2 % 2) + 1) * eq.R_lmn[k]
raxis_cc = file.createVariable("raxis_cc", np.float64, ("n_tor",))
raxis_cc.long_name = "cos(n*p) component of magnetic axis R coordinate"
raxis_cc.units = "m"
raxis_cc[:] = R0_n[N:]
if not eq.sym:
raxis_cs = file.createVariable("raxis_cs", np.float64, ("n_tor",))
raxis_cs.long_name = "sin(n*p) component of magnetic axis R coordinate"
raxis_cs.units = "m"
raxis_cs[0] = 0.0
raxis_cs[1:] = -R0_n[0:N][::-1]
# Z axis
idx = np.where(eq.Z_basis.modes[:, 1] == 0)[0]
Z0_n = np.zeros((2 * N + 1,))
for k in idx:
(l, m, n) = eq.Z_basis.modes[k, :]
Z0_n[n + N] += (-2 * (l // 2 % 2) + 1) * eq.Z_lmn[k]
zaxis_cs = file.createVariable("zaxis_cs", np.float64, ("n_tor",))
zaxis_cs.long_name = "sin(n*p) component of magnetic axis Z coordinate"
zaxis_cs.units = "m"
zaxis_cs[0] = 0.0
zaxis_cs[1:] = -Z0_n[0:N][::-1]
if not eq.sym:
zaxis_cc = file.createVariable("zaxis_cc", np.float64, ("n_tor",))
zaxis_cc.long_name = "cos(n*p) component of magnetic axis Z coordinate"
zaxis_cc.units = "m"
zaxis_cc[:] = Z0_n[N:]
# R
timer.start("R")
if verbose > 0:
print("Saving R")
rmnc = file.createVariable("rmnc", np.float64, ("radius", "mn_mode"))
rmnc.long_name = "cos(m*t-n*p) component of cylindrical R, on full mesh"
rmnc.units = "m"
if not eq.sym:
rmns = file.createVariable("rmns", np.float64, ("radius", "mn_mode"))
rmns.long_name = "sin(m*t-n*p) component of cylindrical R, on full mesh"
rmns.units = "m"
r1 = np.ones_like(eq.R_lmn)
r1[eq.R_basis.modes[:, 1] < 0] *= -1
m, n, x_mn = zernike_to_fourier(r1 * eq.R_lmn, basis=eq.R_basis, rho=r_full)
xm, xn, s, c = ptolemy_identity_rev(m, n, x_mn)
rmnc[:] = c
if not eq.sym:
rmns[:] = s
timer.stop("R")
if verbose > 1:
timer.disp("R")
# Z
timer.start("Z")
if verbose > 0:
print("Saving Z")
zmns = file.createVariable("zmns", np.float64, ("radius", "mn_mode"))
zmns.long_name = "sin(m*t-n*p) component of cylindrical Z, on full mesh"
zmns.units = "m"
if not eq.sym:
zmnc = file.createVariable("zmnc", np.float64, ("radius", "mn_mode"))
zmnc.long_name = "cos(m*t-n*p) component of cylindrical Z, on full mesh"
zmnc.units = "m"
z1 = np.ones_like(eq.Z_lmn)
z1[eq.Z_basis.modes[:, 1] < 0] *= -1
m, n, x_mn = zernike_to_fourier(z1 * eq.Z_lmn, basis=eq.Z_basis, rho=r_full)
xm, xn, s, c = ptolemy_identity_rev(m, n, x_mn)
zmns[:] = s
if not eq.sym:
zmnc[:] = c
timer.stop("Z")
if verbose > 1:
timer.disp("Z")
# lambda
timer.start("lambda")
if verbose > 0:
print("Saving lambda")
lmns = file.createVariable("lmns", np.float64, ("radius", "mn_mode"))
lmns.long_name = "sin(m*t-n*p) component of lambda, on half mesh"
lmns.units = "rad"
if not eq.sym:
lmnc = file.createVariable("lmnc", np.float64, ("radius", "mn_mode"))
lmnc.long_name = "cos(m*t-n*p) component of lambda, on half mesh"
lmnc.units = "rad"
l1 = np.ones_like(eq.L_lmn)
l1[eq.L_basis.modes[:, 1] >= 0] *= -1
m, n, x_mn = zernike_to_fourier(l1 * eq.L_lmn, basis=eq.L_basis, rho=r_half)
xm, xn, s, c = ptolemy_identity_rev(m, n, x_mn)
lmns[0, :] = 0
lmns[1:, :] = s
if not eq.sym:
lmnc[0, :] = 0
lmnc[1:, :] = c
timer.stop("lambda")
if verbose > 1:
timer.disp("lambda")
# derived quantities (approximate conversion)
sin_basis = DoubleFourierSeries(M=M_nyq, N=N_nyq, NFP=NFP, sym="sin")
cos_basis = DoubleFourierSeries(M=M_nyq, N=N_nyq, NFP=NFP, sym="cos")
full_basis = DoubleFourierSeries(M=M_nyq, N=N_nyq, NFP=NFP, sym=None)
if eq.sym:
sin_transform = Transform(
grid=grid_lcfs, basis=sin_basis, build=False, build_pinv=True
)
cos_transform = Transform(
grid=grid_lcfs, basis=cos_basis, build=False, build_pinv=True
)
def cosfit(x):
y = cos_transform.fit(x)
return np.where(cos_transform.basis.modes[:, 1] < 0, -y, y)
def sinfit(x):
y = sin_transform.fit(x)
return np.where(sin_transform.basis.modes[:, 1] < 0, -y, y)
else:
full_transform = Transform(
grid=grid_lcfs, basis=full_basis, build=False, build_pinv=True
)
def fullfit(x):
y = full_transform.fit(x)
return np.where(full_transform.basis.modes[:, 1] < 0, -y, y)
rmin_surf = file.createVariable("rmin_surf", np.float64)
rmin_surf.long_name = "minimum R coordinate range"
rmin_surf.units = "m"
rmin_surf[:] = np.amin(data_lcfs["R"])
rmax_surf = file.createVariable("rmax_surf", np.float64)
rmax_surf.long_name = "maximum R coordinate range"
rmax_surf.units = "m"
rmax_surf[:] = np.amax(data_lcfs["R"])
zmax_surf = file.createVariable("zmax_surf", np.float64)
zmax_surf.long_name = "maximum Z coordinate range"
zmax_surf.units = "m"
zmax_surf[:] = np.amax(np.abs(data_lcfs["Z"]))
# Jacobian
timer.start("Jacobian")
if verbose > 0:
print("Saving Jacobian")
gmnc = file.createVariable("gmnc", np.float64, ("radius", "mn_mode_nyq"))
gmnc.long_name = "cos(m*t-n*p) component of Jacobian, on half mesh"
gmnc.units = "m"
m = cos_basis.modes[:, 1]
n = cos_basis.modes[:, 2]
if not eq.sym:
gmns = file.createVariable("gmns", np.float64, ("radius", "mn_mode_nyq"))
gmns.long_name = "sin(m*t-n*p) component of Jacobian, on half mesh"
gmns.units = "m"
m = full_basis.modes[:, 1]
n = full_basis.modes[:, 2]
# d(rho)/d(s) = 1/(2*rho)
data = (
(data_half["sqrt(g)"] / (2 * data_half["rho"]))
.reshape(
(grid_half.num_theta, grid_half.num_rho, grid_half.num_zeta), order="F"
)
.transpose((1, 0, 2))
.reshape((grid_half.num_rho, -1), order="F")
)
x_mn = np.zeros((surfs - 1, m.size))
for i in range(surfs - 1):
if eq.sym:
x_mn[i, :] = cosfit(data[i, :])
else:
x_mn[i, :] = fullfit(data[i, :])
xm, xn, s, c = ptolemy_identity_rev(m, n, x_mn)
gmnc[0, :] = 0
gmnc[1:, :] = -c # negative sign for negative Jacobian
if not eq.sym:
gmns[0, :] = 0
gmns[1:, :] = -s
timer.stop("Jacobian")
if verbose > 1:
timer.disp("Jacobian")
# |B|
timer.start("|B|")
if verbose > 0:
print("Saving |B|")
bmnc = file.createVariable("bmnc", np.float64, ("radius", "mn_mode_nyq"))
bmnc.long_name = "cos(m*t-n*p) component of |B|, on half mesh"
bmnc.units = "T"
m = cos_basis.modes[:, 1]
n = cos_basis.modes[:, 2]
if not eq.sym:
bmns = file.createVariable("bmns", np.float64, ("radius", "mn_mode_nyq"))
bmns.long_name = "sin(m*t-n*p) component of |B|, on half mesh"
bmns.units = "T"
m = full_basis.modes[:, 1]
n = full_basis.modes[:, 2]
data = (
data_half["|B|"]
.reshape(
(grid_half.num_theta, grid_half.num_rho, grid_half.num_zeta), order="F"
)
.transpose((1, 0, 2))
.reshape((grid_half.num_rho, -1), order="F")
)
x_mn = np.zeros((surfs - 1, m.size))
for i in range(surfs - 1):
if eq.sym:
x_mn[i, :] = cosfit(data[i, :])
else:
x_mn[i, :] = full_transform.fit(data[i, :])
xm, xn, s, c = ptolemy_identity_rev(m, n, x_mn)
bmnc[0, :] = 0
bmnc[1:, :] = c
if not eq.sym:
bmns[0, :] = 0
bmns[1:, :] = s
timer.stop("|B|")
if verbose > 1:
timer.disp("|B|")
# B^theta
timer.start("B^theta")
if verbose > 0:
print("Saving B^theta")
bsupumnc = file.createVariable(
"bsupumnc", np.float64, ("radius", "mn_mode_nyq")
)
bsupumnc.long_name = "cos(m*t-n*p) component of B^theta, on half mesh"
bsupumnc.units = "T/m"
m = cos_basis.modes[:, 1]
n = cos_basis.modes[:, 2]
if not eq.sym:
bsupumns = file.createVariable(
"bsupumns", np.float64, ("radius", "mn_mode_nyq")
)
bsupumns.long_name = "sin(m*t-n*p) component of B^theta, on half mesh"
bsupumns.units = "T/m"
m = full_basis.modes[:, 1]
n = full_basis.modes[:, 2]
data = (
data_half["B^theta"]
.reshape(
(grid_half.num_theta, grid_half.num_rho, grid_half.num_zeta), order="F"
)
.transpose((1, 0, 2))
.reshape((grid_half.num_rho, -1), order="F")
)
x_mn = np.zeros((surfs - 1, m.size))
for i in range(surfs - 1):
if eq.sym:
x_mn[i, :] = cosfit(data[i, :])
else:
x_mn[i, :] = full_transform.fit(data[i, :])
xm, xn, s, c = ptolemy_identity_rev(m, n, x_mn)
bsupumnc[0, :] = 0
bsupumnc[1:, :] = -c # negative sign for negative Jacobian
if not eq.sym:
bsupumns[0, :] = 0
bsupumns[1:, :] = -s
timer.stop("B^theta")
if verbose > 1:
timer.disp("B^theta")
# B^zeta
timer.start("B^zeta")
if verbose > 0:
print("Saving B^zeta")
bsupvmnc = file.createVariable(
"bsupvmnc", np.float64, ("radius", "mn_mode_nyq")
)
bsupvmnc.long_name = "cos(m*t-n*p) component of B^zeta, on half mesh"
bsupvmnc.units = "T/m"
m = cos_basis.modes[:, 1]
n = cos_basis.modes[:, 2]
if not eq.sym:
bsupvmns = file.createVariable(
"bsupvmns", np.float64, ("radius", "mn_mode_nyq")
)
bsupvmns.long_name = "sin(m*t-n*p) component of B^zeta, on half mesh"
bsupvmns.units = "T/m"
m = full_basis.modes[:, 1]
n = full_basis.modes[:, 2]
data = (
data_half["B^zeta"]
.reshape(
(grid_half.num_theta, grid_half.num_rho, grid_half.num_zeta), order="F"
)
.transpose((1, 0, 2))