VMEC 8.52
3D Equilibrium solver with nested flux surfaces.
Loading...
Searching...
No Matches
fileout.f90
Go to the documentation of this file.
1
3
7SUBROUTINE fileout(ier_flag)
8 USE vmec_main
10 USE realspace
12 USE vforces
13 USE xstuff, ONLY: xc, gc, xsave
14 IMPLICIT NONE
15
16 INTEGER, INTENT(inout) :: ier_flag
17
18 INTEGER :: istat
19 INTEGER :: loc_ier_flag
20 INTEGER :: js
21 INTEGER :: first0
22 LOGICAL :: lterm
23
24 REAL(rprec), ALLOCATABLE :: br_out(:), bz_out(:)
25
26 CHARACTER(LEN=*), PARAMETER, DIMENSION(0:10) :: werror = (/ &
27 'EXECUTION TERMINATED NORMALLY ', &
28 'INITIAL JACOBIAN CHANGED SIGN (IMPROVE INITIAL GUESS) ', &
29 'FORCE RESIDUALS EXCEED FTOL: MORE ITERATIONS REQUIRED ', &
30 .ne..ne.'VMEC INDATA ERROR: NCURR1 but BLOAT1. ', &
31 'MORE THAN 75 JACOBIAN ITERATIONS (DECREASE DELT) ', &
32 'ERROR READING INPUT FILE OR NAMELIST ', &
33 'NEW AXIS GUESS STILL FAILED TO GIVE GOOD JACOBIAN ', &
34 'PHIEDGE HAS WRONG SIGN IN VACUUM SUBROUTINE ', &
35 'NS ARRAY MUST NOT BE ALL ZEROES ', &
36 'ERROR READING MGRID FILE ', &
37 'VAC-VMEC I_TOR MISMATCH : BOUNDARY MAY ENCLOSE EXT. COIL ' /)
38
39 ! COMPUTE REMAINING COVARIANT COMPONENT OF B (BSUBS),
40 ! CYLINDRICAL COMPONENTS OF B (BR, BPHI, BZ), AND
41 ! AVERAGE EQUILIBRIUM PROPERTIES AT END OF RUN
42 iequi = 1
43 lterm = (ier_flag.eq.norm_term_flag .or. ier_flag.eq.successful_term_flag)
44
45 loc_ier_flag = ier_flag
46 if (ier_flag .eq. successful_term_flag) then
47 loc_ier_flag = norm_term_flag
48 end if
49
50 IF (lterm) THEN
51 ! Must save first value if in "restart" mode
52 first0 = first
53 CALL funct3d (istat)
54
55 ! The sign of the jacobian MUST multiply phi to get the physically correct toroidal flux
56 ! Below few lines compute the toroidal flux profile from phip by quadrature
57 phi(1) = zero
58 DO js = 2, ns
59 phi(js) = phi(js-1) + phip(js)
60 END DO
61 phi = (signgs*twopi*hs)*phi
62
63 first = first0
64
65 ALLOCATE(br_out(nrzt), bz_out(nrzt), stat=istat)
66
67 gc = xc
68 ! br, bz, bsubu, bsubv, tau, rzl_array, ier_flag
69 CALL eqfor(br_out, bz_out, clmn, blmn, rcon(1,1), gc, ier_flag)
70
71 END IF
72
73 IF (ASSOCIATED(bzmn_o)) THEN
74
75 ! Call WROUT to write output
76 ! or
77 ! error message if lterm = false
78 ! bsq, gsqrt, bsubu, bsubv, bsubs, bsupv,
80 crmn_e, xsave, gc, loc_ier_flag, lterm)
81 ! bsupu, rzl_array, gc_array, ier_flag TODO: lterm omitted in argument list ?
82
83 ! These are the last few lines that appear on screen / in the threed1 file
84 print 120, trim(werror(loc_ier_flag))
85 IF (lterm) print 10, trim(input_extension), ijacob
86 IF (nthreed .gt. 0) THEN
87 WRITE (nthreed,120) trim(werror(loc_ier_flag))
88 IF (lterm) then
89 WRITE (nthreed, 10) trim(input_extension), ijacob
90 end if
91 END IF
92
93 END IF
94
95 10 FORMAT(' FILE : ',a,/,' NUMBER OF JACOBIAN RESETS = ',i4,/)
96120 FORMAT(/1x,a,/)
97
98 IF (ALLOCATED(br_out)) THEN
99 DEALLOCATE (br_out, bz_out)
100 END IF
101
102END SUBROUTINE fileout
subroutine eqfor(br, bz, bsubu, bsubv, tau, rzl_array, ier_flag)
Basis physics analysis and evaluaton of force balance. This is where most of the contents of the thre...
Definition eqfor.f90:16
subroutine fileout(ier_flag)
Write the output files.
Definition fileout.f90:8
subroutine funct3d(ier_flag)
Evaluate the three-dimensional MHD energy functional. Think of this as the "forward model" that tells...
Definition funct3d.f90:14
real(rprec), dimension(:,:), allocatable rcon
spectral condensation term in
Definition realspace.f90:15
real(rprec), dimension(:), allocatable phip
radial derivative of phi/(2*pi) on half-grid
Definition realspace.f90:29
real(rprec), dimension(:), allocatable, target clmn
Definition vforces.f90:17
real(rprec), dimension(:), pointer azmn_o
Definition vforces.f90:27
real(rprec), dimension(:), allocatable, target blmn
Definition vforces.f90:16
real(rprec), dimension(:), pointer czmn_e
Definition vforces.f90:30
real(rprec), dimension(:), pointer bzmn_o
Definition vforces.f90:29
real(rprec), dimension(:), pointer crmn_e
Definition vforces.f90:23
real(rprec), dimension(:), pointer crmn_o
Definition vforces.f90:24
real(rprec) hs
radial mesh size increment
Definition vmec_main.f90:84
real(rprec), dimension(:), allocatable phi
toroidal magnetic flux
Definition vmec_main.f90:37
integer iequi
counter used to call -EQFOR- at end of run
integer first
"counter" monitoring sign of jacobian; resets R, Z, and Lambda when jacobian changes sign and decreas...
integer ijacob
counter for number of times jacobian changes sign
real(rprec) signgs
sign of Jacobian : must be =1 (right-handed) or =-1 (left-handed)
integer, parameter successful_term_flag
integer, parameter norm_term_flag
real(rprec), dimension(:), allocatable gc
stacked array of R, Z, Lambda Spectral force coefficients (see above for stack order)
Definition xstuff.f90:37
real(rprec), dimension(:), allocatable, target xc
stacked array of scaled R, Z, Lambda Fourier coefficients (see above for stack order)
Definition xstuff.f90:40
real(rprec), dimension(:), allocatable xsave
Definition xstuff.f90:45
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