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