VMEC 8.52
3D Equilibrium solver with nested flux surfaces.
Loading...
Searching...
No Matches
read_wout_mod.f
Go to the documentation of this file.
1!> \file
2!> \brief Reading of \c wout VMEC output file.
3
4!> \brief Reading of \c wout VMEC output file.
6!
7! USE READ_WOUT_MOD to include variables dynamically allocated
8! in this module
9! Call DEALLOCATE_READ_WOUT to free this memory when it is no longer needed
10!
11! Reads in output from VMEC equilibrium code(s), contained in wout file
12!
13! Contained subroutines:
14!
15! read_wout_file wrapper alias called to read/open wout file
16! read_wout_text called by read_wout_file to read text file wout
17! read_wout_nc called by read_wout_file to read netcdf file wout
18!
19! Post-processing routines
20!
21! mse_pitch user-callable function to compute mse pitch angle
22! for the computed equilibrium
23!
24 USE stel_kinds
25 USE mgrid_mod
26
27 IMPLICIT NONE
28C-----------------------------------------------
29C L O C A L P A R A M E T E R S
30C-----------------------------------------------
31! Variable names (vn_...) : put eventually into library, used by read_wout too...
32 CHARACTER(LEN=*), PARAMETER ::
33 1 vn_version = 'version_',
34 2 vn_extension = 'input_extension', vn_mgrid = 'mgrid_file',
35 3 vn_magen = 'wb', vn_therm = 'wp', vn_gam = 'gamma',
36 4 vn_maxr = 'rmax_surf', vn_minr = 'rmin_surf',
37 5 vn_maxz = 'zmax_surf', vn_fp = 'nfp',
38 6 vn_radnod = 'ns', vn_polmod = 'mpol', vn_tormod = 'ntor',
39 7 vn_maxmod = 'mnmax', vn_maxit = 'niter',
40 8 vn_asym = 'lasym', vn_free = 'lfreeb',
41 9 vn_error = 'ier_flag', vn_aspect = 'aspect',
42 a vn_maxmod_nyq = 'mnmax_nyq',
43 b vn_beta = 'betatotal', vn_pbeta = 'betapol',
44 c vn_tbeta = 'betator', vn_abeta = 'betaxis',
45 d vn_b0 = 'b0', vn_rbt0 = 'rbtor0', vn_rbt1 = 'rbtor',
46 e vn_sgs = 'signgs', vn_lar = 'IonLarmor', vn_modb = 'volavgB',
47 f vn_ctor = 'ctor', vn_amin = 'Aminor_p', vn_rmaj = 'Rmajor_p',
48 g vn_vol = 'volume_p', vn_am = 'am', vn_ai = 'ai', vn_ac = 'ac',
49 g vn_ah = 'hot particle fraction', vn_atuname = 'T-perp/T-par',
50 h vn_pmass_type = 'pmass_type', vn_piota_type = 'piota_type',
51 i vn_pcurr_type = 'pcurr_type',
52 j vn_am_aux_s = 'am_aux_s', vn_am_aux_f = 'am_aux_f',
53 k vn_ai_aux_s = 'ai_aux_s', vn_ai_aux_f = 'ai_aux_f',
54 l vn_ac_aux_s = 'ac_aux_s', vn_ac_aux_f = 'ac_aux_f',
55 m vn_mse = 'imse', vn_thom = 'itse',
56 n vn_pmod = 'xm', vn_tmod = 'xn', vn_pmod_nyq = 'xm_nyq',
57 o vn_tmod_nyq = 'xn_nyq',
58 p vn_racc = 'raxis_cc', vn_zacs = 'zaxis_cs',
59 q vn_racs = 'raxis_cs', vn_zacc = 'zaxis_cc', vn_iotaf = 'iotaf',
60 q vn_qfact='q-factor', vn_chi='chi', vn_chipf='chipf',
61 r vn_presf = 'presf', vn_phi = 'phi', vn_phipf = 'phipf',
62 s vn_jcuru = 'jcuru', vn_jcurv = 'jcurv', vn_iotah = 'iotas',
63 t vn_mass = 'mass', vn_presh = 'pres', vn_betah = 'beta_vol',
64 u vn_buco = 'buco', vn_bvco = 'bvco', vn_vp = 'vp',
65 v vn_specw = 'specw', vn_phip = 'phips', vn_jdotb = 'jdotb',
66 w vn_bdotb = 'bdotb', vn_overr = 'over_r',
67 x vn_bgrv = 'bdotgradv', vn_merc = 'DMerc', vn_mshear = 'DShear',
68 y vn_mwell = 'DWell', vn_mcurr = 'DCurr', vn_mgeo = 'DGeod',
69 z vn_equif = 'equif', vn_fsq = 'fsqt', vn_wdot = 'wdot',
70 1 vn_ftolv = 'ftolv', vn_fsql= 'fsql', vn_fsqr = 'fsqr',
71 2 vn_fsqz = 'fsqz',
72 3 vn_extcur = 'extcur', vn_curlab = 'curlabel', vn_rmnc = 'rmnc',
73 4 vn_zmns = 'zmns', vn_lmns = 'lmns', vn_gmnc = 'gmnc',
74 5 vn_bmnc = 'bmnc', vn_bsubumnc = 'bsubumnc',
75 6 vn_bsubvmnc = 'bsubvmnc', vn_bsubsmns = 'bsubsmns',
76 7 vn_bsupumnc = 'bsupumnc', vn_bsupvmnc = 'bsupvmnc',
77 8 vn_rmns = 'rmns', vn_zmnc = 'zmnc',
78 9 vn_lmnc = 'lmnc', vn_gmns = 'gmns', vn_bmns = 'bmns',
79 a vn_bsubumns = 'bsubumns', vn_bsubvmns = 'bsubvmns',
80 b vn_bsubsmnc = 'bsubsmnc', vn_bsupumns = 'bsupumns',
81 c vn_bsupvmns = 'bsupvmns',
82 d vn_bsubumnc_sur = 'bsubumnc_sur',
83 e vn_bsubvmnc_sur = 'bsubvmnc_sur',
84 f vn_bsupumnc_sur = 'bsupumnc_sur',
85 g vn_bsupvmnc_sur = 'bsupvmnc_sur',
86 h vn_bsubumns_sur = 'bsubumns_sur',
87 i vn_bsubvmns_sur = 'bsubvmns_sur',
88 j vn_bsupumns_sur = 'bsupumns_sur',
89 k vn_bsupvmns_sur = 'bsupvmns_sur',
90 d vn_rbc = 'rbc', vn_zbs = 'zbs', vn_rbs = 'rbs', vn_zbc = 'zbc',
91 e vn_potvac = 'potvac'
92
93! Long names (ln_...)
94 CHARACTER(LEN=*), PARAMETER ::
95 1 ln_version = 'VMEC Version',
96 2 ln_extension = 'Input file extension',
97 3 ln_mgrid = 'MGRID file',
98 4 ln_magen = 'Magnetic Energy', ln_therm = 'Thermal Energy',
99 5 ln_gam = 'Gamma', ln_maxr = 'Maximum R', ln_minr = 'Minimum R',
100 6 ln_maxz = 'Maximum Z', ln_fp = 'Field Periods',
101 7 ln_radnod = 'Radial nodes', ln_polmod = 'Poloidal modes',
102 8 ln_tormod = 'Toroidal modes', ln_maxmod = 'Fourier modes',
103 8 ln_maxmod_nyq = 'Fourier modes (Nyquist)',
104 9 ln_maxit = 'Max iterations',
105 1 ln_asym = 'Asymmetry', ln_recon = 'Reconstruction',
106 1 ln_free = 'Free boundary',
107 2 ln_error = 'Error flag', ln_aspect = 'Aspect ratio',
108 3 ln_beta = 'Total beta', ln_pbeta = 'Poloidal beta',
109 4 ln_tbeta = 'Toroidal beta', ln_abeta = 'Beta axis',
110 5 ln_b0 = 'RB-t over R axis', ln_rbt0 = 'RB-t axis',
111 6 ln_rbt1 = 'RB-t edge', ln_sgs = 'Sign jacobian',
112 7 ln_lar = 'Ion Larmor radius', ln_modb = 'avg mod B',
113 8 ln_ctor = 'Toroidal current', ln_amin = 'minor radius',
114 9 ln_rmaj = 'major radius', ln_vol = 'Plasma volume',
115 1 ln_mse = 'Number of MSE points',
116 1 ln_thom = 'Number of Thompson scattering points',
117 1 ln_am = 'Specification parameters for mass(s)',
118 1 ln_ac = 'Specification parameters for <J>(s)',
119 1 ln_ai = 'Specification parameters for iota(s)',
120 1 ln_pmass_type = 'Profile type specifier for mass(s)',
121 1 ln_pcurr_type = 'Profile type specifier for <J>(s)',
122 1 ln_piota_type = 'Profile type specifier for iota(s)',
123 1 ln_am_aux_s = 'Auxiliary-s parameters for mass(s)',
124 1 ln_am_aux_f = 'Auxiliary-f parameters for mass(s)',
125 1 ln_ac_aux_s = 'Auxiliary-s parameters for <J>(s)',
126 1 ln_ac_aux_f = 'Auxiliary-f parameters for <J>(s)',
127 1 ln_ai_aux_s = 'Auxiliary-s parameters for iota(s)',
128 1 ln_ai_aux_f = 'Auxiliary-f parameters for iota(s)',
129 4 ln_pmod = 'Poloidal mode numbers',
130 5 ln_tmod = 'Toroidal mode numbers',
131 4 ln_pmod_nyq = 'Poloidal mode numbers (Nyquist)',
132 5 ln_tmod_nyq = 'Toroidal mode numbers (Nyquist)',
133 5 ln_racc = 'raxis (cosnv)', ln_racs = 'raxis (sinnv)',
134 6 ln_zacs = 'zaxis (sinnv)', ln_zacc = 'zaxis (cosnv)',
135 7 ln_iotaf = 'iota on full mesh',
136 7 ln_qfact = 'q-factor on full mesh',
137 8 ln_presf = 'pressure on full mesh',
138 8 ln_phi = 'Toroidal flux on full mesh',
139 9 ln_phipf = 'd(phi)/ds: Toroidal flux deriv on full mesh',
140 9 ln_chi = 'Poloidal flux on full mesh',
141 9 ln_chipf = 'd(chi)/ds: Poroidal flux deriv on full mesh',
142 9 ln_jcuru = 'j dot gradu full',
143 1 ln_jcurv = 'j dot gradv full', ln_iotah = 'iota half',
144 2 ln_mass = 'mass half', ln_presh = 'pressure half',
145 3 ln_betah = 'beta half', ln_buco = 'bsubu half',
146 4 ln_bvco = 'bsubv half', ln_vp = 'volume deriv half',
147 5 ln_specw = 'Spectral width half',
148 6 ln_phip = 'tor flux deriv over 2pi half',
149 7 ln_jdotb = 'J dot B', ln_bdotb = 'B dot B',
150 7 ln_bgrv = 'B dot grad v',
151 8 ln_merc = 'Mercier criterion', ln_mshear = 'Shear Mercier',
152 9 ln_mwell = 'Well Mercier', ln_mcurr = 'Current Mercier',
153 1 ln_mgeo = 'Geodesic Mercier', ln_equif='Average force balance',
154 1 ln_fsq = 'Residual decay',
155 2 ln_wdot = 'Wdot decay', ln_extcur = 'External coil currents',
156 2 ln_fsqr = 'Residual decay - radial',
157 2 ln_fsqz = 'Residual decay - vertical',
158 2 ln_fsql = 'Residual decay - hoop',
159 2 ln_ftolv = 'Residual decay - requested',
160 3 ln_curlab = 'External current names',
161
162 3 ln_rmnc = 'cosmn component of cylindrical R, full mesh',
163 4 ln_zmns = 'sinmn component of cylindrical Z, full mesh',
164 4 ln_lmns = 'sinmn component of lambda, half mesh',
165 5 ln_gmnc = 'cosmn component of jacobian, half mesh',
166 6 ln_bmnc = 'cosmn component of mod-B, half mesh',
167 6 ln_bsubumnc = 'cosmn covariant u-component of B, half mesh',
168 6 ln_bsubvmnc = 'cosmn covariant v-component of B, half mesh',
169 7 ln_bsubsmns = 'sinmn covariant s-component of B, full mesh',
170
171 8 ln_bsubumnc_sur = 'cosmn bsubu of B, surface',
172 9 ln_bsubvmnc_sur = 'cosmn bsubv of B, surface',
173 a ln_bsupumnc_sur = 'cosmn bsupu of B, surface',
174 b ln_bsupvmnc_sur = 'cosmn bsupv of B, surface',
175
176 7 ln_bsupumnc = 'BSUPUmnc half',
177 8 ln_bsupvmnc = 'BSUPVmnc half',
178
179 3 ln_rmns = 'sinmn component of cylindrical R, full mesh',
180 4 ln_zmnc = 'cosmn component of cylindrical Z, full mesh',
181 4 ln_lmnc = 'cosmn component of lambda, half mesh',
182 5 ln_gmns = 'sinmn component of jacobian, half mesh',
183 6 ln_bmns = 'sinmn component of mod-B, half mesh',
184 6 ln_bsubumns = 'sinmn covariant u-component of B, half mesh',
185 6 ln_bsubvmns = 'sinmn covariant v-component of B, half mesh',
186 7 ln_bsubsmnc = 'cosmn covariant s-component of B, full mesh',
187
188 8 ln_bsubumns_sur = 'sinmn bsubu of B, surface',
189 9 ln_bsubvmns_sur = 'sinmn bsubv of B, surface',
190 a ln_bsupumns_sur = 'sinmn bsupu of B, surface',
191 b ln_bsupvmns_sur = 'sinmn bsupv of B, surface',
192
193 4 ln_bsupumns = 'BSUPUmns half', ln_bsupvmns = 'BSUPVmns half',
194 6 ln_rbc = 'Initial boundary R cos(mu-nv) coefficients',
195 7 ln_zbs = 'Initial boundary Z sin(mu-nv) coefficients',
196 8 ln_rbs = 'Initial boundary R sin(mu-nv) coefficients',
197 9 ln_zbc = 'Initial boundary Z cos(mu-nv) coefficients',
198 1 ln_potvac = 'Vacuum Potential on Boundary'
199!-----------------------------------------------
200! L o c a l V a r i a b l e s
201!-----------------------------------------------
202 INTEGER :: nfp, ns, mpol, ntor, mnmax, mnmax_nyq, niter,
209 3 fsql, fsqr, fsqz, ftolv,
211 5 machsq !SAL
212 REAL(rprec), ALLOCATABLE :: rzl_local(:,:,:,:)
213 REAL(rprec), DIMENSION(:,:), ALLOCATABLE ::
217 REAL(rprec), DIMENSION(:,:), ALLOCATABLE ::
220 REAL(rprec), DIMENSION(:), ALLOCATABLE ::
222 1 qfact, chipf, phi, chi,
230
231 LOGICAL :: lasym, lthreed, lwout_opened=.false.
232 CHARACTER :: mgrid_file*200, input_extension*100
233 CHARACTER :: pmass_type*20, pcurr_type*20, piota_type*20
234
235 INTEGER, PARAMETER :: norm_term_flag=0,
236 1 bad_jacobian_flag=1, more_iter_flag=2, jac75_flag=4
237
238! OVERLOAD SUBROUTINE READ_WOUT_FILE TO ACCEPT BOTH UNIT NO. (OPENED EXTERNALLY)
239! OR FILENAME (HANDLE OPEN/CLOSE HERE)
241 MODULE PROCEDURE readw_and_open
242 END INTERFACE
243
244 PRIVATE :: read_wout_nc
245 PRIVATE :: norm_term_flag, bad_jacobian_flag,
246 1 more_iter_flag, jac75_flag
247
248 CONTAINS
249
250 SUBROUTINE readw_and_open(file_or_extension, ierr, iopen)
251 USE safe_open_mod
252 IMPLICIT NONE
253C-----------------------------------------------
254C D u m m y A r g u m e n t s
255C-----------------------------------------------
256 INTEGER, INTENT(out) :: ierr
257 INTEGER, OPTIONAL :: iopen
258 CHARACTER(LEN=*), INTENT(in) :: file_or_extension
259C-----------------------------------------------
260C L o c a l V a r i a b l e s
261C-----------------------------------------------
262 INTEGER, PARAMETER :: iunit_init = 10
263 LOGICAL :: isnc
264 CHARACTER(len=LEN_TRIM(file_or_extension)+10) :: filename
265C-----------------------------------------------
266!
267! THIS SUBROUTINE READS THE WOUT FILE CREATED BY THE VMEC CODE
268! AND STORES THE DATA IN THE READ_WOUT MODULE
269!
270! FIRST, CHECK IF THIS IS A FULLY-QUALIFIED PATH NAME
271! MAKE SURE wout IS NOT EMBEDDED IN THE NAME (PERVERSE USER...)
272!
273 filename = 'wout'
274 CALL parse_extension(filename, file_or_extension, isnc)
275 CALL flush(6)
276 IF (isnc) THEN
277 CALL read_wout_nc(filename, ierr)
278 END IF
279
280 IF (PRESENT(iopen)) iopen = ierr
281 lwout_opened = (ierr .eq. 0)
282 ! WHEN READING A NETCDF FILE, A BAD RUN MAY PREVENT XN FROM BEING
283 ! READ, SUBSEQUENTLY WE MUST CHECK TO SEE IF XN HAS BEEN ALLOCATED
284 ! BEFORE DOING ANYTHING WITH IT OTHERWISE WE DEFAULT LTHREED TO
285 ! FALSE. - SAL 09/07/11
286 IF (ALLOCATED(xn)) THEN
287 lthreed = any(nint(xn) .ne. 0)
288 ELSE
289 lthreed = .false.
290 END IF
291
292 END SUBROUTINE readw_and_open
293
294 SUBROUTINE read_wout_nc(filename, ierr)
295 USE ezcdf
296 USE stel_constants, ONLY: mu0
297 IMPLICIT NONE
298C-----------------------------------------------
299C D u m m y A r g u m e n t s
300C-----------------------------------------------
301 INTEGER, INTENT(out) :: ierr
302 CHARACTER(LEN=*), INTENT(in) :: filename
303C-----------------------------------------------
304C L o c a l V a r i a b l e s
305C-----------------------------------------------
306 INTEGER :: nwout, ierror
307 INTEGER, DIMENSION(3) :: dimlens
308 REAL(rprec) :: ohs
309 REAL(rprec), DIMENSION(:), ALLOCATABLE :: raxis_cc, raxis_cs,
310 1 zaxis_cs, zaxis_cc
311C-----------------------------------------------
312! Open cdf File
313 CALL cdf_open(nwout,filename,'r', ierr)
314 IF (ierr .ne. 0) THEN
315 print *,' Error opening wout .nc file'
316 RETURN
317 END IF
318
319! Be sure all arrays are deallocated
321
322 ! reset iasym, needed when using the same instance multiple times -jons
323 iasym = 0
324
325! Read in scalar variables
326 CALL cdf_read(nwout, vn_error, ierr_vmec)
327
328 IF (ierr_vmec.ne.norm_term_flag .and. ierr_vmec.ne.more_iter_flag)
329 1 GOTO 1000
330
331 CALL cdf_read(nwout, vn_version, version_)
332 CALL cdf_read(nwout, vn_extension, input_extension)
333 CALL cdf_read(nwout, vn_mgrid, mgrid_file)
334 CALL cdf_read(nwout, vn_magen, wb)
335 CALL cdf_read(nwout, vn_therm, wp)
336 CALL cdf_read(nwout, vn_gam, gamma)
337 CALL cdf_read(nwout, vn_maxr, rmax_surf)
338 CALL cdf_read(nwout, vn_minr, rmin_surf)
339 CALL cdf_read(nwout, vn_maxz, zmax_surf)
340 CALL cdf_read(nwout, vn_fp, nfp)
341 CALL cdf_read(nwout, vn_radnod, ns)
342 CALL cdf_read(nwout, vn_polmod, mpol)
343 CALL cdf_read(nwout, vn_tormod, ntor)
344 CALL cdf_read(nwout, vn_maxmod, mnmax)
345 mnmax_nyq = -1
346 CALL cdf_read(nwout, vn_maxmod_nyq, mnmax_nyq)
347 CALL cdf_read(nwout, vn_maxit, niter)
348 CALL cdf_read(nwout, vn_asym, lasym)
349 IF (lasym) iasym = 1
350 CALL cdf_read(nwout, vn_free, lfreeb)
351 CALL cdf_read(nwout, vn_aspect, aspect)
352 CALL cdf_read(nwout, vn_beta, betatot)
353 CALL cdf_read(nwout, vn_pbeta, betapol)
354 CALL cdf_read(nwout, vn_tbeta, betator)
355 CALL cdf_read(nwout, vn_abeta, betaxis)
356 CALL cdf_read(nwout, vn_b0, b0)
357 CALL cdf_read(nwout, vn_rbt0, rbtor0)
358 CALL cdf_read(nwout, vn_rbt1, rbtor)
359 CALL cdf_read(nwout, vn_sgs, isigng)
360 CALL cdf_read(nwout, vn_lar, ionlarmor)
361 CALL cdf_read(nwout, vn_modb, volavgb)
362 CALL cdf_read(nwout, vn_ctor, itor)
363 CALL cdf_read(nwout, vn_amin, aminor)
364 CALL cdf_read(nwout, vn_rmaj, rmajor)
365 CALL cdf_read(nwout, vn_vol, volume)
366 CALL cdf_read(nwout, vn_ftolv, ftolv)
367 CALL cdf_read(nwout, vn_fsqr, fsqr)
368 CALL cdf_read(nwout, vn_fsqz, fsqz)
369 CALL cdf_read(nwout, vn_fsql, fsql)
370 CALL cdf_read(nwout, vn_pcurr_type, pcurr_type)
371 CALL cdf_read(nwout, vn_piota_type, piota_type)
372 CALL cdf_read(nwout, vn_pmass_type, pmass_type)
373 imse = -1
374 CALL cdf_read(nwout, vn_nextcur, nextcur)
375
376 mgrid_mode = 'N'
377 CALL cdf_inquire(nwout, vn_mgmode, dimlens, ier=ierror)
378 IF (ierror.eq.0) CALL cdf_read(nwout, vn_mgmode, mgrid_mode)
379
380! Inquire existence, dimensions of arrays for allocation
381! 1D Arrays
382 CALL cdf_inquire(nwout, vn_pmod, dimlens)
383 ALLOCATE (xm(dimlens(1)), stat = ierror)
384 CALL cdf_inquire(nwout, vn_tmod, dimlens)
385 ALLOCATE (xn(dimlens(1)), stat = ierror)
386 IF (mnmax_nyq .gt. 0) THEN
387 CALL cdf_inquire(nwout, vn_pmod_nyq, dimlens)
388 ALLOCATE (xm_nyq(dimlens(1)), stat = ierror)
389 CALL cdf_inquire(nwout, vn_tmod_nyq, dimlens)
390 ALLOCATE (xn_nyq(dimlens(1)), stat = ierror)
391 END IF
392
393 CALL cdf_inquire(nwout, vn_racc, dimlens)
394 ALLOCATE (raxis_cc(0:dimlens(1)-1), stat = ierror)
395 CALL cdf_inquire(nwout, vn_zacs, dimlens)
396 ALLOCATE (zaxis_cs(0:dimlens(1)-1), stat = ierror)
397 IF (lasym) THEN
398 CALL cdf_inquire(nwout, vn_racs, dimlens)
399 ALLOCATE (raxis_cs(0:dimlens(1)-1), stat = ierror)
400 CALL cdf_inquire(nwout, vn_zacc, dimlens)
401 ALLOCATE (zaxis_cc(0:dimlens(1)-1), stat = ierror)
402 END IF
403
404! Profile coefficients, dimensioned from 0
405 CALL cdf_inquire(nwout, vn_am, dimlens)
406 ALLOCATE (am(0:dimlens(1)-1), stat = ierror)
407 CALL cdf_inquire(nwout, vn_ac, dimlens)
408 ALLOCATE (ac(0:dimlens(1)-1), stat = ierror)
409 CALL cdf_inquire(nwout, vn_ai, dimlens)
410 ALLOCATE (ai(0:dimlens(1)-1), stat = ierror)
411
412 CALL cdf_inquire(nwout, vn_ac_aux_s, dimlens)
413 ALLOCATE (ac_aux_s(dimlens(1)), stat = ierror)
414 CALL cdf_inquire(nwout, vn_ac_aux_f, dimlens)
415 ALLOCATE (ac_aux_f(dimlens(1)), stat = ierror)
416 CALL cdf_inquire(nwout, vn_ai_aux_s, dimlens)
417 ALLOCATE (ai_aux_s(dimlens(1)), stat = ierror)
418 CALL cdf_inquire(nwout, vn_ai_aux_f, dimlens)
419 ALLOCATE (ai_aux_f(dimlens(1)), stat = ierror)
420 CALL cdf_inquire(nwout, vn_am_aux_s, dimlens)
421 ALLOCATE (am_aux_s(dimlens(1)), stat = ierror)
422 CALL cdf_inquire(nwout, vn_am_aux_f, dimlens)
423 ALLOCATE (am_aux_f(dimlens(1)), stat = ierror)
424
425 CALL cdf_inquire(nwout, vn_iotaf, dimlens)
426 ALLOCATE (iotaf(dimlens(1)), stat = ierror)
427 CALL cdf_inquire(nwout, vn_qfact, dimlens)
428 ALLOCATE (qfact(dimlens(1)), stat = ierror)
429 CALL cdf_inquire(nwout, vn_presf, dimlens)
430 ALLOCATE (presf(dimlens(1)), stat = ierror)
431 CALL cdf_inquire(nwout, vn_phi, dimlens)
432 ALLOCATE (phi(dimlens(1)), stat = ierror)
433 CALL cdf_inquire(nwout, vn_chi, dimlens)
434 ALLOCATE (chi(dimlens(1)), stat = ierror)
435 CALL cdf_inquire(nwout, vn_phipf, dimlens)
436 ALLOCATE (phipf(dimlens(1)), stat = ierror)
437 CALL cdf_inquire(nwout, vn_chipf, dimlens)
438 ALLOCATE (chipf(dimlens(1)), stat = ierror)
439 CALL cdf_inquire(nwout, vn_jcuru, dimlens)
440 ALLOCATE (jcuru(dimlens(1)), stat = ierror)
441 CALL cdf_inquire(nwout, vn_jcurv, dimlens)
442 ALLOCATE (jcurv(dimlens(1)), stat = ierror)
443 CALL cdf_inquire(nwout, vn_iotah, dimlens)
444 ALLOCATE (iotas(dimlens(1)), stat = ierror)
445 CALL cdf_inquire(nwout, vn_mass, dimlens)
446 ALLOCATE (mass(dimlens(1)), stat = ierror)
447 CALL cdf_inquire(nwout, vn_presh, dimlens)
448 ALLOCATE (pres(dimlens(1)), stat = ierror)
449 CALL cdf_inquire(nwout, vn_betah, dimlens)
450 ALLOCATE (beta_vol(dimlens(1)), stat = ierror)
451 CALL cdf_inquire(nwout, vn_buco, dimlens)
452 ALLOCATE (buco(dimlens(1)), stat = ierror)
453 CALL cdf_inquire(nwout, vn_bvco, dimlens)
454 ALLOCATE (bvco(dimlens(1)), stat = ierror)
455 CALL cdf_inquire(nwout, vn_vp, dimlens)
456 ALLOCATE (vp(dimlens(1)), stat = ierror)
457 CALL cdf_inquire(nwout, vn_specw, dimlens)
458 ALLOCATE (specw(dimlens(1)), stat = ierror)
459 CALL cdf_inquire(nwout, vn_phip, dimlens)
460 ALLOCATE (phip(dimlens(1)), stat = ierror)
461 CALL cdf_inquire(nwout, vn_overr, dimlens)
462 ALLOCATE (overr(dimlens(1)), stat = ierror)
463
464 CALL cdf_inquire(nwout, vn_jdotb, dimlens)
465 ALLOCATE (jdotb(dimlens(1)), stat = ierror)
466 CALL cdf_inquire(nwout, vn_bdotb, dimlens)
467 ALLOCATE (bdotb(dimlens(1)), stat = ierror)
468 CALL cdf_inquire(nwout, vn_bgrv, dimlens)
469 ALLOCATE (bdotgradv(dimlens(1)), stat = ierror)
470
471 CALL cdf_inquire(nwout, vn_merc, dimlens)
472 ALLOCATE (dmerc(dimlens(1)), stat = ierror)
473 CALL cdf_inquire(nwout, vn_mshear, dimlens)
474 ALLOCATE (dshear(dimlens(1)), stat = ierror)
475 CALL cdf_inquire(nwout, vn_mwell, dimlens)
476 ALLOCATE (dwell(dimlens(1)), stat = ierror)
477 CALL cdf_inquire(nwout, vn_mcurr, dimlens)
478 ALLOCATE (dcurr(dimlens(1)), stat = ierror)
479 CALL cdf_inquire(nwout, vn_mgeo, dimlens)
480 ALLOCATE (dgeod(dimlens(1)), stat = ierror)
481 CALL cdf_inquire(nwout, vn_equif, dimlens)
482 ALLOCATE (equif(dimlens(1)), stat = ierror)
483
484 CALL cdf_inquire(nwout, vn_fsq, dimlens)
485 ALLOCATE (fsqt(dimlens(1)), stat = ierror)
486 CALL cdf_inquire(nwout, vn_wdot, dimlens)
487 ALLOCATE (wdot(dimlens(1)), stat = ierror)
488
489 IF (nextcur .gt. 0) THEN
490 CALL cdf_inquire(nwout, vn_extcur, dimlens)
491 ALLOCATE (extcur(dimlens(1)), stat = ierror)
492!NOTE: curlabel is an array of CHARACTER(30) strings - defined in mgrid_mod
493! so dimlens(1) == 30 (check this) and dimlens(2) is the number of strings in the array
494 CALL cdf_inquire(nwout, vn_curlab, dimlens)
495 ALLOCATE (curlabel(dimlens(2)), stat = ierror)
496 ! SAL
497 CALL cdf_inquire(nwout, vn_potvac, dimlens, ier = ierror)
498 IF (ierror == 0) ALLOCATE (potvac(1:dimlens(1)), stat = ierror)
499 ENDIF
500
501! 2D Arrays
502 CALL cdf_inquire(nwout, vn_rmnc, dimlens)
503 ALLOCATE (rmnc(dimlens(1),dimlens(2)), stat = ierror)
504 CALL cdf_inquire(nwout, vn_zmns, dimlens)
505 ALLOCATE (zmns(dimlens(1),dimlens(2)), stat = ierror)
506 CALL cdf_inquire(nwout, vn_lmns, dimlens)
507 ALLOCATE (lmns(dimlens(1),dimlens(2)), stat = ierror)
508 CALL cdf_inquire(nwout, vn_gmnc, dimlens)
509 ALLOCATE (gmnc(dimlens(1),dimlens(2)), stat = ierror)
510 CALL cdf_inquire(nwout, vn_bmnc, dimlens)
511 ALLOCATE (bmnc(dimlens(1),dimlens(2)), stat = ierror)
512 CALL cdf_inquire(nwout, vn_bsubumnc, dimlens)
513 ALLOCATE (bsubumnc(dimlens(1),dimlens(2)), stat = ierror)
514 CALL cdf_inquire(nwout, vn_bsubvmnc, dimlens)
515 ALLOCATE (bsubvmnc(dimlens(1),dimlens(2)), stat = ierror)
516 CALL cdf_inquire(nwout, vn_bsubsmns, dimlens)
517 ALLOCATE (bsubsmns(dimlens(1),dimlens(2)), stat = ierror)
518
519! ELIMINATE THESE EVENTUALLY: DON'T NEED THEM
520 CALL cdf_inquire(nwout, vn_bsupumnc, dimlens)
521 ALLOCATE (bsupumnc(dimlens(1),dimlens(2)), stat = ierror)
522 CALL cdf_inquire(nwout, vn_bsupvmnc, dimlens)
523 ALLOCATE (bsupvmnc(dimlens(1),dimlens(2)), stat = ierror)
524
525 IF (.NOT. lasym) GO TO 800
526
527 CALL cdf_inquire(nwout, vn_rmns, dimlens)
528 ALLOCATE (rmns(dimlens(1),dimlens(2)), stat = ierror)
529 CALL cdf_inquire(nwout, vn_zmnc, dimlens)
530 ALLOCATE (zmnc(dimlens(1),dimlens(2)), stat = ierror)
531 CALL cdf_inquire(nwout, vn_lmnc, dimlens)
532 ALLOCATE (lmnc(dimlens(1),dimlens(2)), stat = ierror)
533 CALL cdf_inquire(nwout, vn_gmns, dimlens)
534 ALLOCATE (gmns(dimlens(1),dimlens(2)), stat = ierror)
535 CALL cdf_inquire(nwout, vn_bmns, dimlens)
536 ALLOCATE (bmns(dimlens(1),dimlens(2)), stat = ierror)
537 CALL cdf_inquire(nwout, vn_bsubumns, dimlens)
538 ALLOCATE (bsubumns(dimlens(1),dimlens(2)), stat = ierror)
539 CALL cdf_inquire(nwout, vn_bsubvmns, dimlens)
540 ALLOCATE (bsubvmns(dimlens(1),dimlens(2)), stat = ierror)
541 CALL cdf_inquire(nwout, vn_bsubsmnc, dimlens)
542 ALLOCATE (bsubsmnc(dimlens(1),dimlens(2)), stat = ierror)
543
544! ELIMINATE THESE EVENTUALLY: DO NOT NEED THEM
545 CALL cdf_inquire(nwout, vn_bsupumns, dimlens)
546 ALLOCATE (bsupumns(dimlens(1),dimlens(2)), stat = ierror)
547 CALL cdf_inquire(nwout, vn_bsupvmns, dimlens)
548 ALLOCATE (bsupvmns(dimlens(1),dimlens(2)), stat = ierror)
549
550 800 CONTINUE
551
552! Read Arrays
553 CALL cdf_read(nwout, vn_pmod, xm)
554 CALL cdf_read(nwout, vn_tmod, xn)
555 IF (mnmax_nyq .le. 0) THEN
557 ALLOCATE (xm_nyq(mnmax_nyq), xn_nyq(mnmax_nyq), stat=ierror)
558 xm_nyq = xm; xn_nyq = xn
559 ELSE
560 CALL cdf_read(nwout, vn_pmod_nyq, xm_nyq)
561 CALL cdf_read(nwout, vn_tmod_nyq, xn_nyq)
562 END IF
563
564 mnyq = int(maxval(xm_nyq)); nnyq = int(maxval(abs(xn_nyq)))/nfp
565
566 CALL cdf_read(nwout, vn_racc, raxis_cc)
567 CALL cdf_read(nwout, vn_zacs, zaxis_cs)
568
569 IF (SIZE(raxis_cc) .ne. ntor+1)
570 1 stop 'WRONG SIZE(raxis_cc) in READ_WOUT_NC'
571 ALLOCATE (raxis(0:ntor,2), zaxis(0:ntor,2), stat=ierror)
572 raxis(:,1) = raxis_cc(0:ntor); zaxis(:,1) = zaxis_cs(0:ntor)
573 raxis(:,2) = 0; zaxis(:,2) = 0
574 DEALLOCATE (raxis_cc, zaxis_cs, stat=ierror)
575
576 CALL cdf_read(nwout, vn_rmnc, rmnc)
577 CALL cdf_read(nwout, vn_zmns, zmns)
578 CALL cdf_read(nwout, vn_lmns, lmns)
579 CALL cdf_read(nwout, vn_gmnc, gmnc) !Half mesh
580 CALL cdf_read(nwout, vn_bmnc, bmnc) !Half mesh
581 CALL cdf_read(nwout, vn_bsubumnc, bsubumnc) !Half mesh
582 CALL cdf_read(nwout, vn_bsubvmnc, bsubvmnc) !Half mesh
583 CALL cdf_read(nwout, vn_bsubsmns, bsubsmns) !Full mesh
584! ELIMINATE THESE EVENTUALLY: DON'T NEED THEM (can express in terms of lambdas)
585 CALL cdf_read(nwout, vn_bsupumnc, bsupumnc)
586 CALL cdf_read(nwout, vn_bsupvmnc, bsupvmnc)
587 IF (lasym) THEN
588 CALL cdf_read(nwout, vn_racs, raxis_cs)
589 CALL cdf_read(nwout, vn_zacc, zaxis_cc)
590 raxis(:,2) = raxis_cs; zaxis(:,2) = zaxis_cc
591 DEALLOCATE (raxis_cs, zaxis_cc, stat=ierror)
592 CALL cdf_read(nwout, vn_rmns, rmns)
593 CALL cdf_read(nwout, vn_zmnc, zmnc)
594 CALL cdf_read(nwout, vn_lmnc, lmnc)
595 CALL cdf_read(nwout, vn_gmns, gmns)
596 CALL cdf_read(nwout, vn_bmns, bmns)
597 CALL cdf_read(nwout, vn_bsubumns, bsubumns)
598 CALL cdf_read(nwout, vn_bsubvmns, bsubvmns)
599 CALL cdf_read(nwout, vn_bsubsmnc, bsubsmnc)
600! ELIMINATE THESE EVENTUALLY: DON'T NEED THEM
601 CALL cdf_read(nwout, vn_bsupumns, bsupumns)
602 CALL cdf_read(nwout, vn_bsupvmns, bsupvmns)
603 END IF
604
605 CALL cdf_read(nwout, vn_am, am)
606 CALL cdf_read(nwout, vn_ac, ac)
607 CALL cdf_read(nwout, vn_ai, ai)
608
609 CALL cdf_read(nwout, vn_am_aux_s, am_aux_s)
610 CALL cdf_read(nwout, vn_am_aux_f, am_aux_f)
611 CALL cdf_read(nwout, vn_ac_aux_s, ac_aux_s)
612 CALL cdf_read(nwout, vn_ac_aux_f, ac_aux_f)
613 CALL cdf_read(nwout, vn_ai_aux_s, ai_aux_s)
614 CALL cdf_read(nwout, vn_ai_aux_f, ai_aux_f)
615
616 CALL cdf_read(nwout, vn_iotaf, iotaf)
617 CALL cdf_read(nwout, vn_qfact, qfact)
618 CALL cdf_read(nwout, vn_presf, presf)
619 CALL cdf_read(nwout, vn_phi, phi)
620 CALL cdf_read(nwout, vn_phipf, phipf)
621 CALL cdf_read(nwout, vn_chi, chi)
622 CALL cdf_read(nwout, vn_chipf, chipf)
623 CALL cdf_read(nwout, vn_jcuru, jcuru)
624 CALL cdf_read(nwout, vn_jcurv, jcurv)
625
626
627! HALF-MESH quantities
628! NOTE: jdotb is in units_of_A (1/mu0 incorporated in jxbforce...)
629! prior to version 6.00, this was output in internal VMEC units...
630 CALL cdf_read(nwout, vn_iotah, iotas)
631 CALL cdf_read(nwout, vn_mass, mass)
632 CALL cdf_read(nwout, vn_presh, pres)
633 CALL cdf_read(nwout, vn_betah, beta_vol)
634 CALL cdf_read(nwout, vn_buco, buco)
635 CALL cdf_read(nwout, vn_bvco, bvco)
636 CALL cdf_read(nwout, vn_vp, vp)
637 CALL cdf_read(nwout, vn_specw, specw)
638 CALL cdf_read(nwout, vn_phip, phip)
639 CALL cdf_read(nwout, vn_jdotb, jdotb)
640 CALL cdf_read(nwout, vn_bdotb, bdotb)
641 CALL cdf_read(nwout, vn_bgrv, bdotgradv)
642
643! MERCIER_CRITERION
644 CALL cdf_read(nwout, vn_merc, dmerc)
645 CALL cdf_read(nwout, vn_mshear, dshear)
646 CALL cdf_read(nwout, vn_mwell, dwell)
647 CALL cdf_read(nwout, vn_mcurr, dcurr)
648 CALL cdf_read(nwout, vn_mgeo, dgeod)
649 CALL cdf_read(nwout, vn_equif, equif)
650
651 CALL cdf_read(nwout, vn_fsq, fsqt)
652 CALL cdf_read(nwout, vn_wdot, wdot)
653
654 IF (nextcur .gt. 0) THEN
655 CALL cdf_read(nwout, vn_extcur, extcur)
656 CALL cdf_read(nwout, vn_curlab, curlabel)
657 ENDIF
658
659 !SAL Addition
660 IF (ALLOCATED(potvac)) CALL cdf_read(nwout,vn_potvac, potvac)
661
662 1000 CONTINUE
663
664 CALL cdf_close(nwout, ierr)
665
666 IF (.not.ALLOCATED(bsubumnc)) RETURN !Moved this here because ns may not be set. SAL -09/07/11
667!
668! COMPUTE CONTRAVARIANT CURRENT COMPONENTS IN AMPS
669! ON THE FULL RADIAL MESH, WHERE JACOBIAN = SQRT(G)
670!
671! CURRU = SQRT(G) * J dot grad(u)
672! CURRV = SQRT(G) * J dot grad(v)
673!
674 ohs = (ns-1)
675
676
677 IF (ierror .eq. 0) CALL compute_currents(ierror)
678
679 IF (ierr. ne. 0) print *,"in read_wout_nc ierr=",ierr
680 IF (ierror. ne. 0) print *,"in read_wout_nc ierror=",ierror
681
682 END SUBROUTINE read_wout_nc
683
684 SUBROUTINE compute_currents(ierror)
685 USE stel_constants, ONLY: mu0
686 IMPLICIT NONE
687 INTEGER, INTENT(out) :: ierror
688!-----------------------------------------------
689! L o c a l V a r i a b l e s
690!-----------------------------------------------
691 INTEGER :: js
692 REAL(rprec) :: ohs, hs, shalf(ns), sfull(ns)
693 REAL(rprec), DIMENSION(mnmax_nyq) :: bu1, bu0, bv1, bv0, t1, t2,
694 & t3
695!-----------------------------------------------
696!
697! Computes current harmonics for currXmn == sqrt(g)*JsupX, X = u,v
698! [Corrected above "JsubX" to "JsupX", JDH 2010-08-16]
699
700! NOTE: bsub(s,u,v)mn are on HALF radial grid
701! (in earlier versions, bsubsmn was on FULL radial grid)
702
703!
704 ohs = (ns-1)
705 hs = 1._dp/ohs
706
707 DO js = 2, ns
708 shalf(js) = sqrt(hs*(js-1.5_dp))
709 sfull(js) = sqrt(hs*(js-1))
710 END DO
711
712 ALLOCATE (currumnc(mnmax_nyq,ns), currvmnc(mnmax_nyq,ns), &
713 & stat=ierror)
714 IF (ierror .ne. 0) RETURN
715
716 DO js = 2, ns-1
717 WHERE (mod(int(xm_nyq),2) .EQ. 1)
718 t1 = 0.5_dp*(shalf(js+1)*bsubsmns(:,js+1) + &
719 & shalf(js) *bsubsmns(:,js)) /sfull(js)
720 bu0 = bsubumnc(:,js )/shalf(js)
721 bu1 = bsubumnc(:,js+1)/shalf(js+1)
722 t2 = ohs*(bu1-bu0)*sfull(js)+0.25_dp*(bu0+bu1)/sfull(js)
723 bv0 = bsubvmnc(:,js )/shalf(js)
724 bv1 = bsubvmnc(:,js+1)/shalf(js+1)
725 t3 = ohs*(bv1-bv0)*sfull(js)+0.25_dp*(bv0+bv1)/sfull(js)
726 ELSEWHERE
727 t1 = 0.5_dp*(bsubsmns(:,js+1)+bsubsmns(:,js))
728 t2 = ohs*(bsubumnc(:,js+1)-bsubumnc(:,js))
729 t3 = ohs*(bsubvmnc(:,js+1)-bsubvmnc(:,js))
730 ENDWHERE
731 currumnc(:,js) = -xn_nyq(:)*t1 - t3
732 currvmnc(:,js) = -xm_nyq(:)*t1 + t2
733 END DO
734
735 WHERE (xm_nyq .LE. 1)
736 currvmnc(:,1) = 2*currvmnc(:,2) - currvmnc(:,3)
737 currumnc(:,1) = 2*currumnc(:,2) - currumnc(:,3)
738 ELSEWHERE
739 currvmnc(:,1) = 0
740 currumnc(:,1) = 0
741 ENDWHERE
742
743 currumnc(:,ns) = 2*currumnc(:,ns-1) - currumnc(:,ns-2)
744 currvmnc(:,ns) = 2*currvmnc(:,ns-1) - currvmnc(:,ns-2)
746
747 IF (.NOT.lasym) RETURN
748
749 ALLOCATE (currumns(mnmax_nyq,ns), currvmns(mnmax_nyq,ns), &
750 & stat=ierror)
751
752 DO js = 2, ns-1
753 WHERE (mod(int(xm_nyq),2) .EQ. 1)
754 t1 = 0.5_dp*(shalf(js+1)*bsubsmnc(:,js+1) &
755 & + shalf(js) *bsubsmnc(:,js)) / sfull(js)
756 bu0 = bsubumns(:,js )/shalf(js+1)
757 bu1 = bsubumns(:,js+1)/shalf(js+1)
758 t2 = ohs*(bu1-bu0)*sfull(js) + 0.25_dp*(bu0+bu1)/sfull(js)
759 bv0 = bsubvmns(:,js )/shalf(js)
760 bv1 = bsubvmns(:,js+1)/shalf(js+1)
761 t3 = ohs*(bv1-bv0)*sfull(js)+0.25_dp*(bv0+bv1)/sfull(js)
762 ELSEWHERE
763 t1 = 0.5_dp*(bsubsmnc(:,js+1) + bsubsmnc(:,js))
764 t2 = ohs*(bsubumns(:,js+1)-bsubumns(:,js))
765 t3 = ohs*(bsubvmns(:,js+1)-bsubvmns(:,js))
766 END WHERE
767 currumns(:,js) = xn_nyq(:)*t1 - t3
768 currvmns(:,js) = xm_nyq(:)*t1 + t2
769 END DO
770
771 WHERE (xm_nyq .LE. 1)
772 currvmns(:,1) = 2*currvmns(:,2) - currvmns(:,3)
773 currumns(:,1) = 2*currumns(:,2) - currumns(:,3)
774 ELSEWHERE
775 currvmns(:,1) = 0
776 currumns(:,1) = 0
777 END WHERE
778 currumns(:,ns) = 2*currumns(:,ns-1) - currumns(:,ns-2)
779 currvmns(:,ns) = 2*currvmns(:,ns-1) - currvmns(:,ns-2)
781
782 END SUBROUTINE compute_currents
783
785 IMPLICIT NONE
786!-----------------------------------------------
787! L o c a l V a r i a b l e s
788!-----------------------------------------------
789 INTEGER :: istat(10)
790!-----------------------------------------------
791 istat = 0
792 lwout_opened = .false.
793
794 IF (ALLOCATED(extcur)) DEALLOCATE (extcur,
795 1 stat = istat(1))
796 IF (ALLOCATED(curlabel)) DEALLOCATE (curlabel,
797 1 stat = istat(1))
798 IF (ALLOCATED(overr)) DEALLOCATE (overr, stat = istat(2))
799
800 IF (ALLOCATED(xm)) DEALLOCATE (xm, xn, xm_nyq, xn_nyq,
803 3 pres, beta_vol, phip, buco, bvco, phi, vp, jcuru, am, ac, ai,
805 5 bdotb, bdotgradv, raxis, zaxis, fsqt, wdot, stat = istat(3))
806
807 IF (ALLOCATED(chipf)) DEALLOCATE (chipf, chi)
808
809 IF (ALLOCATED(am_aux_s)) DEALLOCATE (am_aux_s, am_aux_f, ac_aux_s,
810 1 ac_aux_f, ai_aux_s, ai_aux_f, stat=istat(6))
811
812 IF (ALLOCATED(rmns)) DEALLOCATE (rmns, zmnc, lmnc,
814 2 bsupumns, bsupvmns, stat=istat(5))
815
816 IF (ALLOCATED(currumnc)) DEALLOCATE (currumnc)
817 IF (ALLOCATED(currumns)) DEALLOCATE (currumns, currvmns)
818 IF (ALLOCATED(rzl_local)) DEALLOCATE (rzl_local)
819
820 ! SAL Addition
821 IF (ALLOCATED(potvac)) DEALLOCATE(potvac)
822
823 IF (any(istat .ne. 0)) THEN
824 print *,istat
825 stop 'Deallocation error in read_wout_deallocate'
826 END IF
827
828 END SUBROUTINE read_wout_deallocate
829
830 SUBROUTINE tosuvspace (s_in, u_in, v_in, gsqrt,
831 1 bsupu, bsupv, jsupu, jsupv, lam)
832 USE stel_constants, ONLY: zero, one
833 IMPLICIT NONE
834C-----------------------------------------------
835C D u m m y A r g u m e n t s
836C-----------------------------------------------
837 REAL(rprec), INTENT(in) :: s_in, u_in, v_in
838 REAL(rprec), INTENT(out), OPTIONAL :: gsqrt, bsupu, bsupv,
839 1 jsupu, jsupv, lam
840C-----------------------------------------------
841C L o c a l V a r i a b l e s
842C-----------------------------------------------
843 REAL(rprec), PARAMETER :: c1p5 = 1.5_dp
844 INTEGER :: m, n, n1, mn, ipresent, jslo, jshi
845 REAL(rprec) :: hs1, wlo, whi, wlo_odd, whi_odd
846 REAL(rprec), DIMENSION(mnmax_nyq) :: gmnc1, gmns1, bsupumnc1,
847 1 bsupumns1, bsupvmnc1, bsupvmns1, jsupumnc1, jsupumns1,
848 2 jsupvmnc1, jsupvmns1, wmins, wplus, lammns1, lammnc1
849 REAL(rprec) :: cosu, sinu, cosv, sinv, tcosmn, tsinmn, sgn
850 REAL(rprec) :: cosmu(0:mnyq), sinmu(0:mnyq),
851 1 cosnv(0:nnyq), sinnv(0:nnyq)
852 LOGICAL :: lgsqrt, lbsupu, lbsupv, ljsupu, ljsupv, llam
853C-----------------------------------------------
854!
855! COMPUTE VARIOUS HALF/FULL-RADIAL GRID QUANTITIES AT THE INPUT POINT
856! (S, U, V) , WHERE
857! S = normalized toroidal flux (0 - 1),
858! U = poloidal angle
859! V = N*phi = toroidal angle * no. field periods
860!
861! HALF-RADIAL GRID QUANTITIES
862! gsqrt, bsupu, bsupv
863!
864! FULL-RADIAL GRID QUANTITIES
865! dbsubuds, dbsubvds, dbsubsdu, dbsubsdv
866!
867C-----------------------------------------------
868 IF (s_in.lt.zero .or. s_in.gt.one) THEN
869 WRITE(6, *)
870 1 ' In tosuvspace, s(flux) must be between 0 and 1'
871 RETURN
872 END IF
873
874 IF (.not.lwout_opened) THEN
875 WRITE(6, *)
876 1 ' tosuvspace can only be called AFTER opening wout file!'
877 RETURN
878 END IF
879
880!
881! SETUP TRIG ARRAYS
882!
883 cosu = cos(u_in); sinu = sin(u_in)
884 cosv = cos(v_in); sinv = sin(v_in)
885
886 cosmu(0) = 1; sinmu(0) = 0
887 cosnv(0) = 1; sinnv(0) = 0
888 DO m = 1, mnyq
889 cosmu(m) = cosmu(m-1)*cosu - sinmu(m-1)*sinu
890 sinmu(m) = sinmu(m-1)*cosu + cosmu(m-1)*sinu
891 END DO
892
893 DO n = 1, nnyq
894 cosnv(n) = cosnv(n-1)*cosv - sinnv(n-1)*sinv
895 sinnv(n) = sinnv(n-1)*cosv + cosnv(n-1)*sinv
896 END DO
897
898
899!
900! FIND INTERPOLATED s VALUE AND COMPUTE INTERPOLATION WEIGHTS wlo, whi
901! RECALL THAT THESE QUANTITIES ARE ON THE HALF-RADIAL GRID...
902! s-half(j) = (j-1.5)*hs, for j = 2,...ns
903!
904 hs1 = one/(ns-1)
905 jslo = int(c1p5 + s_in/hs1)
906 jshi = jslo+1
907 wlo = (hs1*(jshi-c1p5) - s_in)/hs1
908 whi = 1 - wlo
909 IF (jslo .eq. ns) THEN
910! USE Xhalf(ns+1) = 2*Xhalf(ns) - Xhalf(ns-1) FOR "GHOST" POINT VALUE 1/2hs OUTSIDE EDGE
911! THEN, X = wlo*Xhalf(ns) + whi*Xhalf(ns+1) == Xhalf(ns) + whi*(Xhalf(ns) - Xhalf(ns-1))
912 jshi = jslo-1
913 wlo = 1+whi; whi = -whi
914 ELSE IF (jslo .eq. 1) THEN
915 jslo = 2
916 END IF
917
918!
919! FOR ODD-m MODES X ~ SQRT(s), SO INTERPOLATE Xmn/SQRT(s)
920!
921 whi_odd = whi*sqrt(s_in/(hs1*(jshi-c1p5)))
922 IF (jslo .ne. 1) THEN
923 wlo_odd = wlo*sqrt(s_in/(hs1*(jslo-c1p5)))
924 ELSE
925 wlo_odd = 0
926 whi_odd = sqrt(s_in/(hs1*(jshi-c1p5)))
927 END IF
928
929 WHERE (mod(nint(xm_nyq(:)),2) .eq. 0)
930 wmins = wlo
931 wplus = whi
932 ELSEWHERE
933 wmins = wlo_odd
934 wplus = whi_odd
935 END WHERE
936
937 ipresent = 0
938 lgsqrt = PRESENT(gsqrt)
939 IF (lgsqrt) THEN
940 gsqrt = 0 ; ipresent = ipresent+1
941 gmnc1 = wmins*gmnc(:,jslo) + wplus*gmnc(:,jshi)
942 IF (lasym)
943 1 gmns1 = wmins*gmns(:,jslo) + wplus*gmns(:,jshi)
944 END IF
945 lbsupu = PRESENT(bsupu)
946 IF (lbsupu) THEN
947 bsupu = 0 ; ipresent = ipresent+1
948 bsupumnc1 = wmins*bsupumnc(:,jslo) + wplus*bsupumnc(:,jshi)
949 IF (lasym)
950 1 bsupumns1 = wmins*bsupumns(:,jslo) + wplus*bsupumns(:,jshi)
951 END IF
952 lbsupv = PRESENT(bsupv)
953 IF (lbsupv) THEN
954 bsupv = 0 ; ipresent = ipresent+1
955 bsupvmnc1 = wmins*bsupvmnc(:,jslo) + wplus*bsupvmnc(:,jshi)
956 IF (lasym)
957 1 bsupvmns1 = wmins*bsupvmns(:,jslo) + wplus*bsupvmns(:,jshi)
958 END IF
959 llam = PRESENT(lam)
960 IF (llam) THEN
961 lam = 0 ; ipresent = ipresent+1
962 lammns1 = wmins*lmns(:,jslo) + wplus*lmns(:,jshi)
963 IF (lasym)
964 1 lammnc1 = wmins*lmnc(:,jslo) + wplus*lmnc(:,jshi)
965 END IF
966
967 IF (ipresent .eq. 0) GOTO 1000
968
969!
970! COMPUTE GSQRT, ... IN REAL SPACE
971! tcosmn = cos(mu - nv); tsinmn = sin(mu - nv)
972!
973 DO mn = 1, mnmax_nyq
974 m = nint(xm_nyq(mn)); n = nint(xn_nyq(mn))/nfp
975 n1 = abs(n); sgn = sign(1,n)
976 tcosmn = cosmu(m)*cosnv(n1) + sgn*sinmu(m)*sinnv(n1)
977 tsinmn = sinmu(m)*cosnv(n1) - sgn*cosmu(m)*sinnv(n1)
978 IF (lgsqrt) gsqrt = gsqrt + gmnc1(mn)*tcosmn
979 IF (lbsupu) bsupu = bsupu + bsupumnc1(mn)*tcosmn
980 IF (lbsupv) bsupv = bsupv + bsupvmnc1(mn)*tcosmn
981 IF (llam) lam = lam + lammns1(mn)*tsinmn
982 END DO
983
984 IF (.not.lasym) GOTO 1000
985
986 DO mn = 1, mnmax_nyq
987 m = nint(xm_nyq(mn)); n = nint(xn_nyq(mn))/nfp
988 n1 = abs(n); sgn = sign(1,n)
989 tcosmn = cosmu(m)*cosnv(n1) + sgn*sinmu(m)*sinnv(n1)
990 tsinmn = sinmu(m)*cosnv(n1) - sgn*cosmu(m)*sinnv(n1)
991 IF (lgsqrt) gsqrt = gsqrt + gmns1(mn)*tsinmn
992 IF (lbsupu) bsupu = bsupu + bsupumns1(mn)*tsinmn
993 IF (lbsupv) bsupv = bsupv + bsupvmns1(mn)*tsinmn
994 IF (llam) lam = lam + lammnc1(mn)*tcosmn
995 END DO
996
997 1000 CONTINUE
998
999! FULL-MESH QUANTITIES
1000!
1001! FIND INTERPOLATED s VALUE AND COMPUTE INTERPOLATION WEIGHTS wlo, whi
1002! RECALL THAT THESE QUANTITIES ARE ON THE FULL-RADIAL GRID...
1003! s-full(j) = (j-1)*hs, for j = 1,...ns
1004!
1005 hs1 = one/(ns-1)
1006 jslo = 1+int(s_in/hs1)
1007 jshi = jslo+1
1008 IF (jslo .eq. ns) jshi = ns
1009 wlo = (hs1*(jshi-1) - s_in)/hs1
1010 whi = 1 - wlo
1011!
1012! FOR ODD-m MODES X ~ SQRT(s), SO INTERPOLATE Xmn/SQRT(s)
1013!
1014 whi_odd = whi*sqrt(s_in/(hs1*(jshi-1)))
1015 IF (jslo .ne. 1) THEN
1016 wlo_odd = wlo*sqrt(s_in/(hs1*(jslo-1)))
1017 ELSE
1018 wlo_odd = 0
1019 whi_odd = sqrt(s_in/(hs1*(jshi-1)))
1020 END IF
1021
1022 WHERE (mod(nint(xm_nyq(:)),2) .eq. 0)
1023 wmins = wlo
1024 wplus = whi
1025 ELSEWHERE
1026 wmins = wlo_odd
1027 wplus = whi_odd
1028 END WHERE
1029
1030 ipresent = 0
1031 ljsupu = PRESENT(jsupu)
1032 IF (ljsupu) THEN
1033 IF (.not.lgsqrt) stop 'MUST compute gsqrt for jsupu'
1034 jsupu = 0 ; ipresent = ipresent+1
1035 jsupumnc1 = wmins*currumnc(:,jslo) + wplus*currumnc(:,jshi)
1036 IF (lasym)
1037 1 jsupumns1 = wmins*currumns(:,jslo) + wplus*currumns(:,jshi)
1038 END IF
1039
1040 ljsupv = PRESENT(jsupv)
1041 IF (ljsupv) THEN
1042 IF (.not.lgsqrt) stop 'MUST compute gsqrt for jsupv'
1043 jsupv = 0 ; ipresent = ipresent+1
1044 jsupvmnc1 = wmins*currvmnc(:,jslo) + wplus*currvmnc(:,jshi)
1045 IF (lasym)
1046 1 jsupvmns1 = wmins*currvmns(:,jslo) + wplus*currvmns(:,jshi)
1047 END IF
1048
1049 IF (ipresent .eq. 0) RETURN
1050
1051 DO mn = 1, mnmax_nyq
1052 m = nint(xm_nyq(mn)); n = nint(xn_nyq(mn))/nfp
1053 n1 = abs(n); sgn = sign(1,n)
1054 tcosmn = cosmu(m)*cosnv(n1) + sgn*sinmu(m)*sinnv(n1)
1055 IF (ljsupu) jsupu = jsupu + jsupumnc1(mn)*tcosmn
1056 IF (ljsupv) jsupv = jsupv + jsupvmnc1(mn)*tcosmn
1057 END DO
1058
1059 IF (.not.lasym) GOTO 2000
1060
1061 DO mn = 1, mnmax_nyq
1062 m = nint(xm_nyq(mn)); n = nint(xn_nyq(mn))/nfp
1063 n1 = abs(n); sgn = sign(1,n)
1064 tsinmn = sinmu(m)*cosnv(n1) - sgn*cosmu(m)*sinnv(n1)
1065 IF (ljsupu) jsupu = jsupu + jsupumns1(mn)*tsinmn
1066 IF (ljsupv) jsupv = jsupv + jsupvmns1(mn)*tsinmn
1067 END DO
1068
1069 2000 CONTINUE
1070
1071 IF (ljsupu) jsupu = jsupu/gsqrt
1072 IF (ljsupv) jsupv = jsupv/gsqrt
1073
1074 END SUBROUTINE tosuvspace
1075
1076 SUBROUTINE loadrzl
1077 IMPLICIT NONE
1078C-----------------------------------------------
1079C L o c a l V a r i a b l e s
1080C-----------------------------------------------
1081 INTEGER :: rcc, rss, zsc, zcs, rsc, rcs, zcc, zss
1082 INTEGER :: mpol1, mn, m, n, n1
1083 REAL(rprec) :: sgn
1084C-----------------------------------------------
1085!
1086! Arrays must be stacked (and ns,ntor,mpol ordering imposed)
1087! as coefficients of cos(mu)*cos(nv), etc
1088! Only need R, Z components(not lambda, for now anyhow)
1089!
1090 IF (ALLOCATED(rzl_local)) RETURN
1091
1092 mpol1 = mpol-1
1093 rcc = 1; zsc = 1
1094 IF (.not.lasym) THEN
1095 IF (lthreed) THEN
1096 ntmax = 2
1097 rss = 2; zcs = 2
1098 ELSE
1099 ntmax = 1
1100 END IF
1101 ELSE
1102 IF (lthreed) THEN
1103 ntmax = 4
1104 rss = 2; rsc = 3; rcs = 4
1105 zcs = 2; zcc = 3; zss = 4
1106 ELSE
1107 ntmax = 2
1108 rsc = 2; zcc = 2
1109 END IF
1110 END IF
1111
1112! only ALLOCATE 2*ntmax, don't need lambdas
1113 zsc = 1+ntmax; zcs = zcs+ntmax; zcc = zcc+ntmax; zss = zss+ntmax
1114 ALLOCATE(rzl_local(ns,0:ntor,0:mpol1,2*ntmax), stat=n)
1115 IF (n .ne. 0) stop 'Allocation error in LoadRZL'
1116 rzl_local = 0
1117
1118 DO mn = 1, mnmax
1119 m = nint(xm(mn)); n = nint(xn(mn))/nfp; n1 = abs(n)
1120 sgn = sign(1, n)
1121 rzl_local(:,n1,m,rcc) = rzl_local(:,n1,m,rcc) + rmnc(mn,:)
1122 rzl_local(:,n1,m,zsc) = rzl_local(:,n1,m,zsc) + zmns(mn,:)
1123 IF (lthreed) THEN
1124 rzl_local(:,n1,m,rss) = rzl_local(:,n1,m,rss)
1125 1 + sgn*rmnc(mn,:)
1126 rzl_local(:,n1,m,zcs) = rzl_local(:,n1,m,zcs)
1127 1 - sgn*zmns(mn,:)
1128 END IF
1129 IF (lasym) THEN
1130 rzl_local(:,n1,m,rsc) = rzl_local(:,n1,m,rsc)
1131 1 + rmns(mn,:)
1132 rzl_local(:,n1,m,zcc) = rzl_local(:,n1,m,zcc)
1133 1 + zmnc(mn,:)
1134 IF (lthreed) THEN
1135 rzl_local(:,n1,m,rcs) = rzl_local(:,n1,m,rcs)
1136 1 - sgn*rmns(mn,:)
1137 rzl_local(:,n1,m,zss) = rzl_local(:,n1,m,zss)
1138 1 + sgn*zmnc(mn,:)
1139 END IF
1140 END IF
1141 END DO
1142
1143! ADDED by SAL for Vecpot calc
1144 IF (.not. ALLOCATED(chi)) ALLOCATE(chi(1:ns))
1145 DO mn = 1, ns
1146 chi(mn) = sum(iotaf(1:mn)*phipf(1:mn))
1147 END DO
1148
1149 END SUBROUTINE loadrzl
1150
1151 END MODULE read_wout_mod
1152
character(len= *), parameter vn_mgmode
Definition mgrid_mod.f:16
character(len=30), dimension(:), allocatable curlabel
Definition mgrid_mod.f:95
character(len=1) mgrid_mode
Definition mgrid_mod.f:102
integer nextcur
Definition mgrid_mod.f:78
character(len= *), parameter vn_nextcur
Definition mgrid_mod.f:16
Reading of wout VMEC output file.
character(len= *), parameter vn_bmns
subroutine readw_and_open(file_or_extension, ierr, iopen)
character(len= *), parameter vn_piota_type
real(rprec) aminor
character(len= *), parameter vn_bsupumns_sur
character pcurr_type
character(len= *), parameter vn_bsubvmns_sur
character(len= *), parameter vn_aspect
real(rprec), dimension(:,:), allocatable lmns
character(len= *), parameter vn_bsupumnc
character(len= *), parameter ln_ctor
character(len= *), parameter ln_tmod
character(len= *), parameter ln_vp
character(len= *), parameter vn_rmns
real(rprec) flmwgt
character(len= *), parameter ln_beta
character(len= *), parameter vn_radnod
character(len= *), parameter ln_bsupumns_sur
character(len= *), parameter ln_curlab
real(rprec), dimension(:), allocatable qfact
character(len= *), parameter ln_pbeta
real(rprec), dimension(:), allocatable y2thom
character(len= *), parameter vn_zacs
character(len= *), parameter ln_am_aux_f
real(rprec), dimension(:), allocatable fsqt
character(len= *), parameter vn_rbt1
character(len= *), parameter ln_maxmod_nyq
character(len= *), parameter vn_bsupumns
character(len= *), parameter vn_bsubsmns
character(len= *), parameter ln_extcur
character(len= *), parameter vn_pmass_type
character(len= *), parameter ln_zbs
character(len= *), parameter ln_fsq
character(len= *), parameter ln_amin
character(len= *), parameter ln_bsubvmns
character(len= *), parameter ln_ftolv
real(rprec) aspect
character(len= *), parameter ln_chipf
subroutine tosuvspace(s_in, u_in, v_in, gsqrt, bsupu, bsupv, jsupu, jsupv, lam)
real(rprec), dimension(:,:), allocatable bsubvmnc
character(len= *), parameter ln_maxit
character(len= *), parameter vn_lmnc
real(rprec), dimension(:,:), allocatable currvmns
real(rprec), dimension(:,:), allocatable zmns
character(len= *), parameter ln_bsubumns
character(len= *), parameter ln_jcurv
character(len= *), parameter vn_tbeta
character(len= *), parameter ln_bsubvmnc
character(len= *), parameter ln_presh
character(len= *), parameter ln_bsupumnc_sur
real(rprec), dimension(:), allocatable presf
real(rprec), dimension(:), allocatable rthom
real(rprec), dimension(:,:), allocatable bsubumnc
character(len= *), parameter ln_merc
real(rprec), dimension(:,:), allocatable bsupumnc
character(len= *), parameter ln_zacc
real(rprec), dimension(:), allocatable extcur
real(rprec), dimension(:,:), allocatable zaxis
character(len= *), parameter ln_thom
character(len= *), parameter ln_rmns
character(len= *), parameter vn_overr
character(len= *), parameter vn_ftolv
character(len= *), parameter ln_therm
character(len= *), parameter ln_bsubumnc_sur
real(rprec), dimension(:), allocatable pres
character(len= *), parameter vn_curlab
real(rprec), dimension(:), allocatable datathom
character(len= *), parameter vn_specw
character(len= *), parameter ln_lmns
real(rprec), dimension(:,:), allocatable rmnc
character(len= *), parameter vn_vp
real(rprec), dimension(:), allocatable am_aux_f
character(len= *), parameter vn_vol
character(len= *), parameter vn_magen
real(rprec), dimension(:), allocatable presmid
real(rprec), dimension(:,:), allocatable rmns
real(rprec) fsql
character(len= *), parameter ln_fsql
character(len= *), parameter vn_maxz
character(len= *), parameter vn_bsupvmns_sur
real(rprec), dimension(:), allocatable xn_nyq
character(len= *), parameter ln_racc
character(len= *), parameter vn_ah
character(len= *), parameter vn_abeta
real(rprec), dimension(:), allocatable buco
real(rprec), dimension(:,:), allocatable currumnc
character pmass_type
character(len= *), parameter vn_ac
character(len= *), parameter vn_presf
real(rprec), dimension(:,:), allocatable gmnc
character(len= *), parameter vn_chipf
real(rprec), dimension(:,:), allocatable bsupvmnc
real(rprec) fsqr
character(len= *), parameter vn_jcuru
character(len= *), parameter ln_recon
character(len= *), parameter ln_ai_aux_s
real(rprec), dimension(:), allocatable xn
character(len= *), parameter vn_free
character(len= *), parameter vn_version
character(len= *), parameter vn_buco
real(rprec) betator
real(rprec) fsqz
character(len= *), parameter vn_mse
character(len= *), parameter vn_phipf
character(len= *), parameter ln_bsubsmns
character(len= *), parameter ln_pmod
real(rprec), dimension(:), allocatable dcurr
real(rprec), dimension(:), allocatable equif
real(rprec), dimension(:), allocatable phipf
real(rprec) volume
character(len= *), parameter vn_ac_aux_f
character(len= *), parameter ln_radnod
character(len= *), parameter vn_bgrv
real(rprec), dimension(:,:), allocatable bmns
character(len= *), parameter vn_mgrid
character(len= *), parameter ln_extension
character(len= *), parameter vn_zmnc
real(rprec), dimension(:), allocatable pknots
character(len= *), parameter ln_chi
character(len= *), parameter vn_iotah
character(len= *), parameter vn_lmns
real(rprec), dimension(:), allocatable dsiobt
real(rprec), dimension(:), allocatable rstark
character(len= *), parameter ln_rbs
real(rprec), dimension(:), allocatable phip
character(len= *), parameter ln_ai
real(rprec), dimension(:,:), allocatable currumns
real(rprec), dimension(:), allocatable mass
character(len= *), parameter vn_bvco
character(len= *), parameter vn_fsqr
character(len= *), parameter vn_zacc
character(len= *), parameter ln_mse
character(len= *), parameter vn_bsubvmnc_sur
character(len= *), parameter vn_gam
real(rprec), dimension(:), allocatable ai
character(len= *), parameter ln_ac_aux_f
real(rprec), dimension(:), allocatable xm_nyq
real(rprec), dimension(:), allocatable bvco
character(len= *), parameter vn_betah
character(len= *), parameter ln_fsqr
real(rprec), dimension(:), allocatable jcurv
character(len= *), parameter ln_bsupvmns_sur
character(len= *), parameter ln_maxr
character(len= *), parameter vn_bsupvmnc
real(rprec), dimension(:), allocatable ystark
character(len= *), parameter ln_mgeo
real(rprec), dimension(:), allocatable y2stark
real(rprec), dimension(:), allocatable jdotb
character(len= *), parameter vn_mshear
character(len= *), parameter ln_bgrv
character(len= *), parameter vn_merc
real(rprec) volavgb
character(len= *), parameter vn_b0
real(rprec), dimension(:), allocatable xm
real(rprec), dimension(:), allocatable vp
character(len= *), parameter ln_zmnc
real(rprec), dimension(:), allocatable ythom
character(len= *), parameter ln_ac_aux_s
character(len= *), parameter vn_ai_aux_s
real(rprec), dimension(:), allocatable iotaf
character(len= *), parameter ln_free
subroutine compute_currents(ierror)
character(len= *), parameter ln_zbc
character(len= *), parameter vn_rbc
real(rprec) rmajor
character(len= *), parameter vn_zbs
character(len= *), parameter ln_phip
character(len= *), parameter ln_bsupvmns
character(len= *), parameter ln_buco
character(len= *), parameter vn_gmnc
character(len= *), parameter ln_potvac
real(rprec), dimension(:), allocatable chipf
real(rprec), dimension(:), allocatable beta_vol
real(rprec), dimension(:), allocatable specw
character(len= *), parameter vn_bsupumnc_sur
character(len= *), parameter vn_qfact
real(rprec), dimension(:), allocatable dwell
character(len= *), parameter ln_bmns
character(len= *), parameter vn_maxr
real(rprec) msewgt
character(len= *), parameter vn_bsubsmnc
real(rprec), dimension(:), allocatable jcuru
character(len= *), parameter vn_jcurv
real(rprec), dimension(:,:), allocatable bmnc
real(rprec), dimension(:), allocatable bdotb
character(len= *), parameter vn_zmns
character(len= *), parameter vn_jdotb
character(len= *), parameter ln_bsupvmnc
real(rprec), dimension(:), allocatable ac
character(len= *), parameter ln_equif
character(len= *), parameter ln_am
real(rprec), dimension(:), allocatable ac_aux_f
character(len= *), parameter vn_am_aux_s
real(rprec), dimension(:,:), allocatable currvmnc
character(len= *), parameter vn_tmod_nyq
character(len= *), parameter ln_gmns
real(rprec), dimension(:), allocatable overr
character(len= *), parameter ln_version
real(rprec) pfac
character(len= *), parameter vn_mwell
character(len= *), parameter ln_jcuru
character(len= *), parameter ln_bvco
character(len= *), parameter vn_fsql
character mgrid_file
character(len= *), parameter ln_ac
character(len= *), parameter vn_amin
character(len= *), parameter vn_rmnc
real(rprec), dimension(:), allocatable rmid
character(len= *), parameter ln_bsupvmnc_sur
real(rprec), dimension(:), allocatable chi
character(len= *), parameter vn_mgeo
character(len= *), parameter vn_bsubumns
character(len= *), parameter vn_bsupvmns
character(len= *), parameter vn_mass
character(len= *), parameter vn_modb
real(rprec) zmax_surf
character(len= *), parameter ln_vol
real(rprec) b0
character(len= *), parameter ln_racs
character(len= *), parameter vn_gmns
real(rprec), dimension(:,:), allocatable zmnc
character(len= *), parameter vn_atuname
character(len= *), parameter vn_am
character(len= *), parameter vn_fsq
character(len= *), parameter ln_rbt1
character(len= *), parameter vn_am_aux_f
real(rprec) betaxis
character(len= *), parameter vn_lar
character(len= *), parameter vn_rmaj
character(len= *), parameter ln_specw
real(rprec), dimension(:), allocatable ac_aux_s
logical lwout_opened
character(len= *), parameter vn_bmnc
character(len= *), parameter ln_bdotb
real(rprec) wb
character(len= *), parameter vn_ai
character(len= *), parameter ln_phipf
character(len= *), parameter ln_polmod
character(len= *), parameter ln_bmnc
character(len= *), parameter ln_mcurr
character(len= *), parameter vn_racc
real(rprec), dimension(:,:), allocatable bsupvmns
character(len= *), parameter ln_bsupumns
character(len= *), parameter vn_pmod
character(len= *), parameter ln_iotah
character(len= *), parameter ln_error
character(len= *), parameter ln_piota_type
real(rprec), dimension(:), allocatable dgeod
character(len= *), parameter vn_wdot
character(len= *), parameter vn_bsubumns_sur
character(len= *), parameter vn_asym
character(len= *), parameter vn_extension
character(len= *), parameter ln_maxz
real(rprec), dimension(:), allocatable phi
real(rprec) tswgt
character(len= *), parameter vn_tmod
character(len= *), parameter ln_mwell
real(rprec) rbtor0
real(rprec), dimension(:), allocatable am
character(len= *), parameter ln_mshear
real(rprec), dimension(:,:), allocatable bbc
real(rprec) rmin_surf
real(rprec) betatot
character(len= *), parameter ln_fsqz
real(rprec), dimension(:), allocatable wdot
real(rprec), dimension(:), allocatable am_aux_s
character(len= *), parameter vn_beta
character(len= *), parameter vn_pcurr_type
character(len= *), parameter ln_rmaj
character(len= *), parameter ln_rmnc
subroutine loadrzl
real(rprec), dimension(:), allocatable shear
character(len= *), parameter ln_lar
character(len= *), parameter vn_bsubvmns
character(len= *), parameter ln_pmass_type
real(rprec), dimension(:), allocatable curmid
character(len= *), parameter vn_tormod
real(rprec) rmax_surf
real(rprec), dimension(:), allocatable anglemse
real(rprec), dimension(:), allocatable iotas
character(len= *), parameter ln_magen
character(len= *), parameter ln_minr
character(len= *), parameter ln_presf
character(len= *), parameter ln_b0
character(len= *), parameter vn_equif
real(rprec), dimension(:), allocatable potvac
integer imatch_phiedge
real(rprec), dimension(:), allocatable datastark
real(rprec) version_
character(len= *), parameter ln_rbt0
real(rprec), dimension(:), allocatable ai_aux_f
character(len= *), parameter ln_fp
real(rprec), dimension(:,:), allocatable bsubsmnc
character input_extension
character(len= *), parameter vn_error
real(rprec) delphid
character(len= *), parameter vn_minr
character(len= *), parameter vn_extcur
character(len= *), parameter vn_pmod_nyq
character(len= *), parameter ln_mass
character(len= *), parameter vn_bsubumnc
real(rprec), dimension(:), allocatable qmid
character(len= *), parameter vn_presh
real(rprec), dimension(:,:), allocatable bsubumns
character(len= *), parameter vn_sgs
character(len= *), parameter ln_pmod_nyq
character(len= *), parameter vn_phi
character(len= *), parameter ln_rbc
character(len= *), parameter ln_aspect
character(len= *), parameter ln_bsubumns_sur
subroutine read_wout_deallocate
real(rprec) bcwgt
character(len= *), parameter vn_ctor
character(len= *), parameter ln_ai_aux_f
character(len= *), parameter vn_bdotb
character(len= *), parameter vn_racs
character(len= *), parameter vn_pbeta
character(len= *), parameter ln_bsubvmnc_sur
character(len= *), parameter vn_chi
real(rprec) machsq
character(len= *), parameter vn_polmod
real(rprec), dimension(:,:), allocatable bsubvmns
real(rprec), dimension(:), allocatable sknots
character(len= *), parameter vn_mcurr
character(len= *), parameter ln_zmns
real(rprec), dimension(:,:), allocatable bsupumns
real(rprec), dimension(:), allocatable dshear
character(len= *), parameter ln_bsubumnc
character(len= *), parameter vn_maxmod_nyq
character(len= *), parameter vn_bsubvmnc
character(len= *), parameter ln_bsupumnc
character(len= *), parameter ln_sgs
character(len= *), parameter ln_wdot
character(len= *), parameter ln_pcurr_type
character(len= *), parameter vn_thom
character(len= *), parameter ln_tbeta
real(rprec), dimension(:), allocatable qmeas
real(rprec) phidiam
character(len= *), parameter ln_maxmod
character(len= *), parameter vn_maxit
character(len= *), parameter vn_ac_aux_s
character(len= *), parameter vn_rbt0
real(rprec) betapol
character(len= *), parameter ln_lmnc
real(rprec), dimension(:), allocatable alfa
character(len= *), parameter ln_tmod_nyq
character(len= *), parameter vn_fp
character(len= *), parameter ln_bsubsmnc
character(len= *), parameter ln_mgrid
real(rprec), dimension(:), allocatable dmerc
real(rprec) wp
real(rprec) gamma
character(len= *), parameter ln_iotaf
real(rprec), dimension(:,:), allocatable lmnc
character(len= *), parameter vn_zbc
character(len= *), parameter vn_potvac
character(len= *), parameter ln_am_aux_s
character(len= *), parameter ln_phi
character piota_type
character(len= *), parameter ln_asym
character(len= *), parameter ln_tormod
character(len= *), parameter ln_abeta
character(len= *), parameter ln_qfact
character(len= *), parameter vn_maxmod
real(rprec), dimension(:), allocatable ai_aux_s
character(len= *), parameter vn_iotaf
real(rprec) rbtor
character(len= *), parameter ln_betah
real(rprec), dimension(:,:,:,:), allocatable rzl_local
real(rprec), dimension(:,:), allocatable raxis
real(rprec), dimension(:), allocatable bdotgradv
character(len= *), parameter vn_therm
character(len= *), parameter vn_rbs
character(len= *), parameter ln_gmnc
character(len= *), parameter vn_phip
character(len= *), parameter vn_fsqz
character(len= *), parameter vn_bsupvmnc_sur
character(len= *), parameter ln_modb
character(len= *), parameter ln_bsubvmns_sur
real(rprec) itor
character(len= *), parameter ln_gam
real(rprec), dimension(:,:), allocatable gmns
real(rprec) ionlarmor
real(rprec), dimension(:,:), allocatable bsubsmns
character(len= *), parameter ln_jdotb
character(len= *), parameter ln_zacs
character(len= *), parameter vn_ai_aux_f
character(len= *), parameter vn_bsubumnc_sur
real(rprec) ftolv
fault-tolerant file opening routines
real(dp), parameter mu0
real(dp), parameter one
real(dp), parameter zero
integer, parameter rprec
subroutine parse_extension(file_to_parse, file_or_extension, lnc)
Parse the first command-line argument into a filename.