VMEC 8.52
3D Equilibrium solver with nested flux surfaces.
Loading...
Searching...
No Matches
wrout.f90
Go to the documentation of this file.
1
3
16SUBROUTINE wrout(bsq, gsqrt, bsubu, bsubv, bsubs, bsupv, bsupu, rzl_array, gc_array, ier_flag)
17 USE vmec_main
18 USE vparams, p5 => cp5, two => c2p0
20 USE vmec_params
21 USE vmercier
23 USE xstuff
24 USE vmec_io
25 USE realspace, ONLY: phip, gsqrta=>z1, z1=>z1
26 USE vforces, ONLY: bsupua=>brmn_e, bsupva=>czmn_o, bsqa=>bzmn_e, &
27 bsubsa=>armn_e, bsubua=>azmn_e, bsubva=>armn_o
28 USE vacmod, ONLY: potvac
29 USE ezcdf
59 ln_mse, ln_thom, ln_next, &
74
76 USE mgrid_mod
77
78 IMPLICIT NONE
79
80 INTEGER, INTENT(in) :: ier_flag
81
82 ! reverse ns, mnmax for backwards compatibility
83 REAL(rprec), DIMENSION(mnmax,ns,3*MAX(ntmax/2,1)), INTENT(inout), TARGET :: rzl_array, gc_array ! only used as temporary storage !
84 REAL(rprec), DIMENSION(ns,nznt), INTENT(inout) :: bsq, gsqrt, bsubu, bsubv, bsubs, bsupv, bsupu
85
86 REAL(rprec) :: qfact(ns)
87
88 CHARACTER(LEN=*), PARAMETER, DIMENSION(1) :: &
89 r1dim = (/'radius'/), &
90 mn1dim = (/'mn_mode'/), &
91 mn2dim = (/'mn_mode_nyq'/), &
92 mnpddim = (/'mnpd'/), &
93 currg = (/'ext_current'/), &
94 currl = (/'current_label'/)
95 CHARACTER(LEN=*), PARAMETER, DIMENSION(2) :: &
96 r2dim = (/'mn_mode', 'radius ' /), &
97 r3dim = (/'mn_mode_nyq','radius '/)
98
99 INTEGER :: j, js, jlk, mn, lk, &
100 m, n, k, iwout0, n1, nwout, istat, indx1(1)
101 REAL(rprec) :: dmult, tcosi, tsini, vversion, sgn, tmult, &
102 ftolx1
103 REAL(rprec), POINTER, DIMENSION(:,:) :: rmnc, rmns, zmns, zmnc, lmns, lmnc
104 REAL(rprec), ALLOCATABLE, DIMENSION(:,:) :: gmnc, bmnc, gmns, bmns, &
105 bsubumnc, bsubvmnc, bsubsmns, bsubumns, bsubvmns, bsubsmnc
106 REAL(rprec), DIMENSION(mnmax) :: rmnc1, zmns1, lmns1, &
107 rmns1, zmnc1, lmnc1, bmodmn, bmodmn1
108 REAL(rprec), DIMENSION(:), ALLOCATABLE :: gmn, bmn, &
109 bsubumn, bsubvmn, bsubsmn, bsupumn, bsupvmn
110 CHARACTER(LEN=120) :: wout_file
111 REAL(rprec), DIMENSION(:), ALLOCATABLE :: xfinal
112
113 ! ELIMINATE THESE EVENTUALLY
114 REAL(rprec), ALLOCATABLE, DIMENSION(:,:) :: bsupumnc, bsupumns, bsupvmnc, bsupvmns
115
116
117 ! THIS SUBROUTINE CREATES THE FILE WOUT.
118 ! IT CONTAINS THE CYLINDRICAL COORDINATE SPECTRAL COEFFICIENTS
119 ! RMN,ZMN (full), LMN (half_mesh - CONVERTED FROM INTERNAL full REPRESENTATION),
120 ! AS WELL AS COEFFICIENTS (ON NYQ MESH) FOR COMPUTED QUANTITIES:
121 ! BSQ, BSUPU,V, BSUBU,V, GSQRT (HALF); BSUBS (FULL-CONVERTED IN JXBFORCE)
122
123
124 ! Pointer assignments for storage arrays
125 n1 = max(1,ntmax/2)
126 rmnc => rzl_array(:,:,1) ! store COS(mu-nv) components
127 zmns => rzl_array(:,:,1+n1) ! store SIN(mu-nv)
128 lmns => rzl_array(:,:,1+2*n1) ! store SIN(mu-nv)
129
130 IF (lasym) THEN
131 rmns => gc_array(:,:,1) ! store SIN(mu-nv)
132 zmnc => gc_array(:,:,1+n1) ! store COS(mu-nv)
133 lmnc => gc_array(:,:,1+2*n1) ! store COS(mu-nv)
134 END IF
135
136 ALLOCATE (gmn(mnmax_nyq), bmn(mnmax_nyq), &
137 bsubumn(mnmax_nyq), bsubvmn(mnmax_nyq), bsubsmn(mnmax_nyq), &
138 bsupumn(mnmax_nyq), bsupvmn(mnmax_nyq), stat=istat)
139
140 ALLOCATE (gmnc(mnmax_nyq,ns), bmnc(mnmax_nyq,ns), &
141 bsubumnc(mnmax_nyq,ns), bsubvmnc(mnmax_nyq,ns), &
142 bsubsmns(mnmax_nyq,ns), bsupumnc(mnmax_nyq,ns), &
143 bsupvmnc(mnmax_nyq,ns), stat=istat)
144 IF (lasym) THEN
145 ALLOCATE (gmns(mnmax_nyq,ns), bmns(mnmax_nyq,ns), &
146 bsubumns(mnmax_nyq,ns), bsubvmns(mnmax_nyq,ns), &
147 bsubsmnc(mnmax_nyq,ns), bsupumns(mnmax_nyq,ns), &
148 bsupvmns(mnmax_nyq,ns), stat=istat)
149 END IF
150 IF (istat .ne. 0) stop 'Error allocating arrays in VMEC WROUT'
151
152 ! IF (nextcur .eq. 0) THEN
153 ! DO j = SIZE(extcur), 1, -1
154 ! IF (extcur(j) .ne. zero) THEN
155 ! nextcur = j
156 ! EXIT
157 ! END IF
158 ! END DO
159 ! END IF
160
161 ! ftol info evaluated here!
162 indx1=maxloc(ns_array)
163 ftolx1=ftol_array(indx1(1))
164
165 ! NYQUIST FREQUENCY REQUIRES FACTOR OF 1/2
166 IF (mnyq .ne. 0) cosmui(:,mnyq) = p5*cosmui(:,mnyq)
167 IF (nnyq .ne. 0) cosnv(:,nnyq) = p5*cosnv(:,nnyq)
168
169 ! use wout_file as temporary storage for parsing version number to REAL
170 wout_file = version_
171 READ (wout_file, *) vversion
172
173 wout_file = 'wout_' // trim(input_extension) // '.nc'
174 CALL cdf_open(nwout, wout_file, 'w', iwout0)
175 IF (iwout0 .ne. 0) stop 'Error opening wout.nc file VMEC WROUT'
176
177 !================================
178 ! Define Variables
179 !================================
180 ! Scalars
181 CALL cdf_define(nwout, vn_version, vversion)
182 CALL cdf_define(nwout, vn_extension, input_extension)
183 CALL cdf_define(nwout, vn_mgrid, mgrid_file)
184 CALL cdf_define(nwout, vn_pcurr_type, pcurr_type)
185 CALL cdf_define(nwout, vn_pmass_type, pmass_type)
186 CALL cdf_define(nwout, vn_piota_type, piota_type)
187 CALL cdf_define(nwout, vn_magen, wb)
188 CALL cdf_define(nwout, vn_therm, wp)
189 CALL cdf_define(nwout, vn_gam, gamma)
190 CALL cdf_define(nwout, vn_maxr, rmax_surf)
191 CALL cdf_define(nwout, vn_minr, rmin_surf)
192 CALL cdf_define(nwout, vn_maxz, zmax_surf)
193 CALL cdf_define(nwout, vn_fp, nfp)
194 CALL cdf_define(nwout, vn_radnod, ns)
195 CALL cdf_define(nwout, vn_polmod, mpol)
196 CALL cdf_define(nwout, vn_tormod, ntor)
197 CALL cdf_define(nwout, vn_maxmod, mnmax)
198 CALL cdf_define(nwout, vn_maxmod_nyq, mnmax_nyq)
199 CALL cdf_define(nwout, vn_maxit, iter2)
200 CALL cdf_define(nwout, vn_asym, lasym)
201 CALL cdf_define(nwout, vn_free, lfreeb)
202 CALL cdf_define(nwout, vn_error, ier_flag)
203 CALL cdf_define(nwout, vn_aspect, aspect)
204 CALL cdf_define(nwout, vn_beta, betatot)
205 CALL cdf_define(nwout, vn_pbeta, betapol)
206 CALL cdf_define(nwout, vn_tbeta, betator)
207 CALL cdf_define(nwout, vn_abeta, betaxis)
208 CALL cdf_define(nwout, vn_b0, b0)
209 CALL cdf_define(nwout, vn_rbt0, rbtor0)
210 CALL cdf_define(nwout, vn_rbt1, rbtor)
211 CALL cdf_define(nwout, vn_sgs, nint(signgs))
212 CALL cdf_define(nwout, vn_lar, ionlarmor)
213 CALL cdf_define(nwout, vn_modb, volavgb)
214 CALL cdf_define(nwout, vn_ctor, ctor)
215 CALL cdf_define(nwout, vn_amin, aminor_p)
216 CALL cdf_define(nwout, vn_rmaj, rmajor_p)
217 CALL cdf_define(nwout, vn_vol, volume_p)
218 CALL cdf_define(nwout, vn_ftolv, ftolx1)
219 CALL cdf_define(nwout, vn_fsqr, fsqr)
220 CALL cdf_define(nwout, vn_fsqz, fsqz)
221 CALL cdf_define(nwout, vn_fsql, fsql)
222
223 CALL cdf_define(nwout, vn_nextcur, nextcur)
224 CALL cdf_define(nwout, vn_extcur, extcur(1:nextcur), dimname=currg)
225 CALL cdf_define(nwout, vn_mgmode, mgrid_mode)
226
227 ! 1D Arrays
228 CALL cdf_define(nwout, vn_pmod, xm, dimname=mn1dim)
229 CALL cdf_setatt(nwout, vn_pmod, ln_pmod)
230 CALL cdf_define(nwout, vn_tmod, xn, dimname=mn1dim)
231 CALL cdf_setatt(nwout, vn_tmod, ln_tmod)
232 CALL cdf_define(nwout, vn_pmod_nyq, xm_nyq, dimname=mn2dim)
233 CALL cdf_setatt(nwout, vn_pmod_nyq, ln_pmod_nyq)
234 CALL cdf_define(nwout, vn_tmod_nyq, xn_nyq, dimname=mn2dim)
235 CALL cdf_setatt(nwout, vn_tmod_nyq, ln_tmod_nyq)
236
237 CALL cdf_define(nwout, vn_racc, raxis_cc(0:ntor), dimname=(/'n_tor'/))
238 CALL cdf_setatt(nwout, vn_racc, ln_racc)
239 CALL cdf_define(nwout, vn_zacs, zaxis_cs(0:ntor), dimname=(/'n_tor'/))
240 CALL cdf_setatt(nwout, vn_zacs, ln_zacs)
241 IF (lasym) THEN
242 CALL cdf_define(nwout, vn_racs, raxis_cs(0:ntor), dimname=(/'n_tor'/))
243 CALL cdf_setatt(nwout, vn_racs, ln_racs)
244 CALL cdf_define(nwout, vn_zacc, zaxis_cc(0:ntor), dimname=(/'n_tor'/))
245 CALL cdf_setatt(nwout, vn_zacc, ln_zacc)
246 END IF
247
248 j = SIZE(am)-1
249 CALL cdf_define(nwout, vn_am, am(0:j), dimname=(/'preset'/))
250 j = SIZE(ac)-1
251 CALL cdf_define(nwout, vn_ac, ac(0:j), dimname=(/'preset'/))
252 j = SIZE(ai)-1
253 CALL cdf_define(nwout, vn_ai, ai(0:j), dimname=(/'preset'/))
254
255 j = SIZE(am_aux_s)
256 CALL cdf_define(nwout, vn_am_aux_s, am_aux_s(1:j), dimname=(/'ndfmax'/))
257 j = SIZE(am_aux_f)
258 CALL cdf_define(nwout, vn_am_aux_f, am_aux_f(1:j), dimname=(/'ndfmax'/))
259 j = SIZE(ai_aux_s)
260 CALL cdf_define(nwout, vn_ai_aux_s, ai_aux_s(1:j), dimname=(/'ndfmax'/))
261 j = SIZE(ai_aux_f)
262 CALL cdf_define(nwout, vn_ai_aux_f, ai_aux_f(1:j), dimname=(/'ndfmax'/))
263 j = SIZE(ac_aux_s)
264 CALL cdf_define(nwout, vn_ac_aux_s, ac_aux_s(1:j), dimname=(/'ndfmax'/))
265 j = SIZE(ac_aux_f)
266 CALL cdf_define(nwout, vn_ac_aux_f, ac_aux_f(1:j), dimname=(/'ndfmax'/))
267
268
269 CALL cdf_define(nwout, vn_iotaf, iotaf(1:ns), dimname=r1dim)
270 CALL cdf_setatt(nwout, vn_iotaf, ln_iotaf)
271
272 qfact=huge(qfact)
273 WHERE (iotaf(1:ns) .NE. zero) qfact=one/iotaf(1:ns)
274
275 CALL cdf_define(nwout, vn_qfact, qfact(1:ns), dimname=r1dim)
276 CALL cdf_setatt(nwout, vn_iotaf, ln_qfact)
277
278 CALL cdf_define(nwout, vn_presf, presf, dimname=r1dim)
279 CALL cdf_setatt(nwout, vn_presf, ln_presf, units='Pa')
280 CALL cdf_define(nwout, vn_phi, phi, dimname=r1dim)
281 CALL cdf_setatt(nwout, vn_phi, ln_phi, units='wb')
282 CALL cdf_define(nwout, vn_phipf, phipf, dimname=r1dim)
283 CALL cdf_setatt(nwout, vn_phipf, ln_phipf)
284 CALL cdf_define(nwout, vn_chi, chi, dimname=r1dim)
285 CALL cdf_setatt(nwout, vn_chi, ln_chi, units='wb')
286 !CALL cdf_define(nwout, vn_chipf, phipf, dimname=r1dim) ! TODO: wrong quantity !!!
287 CALL cdf_define(nwout, vn_chipf, chipf, dimname=r1dim)
288 CALL cdf_setatt(nwout, vn_chipf, ln_chipf)
289 CALL cdf_define(nwout, vn_jcuru, jcuru, dimname=r1dim)
290 CALL cdf_define(nwout, vn_jcurv, jcurv, dimname=r1dim)
291
292 CALL cdf_define(nwout, vn_iotah, iotas(1:ns), dimname=r1dim)
293 CALL cdf_setatt(nwout, vn_iotah, ln_iotah)
294 CALL cdf_define(nwout, vn_mass, mass, dimname=r1dim)
295 CALL cdf_setatt(nwout, vn_mass, ln_mass)
296 CALL cdf_define(nwout, vn_presh, pres(1:ns), dimname=r1dim)
297 CALL cdf_setatt(nwout, vn_presh, ln_presh, units='Pa')
298 CALL cdf_define(nwout, vn_betah, beta_vol, dimname=r1dim)
299 CALL cdf_define(nwout, vn_buco, buco, dimname=r1dim)
300 CALL cdf_define(nwout, vn_bvco, bvco, dimname=r1dim)
301 CALL cdf_define(nwout, vn_vp, vp(1:ns), dimname=r1dim)
302 CALL cdf_define(nwout, vn_specw, specw, dimname=r1dim) ! NOTE: specw on full grid
303 CALL cdf_define(nwout, vn_phip, phips(1:ns), dimname=r1dim)
304 CALL cdf_define(nwout, vn_overr, overr(1:ns), dimname=r1dim)
305
306 CALL cdf_define(nwout, vn_jdotb, jdotb, dimname=r1dim)
307 CALL cdf_define(nwout, vn_bgrv, bdotgradv, dimname=r1dim)
308
309 CALL cdf_define(nwout, vn_merc, dmerc, dimname=r1dim)
310 CALL cdf_define(nwout, vn_mshear, dshear, dimname=r1dim)
311 CALL cdf_define(nwout, vn_mwell, dwell, dimname=r1dim)
312 CALL cdf_define(nwout, vn_mcurr, dcurr, dimname=r1dim)
313 CALL cdf_define(nwout, vn_mgeo, dgeod, dimname=r1dim)
314 CALL cdf_define(nwout, vn_equif, equif, dimname=r1dim)
315
316 IF (lfreeb .and. nextcur.gt.0 .and. ALLOCATED(curlabel)) THEN
317 CALL cdf_define(nwout, vn_curlab, curlabel(1:nextcur), dimname=currl)
318 ! SAL for potvac
319 CALL cdf_define(nwout, vn_potvac, potvac, dimname=mnpddim)
320 ENDIF
321
322 ! 2D Arrays
323 CALL cdf_define(nwout, vn_rmnc, rmnc, dimname=r2dim)
324 CALL cdf_setatt(nwout, vn_rmnc, ln_rmnc, units='m')
325 CALL cdf_define(nwout, vn_zmns, zmns, dimname=r2dim)
326 CALL cdf_setatt(nwout, vn_zmns, ln_zmns, units='m')
327 CALL cdf_define(nwout, vn_lmns, lmns, dimname=r2dim)
328 CALL cdf_setatt(nwout, vn_lmns, ln_lmns)
329 CALL cdf_define(nwout, vn_gmnc, gmnc, dimname=r3dim)
330 CALL cdf_setatt(nwout, vn_gmnc, ln_gmnc)
331 CALL cdf_define(nwout, vn_bmnc, bmnc, dimname=r3dim)
332 CALL cdf_setatt(nwout, vn_bmnc, ln_bmnc)
333 CALL cdf_define(nwout, vn_bsubumnc, bsubumnc, dimname=r3dim)
334 CALL cdf_setatt(nwout, vn_bsubumnc, ln_bsubumnc)
335 CALL cdf_define(nwout, vn_bsubvmnc, bsubvmnc, dimname=r3dim)
336 CALL cdf_setatt(nwout, vn_bsubvmnc, ln_bsubvmnc)
337 CALL cdf_define(nwout, vn_bsubsmns, bsubsmns, dimname=r3dim)
338 CALL cdf_setatt(nwout, vn_bsubsmns, ln_bsubsmns)
339
340 ! ELIMINATE THESE EVENTUALLY: DON'T NEED THEM - CAN COMPUTE FROM GSQRT
341 CALL cdf_define(nwout, vn_bsupumnc, bsupumnc, dimname=r3dim)
342 CALL cdf_define(nwout, vn_bsupvmnc, bsupvmnc, dimname=r3dim)
343
344 IF (lasym) then
345 CALL cdf_define(nwout, vn_rmns, rmns, dimname=r2dim)
346 CALL cdf_setatt(nwout, vn_rmns, ln_rmns, units='m')
347 CALL cdf_define(nwout, vn_zmnc, zmnc, dimname=r2dim)
348 CALL cdf_setatt(nwout, vn_zmnc, ln_zmnc, units='m')
349 CALL cdf_define(nwout, vn_lmnc, lmnc, dimname=r2dim)
350 CALL cdf_setatt(nwout, vn_lmnc, ln_lmnc)
351 CALL cdf_define(nwout, vn_gmns, gmns, dimname=r3dim)
352 CALL cdf_setatt(nwout, vn_gmns, ln_gmns)
353 CALL cdf_define(nwout, vn_bmns, bmns, dimname=r3dim)
354 CALL cdf_setatt(nwout, vn_bmns, ln_bmns)
355 CALL cdf_define(nwout, vn_bsubumns, bsubumns, dimname=r3dim)
356 CALL cdf_setatt(nwout, vn_bsubumns, ln_bsubumns)
357 CALL cdf_define(nwout, vn_bsubvmns, bsubvmns, dimname=r3dim)
358 CALL cdf_setatt(nwout, vn_bsubvmns, ln_bsubvmns)
359 CALL cdf_define(nwout, vn_bsubsmnc, bsubsmnc, dimname=r3dim)
360 CALL cdf_setatt(nwout, vn_bsubsmnc, ln_bsubsmnc)
361
362 ! ELIMINATE THESE EVENTUALLY: DON'T NEED THEM
363 CALL cdf_define(nwout, vn_bsupumns, bsupumns, dimname=r3dim)
364 CALL cdf_define(nwout, vn_bsupvmns, bsupvmns, dimname=r3dim)
365 end if
366
367 !================================
368 ! Write Variables
369 !================================
370 ! Scalars
371 CALL cdf_write(nwout, vn_version, vversion)
372 CALL cdf_write(nwout, vn_extension, input_extension)
373 CALL cdf_write(nwout, vn_mgrid, mgrid_file)
374 CALL cdf_write(nwout, vn_pcurr_type, pcurr_type)
375 CALL cdf_write(nwout, vn_piota_type, piota_type)
376 CALL cdf_write(nwout, vn_pmass_type, pmass_type)
377 CALL cdf_write(nwout, vn_magen, wb)
378 CALL cdf_write(nwout, vn_therm, wp)
379 CALL cdf_write(nwout, vn_gam, gamma)
380 CALL cdf_write(nwout, vn_maxr, rmax_surf)
381 CALL cdf_write(nwout, vn_minr, rmin_surf)
382 CALL cdf_write(nwout, vn_maxz, zmax_surf)
383 CALL cdf_write(nwout, vn_fp, nfp)
384 CALL cdf_write(nwout, vn_radnod, ns)
385 CALL cdf_write(nwout, vn_polmod, mpol)
386 CALL cdf_write(nwout, vn_tormod, ntor)
387 CALL cdf_write(nwout, vn_maxmod, mnmax)
388 CALL cdf_write(nwout, vn_maxmod_nyq, mnmax_nyq)
389 CALL cdf_write(nwout, vn_maxit, iter2)
390 CALL cdf_write(nwout, vn_asym, lasym)
391 CALL cdf_write(nwout, vn_free, lfreeb)
392 CALL cdf_write(nwout, vn_error, ier_flag)
393
394 CALL cdf_write(nwout, vn_aspect, aspect)
395 CALL cdf_write(nwout, vn_beta, betatot)
396 CALL cdf_write(nwout, vn_pbeta, betapol)
397 CALL cdf_write(nwout, vn_tbeta, betator)
398 CALL cdf_write(nwout, vn_abeta, betaxis)
399 CALL cdf_write(nwout, vn_b0, b0)
400 CALL cdf_write(nwout, vn_rbt0, rbtor0)
401 CALL cdf_write(nwout, vn_rbt1, rbtor)
402 CALL cdf_write(nwout, vn_sgs, nint(signgs))
403 CALL cdf_write(nwout, vn_lar, ionlarmor)
404 CALL cdf_write(nwout, vn_modb, volavgb)
405 CALL cdf_write(nwout, vn_ctor, ctor/mu0)
406 CALL cdf_write(nwout, vn_amin, aminor_p)
407 CALL cdf_write(nwout, vn_rmaj, rmajor_p)
408 CALL cdf_write(nwout, vn_vol, volume_p)
409 CALL cdf_write(nwout, vn_ftolv, ftolx1)
410 CALL cdf_write(nwout, vn_fsql, fsql)
411 CALL cdf_write(nwout, vn_fsqr, fsqr)
412 CALL cdf_write(nwout, vn_fsqz, fsqz)
413
414 CALL cdf_write(nwout, vn_nextcur, nextcur)
415 IF (nextcur .gt. 0) THEN
416 CALL cdf_write(nwout, vn_extcur, extcur(1:nextcur))
417 CALL cdf_write(nwout, vn_mgmode, mgrid_mode)
418 ENDIF
419 IF (lfreeb) THEN
420 ! TODO: write current labels
421 CALL cdf_write(nwout, vn_potvac, potvac)
422 END IF
423
424 ! 1D Arrays
425 CALL cdf_write(nwout, vn_pmod, xm)
426 CALL cdf_write(nwout, vn_tmod, xn)
427 CALL cdf_write(nwout, vn_pmod_nyq, xm_nyq)
428 CALL cdf_write(nwout, vn_tmod_nyq, xn_nyq)
429
430
431 ! here start the convert-to-regular-Fourier-output part
432
433
434 ALLOCATE (xfinal(neqs), stat=js) ! re-use js for allocation return code
435 IF (js .ne. 0) stop 'Allocation error for xfinal in WROUT!'
436 xfinal = xc ! ignore passed rzl_array !!!
437 ! TODO: scalxc ???
438 ! --> rzl_array is re-used to store temporary symmetric rmnc, ... !
439 ! --> gc_array is re-used to store temporary asymmetric rmns, ... ! (for lasym)
440
441 ! MUST CONVERT m=1 MODES... FROM INTERNAL TO PHYSICAL FORM
442 ! Extrapolation of m=0 Lambda (cs) modes, which are not evolved at j=1, done in CONVERT
443
444 IF (lthreed) CALL convert_sym (xfinal(1+mns*(rss-1)), xfinal(1+irzloff+mns*(zcs-1)))
445 IF (lasym) CALL convert_asym (xfinal(1+mns*(rsc-1)), xfinal(1+irzloff+mns*(zcc-1)))
446
447 ! CONVERT TO rmnc, zmns, lmns, etc EXTERNAL representation (without internal mscale, nscale)
448 ! IF B^v ~ phip + lamu, MUST DIVIDE BY phipf(js) below to maintain old-style format
449 radius1: DO js = 1, ns
450
451 CALL convert (rmnc1, zmns1, lmns1, rmns1, zmnc1, lmnc1, xfinal, js)
452
453 rmnc(:,js) = rmnc1(:)
454 zmns(:,js) = zmns1(:)
455 lmns(:,js) = (lmns1(:)/phipf(js)) * lamscale
456 IF (lasym) THEN
457 rmns(:,js) = rmns1(:)
458 zmnc(:,js) = zmnc1(:)
459 lmnc(:,js) = (lmnc1(:)/phipf(js)) * lamscale
460 END IF
461
462 END DO radius1
463
464 DEALLOCATE (xfinal)
465
466 ! INTERPOLATE LAMBDA ONTO HALF-MESH FOR BACKWARDS CONSISTENCY WITH EARLIER VERSIONS OF VMEC
467 ! AND SMOOTHS POSSIBLE UNPHYSICAL "WIGGLE" ON RADIAL MESH
468 WHERE (nint(xm) .le. 1) lmns(:,1) = lmns(:,2) ! constant extrapolation of lambda
469 DO js = ns,2,-1
470 WHERE (mod(nint(xm),2) .eq. 0)
471 lmns(:,js) = p5*(lmns(:,js) + lmns(:,js-1))
472 ELSEWHERE
473 lmns(:,js) = p5*(sm(js)*lmns(:,js) + sp(js-1)*lmns(:,js-1))
474 END WHERE
475 END DO
476
477 lmns(:,1) = 0 ! NOTE: needed constant extrapolation above only in averaging step!
478 raxis_cc(0:ntor) = rmnc(1:ntor+1,1)
479 zaxis_cs(0:ntor) = zmns(1:ntor+1,1)
480
481 IF (lasym) then
482 WHERE (nint(xm) .le. 1) lmnc(:,1) = lmnc(:,2)
483 DO js = ns,2,-1
484 WHERE (mod(nint(xm),2) .eq. 0)
485 lmnc(:,js) = p5*(lmnc(:,js) + lmnc(:,js-1))
486 ELSEWHERE
487 lmnc(:,js) = p5*(sm(js)*lmnc(:,js) + sp(js-1)*lmnc(:,js-1))
488 END WHERE
489 END DO
490
491 lmnc(:,1) = 0 ! NOTE: needed constant extrapolation above only in averaging step!
492 raxis_cs(0:ntor) = rmns(1:ntor+1,1)
493 zaxis_cc(0:ntor) = zmnc(1:ntor+1,1)
494 end if
495
496
497
498
499
500 ! here comes the realspace-to-Fourier transform part...
501
502
503 ! COMPUTE |B| = SQRT(|B|**2) and store in bsq, bsqa
504 DO js = 2, ns
505 bsq(js,:nznt) = sqrt(2*abs(bsq(js,:nznt)-pres(js)))
506 END DO
507
508 tmult = p5/r0scale**2
509 !SPH: FIXED THIS 03-05-07 TO CALL symmetrization routine
510 IF (lasym) THEN
511 ! Changed integration norm in fixaray, SPH012314
512 tmult = 2*tmult
513 bsubs(1,:) = 0
514 CALL symoutput (bsq, gsqrt, bsubu, bsubv, bsupu, &
515 bsupv, bsubs, &
516 bsqa, gsqrta, bsubua, bsubva, bsupua, &
517 bsupva, bsubsa )
518 END IF
519
520 ! Fourier-transform derived quantities for each surface individually
521 radius2: DO js = 2, ns
522 gmn = 0
523 bmn = 0
524
525 bsubumn = 0
526 bsubvmn = 0
527 bsubsmn = 0
528
529 bsupumn = 0
530 bsupvmn = 0
531
532 mn2: DO mn = 1, mnmax_nyq
533 n = nint(xn_nyq(mn))/nfp
534 m = nint(xm_nyq(mn))
535 n1 = abs(n)
536
537 dmult = mscale(m)*nscale(n1)*tmult
538 IF (m.eq.0 .or. n.eq.0) dmult = 2*dmult
539
540 sgn = sign(1, n)
541
542 lk = 0
543 DO j = 1, ntheta2
544 DO k = 1, nzeta
545 lk = lk + 1
546
547 ! cos(mu - nv)
548 tcosi = dmult*(cosmui(j,m)*cosnv(k,n1) + sgn*sinmui(j,m)*sinnv(k,n1))
549
550 ! sin(mu - nv)
551 tsini = dmult*(sinmui(j,m)*cosnv(k,n1) - sgn*cosmui(j,m)*sinnv(k,n1))
552
553 bmn(mn) = bmn(mn) + tcosi*bsq(js,lk)
554 gmn(mn) = gmn(mn) + tcosi*gsqrt(js,lk)
555 bsubumn(mn) = bsubumn(mn) + tcosi*bsubu(js,lk)
556 bsubvmn(mn) = bsubvmn(mn) + tcosi*bsubv(js,lk)
557 bsubsmn(mn) = bsubsmn(mn) + tsini*bsubs(js,lk)
558 bsupumn(mn) = bsupumn(mn) + tcosi*bsupu(js,lk)
559 bsupvmn(mn) = bsupvmn(mn) + tcosi*bsupv(js,lk)
560 END DO
561 END DO
562 END DO mn2
563
564 ! needed for freeb_data below
565 IF (js .eq. ns/2) bmodmn = bmn(1:mnmax)
566 IF (js .eq. ns ) bmodmn1 = bmn(1:mnmax)
567
568 gmnc(:,js) = gmn(:)
569 bmnc(:,js) = bmn(:)
570 bsubumnc(:,js) = bsubumn(:)
571 bsubvmnc(:,js) = bsubvmn(:)
572 bsubsmns(:,js) = bsubsmn(:)
573 bsupumnc(:,js) = bsupumn(:)
574 bsupvmnc(:,js) = bsupvmn(:)
575 END DO radius2
576
577 ! endpoint values at magnetic axis
578 gmnc(:,1) = 0
579 bmnc(:,1) = 0
580
581 bsubumnc(:,1) = 0
582 bsubvmnc(:,1) = 0
583
584 ! NOTE: This assumes that bsubs (from which bsubsmns is computed) is on the full grid.
585 ! However, since crmn_o is passed as bsubs in fileout(),
586 ! HERE, bsubs is ACTUALLY on the HALF-grid !!! WTF ???
587 ! (The full-grid bsubs array is a LOCAL variable in jxbforce().)
588 bsubsmns(:,1) = 2*bsubsmns(:,2) - bsubsmns(:,3) ! extrapolation on full grid
589
590 bsupumnc(:,1) = 0
591 bsupvmnc(:,1) = 0
592
593 IF (lasym) then
594 radius3: DO js = 2, ns
595 gmn = 0
596 bmn = 0
597 bsubumn = 0
598 bsubvmn = 0
599 bsubsmn = 0
600 bsupumn = 0
601 bsupvmn = 0
602 mn3: DO mn = 1, mnmax_nyq
603 n = nint(xn_nyq(mn))/nfp
604 m = nint(xm_nyq(mn))
605 n1 = abs(n)
606
607 dmult = mscale(m)*nscale(n1)*tmult
608 IF (m.eq.0 .or. n.eq.0) dmult = 2*dmult
609
610 sgn = sign(1, n)
611
612 lk = 0
613 jlk = js
614 DO j = 1, ntheta2
615 DO k = 1, nzeta
616 lk = lk + 1
617
618 tcosi = dmult*(cosmui(j,m)*cosnv(k,n1) + sgn*sinmui(j,m)*sinnv(k,n1))
619 tsini = dmult*(sinmui(j,m)*cosnv(k,n1) - sgn*cosmui(j,m)*sinnv(k,n1))
620
621 bmn(mn) = bmn(mn) + tsini*bsqa(jlk)
622 gmn(mn) = gmn(mn) + tsini*gsqrta(jlk,0) ! TODO: what is this zero index?
623
624 bsubumn(mn) = bsubumn(mn) + tsini*bsubua(jlk)
625 bsubvmn(mn) = bsubvmn(mn) + tsini*bsubva(jlk)
626 bsubsmn(mn) = bsubsmn(mn) + tcosi*bsubsa(jlk)
627
628 bsupumn(mn) = bsupumn(mn) + tsini*bsupua(jlk)
629 bsupvmn(mn) = bsupvmn(mn) + tsini*bsupva(jlk)
630 jlk = jlk+ns
631 END DO
632 END DO
633 END DO mn3
634
635 gmns(:,js) = gmn(:)
636 bmns(:,js) = bmn(:)
637 bsubumns(:,js) = bsubumn(:)
638 bsubvmns(:,js) = bsubvmn(:)
639 bsubsmnc(:,js) = bsubsmn(:)
640 bsupumns(:,js) = bsupumn(:)
641 bsupvmns(:,js) = bsupvmn(:)
642 END DO radius3
643
644 gmns(:,1) = 0
645 bmns(:,1) = 0
646 bsubumns(:,1) = 0
647 bsubvmns(:,1) = 0
648 bsubsmnc(:,1) = 2*bsubsmnc(:,2) - bsubsmnc(:,3)
649 bsupumns(:,1) = 0
650 bsupvmns(:,1) = 0
651 end if
652
653 ! WRITE OUT ARRAYS
654 CALL cdf_write(nwout, vn_racc, raxis_cc(0:ntor))
655 CALL cdf_write(nwout, vn_zacs, zaxis_cs(0:ntor))
656 CALL cdf_write(nwout, vn_rmnc, rmnc)
657 CALL cdf_write(nwout, vn_zmns, zmns)
658 CALL cdf_write(nwout, vn_lmns, lmns)
659 CALL cdf_write(nwout, vn_gmnc, gmnc) !Half mesh
660 CALL cdf_write(nwout, vn_bmnc, bmnc) !Half mesh
661 CALL cdf_write(nwout, vn_bsubumnc, bsubumnc) !Half mesh
662 CALL cdf_write(nwout, vn_bsubvmnc, bsubvmnc) !Half mesh
663 CALL cdf_write(nwout, vn_bsubsmns, bsubsmns) !Full mesh
664
665 ! GET RID OF THESE EVENTUALLY: DON'T NEED THEM (can express in terms of lambdas)
666 CALL cdf_write(nwout, vn_bsupumnc, bsupumnc)
667 CALL cdf_write(nwout, vn_bsupvmnc, bsupvmnc)
668
669 ! FULL-MESH quantities
670 ! NOTE: jdotb is in units_of_A (1/mu0 incorporated in jxbforce...)
671 ! prior to version 6.00, this was output in internal VMEC units...
672 j = SIZE(am)-1
673 CALL cdf_write(nwout, vn_am, am(0:j))
674 j = SIZE(ac)-1
675 CALL cdf_write(nwout, vn_ac, ac(0:j))
676 j = SIZE(ai)-1
677 CALL cdf_write(nwout, vn_ai, ai(0:j))
678
679 j = SIZE(am_aux_s)
680 CALL cdf_write(nwout, vn_am_aux_s, am_aux_s(1:j))
681 j = SIZE(am_aux_f)
682 CALL cdf_write(nwout, vn_am_aux_f, am_aux_f(1:j))
683 j = SIZE(ac_aux_s)
684 CALL cdf_write(nwout, vn_ac_aux_s, ac_aux_s(1:j))
685 j = SIZE(ac_aux_f)
686 CALL cdf_write(nwout, vn_ac_aux_f, ac_aux_f(1:j))
687 j = SIZE(ai_aux_s)
688 CALL cdf_write(nwout, vn_ai_aux_s, ai_aux_s(1:j))
689 j = SIZE(ai_aux_f)
690 CALL cdf_write(nwout, vn_ai_aux_f, ai_aux_f(1:j))
691
692 CALL cdf_write(nwout, vn_iotaf, iotaf(1:ns))
693 CALL cdf_write(nwout, vn_qfact, qfact(1:ns))
694 CALL cdf_write(nwout, vn_presf, presf/mu0) ! NOTE: scaling !!!
695 CALL cdf_write(nwout, vn_phi, phi)
696 CALL cdf_write(nwout, vn_phipf, twopi*signgs*phipf) ! NOTE: scaling !!!
697 CALL cdf_write(nwout, vn_chi, chi)
698 CALL cdf_write(nwout, vn_chipf, twopi*signgs*chipf) ! NOTE: scaling !!!
699 CALL cdf_write(nwout, vn_jcuru, jcuru/mu0) ! NOTE: scaling !!!
700 CALL cdf_write(nwout, vn_jcurv, jcurv/mu0) ! NOTE: scaling !!!
701 CALL cdf_write(nwout, vn_jdotb, jdotb)
702 CALL cdf_write(nwout, vn_bgrv, bdotgradv)
703
704 ! HALF-MESH quantities
705 iotas(1) = 0
706 mass(1) = 0
707 pres(1) = 0
708 phip(1) = 0
709 buco(1) = 0
710 bvco(1) = 0
711 vp(1) = 0
712 overr(1) = 0
713 specw(1) = 1
714 beta_vol(1) = 0
715
716 CALL cdf_write(nwout, vn_iotah, iotas(1:ns))
717 CALL cdf_write(nwout, vn_mass, mass/mu0) ! NOTE: scaling !!!
718 CALL cdf_write(nwout, vn_presh, pres(1:ns)/mu0) ! NOTE: scaling !!!
719 CALL cdf_write(nwout, vn_betah, beta_vol)
720 CALL cdf_write(nwout, vn_buco, buco)
721 CALL cdf_write(nwout, vn_bvco, bvco)
722 CALL cdf_write(nwout, vn_vp, vp(1:ns))
723 CALL cdf_write(nwout, vn_specw, specw)
724 CALL cdf_write(nwout, vn_phip, phips(1:ns))
725 CALL cdf_write(nwout, vn_overr, overr(1:ns))
726
727 ! MERCIER_CRITERION
728 CALL cdf_write(nwout, vn_merc, dmerc(1:ns))
729 CALL cdf_write(nwout, vn_mshear, dshear(1:ns))
730 CALL cdf_write(nwout, vn_mwell, dwell(1:ns))
731 CALL cdf_write(nwout, vn_mcurr, dcurr(1:ns))
732 CALL cdf_write(nwout, vn_mgeo, dgeod(1:ns))
733
734 CALL cdf_write(nwout, vn_equif, equif(1:ns))
735
736 IF (lasym) THEN
737 CALL cdf_write(nwout, vn_racs, raxis_cs(0:ntor))
738 CALL cdf_write(nwout, vn_zacc, zaxis_cc(0:ntor))
739
740 CALL cdf_write(nwout, vn_rmns, rmns)
741 CALL cdf_write(nwout, vn_zmnc, zmnc)
742 CALL cdf_write(nwout, vn_lmnc, lmnc)
743
744 CALL cdf_write(nwout, vn_gmns, gmns)
745 CALL cdf_write(nwout, vn_bmns, bmns)
746
747 CALL cdf_write(nwout, vn_bsubumns, bsubumns)
748 CALL cdf_write(nwout, vn_bsubvmns, bsubvmns)
749 CALL cdf_write(nwout, vn_bsubsmnc, bsubsmnc)
750
751 ! GET RID OF THESE EVENTUALLY: DON'T NEED THEM
752 CALL cdf_write(nwout, vn_bsupumns, bsupumns)
753 CALL cdf_write(nwout, vn_bsupvmns, bsupvmns)
754 END IF
755
756 CALL cdf_close(nwout)
757
758 ! RESTORE nyq ENDPOINT VALUES
759 IF (mnyq .ne. 0) cosmui(:,mnyq) = 2*cosmui(:,mnyq)
760 IF (nnyq .ne. 0) cosnv(:,nnyq) = 2*cosnv(:,nnyq)
761
762 ! DEALLOCATIONS
763 IF (ALLOCATED(gmnc)) DEALLOCATE(gmnc, bmnc, bsubumnc, bsubvmnc, &
764 bsubsmns, bsupumnc, bsupvmnc )
765 IF (ALLOCATED(gmns)) DEALLOCATE(gmns, bmns, bsubumns, bsubvmns, &
766 bsubsmnc, bsupumns, bsupvmns )
767 IF (ALLOCATED(gmn)) DEALLOCATE (gmn, bmn, bsubumn, bsubvmn, &
768 bsubsmn, bsupumn, bsupvmn, stat=istat)
769
770 ! FREE BOUNDARY DATA
771 ! rmnc1, ... contain the LCFS values at this point
772 CALL freeb_data(rmnc1, zmns1, rmns1, zmnc1, bmodmn, bmodmn1)
773
774 rzl_array = 0
775
776END SUBROUTINE wrout
subroutine convert(rmnc, zmns, lmns, rmns, zmnc, lmnc, rzl_array, js)
Convert internal mode representation to standard form for output (coefficients of cos(mu-nv),...
Definition convert.f90:17
subroutine freeb_data(rmnc, zmns, rmns, zmnc, bmodmn, bmodmn1)
Write out edge values of fields.
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
character(len= *), parameter vn_piota_type
character(len= *), parameter vn_aspect
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
character(len= *), parameter ln_beta
character(len= *), parameter vn_radnod
character(len= *), parameter ln_curlab
character(len= *), parameter ln_pbeta
character(len= *), parameter vn_zacs
character(len= *), parameter ln_am_aux_f
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_chipf
character(len= *), parameter ln_maxit
character(len= *), parameter vn_lmnc
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_merc
character(len= *), parameter ln_zacc
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 vn_curlab
character(len= *), parameter vn_specw
character(len= *), parameter ln_lmns
character(len= *), parameter vn_vp
character(len= *), parameter vn_vol
character(len= *), parameter vn_magen
character(len= *), parameter vn_maxz
character(len= *), parameter ln_racc
character(len= *), parameter vn_abeta
character(len= *), parameter vn_ac
character(len= *), parameter vn_presf
character(len= *), parameter vn_chipf
character(len= *), parameter vn_jcuru
character(len= *), parameter ln_ai_aux_s
character(len= *), parameter vn_free
character(len= *), parameter vn_version
character(len= *), parameter vn_buco
character(len= *), parameter vn_phipf
character(len= *), parameter ln_bsubsmns
character(len= *), parameter ln_pmod
character(len= *), parameter vn_ac_aux_f
character(len= *), parameter ln_radnod
character(len= *), parameter vn_bgrv
character(len= *), parameter vn_mgrid
character(len= *), parameter ln_extension
character(len= *), parameter vn_zmnc
character(len= *), parameter ln_chi
character(len= *), parameter vn_iotah
character(len= *), parameter vn_lmns
character(len= *), parameter ln_rbs
character(len= *), parameter ln_ai
character(len= *), parameter vn_bvco
character(len= *), parameter vn_fsqr
character(len= *), parameter vn_zacc
character(len= *), parameter ln_mse
character(len= *), parameter vn_gam
character(len= *), parameter ln_ac_aux_f
character(len= *), parameter vn_betah
character(len= *), parameter ln_maxr
character(len= *), parameter vn_bsupvmnc
character(len= *), parameter ln_mgeo
character(len= *), parameter vn_mshear
character(len= *), parameter ln_bgrv
character(len= *), parameter vn_merc
character(len= *), parameter vn_b0
character(len= *), parameter ln_zmnc
character(len= *), parameter ln_ac_aux_s
character(len= *), parameter vn_ai_aux_s
character(len= *), parameter ln_free
character(len= *), parameter ln_zbc
character(len= *), parameter vn_rbc
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
character(len= *), parameter vn_qfact
character(len= *), parameter ln_bmns
character(len= *), parameter vn_maxr
character(len= *), parameter vn_bsubsmnc
character(len= *), parameter vn_jcurv
character(len= *), parameter vn_zmns
character(len= *), parameter vn_jdotb
character(len= *), parameter ln_bsupvmnc
character(len= *), parameter ln_equif
character(len= *), parameter ln_am
character(len= *), parameter vn_am_aux_s
character(len= *), parameter vn_tmod_nyq
character(len= *), parameter ln_gmns
character(len= *), parameter ln_version
character(len= *), parameter vn_mwell
character(len= *), parameter ln_jcuru
character(len= *), parameter ln_bvco
character(len= *), parameter vn_fsql
character(len= *), parameter ln_ac
character(len= *), parameter vn_amin
character(len= *), parameter vn_rmnc
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
character(len= *), parameter ln_vol
character(len= *), parameter ln_racs
character(len= *), parameter vn_gmns
character(len= *), parameter vn_am
character(len= *), parameter vn_fsq
character(len= *), parameter ln_rbt1
character(len= *), parameter vn_am_aux_f
character(len= *), parameter vn_lar
character(len= *), parameter vn_rmaj
character(len= *), parameter ln_specw
character(len= *), parameter vn_bmnc
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
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
character(len= *), parameter vn_wdot
character(len= *), parameter vn_asym
character(len= *), parameter vn_extension
character(len= *), parameter ln_maxz
character(len= *), parameter vn_tmod
character(len= *), parameter ln_mwell
character(len= *), parameter ln_mshear
character(len= *), parameter vn_beta
character(len= *), parameter vn_pcurr_type
character(len= *), parameter ln_rmaj
character(len= *), parameter ln_rmnc
character(len= *), parameter ln_lar
character(len= *), parameter vn_bsubvmns
character(len= *), parameter ln_pmass_type
character(len= *), parameter vn_tormod
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
character(len= *), parameter ln_rbt0
character(len= *), parameter ln_fp
character(len= *), parameter vn_error
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
character(len= *), parameter vn_presh
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 vn_ctor
character(len= *), parameter ln_ai_aux_f
character(len= *), parameter vn_racs
character(len= *), parameter vn_pbeta
character(len= *), parameter vn_chi
character(len= *), parameter vn_polmod
character(len= *), parameter vn_mcurr
character(len= *), parameter ln_zmns
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 ln_tbeta
character(len= *), parameter ln_maxmod
character(len= *), parameter vn_maxit
character(len= *), parameter vn_ac_aux_s
character(len= *), parameter vn_rbt0
character(len= *), parameter ln_lmnc
character(len= *), parameter ln_tmod_nyq
character(len= *), parameter vn_fp
character(len= *), parameter ln_bsubsmnc
character(len= *), parameter ln_mgrid
character(len= *), parameter ln_iotaf
character(len= *), parameter vn_zbc
character(len= *), parameter vn_potvac
character(len= *), parameter ln_am_aux_s
character(len= *), parameter ln_phi
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
character(len= *), parameter vn_iotaf
character(len= *), parameter ln_betah
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 ln_modb
character(len= *), parameter ln_gam
character(len= *), parameter ln_jdotb
character(len= *), parameter ln_zacs
character(len= *), parameter vn_ai_aux_f
real(rprec), dimension(:,:), allocatable, target z1
Definition realspace.f90:11
real(rprec), dimension(:), allocatable phip
radial derivative of phi/(2*pi) on half-grid
Definition realspace.f90:29
fault-tolerant file opening routines
real(rprec), dimension(:), allocatable, target potvac
Definition vacmod.f90:29
real(rprec), parameter p5
Definition vacmod.f90:13
real(rprec), dimension(:), pointer bzmn_e
Definition vforces.f90:28
real(rprec), dimension(:), pointer brmn_e
Definition vforces.f90:21
real(rprec), dimension(:), pointer azmn_e
Definition vforces.f90:26
real(rprec), dimension(:), pointer czmn_o
Definition vforces.f90:31
real(rprec), dimension(:), pointer armn_o
Definition vforces.f90:20
real(rprec), dimension(:), pointer armn_e
Definition vforces.f90:19
real(rprec), dimension(0:20) ac
array of coefficients in phi-series for the quantity d(Icurv)/ds = toroidal current density * Vprime,...
logical lfreeb
character(len=20) pcurr_type
real(rprec), dimension(ndatafmax) ac_aux_f
real(rprec), dimension(ndatafmax) ai_aux_s
character(len=20) pmass_type
real(rprec) gamma
real(rprec), dimension(ndatafmax) am_aux_s
character(len=20) piota_type
real(rprec), dimension(ndatafmax) ai_aux_f
real(rprec), dimension(0:20) am
array of coefficients in phi-series for mass (NWT/m**2)
real(rprec), dimension(ndatafmax) am_aux_f
real(rprec), dimension(0:ntord) raxis_cs
character(len=200) mgrid_file
real(rprec), dimension(ndatafmax) ac_aux_s
character(len=100) input_extension
real(rprec), dimension(0:ntord) raxis_cc
logical lasym
real(rprec), dimension(nigroup) extcur
real(rprec), dimension(0:20) ai
array of coefficients in phi-series for iota (ncurr=0)
integer nfp
integer ntor
integer nzeta
integer mpol
real(rprec), dimension(0:ntord) zaxis_cs
real(rprec), dimension(100) ftol_array
real(rprec), dimension(0:ntord) zaxis_cc
integer, dimension(100) ns_array
real(rprec) ionlarmor
Definition vmec_io.f90:9
real(rprec) rmax_surf
Definition vmec_io.f90:22
real(rprec) aminor_p
Definition vmec_io.f90:10
real(rprec) betapol
Definition vmec_io.f90:13
real(rprec) zmax_surf
Definition vmec_io.f90:24
real(rprec) betator
Definition vmec_io.f90:14
real(rprec) betaxis
Definition vmec_io.f90:15
real(rprec) rmin_surf
Definition vmec_io.f90:23
real(rprec) betatot
Definition vmec_io.f90:12
real(rprec) b0
Definition vmec_io.f90:16
real(rprec) volavgb
Definition vmec_io.f90:8
real(rprec) volume_p
Definition vmec_io.f90:17
real(rprec) rmajor_p
Definition vmec_io.f90:11
real(rprec), dimension(:), allocatable equif
radial force balance error: grad(p) - <j x B>
Definition vmec_main.f90:45
integer irzloff
offset in xc array between R,Z,L components
real(rprec) rbtor
poloidal current at LCFS
real(rprec), dimension(:), allocatable vp
radial derivative of enclosed volume
Definition vmec_main.f90:56
real(rprec), dimension(:), allocatable bvco
enclosed poloidal current profile
Definition vmec_main.f90:43
real(rprec), dimension(:), allocatable buco
enclosed toroidal current profile
Definition vmec_main.f90:42
real(rprec), dimension(:), allocatable sp
shalf(i+1)/sfull(i)
Definition vmec_main.f90:26
real(rprec) fsql
Definition vmec_main.f90:96
real(rprec), dimension(:), allocatable jcuru
poloidal current density
Definition vmec_main.f90:39
real(rprec) wp
kinetic/thermal energy (from pressure)
real(rprec) fsqr
Definition vmec_main.f90:94
real(rprec) fsqz
Definition vmec_main.f90:95
real(rprec) ctor
toroidal current (?)
real(rprec), dimension(:), allocatable presf
pressure profile on full-grid, mass/phip**gamma
Definition vmec_main.f90:66
real(rprec), dimension(:), allocatable phi
toroidal magnetic flux
Definition vmec_main.f90:37
integer neqs
total number of equations to evolve (size of xc)
real(rprec), dimension(:), allocatable overr
Definition vmec_main.f90:54
real(rprec), dimension(:), allocatable sm
shalf(i)/sfull(i)
Definition vmec_main.f90:25
real(rprec), dimension(:), allocatable jdotb
Definition vmec_main.f90:41
real(rprec) aspect
Definition vmec_main.f90:86
real(rprec), dimension(:), allocatable bdotgradv
Definition vmec_main.f90:44
real(rprec), dimension(:), allocatable specw
spectral width (diagnostic)
Definition vmec_main.f90:46
real(rprec), dimension(:), allocatable iotaf
rotational transform (full grid)
Definition vmec_main.f90:34
real(rprec) rbtor0
poloidal current at magnetic axis
real(rprec), dimension(:), allocatable pres
pressure profile
Definition vmec_main.f90:55
real(rprec), dimension(:), allocatable chipf
radial derivative of poloidal magnetic flux (full grid)
Definition vmec_main.f90:36
real(rprec) wb
magnetic energy: volume integral over B^2/2
logical lthreed
real(rprec), dimension(:), allocatable chi
poloidal magnetic flux
Definition vmec_main.f90:64
real(rprec) r0scale
Definition vmec_main.f90:90
real(rprec), dimension(:), allocatable phips
toroidal flux (same as phip), one-dimensional array
Definition vmec_main.f90:68
real(rprec), dimension(:), allocatable mass
mass profile on half-grid
Definition vmec_main.f90:71
real(rprec), dimension(:), allocatable beta_vol
Definition vmec_main.f90:38
integer iter2
total number of iterations
real(rprec), dimension(:), allocatable iotas
rotational transform , on half radial mesh
Definition vmec_main.f90:69
real(rprec), dimension(:), allocatable phipf
radial derivative of toroidal magnetic flux (full grid)
Definition vmec_main.f90:35
real(rprec), dimension(:), allocatable jcurv
toroidal current density
Definition vmec_main.f90:40
real(rprec) signgs
sign of Jacobian : must be =1 (right-handed) or =-1 (left-handed)
real(rprec), dimension(:), allocatable mscale
array for norming theta-trig functions (internal use only) so that the discrete SUM[cos(mu)*cos(m'u)]...
real(rprec) lamscale
real(rprec), dimension(:), allocatable nscale
array for norming zeta -trig functions (internal use only)
character(len= *), parameter version_
integer ntmax
number of contributing Fourier basis function (can be 1, 2 or 4); assigned in read_indata()
real(rprec), dimension(:,:), allocatable sinmui
real(rprec), dimension(:), allocatable, target xm_nyq
real(rprec), dimension(:), allocatable, target xm
real(rprec), dimension(:,:), allocatable cosmui
real(rprec), dimension(:,:), allocatable sinnv
real(rprec), dimension(:), allocatable, target xn
real(rprec), dimension(:), allocatable, target xn_nyq
real(rprec), dimension(:,:), allocatable cosnv
real(rprec), dimension(nsd) dgeod
Definition vmercier.f90:12
real(rprec), dimension(nsd) dmerc
Definition vmercier.f90:11
real(rprec), dimension(nsd) dwell
Definition vmercier.f90:9
real(rprec), dimension(nsd) dcurr
Definition vmercier.f90:10
real(rprec), dimension(nsd) dshear
Definition vmercier.f90:8
real(rprec), parameter cp5
Definition vparams.f90:37
real(rprec), parameter c2p0
Definition vparams.f90:40
real(rprec), dimension(:), allocatable, target xc
stacked array of scaled R, Z, Lambda Fourier coefficients (see above for stack order)
Definition xstuff.f90:40
subroutine symoutput(bsq, gsqrt, bsubu, bsubv, bsupu, bsupv, bsubs, bsqa, gsqrta, bsubua, bsubva, bsupua, bsupva, bsubsa)
Symmetrize some quantities so that they can be output (?)
Definition symoutput.f90:22
subroutine convert_sym(rmnss, zmncs)
Convert from internal representation to "physical" rmnss, zmncs Fourier form.
Definition totzsp.f90:217
subroutine convert_asym(rmnsc, zmncc)
Convert from internal representation to "physical" rmnsc, zmncc Fourier form.
Definition totzsp.f90:380
subroutine wrout(bsq, gsqrt, bsubu, bsubv, bsubs, bsupv, bsupu, rzl_array, gc_array, ier_flag)
Write the output files of VMEC.
Definition wrout.f90:17