16SUBROUTINE wrout(bsq, gsqrt, bsubu, bsubv, bsubs, bsupv, bsupu, rzl_array, gc_array, ier_flag)
80 INTEGER,
INTENT(in) :: ier_flag
83 REAL(rprec),
DIMENSION(mnmax,ns,3*MAX(ntmax/2,1)),
INTENT(inout),
TARGET :: rzl_array, gc_array
84 REAL(rprec),
DIMENSION(ns,nznt),
INTENT(inout) :: bsq, gsqrt, bsubu, bsubv, bsubs, bsupv, bsupu
86 REAL(rprec) :: qfact(ns)
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 '/)
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, &
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
114 REAL(rprec),
ALLOCATABLE,
DIMENSION(:,:) :: bsupumnc, bsupumns, bsupvmnc, bsupvmns
126 rmnc => rzl_array(:,:,1)
127 zmns => rzl_array(:,:,1+n1)
128 lmns => rzl_array(:,:,1+2*n1)
131 rmns => gc_array(:,:,1)
132 zmnc => gc_array(:,:,1+n1)
133 lmnc => gc_array(:,:,1+2*n1)
150 IF (istat .ne. 0) stop
'Error allocating arrays in VMEC WROUT'
171 READ (wout_file, *) vversion
174 CALL cdf_open(nwout, wout_file,
'w', iwout0)
175 IF (iwout0 .ne. 0) stop
'Error opening wout.nc file VMEC WROUT'
202 CALL cdf_define(nwout,
vn_error, ier_flag)
208 CALL cdf_define(nwout,
vn_b0,
b0)
218 CALL cdf_define(nwout,
vn_ftolv, ftolx1)
228 CALL cdf_define(nwout,
vn_pmod,
xm, dimname=mn1dim)
230 CALL cdf_define(nwout,
vn_tmod,
xn, dimname=mn1dim)
249 CALL cdf_define(nwout,
vn_am,
am(0:j), dimname=(/
'preset'/))
251 CALL cdf_define(nwout,
vn_ac,
ac(0:j), dimname=(/
'preset'/))
253 CALL cdf_define(nwout,
vn_ai,
ai(0:j), dimname=(/
'preset'/))
273 WHERE (
iotaf(1:ns) .NE. zero) qfact=one/
iotaf(1:ns)
275 CALL cdf_define(nwout,
vn_qfact, qfact(1:ns), dimname=r1dim)
280 CALL cdf_define(nwout,
vn_phi,
phi, dimname=r1dim)
284 CALL cdf_define(nwout,
vn_chi,
chi, dimname=r1dim)
294 CALL cdf_define(nwout,
vn_mass,
mass, dimname=r1dim)
296 CALL cdf_define(nwout,
vn_presh,
pres(1:ns), 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)
303 CALL cdf_define(nwout,
vn_phip,
phips(1:ns), dimname=r1dim)
324 CALL cdf_define(nwout,
vn_rmnc, rmnc, dimname=r2dim)
326 CALL cdf_define(nwout,
vn_zmns, zmns, dimname=r2dim)
328 CALL cdf_define(nwout,
vn_lmns, lmns, dimname=r2dim)
330 CALL cdf_define(nwout,
vn_gmnc, gmnc, dimname=r3dim)
332 CALL cdf_define(nwout,
vn_bmnc, bmnc, dimname=r3dim)
334 CALL cdf_define(nwout,
vn_bsubumnc, bsubumnc, dimname=r3dim)
336 CALL cdf_define(nwout,
vn_bsubvmnc, bsubvmnc, dimname=r3dim)
338 CALL cdf_define(nwout,
vn_bsubsmns, bsubsmns, dimname=r3dim)
342 CALL cdf_define(nwout,
vn_bsupumnc, bsupumnc, dimname=r3dim)
343 CALL cdf_define(nwout,
vn_bsupvmnc, bsupvmnc, dimname=r3dim)
346 CALL cdf_define(nwout,
vn_rmns, rmns, dimname=r2dim)
348 CALL cdf_define(nwout,
vn_zmnc, zmnc, dimname=r2dim)
350 CALL cdf_define(nwout,
vn_lmnc, lmnc, dimname=r2dim)
352 CALL cdf_define(nwout,
vn_gmns, gmns, dimname=r3dim)
354 CALL cdf_define(nwout,
vn_bmns, bmns, dimname=r3dim)
356 CALL cdf_define(nwout,
vn_bsubumns, bsubumns, dimname=r3dim)
358 CALL cdf_define(nwout,
vn_bsubvmns, bsubvmns, dimname=r3dim)
360 CALL cdf_define(nwout,
vn_bsubsmnc, bsubsmnc, dimname=r3dim)
364 CALL cdf_define(nwout,
vn_bsupumns, bsupumns, dimname=r3dim)
365 CALL cdf_define(nwout,
vn_bsupvmns, bsupvmns, dimname=r3dim)
393 CALL cdf_write(nwout,
vn_error, ier_flag)
400 CALL cdf_write(nwout,
vn_b0,
b0)
410 CALL cdf_write(nwout,
vn_ftolv, ftolx1)
435 ALLOCATE (xfinal(
neqs), stat=js)
436 IF (js .ne. 0) stop
'Allocation error for xfinal in WROUT!'
450 radius1:
DO js = 1, ns
452 CALL convert (rmnc1, zmns1, lmns1, rmns1, zmnc1, lmnc1, xfinal, js)
454 rmnc(:,js) = rmnc1(:)
455 zmns(:,js) = zmns1(:)
458 rmns(:,js) = rmns1(:)
459 zmnc(:,js) = zmnc1(:)
469 WHERE (nint(
xm) .le. 1) lmns(:,1) = lmns(:,2)
471 WHERE (mod(nint(
xm),2) .eq. 0)
472 lmns(:,js) =
p5*(lmns(:,js) + lmns(:,js-1))
474 lmns(:,js) =
p5*(
sm(js)*lmns(:,js) +
sp(js-1)*lmns(:,js-1))
483 WHERE (nint(
xm) .le. 1) lmnc(:,1) = lmnc(:,2)
485 WHERE (mod(nint(
xm),2) .eq. 0)
486 lmnc(:,js) =
p5*(lmnc(:,js) + lmnc(:,js-1))
488 lmnc(:,js) =
p5*(
sm(js)*lmnc(:,js) +
sp(js-1)*lmnc(:,js-1))
506 bsq(js,:nznt) = sqrt(2*abs(bsq(js,:nznt)-
pres(js)))
515 CALL symoutput (bsq, gsqrt, bsubu, bsubv, bsupu, &
517 bsqa, gsqrta, bsubua, bsubva, bsupua, &
522 radius2:
DO js = 2, ns
539 IF (m.eq.0 .or. n.eq.0) dmult = 2*dmult
554 bmn(mn) = bmn(mn) + tcosi*bsq(js,lk)
555 gmn(mn) = gmn(mn) + tcosi*gsqrt(js,lk)
556 bsubumn(mn) = bsubumn(mn) + tcosi*bsubu(js,lk)
557 bsubvmn(mn) = bsubvmn(mn) + tcosi*bsubv(js,lk)
558 bsubsmn(mn) = bsubsmn(mn) + tsini*bsubs(js,lk)
559 bsupumn(mn) = bsupumn(mn) + tcosi*bsupu(js,lk)
560 bsupvmn(mn) = bsupvmn(mn) + tcosi*bsupv(js,lk)
566 IF (js .eq. ns/2) bmodmn = bmn(1:
mnmax)
567 IF (js .eq. ns ) bmodmn1 = bmn(1:
mnmax)
571 bsubumnc(:,js) = bsubumn(:)
572 bsubvmnc(:,js) = bsubvmn(:)
573 bsubsmns(:,js) = bsubsmn(:)
574 bsupumnc(:,js) = bsupumn(:)
575 bsupvmnc(:,js) = bsupvmn(:)
589 bsubsmns(:,1) = 2*bsubsmns(:,2) - bsubsmns(:,3)
595 radius3:
DO js = 2, ns
609 IF (m.eq.0 .or. n.eq.0) dmult = 2*dmult
622 bmn(mn) = bmn(mn) + tsini*bsqa(jlk)
623 gmn(mn) = gmn(mn) + tsini*gsqrta(jlk,0)
625 bsubumn(mn) = bsubumn(mn) + tsini*bsubua(jlk)
626 bsubvmn(mn) = bsubvmn(mn) + tsini*bsubva(jlk)
627 bsubsmn(mn) = bsubsmn(mn) + tcosi*bsubsa(jlk)
629 bsupumn(mn) = bsupumn(mn) + tsini*bsupua(jlk)
630 bsupvmn(mn) = bsupvmn(mn) + tsini*bsupva(jlk)
638 bsubumns(:,js) = bsubumn(:)
639 bsubvmns(:,js) = bsubvmn(:)
640 bsubsmnc(:,js) = bsubsmn(:)
641 bsupumns(:,js) = bsupumn(:)
642 bsupvmns(:,js) = bsupvmn(:)
649 bsubsmnc(:,1) = 2*bsubsmnc(:,2) - bsubsmnc(:,3)
657 CALL cdf_write(nwout,
vn_rmnc, rmnc)
658 CALL cdf_write(nwout,
vn_zmns, zmns)
659 CALL cdf_write(nwout,
vn_lmns, lmns)
660 CALL cdf_write(nwout,
vn_gmnc, gmnc)
661 CALL cdf_write(nwout,
vn_bmnc, bmnc)
674 CALL cdf_write(nwout,
vn_am,
am(0:j))
676 CALL cdf_write(nwout,
vn_ac,
ac(0:j))
678 CALL cdf_write(nwout,
vn_ai,
ai(0:j))
694 CALL cdf_write(nwout,
vn_qfact, qfact(1:ns))
724 CALL cdf_write(nwout,
vn_vp,
vp(1:ns))
742 CALL cdf_write(nwout,
vn_rmns, rmns)
743 CALL cdf_write(nwout,
vn_zmnc, zmnc)
744 CALL cdf_write(nwout,
vn_lmnc, lmnc)
746 CALL cdf_write(nwout,
vn_gmns, gmns)
747 CALL cdf_write(nwout,
vn_bmns, bmns)
758 CALL cdf_close(nwout)
765 IF (
ALLOCATED(gmnc))
DEALLOCATE(gmnc, bmnc, bsubumnc, bsubvmnc, &
766 bsubsmns, bsupumnc, bsupvmnc )
767 IF (
ALLOCATED(gmns))
DEALLOCATE(gmns, bmns, bsubumns, bsubvmns, &
768 bsubsmnc, bsupumns, bsupvmns )
769 IF (
ALLOCATED(gmn))
DEALLOCATE (gmn, bmn, bsubumn, bsubvmn, &
770 bsubsmn, bsupumn, bsupvmn, stat=istat)
774 CALL freeb_data(rmnc1, zmns1, rmns1, zmnc1, bmodmn, bmodmn1)
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),...
subroutine freeb_data(rmnc, zmns, rmns, zmnc, bmodmn, bmodmn1)
Write out edge values of fields.
character(len= *), parameter vn_mgmode
character(len=30), dimension(:), allocatable curlabel
character(len=1) mgrid_mode
character(len= *), parameter vn_nextcur
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 ln_bdotb
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_bdotb
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
real(rprec), dimension(:), allocatable phip
radial derivative of phi/(2*pi) on half-grid
fault-tolerant file opening routines
real(rprec), dimension(:), allocatable, target potvac
real(rprec), parameter p5
real(rprec), dimension(:), pointer bzmn_e
real(rprec), dimension(:), pointer brmn_e
real(rprec), dimension(:), pointer azmn_e
real(rprec), dimension(:), pointer czmn_o
real(rprec), dimension(:), pointer armn_o
real(rprec), dimension(:), pointer armn_e
real(rprec), dimension(:), allocatable equif
radial force balance error: grad(p) - <j x B>
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
real(rprec), dimension(:), allocatable bvco
enclosed poloidal current profile
real(rprec), dimension(:), allocatable buco
enclosed toroidal current profile
real(rprec), dimension(:), allocatable sp
shalf(i+1)/sfull(i)
real(rprec), dimension(:), allocatable jcuru
poloidal current density
real(rprec) wp
kinetic/thermal energy (from pressure)
real(rprec) ctor
toroidal current (?)
real(rprec), dimension(:), allocatable presf
pressure profile on full-grid, mass/phip**gamma
real(rprec), dimension(:), allocatable phi
toroidal magnetic flux
integer neqs
total number of equations to evolve (size of xc)
real(rprec), dimension(:), allocatable overr
real(rprec), dimension(:), allocatable sm
shalf(i)/sfull(i)
real(rprec), dimension(:), allocatable jdotb
real(rprec), dimension(:), allocatable bdotgradv
real(rprec), dimension(:), allocatable specw
spectral width (diagnostic)
real(rprec), dimension(:), allocatable iotaf
rotational transform (full grid)
real(rprec) rbtor0
poloidal current at magnetic axis
real(rprec), dimension(:), allocatable pres
pressure profile
real(rprec), dimension(:), allocatable chipf
radial derivative of poloidal magnetic flux (full grid)
real(rprec) wb
magnetic energy: volume integral over B^2/2
real(rprec), dimension(:), allocatable chi
poloidal magnetic flux
real(rprec), dimension(:), allocatable phips
toroidal flux (same as phip), one-dimensional array
real(rprec), dimension(:), allocatable mass
mass profile on half-grid
real(rprec), dimension(:), allocatable bdotb
real(rprec), dimension(:), allocatable beta_vol
integer iter2
total number of iterations
real(rprec), dimension(:), allocatable iotas
rotational transform , on half radial mesh
real(rprec), dimension(:), allocatable phipf
radial derivative of toroidal magnetic flux (full grid)
real(rprec), dimension(:), allocatable jcurv
toroidal current density
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), 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
real(rprec), dimension(nsd) dmerc
real(rprec), dimension(nsd) dwell
real(rprec), dimension(nsd) dcurr
real(rprec), dimension(nsd) dshear
real(rprec), parameter cp5
real(rprec), parameter c2p0
real(rprec), dimension(:), allocatable, target xc
stacked array of scaled R, Z, Lambda Fourier coefficients (see above for stack order)
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 (?)
subroutine convert_sym(rmnss, zmncs)
Convert from internal representation to "physical" rmnss, zmncs Fourier form.
subroutine convert_asym(rmnsc, zmncc)
Convert from internal representation to "physical" rmnsc, zmncc Fourier form.
subroutine wrout(bsq, gsqrt, bsubu, bsubv, bsubs, bsupv, bsupu, rzl_array, gc_array, ier_flag)
Write the output files of VMEC.