VMEC 8.52
3D Equilibrium solver with nested flux surfaces.
Loading...
Searching...
No Matches
readin.f90
Go to the documentation of this file.
1
3
8 SUBROUTINE readin(input_file, ier_flag)
9 USE vmec_main
10 USE vmec_params
11 USE mgrid_mod, ONLY: nextcur, curlabel, read_mgrid, &
13
14 use dbgout
15 use theta_flip
16
17 IMPLICIT NONE
18
19 INTEGER, INTENT(inout) :: ier_flag
20 CHARACTER(LEN=*) :: input_file
21
22 INTEGER :: n, iunit, ier_flag_init, i, ni, m, nsmin, mj, isgn, ioff, joff
23 REAL(rprec), DIMENSION(:,:), POINTER :: rbcc, rbss, rbcs, rbsc, zbcs, zbsc, zbcc, zbss
24 REAL(rprec) :: rtest, ztest, tzc, trc, delta
25 REAL(rprec), ALLOCATABLE :: temp(:)
26 CHARACTER(LEN=100) :: line, line2
27
28 ier_flag_init = ier_flag
29 ier_flag = norm_term_flag
30
31 ! READ IN DATA FROM INDATA FILE
32 CALL read_indata(input_file, iunit, ier_flag)
33 IF (ier_flag .ne. norm_term_flag) then
34 RETURN
35 end if
36
37 ! Open output files here, print out heading to threed1 file
38 CALL heading(input_extension)
39
40 IF (lfreeb) THEN
41 ! READ IN AND STORE (FOR SEQUENTIAL RUNNING) MAGNETIC FIELD DATA FROM MGRID_FILE
42 CALL read_mgrid (mgrid_file, extcur, nzeta, nfp, .true., ier_flag)
43
44 ! check again for lfreeb
45 ! --> might have been reset to .false. by read_mgrid
46 ! if mgrid file was not found / could not be opened
47 IF (lfreeb) THEN
48 IF (ier_flag .ne. norm_term_flag) then
49 RETURN
50 end if
51
52 WRITE (nthreed,20) nr0b, nz0b, np0b, rminb, rmaxb, zminb, zmaxb, trim(mgrid_file)
5320 FORMAT(//,' VACUUM FIELD PARAMETERS:',/,1x,24('-'),/, &
54 ' nr-grid nz-grid np-grid rmin rmax zmin', &
55 ' zmax input-file',/,3i9,4f10.3,5x,a)
56 END IF
57 END IF
58
59 ! PARSE NS_ARRAY
61 nsmin = 1
62 DO WHILE (ns_array(multi_ns_grid) .ge. nsmin .and. & ! increasing or constant ns
63 multi_ns_grid .lt. 100) ! max 100 multi-grid iterations
64 nsmin = max(nsmin, ns_array(multi_ns_grid))
65 IF (nsmin .le. nsd) THEN ! max ns is given by nsd
67 ELSE
68 ! Optimizer, Boozer code overflows otherwise
69 ns_array(multi_ns_grid) = nsd
70 nsmin = nsd
71 print *,' NS_ARRAY ELEMENTS CANNOT EXCEED ',nsd
72 print *,' CHANGING NS_ARRAY(',multi_ns_grid,') to ', nsd
73 END IF
74 END DO
76
77 ! consistency check on requested number of flux surfaces
78 IF (all(ns_array.eq.0)) THEN
79 stop 'NS ARRAY MUST NOT BE ALL ZEROES' ! was ier_flag=ns_error_flag
80 END IF
81
82 ns_maxval = nsmin
83
84 ! Convert to Internal units
85 currv = mu0*curtor
86
87 ! WRITE OUT DATA TO THREED1 FILE
88
89 !SPH121912 - SCALING TO RENDER LAMSCALE=1
90 ! delta = twopi/phiedge !phiedge=>twopi
91 ! phiedge = phiedge*delta
92 ! bcrit = bcrit*delta
93 ! curtor = curtor*delta
94 ! extcur = extcur*delta
95 ! am = am*delta**2
96
97 WRITE (nthreed,100) &
98 ns_array(multi_ns_grid),ntheta1,nzeta,mpol,ntor, &
99 nfp,gamma,spres_ped,phiedge,curtor
100100 FORMAT(/,' COMPUTATION PARAMETERS: (u = theta, v = zeta)'/, &
101 1x,45('-'),/, &
102 ' ns nu nv mu mv',/,5i7//, &
103 ' CONFIGURATION PARAMETERS:',/,1x,39('-'),/, &
104 ' nfp gamma spres_ped phiedge(wb) curtor(A)',/, &
105 i7,1p,e11.3,2e15.3,e14.3/)
106
107 IF (nvacskip.le.0) then
108 ! default value for nvacskip: number of field periods ... ?
109 ! how does that work for e.g. LHD with nfp=10 ???
110 nvacskip = nfp
111 end if
112
113 WRITE (nthreed,110) ncurr,niter_array(multi_ns_grid),ns_array(1), &
114 nstep,nvacskip,ftol_array(multi_ns_grid),tcon0,lasym,lconm1, &
115 mfilter_fbdy,nfilter_fbdy
116110 FORMAT(' RUN CONTROL PARAMETERS:',/,1x,23('-'),/, &
117 ' ncurr niter nsin nstep nvacskip ftol tcon0',&
118 ' lasym lconm1',/, &
119 4i7,i10,1p,2e10.2,2l9,/, &
120 ' mfilter_fbdy nfilter_fbdy',/, &
121 2(6x,i7),/)
122
123 IF (nextcur .gt. 0) THEN
124 WRITE(nthreed, "(' EXTERNAL CURRENTS',/,1x,17('-'))")
125 ni = 0
126 IF (ALLOCATED(curlabel)) ni = maxval(len_trim(curlabel(1:nextcur)))
127 ni = max(ni+4, 14)
128 WRITE (line, '(a,i2.2,a)') "(5a",ni,")"
129 WRITE (line2, '(a,i2.2,a)') "(5(",ni-12,"x,1p,e12.4))"
130 DO i = 1,nextcur,5
131 ni = min(i+4, nextcur)
132 IF (ALLOCATED(curlabel)) then
133 WRITE (nthreed, line, iostat=mj) (trim(curlabel(n)),n=i,ni)
134 end if
135 WRITE (nthreed, line2,iostat=mj) (extcur(n), n=i,ni)
136 ENDDO
137 WRITE (nthreed, *)
138 ENDIF
139
140 IF (bloat .ne. one) THEN
141 WRITE (nthreed,'(" Profile Bloat Factor: ",1pe11.4)') bloat
142 phiedge = phiedge*bloat
143 ENDIF
144
145 IF (pres_scale .ne. one) THEN
146 WRITE (nthreed,'(" Pressure profile factor: ",1pe11.4,'// &
147 '" (multiplier for pressure)")') pres_scale
148 END IF
149
150 ! Print out am array
151 WRITE(nthreed,130)
152 WRITE(nthreed,131) trim(pmass_type)
153 WRITE(nthreed,132)
154130 FORMAT(' MASS PROFILE COEFFICIENTS am - newton/m**2 (EXPANSION IN NORMALIZED RADIUS):')
155131 FORMAT(' PMASS parameterization type is ''', a,'''')
156132 FORMAT(1x,35('-'))
157 WRITE(nthreed,135)(am(i-1),i=1, SIZE(am))
158135 FORMAT(1p,6e12.3)
159 SELECT CASE(trim(pmass_type))
160 CASE ('Akima_spline','cubic_spline')
161 WRITE(nthreed,"(' am_aux_s is' )")
162 WRITE(nthreed,135)(am_aux_s(i),i=1, SIZE(am_aux_s))
163 WRITE(nthreed,"(' am_aux_f is' )")
164 WRITE(nthreed,135)(am_aux_f(i),i=1, SIZE(am_aux_f))
165 END SELECT
166
167 IF (ncurr.eq.0) THEN
168 WRITE(nthreed,140)
169 WRITE(nthreed,143) trim(piota_type)
170 WRITE(nthreed,132)
171 ! Print out ai array
172 WRITE(nthreed,135)(ai(i-1),i=1, SIZE(ai))
173 SELECT CASE(trim(piota_type))
174 CASE ('Akima_spline','cubic_spline')
175 WRITE(nthreed,"(' ai_aux_s is' )")
176 WRITE(nthreed,135)(ai_aux_s(i),i=1, SIZE(ai_aux_s))
177 WRITE(nthreed,"(' ai_aux_f is' )")
178 WRITE(nthreed,135)(ai_aux_f(i),i=1, SIZE(ai_aux_f))
179 END SELECT
180 ELSE
181 ! Print out ac array
182 WRITE(nthreed,145)
183 WRITE(nthreed,146) trim(pcurr_type)
184 WRITE(nthreed,147)
185 WRITE(nthreed,135)(ac(i-1),i=1, SIZE(ac))
186 SELECT CASE(trim(pcurr_type))
187 CASE ('Akima_spline_Ip','Akima_spline_I', &
188 'cubic_spline_Ip','cubic_spline_I')
189 WRITE(nthreed,"(' ac_aux_s is' )")
190 WRITE(nthreed,135)(ac_aux_s(i),i=1, SIZE(ac_aux_s))
191 WRITE(nthreed,"(' ac_aux_f is' )")
192 WRITE(nthreed,135)(ac_aux_f(i),i=1, SIZE(ac_aux_f))
193 END SELECT
194 ENDIF
195
196140 FORMAT(/' IOTA PROFILE COEFFICIENTS ai (EXPANSION IN NORMALIZED RADIUS):')
197143 FORMAT(' PIOTA parameterization type is ''', a,'''')
198145 FORMAT(/' TOROIDAL CURRENT DENSITY (*V'') COEFFICIENTS ac (EXPANSION IN NORMALIZED RADIUS):')
199146 FORMAT(' PCURR parameterization type is ''', a,'''')
200147 FORMAT(1x,38('-'))
201
202 WRITE(nthreed,150)
203 WRITE(nthreed,135)(aphi(i),i=1, SIZE(aphi))
204150 FORMAT(/' NORMALIZED TOROIDAL FLUX COEFFICIENTS aphi (EXPANSION IN S):',/,1x,35('-'))
205
206! Fourier Boundary Coefficients (header written always, but data only if in free-boundary mode)
207 WRITE(nthreed,180)
208180 FORMAT(/,' R-Z FOURIER BOUNDARY COEFFICIENTS AND', &
209 ' MAGNETIC AXIS INITIAL GUESS',/, &
210 ' R = RBC*COS(m*u - n*v) + RBS*SIN(m*u - n*v),', &
211 ' Z = ZBC*COS(m*u - n*v) + ZBS*SIN(m*u-n*v)'/1x,86('-'), &
212 /,' nb mb rbc rbs zbc zbs ',&
213 ' raxis(c) raxis(s) zaxis(c) zaxis(s)')
214
215 IF (lasym) THEN
216 ! CONVERT TO REPRESENTATION WITH RBS(m=1) = ZBC(m=1)
217
218 delta = atan2(rbs(0,1) - zbc(0,1), rbc(0,1) + zbs(0,1))
219 IF (delta .ne. zero) THEN
220 DO m = 0,mpol1
221 DO n = -ntor,ntor
222
223 trc = rbc(n,m)*cos(m*delta) + rbs(n,m)*sin(m*delta)
224 rbs(n,m) = rbs(n,m)*cos(m*delta) - rbc(n,m)*sin(m*delta)
225 rbc(n,m) = trc
226
227 tzc = zbc(n,m)*cos(m*delta) + zbs(n,m)*sin(m*delta)
228 zbs(n,m) = zbs(n,m)*cos(m*delta) - zbc(n,m)*sin(m*delta)
229 zbc(n,m) = tzc
230 ENDDO
231 ENDDO
232 ENDIF
233 ENDIF
234
235 ! ALLOCATE MEMORY FOR NU, NV, MPOL, NTOR SIZED ARRAYS
236 CALL allocate_nunv
237
238 ! CONVERT TO INTERNAL REPRESENTATION OF MODES
239 !
240 ! R = RBCC*COS(M*U)*COS(N*V) + RBSS*SIN(M*U)*SIN(N*V)
241 ! + RBCS*COS(M*U)*SIN(N*V) + RBSC*SIN(M*U)*COS(N*V)
242 ! Z = ZBCS*COS(M*U)*SIN(N*V) + ZBSC*SIN(M*U)*COS(N*V)
243 ! + ZBCC*COS(M*U)*COS(N*V) + ZBSS*SIN(M*U)*SIN(N*V)
244 !
245 ! POINTER ASSIGNMENTS
246 ! ARRAY STACKING ORDER DETERMINED HERE
247 rbcc => rmn_bdy(:,:,rcc)
248 zbsc => zmn_bdy(:,:,zsc)
249 IF (lthreed) THEN
250 rbss => rmn_bdy(:,:,rss)
251 zbcs => zmn_bdy(:,:,zcs)
252 END IF
253
254 IF (lasym) THEN
255 rbsc => rmn_bdy(:,:,rsc)
256 zbcc => zmn_bdy(:,:,zcc)
257 IF (lthreed) THEN
258 rbcs => rmn_bdy(:,:,rcs)
259 zbss => zmn_bdy(:,:,zss)
260 END IF
261 ENDIF
262
263 rmn_bdy = 0.0_dp
264 zmn_bdy = 0.0_dp
265
266 ! NOTE: INDICES START AT 1, NOT 0, FOR POINTERS,
267 ! EVEN THOUGH THEY START AT ZERO FOR RMN_BDY
268 ioff = lbound(rbcc,1) ! index offset in dimension 1, which is n
269 joff = lbound(rbcc,2) ! index offset in dimension 2, which is m
270
271 ! go over all mode number combinations that could be specified in input file for given Fourier resolution
272 DO m=0,mpol1
273
274 ! target index in rbcc, ... along m --> m
275 mj = m+joff
276
277 ! in free-boundary, skip initial guess with Fourier mode number m larger than mfilter_fbdy
278 IF (lfreeb .and. (mfilter_fbdy.gt.1 .and. m.gt.mfilter_fbdy)) cycle
279
280 DO n=-ntor,ntor
281
282 ! in free-boundary, skip initial guess with Fourier mode number n larger than nfilter_fbdy
283 IF (lfreeb .and. (nfilter_fbdy.gt.0 .and. abs(n).gt.nfilter_fbdy)) cycle
284
285 ! target index in rbcc, ... along n --> |n|
286 ni = abs(n) + ioff
287
288 ! these two loops over m and n simply go over the full Rbc, ... arrays
289
290 ! rbc, zbs etc are defined for the full 2d Fourier kernel with possibly negative n
291 ! this sign is used to convert it to positive-only Fourier numbers
292 IF (n .eq. 0) THEN
293 isgn = 0
294 ELSE IF (n .gt. 0) THEN
295 isgn = 1
296 ELSE
297 isgn = -1
298 END IF
299
300 rbcc(ni,mj) = rbcc(ni,mj) + rbc(n,m)
301 IF (m .gt. 0) zbsc(ni,mj) = zbsc(ni,mj) + zbs(n,m)
302
303 IF (lthreed) THEN
304 IF (m .gt. 0) then
305 rbss(ni,mj) = rbss(ni,mj) + isgn*rbc(n,m)
306 end if
307 zbcs(ni,mj) = zbcs(ni,mj) - isgn*zbs(n,m)
308 END IF
309
310 IF (lasym) THEN
311 IF (m .gt. 0) rbsc(ni,mj) = rbsc(ni,mj) + rbs(n,m)
312 zbcc(ni,mj) = zbcc(ni,mj) + zbc(n,m)
313
314 IF (lthreed) THEN
315 rbcs(ni,mj) = rbcs(ni,mj) - isgn*rbs(n,m)
316 IF (m .gt. 0) zbss(ni,mj) = zbss(ni,mj) + isgn*zbc(n,m)
317 END IF
318 END IF
319
320
321
322 ! Fourier coefficient assigments are done here,
323 ! only diagnostics below in this loop
324
325
326
327 ! TODO: why exit only here or is this leftover from something else?
328 IF (ier_flag_init .ne. norm_term_flag) then
329 print *, "skip printing of initial guess for geometry, since ier_flag_init = ",ier_flag_init
330 cycle
331 end if
332
333 ! compute sum of quantities to be printed below
334 ! --> check if any of rbc, ..., zbs is != 0 for current m,n values
335 trc = abs(rbc(n,m)) + abs(rbs(n,m)) + abs(zbc(n,m)) + abs(zbs(n,m))
336
337 IF (m .eq. 0) THEN
338 ! m=0: also axis coefficients can be printed
339
340 ! above n-loop is re-used for printing:
341 ! negative n come before positive n, so any |n| will have been completely assigned also its negative counterparts at this point
342 IF (n .lt. 0) cycle
343
344 IF (trc.eq.zero .and. abs(raxis_cc(n)).eq.zero .and. abs(zaxis_cs(n)).eq.zero) then
345 ! all values to be printed are zero, so skip printing them
346 cycle
347 end if
348
349 WRITE (nthreed,195) n, m, rbc(n,m), rbs(n,m), zbc(n,m), zbs(n,m), &
350 raxis_cc(n), raxis_cs(n), zaxis_cc(n), zaxis_cs(n)
351 ELSE
352 ! m > 0, so only print boundary coeffs (axis does not have coeffs with m>1
353
354 IF (trc .eq. zero) then
355 ! all values to be printed are zero, so skip printing them
356 cycle
357 end if
358 WRITE (nthreed,195) n, m, rbc(n,m), rbs(n,m), zbc(n,m), zbs(n,m)
359 END IF
360
361 END DO
362 END DO
363195 FORMAT(i5,i4,1p,8e12.4)
364
365 ! CHECK SIGN OF JACOBIAN (SHOULD BE SAME AS SIGNGS)
366 m = 1
367 mj = m+joff
368 rtest = sum(rbcc(1:ntor1,mj))
369 ztest = sum(zbsc(1:ntor1,mj))
370 lflip=(rtest*ztest .lt. zero)
371 signgs = -1.0_dp
372 IF (lflip) then
373 ! (rtest*ztest < 0) means that rtest and ztest have different signs
375 end if
376
377 rtest = sum(rbcc(1:ntor1,mj))
378 ztest = sum(zbsc(1:ntor1,mj))
379 if (rtest*ztest .lt. zero) then
380 stop "flip_theta did not fix it!"
381 end if
382
383 ! CONVERT TO INTERNAL FORM FOR (CONSTRAINED) m=1 MODES
384 ! INTERNALLY, FOR m=1: XC(rss) = .5(RSS+ZCS), XC(zcs) = .5(RSS-ZCS)
385 ! WITH XC(zcs) -> 0 FOR POLAR CONSTRAINT
386 ! FOR ASYMMETRIC CASE, XC(rsc) = .5(RSC+ZCC), XC(zcc) = .5(RSC-ZCC)
387 ! WITH XC(zss) -> 0 FOR POLAR CONSTRAINT ! TODO: (jons) maybe XC(zcc)-->0 is meant here?
388 ! (see convert_sym, convert_asym in totzsp.f90 file)
389 IF (lconm1 .AND. (lthreed.OR.lasym)) THEN
390 ! the length could be also defined as ntor+1
391 ! (which would be more explicit, but who needs this ...?)
392 ALLOCATE (temp(SIZE(rbcc,1)))
393 IF (lthreed) THEN
394 mj = 1+joff ! index of m=1 modes for all n
395 temp = rbss(:,mj)
396 rbss(:,mj) = cp5*(temp(:) + zbcs(:,mj))
397 zbcs(:,mj) = cp5*(temp(:) - zbcs(:,mj))
398 END IF
399 IF (lasym) THEN
400 mj = 1+joff ! index of m=1 modes for all n
401 temp = rbsc(:,mj)
402 rbsc(:,mj) = cp5*(temp(:) + zbcc(:,mj))
403 zbcc(:,mj) = cp5*(temp(:) - zbcc(:,mj))
404 END IF
405 DEALLOCATE (temp)
406 END IF
407
408 if (open_dbg_context("readin_boundary")) then
409
410 if (.not. lasym) then
411 if (.not. lthreed) then
412 call add_real_2d("rbcc", ntor1, mpol, rbcc)
413 call add_real_2d("zbsc", ntor1, mpol, zbsc)
414 call add_null("rbss")
415 call add_null("zbcs")
416 call add_null("rbsc")
417 call add_null("zbcc")
418 call add_null("rbcs")
419 call add_null("zbss")
420 else
421 call add_real_2d("rbcc", ntor1, mpol, rbcc)
422 call add_real_2d("zbsc", ntor1, mpol, zbsc)
423 call add_real_2d("rbss", ntor1, mpol, rbss)
424 call add_real_2d("zbcs", ntor1, mpol, zbcs)
425 call add_null("rbsc")
426 call add_null("zbcc")
427 call add_null("rbcs")
428 call add_null("zbss")
429 end if
430 else
431 if (.not. lthreed) then
432 call add_real_2d("rbcc", ntor1, mpol, rbcc)
433 call add_real_2d("zbsc", ntor1, mpol, zbsc)
434 call add_null("rbss")
435 call add_null("zbcs")
436 call add_real_2d("rbsc", ntor1, mpol, rbsc)
437 call add_real_2d("zbcc", ntor1, mpol, zbcc)
438 call add_null("rbcs")
439 call add_null("zbss")
440 else
441 call add_real_2d("rbcc", ntor1, mpol, rbcc)
442 call add_real_2d("zbsc", ntor1, mpol, zbsc)
443 call add_real_2d("rbss", ntor1, mpol, rbss)
444 call add_real_2d("zbcs", ntor1, mpol, zbcs)
445 call add_real_2d("rbsc", ntor1, mpol, rbsc)
446 call add_real_2d("zbcc", ntor1, mpol, zbcc)
447 call add_real_2d("rbcs", ntor1, mpol, rbcs)
448 call add_real_2d("zbss", ntor1, mpol, zbss)
449 end if
450 end if
451
452 call close_dbg_out()
453 end if
454
455END SUBROUTINE readin
subroutine allocate_nunv
allocate arrays depending on the number of Fourier coefficients nunv
subroutine heading(extension)
Open output files and print banner message at the top.
Definition heading.f90:8
logical function open_dbg_context(context_name, repetition, id)
check if any output is desired for the current iteration check if the given context should be openend...
Definition dbgout.f90:17
character(len=30), dimension(:), allocatable curlabel
Definition mgrid_mod.f:95
subroutine read_mgrid(mgrid_file, extcur, nv, nfp, lscreen, ier_flag)
Definition mgrid_mod.f:114
integer nr0b
Definition mgrid_mod.f:77
real(rprec) zminb
Definition mgrid_mod.f:84
integer nextcur
Definition mgrid_mod.f:78
real(rprec) zmaxb
Definition mgrid_mod.f:84
real(rprec) rminb
Definition mgrid_mod.f:84
real(rprec) rmaxb
Definition mgrid_mod.f:84
integer nz0b
Definition mgrid_mod.f:77
integer np0b
Definition mgrid_mod.f:77
subroutine flip_theta(rmn, zmn, lmn)
Flip the definition of the poloidal angle in the user-provided initial guess for the LCFS geometry.
real(rprec), dimension(:,:,:), allocatable, target zmn_bdy
logical lconm1
integer multi_ns_grid
real(rprec), dimension(:,:,:), allocatable, target rmn_bdy
real(rprec) currv
toroidal current (?)
Definition vmec_main.f90:85
logical lthreed
logical lflip
from init_geometry
real(rprec) signgs
sign of Jacobian : must be =1 (right-handed) or =-1 (left-handed)
integer, parameter norm_term_flag
subroutine read_indata(in_file, iunit, ier_flag)
Read the INDATA namelist from a given input file.
subroutine readin(input_file, ier_flag)
Read the input file.
Definition readin.f90:9