V3FIT
restart_mod.f90
Go to the documentation of this file.
1 !*******************************************************************************
4 !
5 ! Note separating the Doxygen comment block here so detailed decription is
6 ! found in the Module not the file.
7 !
9 !*******************************************************************************
10  MODULE restart_mod
11  USE ezcdf
12  USE stel_kinds
13  USE utilities, ONLY: gradientfull, to_full_mesh
14  USE metrics, ONLY: tolowerh
15  USE descriptor_mod, ONLY: iam
16  USE shared_data, ONLY: lasym, lrecon, unit_out
17  USE v3_utilities, ONLY: assert_eq
18  USE stel_constants, ONLY: mu0
19 
20  IMPLICIT NONE
21 
22 !*******************************************************************************
23 ! restart module parameters
24 !*******************************************************************************
26  INTEGER, PARAMETER :: restart_version = 0
27 
28 ! Flag will be placed in the last bits. This way the version number and flags
29 ! can fit in the same memory.
31  INTEGER, PARAMETER :: restart_lasym = 31
33  INTEGER, PARAMETER :: restart_lrecon = 30
34 
36  CHARACTER (LEN=*), DIMENSION(1), PARAMETER :: &
37  & radial_dim = (/ 'radius' /)
39  CHARACTER (LEN=*), DIMENSION(3), PARAMETER :: &
40  & restart_dims = (/ 'm-mode', 'n-mode', radial_dim(1) /)
41 
43  CHARACTER (len=*), PARAMETER :: vn_nsin = 'nrad'
45  CHARACTER (len=*), PARAMETER :: vn_mpolin = 'mpol'
47  CHARACTER (len=*), PARAMETER :: vn_ntorin = 'ntor'
49  CHARACTER (len=*), PARAMETER :: vn_nfpin = 'nfp'
51  CHARACTER (len=*), PARAMETER :: vn_wout = 'wout_file'
53  CHARACTER (len=*), PARAMETER :: vn_flags = 'state_flags'
54 
56  CHARACTER (len=*), PARAMETER :: vn_jbsupss = 'JBsupssh(m,n,r)'
58  CHARACTER (len=*), PARAMETER :: vn_jbsupuc = 'JBsupuch(m,n,r)'
60  CHARACTER (len=*), PARAMETER :: vn_jbsupvc = 'JBsupvch(m,n,r)'
62  CHARACTER (len=*), PARAMETER :: vn_jbsupsc = 'JBsupsch(m,n,r)'
64  CHARACTER (len=*), PARAMETER :: vn_jbsupus = 'JBsupush(m,n,r)'
66  CHARACTER (len=*), PARAMETER :: vn_jbsupvs = 'JBsupvsh(m,n,r)'
68  CHARACTER (len=*), PARAMETER :: vn_jpresc = 'jpresch(m,n,r)'
70  CHARACTER (len=*), PARAMETER :: vn_jpress = 'jpressh(m,n,r)'
71 
73  CHARACTER (len=*), PARAMETER :: vn_rmnc = 'rmnc(m,n,r)'
75  CHARACTER (len=*), PARAMETER :: vn_rmns = 'rmns(m,n,r)'
77  CHARACTER (len=*), PARAMETER :: vn_zmnc = 'zmnc(m,n,r)'
79  CHARACTER (len=*), PARAMETER :: vn_zmns = 'zmns(m,n,r)'
80 
82  CHARACTER (len=*), PARAMETER :: vn_chipf = 'chipf(r)'
84  CHARACTER (len=*), PARAMETER :: vn_phipf = 'phipf(r)'
85 
87  CHARACTER (len=*), PARAMETER :: vn_bsupsmns = 'bsupsmnsh(m,n,r)'
89  CHARACTER (len=*), PARAMETER :: vn_bsupsmnc = 'bsupsmnch(m,n,r)'
91  CHARACTER (len=*), PARAMETER :: vn_bsupumns = 'bsupumnsh(m,n,r)'
93  CHARACTER (len=*), PARAMETER :: vn_bsupumnc = 'bsupumnch(m,n,r)'
95  CHARACTER (len=*), PARAMETER :: vn_bsupvmns = 'bsupvmnsh(m,n,r)'
97  CHARACTER (len=*), PARAMETER :: vn_bsupvmnc = 'bsupvmnch(m,n,r)'
99  CHARACTER (len=*), PARAMETER :: vn_bsubsmns = 'bsubsmnsh(m,n,r)'
101  CHARACTER (len=*), PARAMETER :: vn_bsubsmnc = 'bsubsmnch(m,n,r)'
103  CHARACTER (len=*), PARAMETER :: vn_bsubumns = 'bsubumnsh(m,n,r)'
105  CHARACTER (len=*), PARAMETER :: vn_bsubumnc = 'bsubumnch(m,n,r)'
107  CHARACTER (len=*), PARAMETER :: vn_bsubvmns = 'bsubvmnsh(m,n,r)'
109  CHARACTER (len=*), PARAMETER :: vn_bsubvmnc = 'bsubvmnch(m,n,r)'
111  CHARACTER (len=*), PARAMETER :: vn_pmns = 'pmnsh(m,n,r)'
113  CHARACTER (len=*), PARAMETER :: vn_pmnc = 'pmnch(m,n,r)'
115  CHARACTER (len=*), PARAMETER :: vn_jksupsmns = 'jksupsmnsf(m,n,r)'
117  CHARACTER (len=*), PARAMETER :: vn_jksupsmnc = 'jksupsmncf(m,n,r)'
119  CHARACTER (len=*), PARAMETER :: vn_jksupumns = 'jksupumnsf(m,n,r)'
121  CHARACTER (len=*), PARAMETER :: vn_jksupumnc = 'jksupumncf(m,n,r)'
123  CHARACTER (len=*), PARAMETER :: vn_jksupvmns = 'jksupvmnsf(m,n,r)'
125  CHARACTER (len=*), PARAMETER :: vn_jksupvmnc = 'jksupvmncf(m,n,r)'
126 
128  CHARACTER (len=*), PARAMETER :: vn_p_factor = 'p_factor'
130  CHARACTER (len=*), PARAMETER :: vn_b_factor = 'b_factor'
132  CHARACTER (len=*), PARAMETER :: vn_wb = 'wb'
134  CHARACTER (len=*), PARAMETER :: vn_wp = 'wp'
136  CHARACTER (len=*), PARAMETER :: vn_rmajor = 'rmajor'
138  CHARACTER (len=*), PARAMETER :: vn_curtor = 'curtor'
139 
141  CHARACTER (len=*), PARAMETER :: vn_p_max = 'p_max'
143  CHARACTER (len=*), PARAMETER :: vn_p_min = 'p_min'
144 
145  CONTAINS
146 
147 !*******************************************************************************
148 ! UTILITY SUBROUTINES
149 !*******************************************************************************
150 !-------------------------------------------------------------------------------
162 !-------------------------------------------------------------------------------
163  FUNCTION restart_read(restart_ext, wout_file, mpolin, ntorin, nfpin, &
164  nsin)
165  USE quantities, ONLY: jbsupsmnsh, jbsupsmnch, &
166  jbsupumnsh, jbsupumnch, &
167  jbsupvmnsh, jbsupvmnch, &
168  jpmnsh, jpmnch, &
169  b_factor, p_factor, alloc_quantities
170  USE island_params, ONLY: chipf => chipf_i, phipf => phipf_i, &
171  & wb => wb_i, wp => wp_i, nfp_i, gamma, &
172  & gnorm => gnorm_i, rmajor => rmajor_i, &
173  & fourier_context, nu_i, nv_i
174  USE vmec_info, ONLY: rmnc => rmnc_i, zmns => zmns_i, &
175  & rmns => rmns_i, zmnc => zmnc_i, &
176  & vmec_info_construct_island
177  USE metrics, ONLY: loadgrid
178  USE fourier, ONLY: m0, m1, fourier_class
179 
180  IMPLICIT NONE
181 
182 ! Declare Arguments
183  INTEGER :: restart_read
184  CHARACTER (len=*), INTENT(in) :: restart_ext
185  CHARACTER (len=*), INTENT(inout) :: wout_file
186  INTEGER, INTENT(in) :: mpolin
187  INTEGER, INTENT(in) :: ntorin
188  INTEGER, INTENT(in) :: nsin
189  INTEGER, INTENT(in) :: nfpin
190 
191 ! local variables
192  INTEGER :: flags
193  INTEGER :: ncid
194  INTEGER :: ns
195  INTEGER :: mpol
196  INTEGER :: ntor
197  INTEGER :: nfp
198  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: tempmn_r
199  REAL (dp), DIMENSION(:), ALLOCATABLE :: temp_r
200  INTEGER :: nmin
201  INTEGER :: status
202  CHARACTER (LEN=256) :: filename
203  INTEGER :: s
204 
205 ! Start of executable code
206  restart_read = 1
207 
208  CALL alloc_quantities
209 
210  filename = 'siesta_' // trim(restart_ext) // '.nc'
211  CALL cdf_open(ncid, trim(filename), 'r', status)
212  CALL assert_eq(0, status, 'Failed to read restart file: ' // &
213  & trim(filename))
214 
215  CALL cdf_read(ncid, vn_flags, flags)
216 
217 ! The namelist inut file can change the size of the radial grid and the
218 ! poloidal and toroidal modes.
219  CALL cdf_read(ncid, vn_nsin, ns)
220  CALL cdf_read(ncid, vn_mpolin, mpol)
221  CALL cdf_read(ncid, vn_ntorin, ntor)
222  CALL cdf_read(ncid, vn_nfpin, nfp)
223 
224  IF (nfpin .lt. 1) THEN
225  nfp_i = nfp
226  ELSE
227  nfp_i = nfpin
228  END IF
229 
230  nmin = min(ntor, ntorin)
231 
232  CALL cdf_read(ncid, vn_wout, wout_file)
233 
234 ! The namelist input file may turn the asymmetric terms on and off.
235  ALLOCATE(tempmn_r(0:mpol,-ntor:ntor,ns))
236  ALLOCATE(temp_r(ns))
237  CALL vmec_info_construct_island(mpolin, ntorin, nsin, lasym)
238 
239  ALLOCATE(chipf(nsin))
240  CALL cdf_read(ncid, vn_chipf, temp_r)
241  CALL interpit_1d(temp_r, chipf, ns, nsin, .false., 1)
242 
243  ALLOCATE(phipf(nsin))
244  CALL cdf_read(ncid, vn_phipf, temp_r)
245  CALL interpit_1d(temp_r, phipf, ns, nsin, .false., 1)
246 
247  CALL cdf_read(ncid, vn_jbsupss, tempmn_r)
248  jbsupsmnsh(:,:,1) = 0
249  CALL interpit(tempmn_r, jbsupsmnsh, ns, nsin, mpol, mpolin, &
250  & ntor, ntorin, nfp, nfpin, .true.)
251  jbsupsmnsh(m1,-nmin:nmin,1) = jbsupsmnsh(m1,-nmin:nmin,2)
252 
253  CALL cdf_read(ncid, vn_jbsupuc, tempmn_r)
254  jbsupumnch(:,:,1) = 0
255  CALL interpit(tempmn_r, jbsupumnch, ns, nsin, mpol, mpolin, &
256  & ntor, ntorin, nfp, nfpin, .true.)
257  jbsupumnch(m1,-nmin:nmin,1) = jbsupumnch(m1,-nmin:nmin,2)
258 
259  CALL cdf_read(ncid, vn_jbsupvc, tempmn_r)
260  jbsupvmnch(:,:,1) = 0
261  CALL interpit(tempmn_r, jbsupvmnch, ns, nsin, mpol, mpolin, &
262  & ntor, ntorin, nfp, nfpin, .true.)
263  jbsupvmnch(m0,-nmin:nmin,1) = jbsupvmnch(m0,-nmin:nmin,2)
264 
265  CALL cdf_read(ncid, vn_jpresc, tempmn_r)
266  jpmnch(:,:,1) = 0
267  CALL interpit(tempmn_r, jpmnch, ns, nsin, mpol, mpolin, &
268  & ntor, ntorin, nfp, nfpin, .true.)
269  jpmnch(m0,-nmin:nmin,1) = jpmnch(m0,-nmin:nmin,2)
270 
271  CALL cdf_read(ncid, vn_rmnc, tempmn_r)
272  CALL interpit(tempmn_r, rmnc, ns, nsin, mpol, mpolin, &
273  ntor, ntorin, nfp, nfpin, .false.)
274 
275  CALL cdf_read(ncid, vn_zmns, tempmn_r)
276  CALL interpit(tempmn_r, zmns, ns, nsin, mpol, mpolin, &
277  ntor, ntorin, nfp, nfpin, .false.)
278 
279  IF (btest(flags, restart_lasym) .and. lasym) THEN
280  jbsupsmnch(:,:,1) = 0
281  CALL cdf_read(ncid, vn_jbsupsc, tempmn_r)
282  CALL interpit(tempmn_r, jbsupsmnch, ns, nsin, mpol, mpolin, &
283  & ntor, ntorin, nfp, nfpin, .true.)
284  jbsupsmnch(m1,-nmin:nmin,1) = jbsupsmnch(m1,-nmin:nmin,2)
285 
286  CALL cdf_read(ncid, vn_jbsupus, tempmn_r)
287  CALL interpit(tempmn_r, jbsupumnsh, ns, nsin, mpol, mpolin, &
288  & ntor, ntorin, nfp, nfpin, .true.)
289  jbsupumnsh(m1,-nmin:nmin,1) = jbsupumnsh(m1,-nmin:nmin,2)
290 
291  CALL cdf_read(ncid, vn_jbsupvs, tempmn_r)
292  CALL interpit(tempmn_r, jbsupvmnsh, ns, nsin, mpol, mpolin, &
293  & ntor, ntorin, nfp, nfpin, .true.)
294  jbsupvmnsh(m0,-nmin:nmin,1) = jbsupvmnsh(m0,-nmin:nmin,2)
295 
296  CALL cdf_read(ncid, vn_jpress, tempmn_r)
297  CALL interpit(tempmn_r, jpmnsh, ns, nsin, mpol, mpolin, &
298  & ntor, ntorin, nfp, nfpin, .true.)
299  jpmnsh(m0,-nmin:nmin,1) = jpmnsh(m0,-nmin:nmin,2)
300 
301  CALL cdf_read(ncid, vn_rmns, tempmn_r)
302  CALL interpit(tempmn_r, rmns, ns, nsin, mpol, mpolin, &
303  ntor, ntorin, nfp, nfpin, .false.)
304 
305  CALL cdf_read(ncid, vn_zmns, tempmn_r)
306  CALL interpit(tempmn_r, zmnc, ns, nsin, mpol, mpolin, &
307  ntor, ntorin, nfp, nfpin, .false.)
308 
309  ELSE IF (lasym) THEN
310  jbsupsmnch = 0
311  jbsupumnsh = 0
312  jbsupvmnsh = 0
313  jpmnsh = 0
314  rmns = 0
315  rmnc = 0
316  END IF
317 
318 ! Read normalization factors.
319  CALL cdf_read(ncid, vn_wb, wb)
320  CALL cdf_read(ncid, vn_wp, wp)
321  CALL cdf_read(ncid, vn_rmajor, rmajor)
322 
323  CALL cdf_close(ncid)
324 
325  fourier_context => fourier_class(mpolin, ntorin, nu_i, nv_i, nfpin, lasym)
326 
327  DO s = 1, nsin
328  rmnc(:,:,s) = rmnc(:,:,s)/fourier_context%orthonorm
329  zmns(:,:,s) = zmns(:,:,s)/fourier_context%orthonorm
330  IF (lasym) THEN
331  rmns(:,:,s) = rmns(:,:,s)/fourier_context%orthonorm
332  zmnc(:,:,s) = zmnc(:,:,s)/fourier_context%orthonorm
333  END IF
334  END DO
335 
336  CALL loadgrid(status)
337 
338  DEALLOCATE(tempmn_r)
339  DEALLOCATE(temp_r)
340 
341 ! Init quantities.
342  p_factor = 1.0/abs(wb)
343  b_factor = sqrt(p_factor)
344  gnorm = abs(wb)/(wb + wp/(gamma - 1.0))
345 
346  WRITE (unit_out, 1000) mpol, ntor, ns
347 
348  restart_read = 2
349 
350 1000 FORMAT(/,' RESTARTED FROM RUN PARAMETERS M: ',i3,' N: ',i3,' NS: ', i3,/)
351 
352  END FUNCTION
353 
354 !-------------------------------------------------------------------------------
361 !-------------------------------------------------------------------------------
362  SUBROUTINE restart_write(restart_ext, wout_file)
363  USE quantities, ONLY: jbsupsmnsh, jbsupsmnch, &
364  jbsupumnsh, jbsupumnch, &
365  jbsupvmnsh, jbsupvmnch, &
366  jpmnsh, jpmnch, &
367  b_factor, p_factor, jacobh
368  USE fourier, ONLY: f_cos, f_sin, f_sum, f_none
369  USE island_params, ONLY: nfp => nfp_i, chipf => chipf_i, &
370  phipf => phipf_i, wb => wb_i, wp => wp_i, &
371  rmajor => rmajor_i, fourier_context, &
372  mpol=>mpol_i, ntor=>ntor_i, ns=>ns_i, &
373  ntheta=>nu_i, nzeta=>nv_i
374  USE shared_data, ONLY: r1_i, z1_i
375  USE utilities, ONLY: curl_htof
376  USE stel_constants, ONLY: one
377 
378  IMPLICIT NONE
379 
380 ! Declare Arguments
381  CHARACTER (len=*), INTENT(in) :: restart_ext
382  CHARACTER (len=*), INTENT(in) :: wout_file
383 
384 ! Local Variables
385  INTEGER :: flags
386  INTEGER :: ncid
387  INTEGER :: nblock
388  INTEGER :: s
389  INTEGER :: m
390  INTEGER :: n
391  INTEGER :: status
392 
393  REAL (dp) :: r0
394  REAL (dp) :: tempscalar
395  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsupsijh
396  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsupuijh
397  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsupvijh
398  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsubsijh
399  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsubuijh
400  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsubvijh
401  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsubsmnh
402  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsubumnh
403  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: bsubvmnh
404  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: jksupsmnf
405  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: jksupumnf
406  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: jksupvmnf
407  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: tempmn_w
408  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: pijh
409 
410 ! Start of executable code.
411  CALL cdf_open(ncid, 'siesta_' // trim(restart_ext) // '.nc', 'w', status)
412  CALL assert_eq(0, status, 'Failed to write restart file.')
413 
414  flags = restart_version
415  IF (lasym) THEN
416  flags = ibset(flags, restart_lasym)
417  END IF
418  IF (lrecon) THEN
419  flags = ibset(flags, restart_lrecon)
420  END IF
421 
422  ALLOCATE(tempmn_w(0:mpol,-ntor:ntor,ns))
423 
424  CALL cdf_define(ncid, vn_flags, flags)
425  CALL cdf_define(ncid, vn_nsin, ns)
426  CALL cdf_define(ncid, vn_mpolin, mpol)
427  CALL cdf_define(ncid, vn_ntorin, ntor)
428  CALL cdf_define(ncid, vn_nfpin, nfp)
429  CALL cdf_define(ncid, vn_wout, wout_file)
430 
431  CALL cdf_define(ncid, vn_jbsupss, jbsupsmnsh, dimname=restart_dims)
432  CALL cdf_define(ncid, vn_jbsupuc, jbsupumnch, dimname=restart_dims)
433  CALL cdf_define(ncid, vn_jbsupvc, jbsupvmnch, dimname=restart_dims)
434  CALL cdf_define(ncid, vn_jpresc, jpmnch, dimname=restart_dims)
435 
436  CALL cdf_define(ncid, vn_rmnc, tempmn_w, dimname=restart_dims)
437  CALL cdf_define(ncid, vn_zmns, tempmn_w, dimname=restart_dims)
438 
439  CALL cdf_define(ncid, vn_chipf, chipf, dimname=radial_dim)
440  CALL cdf_define(ncid, vn_phipf, phipf, dimname=radial_dim)
441 
442  CALL cdf_define(ncid, vn_wb, wb)
443  CALL cdf_define(ncid, vn_wp, wp)
444 
445  IF (lasym) THEN
446  CALL cdf_define(ncid, vn_jbsupsc, jbsupsmnch, dimname=restart_dims)
447  CALL cdf_define(ncid, vn_jbsupus, jbsupumnsh, dimname=restart_dims)
448  CALL cdf_define(ncid, vn_jbsupvs, jbsupvmnsh, dimname=restart_dims)
449  CALL cdf_define(ncid, vn_jpress, jpmnsh, dimname=restart_dims)
450 
451  CALL cdf_define(ncid, vn_rmns, tempmn_w, dimname=restart_dims)
452  CALL cdf_define(ncid, vn_zmnc, tempmn_w, dimname=restart_dims)
453  END IF
454 
455  recon0: IF (lrecon) THEN
456 ! Using variable with jacobian only for the dimension sizes and types. Before
457 ! writting, the jacobian normalization will be removed.
458  CALL cdf_define(ncid, vn_bsupsmns, jbsupsmnsh, dimname=restart_dims)
459  CALL cdf_define(ncid, vn_bsupumnc, jbsupumnch, dimname=restart_dims)
460  CALL cdf_define(ncid, vn_bsupvmnc, jbsupvmnch, dimname=restart_dims)
461  CALL cdf_define(ncid, vn_bsubsmns, jbsupsmnsh, dimname=restart_dims)
462  CALL cdf_define(ncid, vn_bsubumnc, jbsupumnch, dimname=restart_dims)
463  CALL cdf_define(ncid, vn_bsubvmnc, jbsupvmnch, dimname=restart_dims)
464  CALL cdf_define(ncid, vn_jksupsmns, jbsupsmnsh, dimname=restart_dims)
465  CALL cdf_define(ncid, vn_jksupumnc, jbsupumnch, dimname=restart_dims)
466  CALL cdf_define(ncid, vn_jksupvmnc, jbsupvmnch, dimname=restart_dims)
467  CALL cdf_define(ncid, vn_pmnc, jpmnch, dimname=restart_dims)
468  CALL cdf_define(ncid, vn_b_factor, b_factor)
469  CALL cdf_define(ncid, vn_p_factor, p_factor)
470  CALL cdf_define(ncid, vn_rmajor, rmajor)
471  CALL cdf_define(ncid, vn_p_max, b_factor)
472  CALL cdf_define(ncid, vn_p_min, p_factor)
473  CALL cdf_define(ncid, vn_curtor, tempscalar)
474  IF (lasym) THEN
475  CALL cdf_define(ncid, vn_bsupsmnc, jbsupsmnch, dimname=restart_dims)
476  CALL cdf_define(ncid, vn_bsupumns, jbsupumnsh, dimname=restart_dims)
477  CALL cdf_define(ncid, vn_bsupvmns, jbsupvmnsh, dimname=restart_dims)
478  CALL cdf_define(ncid, vn_bsubsmnc, jbsupsmnch, dimname=restart_dims)
479  CALL cdf_define(ncid, vn_bsubumns, jbsupumnsh, dimname=restart_dims)
480  CALL cdf_define(ncid, vn_bsubvmns, jbsupvmnsh, dimname=restart_dims)
481  CALL cdf_define(ncid, vn_jksupsmnc, jbsupsmnch, &
482  dimname=restart_dims)
483  CALL cdf_define(ncid, vn_jksupumns, jbsupumnsh, &
484  dimname=restart_dims)
485  CALL cdf_define(ncid, vn_jksupvmns, jbsupvmnsh, &
486  dimname=restart_dims)
487  CALL cdf_define(ncid, vn_pmns, jpmnsh, dimname=restart_dims)
488  END IF
489  END IF recon0
490 
491  CALL cdf_write(ncid, vn_flags, flags)
492  CALL cdf_write(ncid, vn_nsin, ns)
493  CALL cdf_write(ncid, vn_mpolin, mpol)
494  CALL cdf_write(ncid, vn_ntorin, ntor)
495  CALL cdf_write(ncid, vn_nfpin, nfp)
496  CALL cdf_write(ncid, vn_wout, wout_file)
497 
498  CALL cdf_write(ncid, vn_jbsupss, jbsupsmnsh)
499  CALL cdf_write(ncid, vn_jbsupuc, jbsupumnch)
500  CALL cdf_write(ncid, vn_jbsupvc, jbsupvmnch)
501  CALL cdf_write(ncid, vn_jpresc, jpmnch)
502 
503  CALL cdf_write(ncid, vn_chipf, chipf)
504  CALL cdf_write(ncid, vn_phipf, phipf)
505 
506  CALL fourier_context%tomnsp(r1_i, tempmn_w, f_cos)
507  CALL restart_denormalize(tempmn_w, one)
508  CALL cdf_write(ncid, vn_rmnc, tempmn_w)
509  CALL fourier_context%tomnsp(z1_i, tempmn_w, f_sin)
510  CALL restart_denormalize(tempmn_w, one)
511  CALL cdf_write(ncid, vn_zmns, tempmn_w)
512 
513  IF (lasym) THEN
514  CALL cdf_write(ncid, vn_jbsupsc, jbsupsmnch)
515  CALL cdf_write(ncid, vn_jbsupus, jbsupumnsh)
516  CALL cdf_write(ncid, vn_jbsupvs, jbsupvmnsh)
517  CALL cdf_write(ncid, vn_jpress, jpmnsh)
518 
519  CALL fourier_context%tomnsp(r1_i, tempmn_w, f_sin)
520  CALL restart_denormalize(tempmn_w, one)
521  CALL cdf_write(ncid, vn_rmns, tempmn_w)
522  CALL fourier_context%tomnsp(z1_i, tempmn_w, f_cos)
523  CALL restart_denormalize(tempmn_w, one)
524  CALL cdf_write(ncid, vn_zmnc, tempmn_w)
525  END IF
526 
527  CALL cdf_write(ncid, vn_wb, wb)
528  CALL cdf_write(ncid, vn_wp, wp)
529  CALL cdf_write(ncid, vn_rmajor, rmajor)
530 
531  recon1: IF (lrecon) THEN
532  ALLOCATE(bsupsijh(ntheta,nzeta,ns))
533  ALLOCATE(bsupuijh(ntheta,nzeta,ns))
534  ALLOCATE(bsupvijh(ntheta,nzeta,ns))
535  ALLOCATE(pijh(ntheta,nzeta,ns))
536 
537  bsupsijh = 0
538  bsupuijh = 0
539  bsupvijh = 0
540  pijh = 0
541 
542  CALL fourier_context%toijsp(jbsupsmnsh, bsupsijh, f_none, f_sin)
543  CALL fourier_context%toijsp(jbsupumnch, bsupuijh, f_none, f_cos)
544  CALL fourier_context%toijsp(jbsupvmnch, bsupvijh, f_none, f_cos)
545  CALL fourier_context%toijsp(jpmnch, pijh, f_none, f_cos)
546 
547  CALL cdf_write(ncid, vn_p_factor, p_factor)
548  CALL cdf_write(ncid, vn_b_factor, b_factor)
549 
550  IF (lasym) THEN
551  CALL fourier_context%toijsp(jbsupsmnch, bsupsijh, f_sum, f_cos)
552  CALL fourier_context%toijsp(jbsupumnsh, bsupuijh, f_sum, f_sin)
553  CALL fourier_context%toijsp(jbsupvmnsh, bsupvijh, f_sum, f_sin)
554  CALL fourier_context%toijsp(jpmnsh, pijh, f_sum, f_sin)
555  END IF
556 
557  bsupsijh(:,:,1) = 0
558  bsupuijh(:,:,1) = 0
559  bsupvijh(:,:,1) = 0
560  pijh(:,:,1) = 0
561 
562 ! Remove the jacobian term.
563  bsupsijh = bsupsijh/jacobh
564  bsupuijh = bsupuijh/jacobh
565  bsupvijh = bsupvijh/jacobh
566  pijh = pijh/jacobh
567 
568  tempscalar = maxval(pijh(:,:,2:))
569  CALL cdf_write(ncid, vn_p_max, tempscalar)
570  tempscalar = minval(pijh(:,:,2:))
571  CALL cdf_write(ncid, vn_p_min, tempscalar)
572 
573  tempmn_w = 0
574 
575 ! Remove the orthonorm and b_factor so these quantities can be directly summed.
576  CALL fourier_context%tomnsp(bsupsijh, tempmn_w, f_sin)
577  CALL restart_denormalize(tempmn_w, b_factor)
578  CALL cdf_write(ncid, vn_bsupsmns, tempmn_w)
579 
580  CALL fourier_context%tomnsp(bsupuijh, tempmn_w, f_cos)
581  CALL restart_denormalize(tempmn_w, b_factor)
582  CALL cdf_write(ncid, vn_bsupumnc, tempmn_w)
583 
584  CALL fourier_context%tomnsp(bsupvijh, tempmn_w, f_cos)
585  CALL restart_denormalize(tempmn_w, b_factor)
586  CALL cdf_write(ncid, vn_bsupvmnc, tempmn_w)
587 
588 ! Remove the orthonorm and p_factor so this quantity can be directly summed.
589  CALL fourier_context%tomnsp(pijh, tempmn_w, f_cos)
590  CALL restart_denormalize(tempmn_w, p_factor*mu0)
591  CALL cdf_write(ncid, vn_pmnc, tempmn_w)
592 
593  ALLOCATE(bsubsijh(ntheta,nzeta,ns))
594  ALLOCATE(bsubuijh(ntheta,nzeta,ns))
595  ALLOCATE(bsubvijh(ntheta,nzeta,ns))
596 
597  CALL tolowerh(bsupsijh, bsupuijh, bsupvijh, &
598  bsubsijh, bsubuijh, bsubvijh, 1, ns)
599 
600 ! Need these to compute the currents later.
601  ALLOCATE(bsubsmnh(0:mpol,-ntor:ntor,ns))
602  ALLOCATE(bsubumnh(0:mpol,-ntor:ntor,ns))
603  ALLOCATE(bsubvmnh(0:mpol,-ntor:ntor,ns))
604 
605 ! Remove the orthonorm and b_factor so these quantities can be directly summed.
606  CALL fourier_context%tomnsp(bsubsijh, bsubsmnh, f_sin)
607  CALL restart_denormalize(bsubsmnh, b_factor)
608  CALL cdf_write(ncid, vn_bsubsmns, bsubsmnh)
609 
610  CALL fourier_context%tomnsp(bsubuijh, bsubumnh, f_cos)
611  CALL restart_denormalize(bsubumnh, b_factor)
612  CALL cdf_write(ncid, vn_bsubumnc, bsubumnh)
613 
614  CALL fourier_context%tomnsp(bsubvijh, bsubvmnh, f_cos)
615  CALL restart_denormalize(bsubvmnh, b_factor)
616  CALL cdf_write(ncid, vn_bsubvmnc, bsubvmnh)
617 
618 ! Compute currents on the full mesh.
619  ALLOCATE(jksupsmnf(0:mpol,-ntor:ntor,ns))
620  ALLOCATE(jksupumnf(0:mpol,-ntor:ntor,ns))
621  ALLOCATE(jksupvmnf(0:mpol,-ntor:ntor,ns))
622 
623  CALL curl_htof(bsubsmnh, bsubumnh, bsubvmnh, &
624  & jksupsmnf, jksupumnf, jksupvmnf, &
625  & f_sin, 1, ns, .false., tempscalar)
626 
627  CALL restart_denormalize(jksupsmnf, b_factor)
628  CALL restart_denormalize(jksupumnf, b_factor)
629  CALL restart_denormalize(jksupvmnf, b_factor)
630 
631  CALL cdf_write(ncid, vn_jksupsmns, jksupsmnf)
632  CALL cdf_write(ncid, vn_jksupumnc, jksupumnf)
633  CALL cdf_write(ncid, vn_jksupvmnc, jksupvmnf)
634 
635  CALL cdf_write(ncid, vn_curtor, tempscalar)
636 
637  IF (lasym) THEN
638 ! Remove the orthonorm and b_factor so these quantities can be directly summed.
639  CALL fourier_context%tomnsp(bsupsijh, tempmn_w, f_cos)
640  CALL restart_denormalize(tempmn_w, b_factor)
641  CALL cdf_write(ncid, vn_bsupsmnc, tempmn_w)
642 
643  CALL fourier_context%tomnsp(bsupuijh, tempmn_w, f_sin)
644  CALL restart_denormalize(tempmn_w, b_factor)
645  CALL cdf_write(ncid, vn_bsupumns, tempmn_w)
646 
647  CALL fourier_context%tomnsp(bsupvijh, tempmn_w, f_sin)
648  CALL restart_denormalize(tempmn_w, b_factor)
649  CALL cdf_write(ncid, vn_bsupvmns, tempmn_w)
650 
651 ! Remove the orthonorm and p_factor so this quantity can be directly summed.
652  CALL fourier_context%tomnsp(pijh, tempmn_w, f_sin)
653  CALL restart_denormalize(tempmn_w, p_factor*mu0)
654  CALL cdf_write(ncid, vn_pmns, tempmn_w)
655 
656 ! Remove the orthonorm and b_factor so these quantities can be directly summed.
657  CALL fourier_context%tomnsp(bsubsijh, bsubsmnh, f_cos)
658  CALL restart_denormalize(tempmn_w, b_factor)
659  CALL cdf_write(ncid, vn_bsubsmnc, tempmn_w)
660 
661  CALL fourier_context%tomnsp(bsubuijh, bsubumnh, f_sin)
662  CALL restart_denormalize(tempmn_w, b_factor)
663  CALL cdf_write(ncid, vn_bsubumns, tempmn_w)
664 
665  CALL fourier_context%tomnsp(bsubvijh, bsubvmnh, f_sin)
666  CALL restart_denormalize(tempmn_w, b_factor)
667  CALL cdf_write(ncid, vn_bsubvmns, tempmn_w)
668 
669  CALL curl_htof(bsubsmnh, bsubumnh, bsubvmnh, &
670  & jksupsmnf, jksupumnf, jksupvmnf, &
671  & f_cos, 1, ns, .false., tempscalar)
672 
673  CALL restart_denormalize(jksupsmnf, b_factor)
674  CALL restart_denormalize(jksupumnf, b_factor)
675  CALL restart_denormalize(jksupvmnf, b_factor)
676 
677  CALL cdf_write(ncid, vn_jksupsmnc, jksupsmnf)
678  CALL cdf_write(ncid, vn_jksupumns, jksupumnf)
679  CALL cdf_write(ncid, vn_jksupvmns, jksupvmnf)
680  END IF
681 
682  DEALLOCATE(bsupsijh)
683  DEALLOCATE(bsupuijh)
684  DEALLOCATE(bsupvijh)
685  DEALLOCATE(pijh)
686 
687  DEALLOCATE(bsubsmnh)
688  DEALLOCATE(bsubumnh)
689  DEALLOCATE(bsubvmnh)
690 
691  DEALLOCATE(jksupsmnf)
692  DEALLOCATE(jksupumnf)
693  DEALLOCATE(jksupvmnf)
694 
695  END IF recon1
696 
697  DEALLOCATE(tempmn_w)
698 
699  CALL cdf_close(ncid)
700 
701  END SUBROUTINE
702 
703 !-------------------------------------------------------------------------------
718 !-------------------------------------------------------------------------------
719  SUBROUTINE interpit(aold, anew, &
720  & ns_old, ns_new, &
721  & mpol_old, mpol_new, &
722  & ntor_old, ntor_new, &
723  & nfp_old, nfp_new, lhalf)
724  USE stel_constants, ONLY: one, zero
725 
726  IMPLICIT NONE
727 
728 ! Declare Arguments
729  REAL (dp), INTENT(in) :: aold(0:mpol_old,-ntor_old:ntor_old,ns_old)
730  REAL (dp), INTENT(out) :: anew(0:mpol_new,-ntor_new:ntor_new,ns_new)
731  INTEGER, INTENT(in) :: ns_old
732  INTEGER, INTENT(in) :: ns_new
733  INTEGER, INTENT(in) :: mpol_old
734  INTEGER, INTENT(in) :: mpol_new
735  INTEGER, INTENT(in) :: ntor_old
736  INTEGER, INTENT(in) :: ntor_new
737  INTEGER, INTENT(in) :: nfp_old
738  INTEGER, INTENT(in) :: nfp_new
739  LOGICAL, INTENT(in) :: lhalf
740 
741 ! local variables
742  INTEGER :: real_ntor_min
743  INTEGER :: i_n
744  INTEGER :: i_m
745  REAL (dp) :: ds
746  INTEGER :: parity
747  INTEGER :: mpol_min
748  INTEGER :: ntor_min
749 
750 ! Start of executable code
751  mpol_min = min(mpol_new, mpol_old)
752  real_ntor_min = min(ntor_new*nfp_new, ntor_old*nfp_old)
753 
754 ! If the Fourier dimensions of the both arrays are the same just copy them.
755 ! Otherwise copy only the matcing modes.
756  IF (ns_old .eq. ns_new .and. &
757  mpol_old .eq. mpol_new .and. &
758  ntor_old .eq. ntor_new .and. &
759  nfp_old .eq. nfp_new) THEN
760  anew = aold
761  ELSE
762 ! Only some of the toroidal modes will match if nfp_new and nfp_old are
763 ! different.
764  DO i_n = -real_ntor_min, real_ntor_min
765  IF (i_n/nfp_new .eq. i_n/nfp_old) THEN
766  DO i_m = 0, mpol_min
767  IF (mod(i_m, 2) .EQ. 0) THEN
768  parity = -1
769  ELSE
770  parity = 1
771  END IF
772 
773  CALL interpit_1d(aold(i_m, i_n, :), anew(i_m, i_n, :), &
774  ns_old, ns_new, lhalf, parity)
775  END DO
776  ELSE
777  anew(:, i_n, :) = 0
778  END IF
779  END DO
780  END IF
781 
782  END SUBROUTINE
783 
784 !-------------------------------------------------------------------------------
793 !-------------------------------------------------------------------------------
794  SUBROUTINE interpit_1d(aold, anew, ns_old, ns_new, lhalf, parity)
795  USE stel_constants, ONLY: one, zero
796 
797  IMPLICIT NONE
798 
799 ! Declare Arguments
800  REAL (dp), INTENT(in) :: aold(ns_old)
801  REAL (dp), INTENT(out) :: anew(ns_new)
802  INTEGER, INTENT(in) :: ns_old
803  INTEGER, INTENT(in) :: ns_new
804  LOGICAL, INTENT(in) :: lhalf
805  INTEGER, INTENT(in) :: parity
806 
807 ! local variables
808  REAL (dp), ALLOCATABLE, DIMENSION(:) :: temp_y
809  REAL (dp), ALLOCATABLE, DIMENSION(:) :: temp_x
810  REAL (dp), ALLOCATABLE, DIMENSION(:) :: temp_new
811  REAL (dp), ALLOCATABLE, DIMENSION(:) :: temp
812  REAL (dp) :: ds
813  INTEGER :: i_s
814 
815 ! local parameters
816  REAL (dp), PARAMETER :: lower_b = -1.e30_dp
817  REAL (dp), PARAMETER :: upper_b = -1.e30_dp
818 
819 ! Start of executable code
820  anew = 0
821  IF (ns_old .EQ. ns_new) THEN
822  anew = aold
823  RETURN
824  END IF
825 
826 ! Need to radially interpolate. Extend the arrays from -ns to ns to get the
827 ! correct behavior at s=0.
828  ds = 1.0/(ns_old - 1)
829  IF (lhalf) THEN
830  ALLOCATE(temp_y(2*(ns_old - 1)))
831  ALLOCATE(temp_x(2*(ns_old - 1)))
832  ALLOCATE(temp(2*(ns_old - 1)))
833 
834  DO i_s = 2, ns_old
835  temp_x(ns_old - i_s + 1) = -ds*(i_s - 1.5)
836  temp_y(ns_old + i_s - 2) = ds*(i_s - 1.5)
837  temp_y(ns_old - i_s + 1) = parity*aold(i_s)
838  temp_y(ns_old + i_s - 2) = aold(i_s)
839  END DO
840  ELSE
841  ALLOCATE(temp_y(2*ns_old - 1))
842  ALLOCATE(temp_x(2*ns_old - 1))
843  ALLOCATE(temp(2*ns_old - 1))
844 
845  DO i_s = 1, ns_old
846  temp_x(ns_old - i_s + 1) = -ds*(i_s - 1)
847  temp_x(ns_old + i_s - 1) = ds*(i_s - 1)
848  temp_y(ns_old - i_s + 1) = parity*aold(i_s)
849  temp_y(ns_old + i_s - 1) = aold(i_s)
850  END DO
851  END IF
852 
853  CALL spline(temp_x, temp_y, SIZE(temp_x), lower_b, upper_b, temp)
854 
855  ALLOCATE(temp_new(ns_new))
856  ds = 1.0/(ns_new - 1)
857  IF (lhalf) THEN
858  DO i_s = 2, ns_new
859  temp_new(i_s) = ds*(i_s - 1.5)
860  END DO
861  ELSE
862  DO i_s = 1, ns_new
863  temp_new(i_s) = ds*(i_s - 1)
864  END DO
865  END IF
866 
867  CALL splint(temp_x, temp_y, temp, SIZE(temp_x), temp_new, anew)
868 
869  DEALLOCATE(temp_y)
870  DEALLOCATE(temp_x)
871  DEALLOCATE(temp)
872 
873  END SUBROUTINE
874 
875 !-------------------------------------------------------------------------------
883 !-------------------------------------------------------------------------------
884  SUBROUTINE restart_denormalize(xmn, factor)
885  USE island_params, ONLY: fourier_context, mpol=>mpol_i, ntor=>ntor_i
886 
887  IMPLICIT NONE
888 
889 ! Declare Arguments
890  REAL (dp), DIMENSION(:,:,:) :: xmn
891  REAL (dp) :: factor
892 
893 ! Local Variables
894  INTEGER :: s
895 
896 ! Start of executable code.
897  DO s = 1, SIZE(xmn,3)
898  xmn(:,:,s) = fourier_context%orthonorm(0:mpol,-ntor:ntor)*xmn(:,:,s) &
899  / factor
900  END DO
901 
902  END SUBROUTINE
903 
904  END MODULE
restart_mod::vn_wout
character(len= *), parameter vn_wout
Name for the restart file number of wout file modes.
Definition: restart_mod.f90:51
shared_data::r1_i
real(dp), dimension(:,:,:), allocatable r1_i
R coordinates of the computational grid.
Definition: shared_data.f90:153
restart_mod::vn_jbsupvc
character(len= *), parameter vn_jbsupvc
Name for the restart file jbsupvc.
Definition: restart_mod.f90:60
restart_mod::vn_zmnc
character(len= *), parameter vn_zmnc
Name for the restart file zmnc.
Definition: restart_mod.f90:77
restart_mod::vn_jbsupsc
character(len= *), parameter vn_jbsupsc
Name for the restart file jbsupsc.
Definition: restart_mod.f90:62
quantities
This file contains subroutines for allocating and initializing curvilinear magnetic covariant and pre...
Definition: quantities.f90:11
restart_mod::vn_bsupvmnc
character(len= *), parameter vn_bsupvmnc
Name for the restart file bsupvmnc.
Definition: restart_mod.f90:97
restart_mod::restart_version
integer, parameter restart_version
Version number.
Definition: restart_mod.f90:26
v3_utilities::assert_eq
Definition: v3_utilities.f:62
restart_mod::restart_dims
character(len= *), dimension(3), parameter restart_dims
Fourier Dimension names.
Definition: restart_mod.f90:39
restart_mod::vn_wp
character(len= *), parameter vn_wp
Name for the restart file wb.
Definition: restart_mod.f90:134
restart_mod::vn_rmnc
character(len= *), parameter vn_rmnc
Name for the restart file rmnc.
Definition: restart_mod.f90:73
restart_mod::vn_bsubumnc
character(len= *), parameter vn_bsubumnc
Name for the restart file bsubumnc.
Definition: restart_mod.f90:105
island_params::phipf_i
real(dp), dimension(:), allocatable phipf_i
Radial toroidal flux derivative.
Definition: island_params.f90:68
island_params::fourier_context
type(fourier_class), pointer fourier_context
Fourier transform object.
Definition: island_params.f90:76
restart_mod::vn_bsupvmns
character(len= *), parameter vn_bsupvmns
Name for the restart file bsupvmns.
Definition: restart_mod.f90:95
restart_mod::vn_jpresc
character(len= *), parameter vn_jpresc
Name for the restart file jpresc.
Definition: restart_mod.f90:68
fourier
Module contains subroutines for computing FFT on parallel radial intervals. Converts quantities betwe...
Definition: fourier.f90:13
restart_mod::vn_jksupvmns
character(len= *), parameter vn_jksupvmns
Name for the restart file jksupvmns.
Definition: restart_mod.f90:123
restart_mod::vn_chipf
character(len= *), parameter vn_chipf
Name for the restart file chipf.
Definition: restart_mod.f90:82
fourier::m0
integer, parameter m0
m = 0 mode.
Definition: fourier.f90:58
shared_data::lrecon
logical lrecon
Output extra information to the restart file that will be used by V3FIT.
Definition: shared_data.f90:232
restart_mod::vn_p_factor
character(len= *), parameter vn_p_factor
Name for the restart file p_factor.
Definition: restart_mod.f90:128
restart_mod::interpit
subroutine interpit(aold, anew, ns_old, ns_new, mpol_old, mpol_new, ntor_old, ntor_new, nfp_old, nfp_new, lhalf)
Interpolate fourier quantites from the restart file.
Definition: restart_mod.f90:724
fourier::f_sin
integer, parameter f_sin
Sine parity.
Definition: fourier.f90:27
restart_mod::vn_phipf
character(len= *), parameter vn_phipf
Name for the restart file phipf.
Definition: restart_mod.f90:84
restart_mod::vn_jbsupss
character(len= *), parameter vn_jbsupss
Name for the restart file jbsupss.
Definition: restart_mod.f90:56
restart_mod::vn_jksupsmnc
character(len= *), parameter vn_jksupsmnc
Name for the restart file jksupsmnc.
Definition: restart_mod.f90:117
restart_mod::restart_denormalize
subroutine restart_denormalize(xmn, factor)
Denormalize a quantity so the value written to the restart file can be summed directly.
Definition: restart_mod.f90:784
island_params::wp_i
real(dp) wp_i
Thermal energy. 2Pi*wp read_wout_mod::wp.
Definition: island_params.f90:60
restart_mod::vn_pmns
character(len= *), parameter vn_pmns
Name for the restart file pmns.
Definition: restart_mod.f90:111
island_params::rmajor_i
real(dp) rmajor_i
Major radius.
Definition: island_params.f90:64
island_params
This file contains fix parameters related to the computational grids.
Definition: island_params.f90:10
restart_mod::vn_jbsupuc
character(len= *), parameter vn_jbsupuc
Name for the restart file jbsupuc.
Definition: restart_mod.f90:58
restart_mod::restart_lrecon
integer, parameter restart_lrecon
Bit position for the lrecon flag.
Definition: restart_mod.f90:33
island_params::wb_i
real(dp) wb_i
Magnetic energy. 2Pi*wb read_wout_mod::wb.
Definition: island_params.f90:58
restart_mod
Contains routines for writting the restart file.
Definition: restart_mod.f90:10
island_params::gamma
real(dp), parameter gamma
Adiabatic constant.
Definition: island_params.f90:66
restart_mod::vn_jksupumnc
character(len= *), parameter vn_jksupumnc
Name for the restart file jksupumnc.
Definition: restart_mod.f90:121
restart_mod::vn_nsin
character(len= *), parameter vn_nsin
Name for the restart file number of radial points.
Definition: restart_mod.f90:43
shared_data::unit_out
integer, parameter unit_out
File output io unit.
Definition: shared_data.f90:39
restart_mod::vn_bsubvmns
character(len= *), parameter vn_bsubvmns
Name for the restart file bsubvmns.
Definition: restart_mod.f90:107
restart_mod::vn_flags
character(len= *), parameter vn_flags
Name for the restart file number of state flags modes.
Definition: restart_mod.f90:53
shared_data::lasym
logical lasym
Use non-stellarator symmetry.
Definition: shared_data.f90:230
restart_mod::vn_jbsupvs
character(len= *), parameter vn_jbsupvs
Name for the restart file jbsupus.
Definition: restart_mod.f90:66
restart_mod::vn_bsupsmnc
character(len= *), parameter vn_bsupsmnc
Name for the restart file bsupsmnc.
Definition: restart_mod.f90:89
restart_mod::restart_lasym
integer, parameter restart_lasym
Bit position for the lasym flag.
Definition: restart_mod.f90:31
metrics::loadgrid
subroutine loadgrid(istat)
Load the R, Z and lambda arrays from VMEC.
Definition: metrics.f90:218
restart_mod::vn_zmns
character(len= *), parameter vn_zmns
Name for the restart file zmns.
Definition: restart_mod.f90:79
island_params::chipf_i
real(dp), dimension(:), allocatable chipf_i
Radial poloidal flux derivative.
Definition: island_params.f90:70
restart_mod::vn_jksupsmns
character(len= *), parameter vn_jksupsmns
Name for the restart file jksupsmns.
Definition: restart_mod.f90:115
restart_mod::vn_rmajor
character(len= *), parameter vn_rmajor
Name for the restart file rmajor.
Definition: restart_mod.f90:136
restart_mod::radial_dim
character(len= *), dimension(1), parameter radial_dim
Radial Dimension names.
Definition: restart_mod.f90:36
restart_mod::restart_read
integer function restart_read(restart_ext, wout_file, mpolin, ntorin, nfpi
Reads the restart file.
Definition: restart_mod.f90:164
restart_mod::vn_jpress
character(len= *), parameter vn_jpress
Name for the restart file jpress.
Definition: restart_mod.f90:70
restart_mod::vn_jksupvmnc
character(len= *), parameter vn_jksupvmnc
Name for the restart file jksupvmnc.
Definition: restart_mod.f90:125
metrics::tolowerh
subroutine tolowerh(xsupsij, xsupuij, xsupvij,
Converts to covariant component half grid.
Definition: metrics.f90:539
restart_mod::vn_pmnc
character(len= *), parameter vn_pmnc
Name for the restart file pmnc.
Definition: restart_mod.f90:113
fourier::f_sum
integer, parameter f_sum
Sum fouier real space transformed quantity. This is used when a non-stellarator symmetric parity is b...
Definition: fourier.f90:43
fourier::m1
integer, parameter m1
m = 1 mode.
Definition: fourier.f90:60
restart_mod::vn_bsupumns
character(len= *), parameter vn_bsupumns
Name for the restart file bsupumns.
Definition: restart_mod.f90:91
restart_mod::vn_p_max
character(len= *), parameter vn_p_max
Name for the restart file p_max.
Definition: restart_mod.f90:141
restart_mod::vn_wb
character(len= *), parameter vn_wb
Name for the restart file wb.
Definition: restart_mod.f90:132
restart_mod::vn_mpolin
character(len= *), parameter vn_mpolin
Name for the restart file number of poloidal modes.
Definition: restart_mod.f90:45
restart_mod::vn_p_min
character(len= *), parameter vn_p_min
Name for the restart file p_min.
Definition: restart_mod.f90:143
restart_mod::vn_bsubumns
character(len= *), parameter vn_bsubumns
Name for the restart file bsubumns.
Definition: restart_mod.f90:103
restart_mod::vn_curtor
character(len= *), parameter vn_curtor
Name for the restart file total toroidal current.
Definition: restart_mod.f90:138
restart_mod::vn_jksupumns
character(len= *), parameter vn_jksupumns
Name for the restart file jksupumns.
Definition: restart_mod.f90:119
restart_mod::vn_bsupsmns
character(len= *), parameter vn_bsupsmns
Name for the restart file bsupsmns.
Definition: restart_mod.f90:87
shared_data::z1_i
real(dp), dimension(:,:,:), allocatable z1_i
Z coordinates of the computational grid.
Definition: shared_data.f90:155
metrics
Module for reading in VMEC input and computing metric elements on the SIESTA (sqrt-flux) radial grid....
Definition: metrics.f90:13
restart_mod::vn_ntorin
character(len= *), parameter vn_ntorin
Name for the restart file number of toroidal modes.
Definition: restart_mod.f90:47
fourier::f_cos
integer, parameter f_cos
Cosine parity.
Definition: fourier.f90:25
restart_mod::vn_nfpin
character(len= *), parameter vn_nfpin
Name for the restart file number of field periods.
Definition: restart_mod.f90:49
restart_mod::restart_write
subroutine restart_write(restart_ext, wout_file)
Write the restart file.
Definition: restart_mod.f90:363
restart_mod::vn_b_factor
character(len= *), parameter vn_b_factor
Name for the restart file b_factor.
Definition: restart_mod.f90:130
restart_mod::vn_bsubsmnc
character(len= *), parameter vn_bsubsmnc
Name for the restart file bsubsmnc.
Definition: restart_mod.f90:101
shared_data
This file contains variables and parameters used by many modules in SIESTA.
Definition: shared_data.f90:10
restart_mod::vn_bsupumnc
character(len= *), parameter vn_bsupumnc
Name for the restart file bsupumnc.
Definition: restart_mod.f90:93
fourier::fourier_class
Base class containing fourier memory.
Definition: fourier.f90:77
fourier::f_none
integer, parameter f_none
Do not sum fouier real space transformed quantity. This is used when a stellarator symmetric parity i...
Definition: fourier.f90:32
restart_mod::vn_bsubvmnc
character(len= *), parameter vn_bsubvmnc
Name for the restart file bsubvmnc.
Definition: restart_mod.f90:109
restart_mod::vn_rmns
character(len= *), parameter vn_rmns
Name for the restart file rmns.
Definition: restart_mod.f90:75
restart_mod::vn_jbsupus
character(len= *), parameter vn_jbsupus
Name for the restart file jbsupus.
Definition: restart_mod.f90:64
restart_mod::vn_bsubsmns
character(len= *), parameter vn_bsubsmns
Name for the restart file bsubsmns.
Definition: restart_mod.f90:99