V3FIT
metrics.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 !
11 !*******************************************************************************
12 
13  MODULE metrics
14  USE stel_kinds
15  USE stel_constants
16  USE island_params
17  USE shared_data, ONLY: lasym, r1_i, z1_i, ru_i, zu_i, rv_i, zv_i, &
19  USE read_wout_mod, ns_vmec=>ns, mpol_vmec=>mpol, ntor_vmec=>ntor, &
20  rmnc_vmec=>rmnc, zmns_vmec=>zmns, lmns_vmec=>lmns, &
21  xm_vmec=>xm, xn_vmec=>xn, chipf_vmec=>chipf, & ! MRC 4/1/2016
22  rmns_vmec=>rmns, zmnc_vmec=>zmnc, lmnc_vmec=>lmnc, & ! MRC 12/1/2016
23  phipf_vmec=>phipf, presf_vmec=>presf, nfp_vmec=>nfp, &
24  wb_vmec=>wb, wp_vmec=>wp, gamma_vmec=>gamma, &
25  volume_vmec=>volume, raxis_vmec=>raxis, lasym_vmec=>lasym, &
26  iasym_vmec=>iasym
27  USE descriptor_mod, ONLY: iam
28  USE timer_mod
29  USE utilities, ONLY: to_full_mesh
30  USE siesta_namelist, ONLY: l_vessel
31 
32  IMPLICIT NONE
33 
34 !*******************************************************************************
35 ! metrics module variables
36 !*******************************************************************************
38  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: sqrtg
40  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gss
42  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gsu
44  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gsv
46  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: guu
48  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: guv
50  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gvv
52  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: hss
54  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: hsu
56  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: hsv
58  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: huu
60  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: huv
62  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: hvv
64  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gssf
66  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gsuf
68  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gsvf
70  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: guuf
72  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: guvf
74  REAL(dp), DIMENSION(:,:), ALLOCATABLE :: gvvf
75 
77  REAL(dp) :: rmax
79  REAL(dp) :: rmin
81  REAL(dp) :: zmax
83  REAL(dp) :: zmin
84 
86  REAL(dp) :: skston
88  REAL(dp) :: skstoff
89 
90  CONTAINS
91 
92 !*******************************************************************************
93 ! UTILITY SUBROUTINES
94 !*******************************************************************************
95 !-------------------------------------------------------------------------------
101 !-------------------------------------------------------------------------------
102  SUBROUTINE init_metric_elements()
103  USE timer_mod, ONLY: init_timers
104  USE shared_data, ONLY: torflux, polflux
105  USE island_params, hs=>hs_i, ns=>ns_i
106 
107  IMPLICIT NONE
108 
109 ! Start of executable code
110 ! Compute half-mesh lower and upper metric elements and jacobian.
111  CALL half_mesh_metrics(r1_i, ru_i, rv_i, z1_i, zu_i, zv_i)
112 ! Convert half grid elements to the full grid.
113  CALL full_mesh_metrics
114 
115  END SUBROUTINE init_metric_elements
116 
117 !-------------------------------------------------------------------------------
125 !-------------------------------------------------------------------------------
126  SUBROUTINE set_grid_sizes(mpol_in, ntor_in)
128  USE shared_data, ONLY: mblk_size, ndims, lasym
129 
130  IMPLICIT NONE
131 
132 ! Declare Arguments
133  INTEGER, INTENT(IN) :: mpol_in
134  INTEGER, INTENT(IN) :: ntor_in
135 
136 ! Start of executable code
137  ndims = 3
138  IF (lasym) THEN
139  ndims = 2*ndims
140  END IF
141 
142  mpol_i = mpol_in
143  ntor_i = ntor_in
144 
145 ! Set number of points == number of modes for now! (May want mid-points for
146 ! flux conservation)
147  nu_i = mpol_i + 2
148  nv_i = 2*ntor_i + 2
149 
150 !USE 3/2 (ORSZAG) RULE FOR ANTI-ALIASING OF EVOLUTION EQNS
151 !Suppresses RADIAL grid separation in pressure
152  nu_i = 3*nu_i/2
153 !SPH051617 nu_i = nu_i+MOD(nu_i,2)
154  IF (lasym) THEN
155  nu_i = 2*nu_i
156  ELSE IF (mod(nu_i, 2) .ne. 1) THEN
157  nu_i = nu_i + 1
158  END IF
159 
160  nv_i = 3*nv_i/2
161  nv_i = nv_i + mod(nv_i, 2)
162 
163  IF (ntor_i .EQ. 0) THEN
164  nv_i = 1
165  END IF
166  nuv_i = nu_i*nv_i
167  mnmax_i = (mpol_i + 1)*(2*ntor_i + 1) ! Added RS. Contains total number of modes.
168 
169  mblk_size = ndims*mnmax_i
170  CALL assert(mblk_size .NE. 0, 'mblk_size = 0 in set_grid_sizes')
171 
172  END SUBROUTINE
173 
174 !-------------------------------------------------------------------------------
181 !-------------------------------------------------------------------------------
182  SUBROUTINE surfavg(average, q3d, nsmin, nsmax)
184 
185  IMPLICIT NONE
186 
187 ! Declare Arguments
188  REAL(dp), DIMENSION(nsmin:nsmax), INTENT(out) :: average
189  REAL(dp), DIMENSION(nu_i,nv_i,nsmin:nsmax), INTENT(IN) :: q3d
190  INTEGER, INTENT(IN) :: nsmin
191  INTEGER, INTENT(IN) :: nsmax
192 
193 ! Local variables
194  INTEGER :: l
195  INTEGER :: js
196  REAL(dp), ALLOCATABLE, DIMENSION(:,:,:) :: wint
197 
198 ! Start of executable code
199  ALLOCATE(wint(nu_i,nv_i,nsmin:nsmax))
200  DO l = 1, nu_i
201  wint(l,:,:) = fourier_context%cosmui(l,0)
202  END DO
203  DO js = nsmin, nsmax
204  average(js) = sum(q3d(:,:,js)*wint(:,:,js))
205  END DO
206  DEALLOCATE(wint)
207 
208  END SUBROUTINE
209 
210 !-------------------------------------------------------------------------------
216 !-------------------------------------------------------------------------------
217  SUBROUTINE loadgrid(istat)
218  USE vmec_info
219 
220  IMPLICIT NONE
221 
222 ! Declare Arguments
223  INTEGER, INTENT(out) :: istat
224 
225 ! Start executable code.
226  ALLOCATE (r1_i(nu_i, nv_i, ns_i), z1_i(nu_i, nv_i, ns_i), &
227  ru_i(nu_i, nv_i, ns_i), zu_i(nu_i, nv_i, ns_i), &
228  rv_i(nu_i, nv_i, ns_i), zv_i(nu_i, nv_i, ns_i), &
229  stat = istat)
230 
231  CALL rz_to_ijsp(rmnc_i, zmns_i, .false.)
232  IF (lasym) THEN
233  CALL rz_to_ijsp(rmns_i, zmnc_i, .true.)
234  END IF
235 
236  rmax = maxval(r1_i)
237  rmin = minval(r1_i)
238  zmax = maxval(z1_i)
239  zmin = minval(z1_i)
240  IF (.NOT.lasym) THEN
241  zmin = min(zmin, -zmax) ! ONLY TOP HALF
242  END IF
243 
244  END SUBROUTINE
245 
246 !-------------------------------------------------------------------------------
255 !-------------------------------------------------------------------------------
256  SUBROUTINE rz_to_ijsp(rmn, zmn, lasym)
258  USE fourier, ONLY: f_cos, f_sin, f_none, f_sum, f_du, f_dv
259 
260  IMPLICIT NONE
261 
262 ! Declare Arguments
263  REAL(dp), DIMENSION(:,:,:), INTENT(in) :: rmn
264  REAL(dp), DIMENSION(:,:,:), INTENT(in) :: zmn
265  LOGICAL, INTENT(in) :: lasym
266 
267 ! Local variables
268  INTEGER :: sum
269  INTEGER :: m
270  INTEGER :: n
271  INTEGER :: rmode
272  INTEGER :: zmode
273 
274 ! Start executable code.
275  IF (lasym) THEN
276  sum = f_sum
277  rmode = f_sin
278  zmode = f_cos
279  ELSE
280  sum = f_none
281  rmode = f_cos
282  zmode = f_sin
283  END IF
284  m = ior(sum, f_du)
285  n = ior(sum, f_dv)
286 
287  CALL fourier_context%toijsp(rmn, r1_i, sum, rmode)
288  CALL fourier_context%toijsp(rmn, ru_i, m, rmode)
289  CALL fourier_context%toijsp(rmn, rv_i, n, rmode)
290 
291  CALL fourier_context%toijsp(zmn, z1_i, sum, zmode)
292  CALL fourier_context%toijsp(zmn, zu_i, m, zmode)
293  CALL fourier_context%toijsp(zmn, zv_i, n, zmode)
294 
295  END SUBROUTINE
296 
297 !-------------------------------------------------------------------------------
306 !-------------------------------------------------------------------------------
307  SUBROUTINE half_mesh_metrics(r1_i, ru_i, rv_i, z1_i, zu_i, zv_i)
308  USE island_params, nuv=>nuv_i, nu=>nu_i, nv=>nv_i!, ns=>ns_i
309 
310  IMPLICIT NONE
311 
312 ! Declare Arguments
313  REAL(dp), DIMENSION(nuv,ns_i), INTENT(in) :: r1_i
314  REAL(dp), DIMENSION(nuv,ns_i), INTENT(in) :: ru_i
315  REAL(dp), DIMENSION(nuv,ns_i), INTENT(in) :: rv_i
316  REAL(dp), DIMENSION(nuv,ns_i), INTENT(in) :: z1_i
317  REAL(dp), DIMENSION(nuv,ns_i), INTENT(in) :: zu_i
318  REAL(dp), DIMENSION(nuv,ns_i), INTENT(in) :: zv_i
319 
320 ! Local variables
321  INTEGER :: js
322  INTEGER :: istat
323  INTEGER :: imin(2)
324  INTEGER :: imax(2)
325  REAL(dp) :: mintest
326  REAL(dp) :: maxtest
327  REAL(dp), DIMENSION(nuv) :: r12
328  REAL(dp), DIMENSION(nuv) :: rs12
329  REAL(dp), DIMENSION(nuv) :: ru12
330  REAL(dp), DIMENSION(nuv) :: rv12
331  REAL(dp), DIMENSION(nuv) :: zs12
332  REAL(dp), DIMENSION(nuv) :: zu12
333  REAL(dp), DIMENSION(nuv) :: zv12
334 
335 ! Local parameters.
336  REAL(dp), PARAMETER :: p5 = 1.d0/2.d0
337  REAL(dp), PARAMETER :: zero = 0
338 
339 ! Start executable code.
340 
341 ! ALLOCATE METRIC ELEMENT ARRAYS
342  ALLOCATE(sqrtg(nuv,ns_i), &
343  gss(nuv,ns_i), gsu(nuv,ns_i), gsv(nuv,ns_i), &
344  guu(nuv,ns_i), guv(nuv,ns_i), gvv(nuv,ns_i), stat=istat)
345 
346 ! COMPUTE ALL ON THE HALF MESH
347  sqrtg(:,1) = 0
348 
349  DO js = 2, ns_i
350  r12 = p5*(r1_i(:,js) + r1_i(:,js - 1))
351  rs12 = (r1_i(:,js) - r1_i(:,js - 1))*ohs_i
352  ru12 = p5*(ru_i(:,js) + ru_i(:,js - 1))
353  rv12 = p5*(rv_i(:,js) + rv_i(:,js - 1))
354  zs12 = (z1_i(:,js) - z1_i(:,js - 1))*ohs_i
355  zu12 = p5*(zu_i(:,js) + zu_i(:,js - 1))
356  zv12 = p5*(zv_i(:,js) + zv_i(:,js - 1))
357  gss(:,js) = rs12*rs12 + zs12*zs12
358  gsu(:,js) = rs12*ru12 + zs12*zu12
359  gsv(:,js) = rs12*rv12 + zs12*zv12
360  guu(:,js) = ru12*ru12 + zu12*zu12
361  guv(:,js) = ru12*rv12 + zu12*zv12
362  gvv(:,js) = rv12*rv12 + r12*r12 + zv12*zv12
363  sqrtg(:,js) = r12*(ru12*zs12 - rs12*zu12)
364  END DO
365 
366  sqrtg(:,1) = sqrtg(:,2)
367 
368 !NEED THESE FOR CURRENT CALCULATION AT ORIGIN
369  gss(:,1) = gss(:,2)
370  gsu(:,1) = 0
371  gsv(:,1) = gsv(:,2)
372  guu(:,1) = 0
373  guv(:,1) = 0
374  gvv(:,1) = gvv(:,2)
375 
376  mintest = minval(sqrtg(:,2:))
377  maxtest = maxval(sqrtg(:,2:))
378 
379  IF (mintest*maxtest .LE. zero) THEN
380  imin = minloc(sqrtg(:,2:))
381  imax = maxloc(sqrtg(:,2:))
382  IF (iam .EQ. 0) THEN
383  WRITE(*,1000) mintest, imin(2), maxtest, imax(2)
384  END IF
385  CALL assert(mintest*maxtest.GT.zero, &
386  ' Jacobian changed sign in half_mesh_metrics')
387  END IF
388 
389 1000 FORMAT(' MIN-G: ',1p,e12.4,' AT JS: ',i5, &
390  ' MAX-G: ',1p,e12.4,' AT JS: ',i5)
391 
392  END SUBROUTINE
393 
394 !-------------------------------------------------------------------------------
399 !-------------------------------------------------------------------------------
400  SUBROUTINE full_mesh_metrics
401  USE island_params, nuv=>nuv_i, nzeta=>nv_i, ntheta=>nu_i!, ns=>ns_i
402 
403  IMPLICIT NONE
404 
405 ! Local variables
406  INTEGER :: istat
407  INTEGER :: js
408  REAL (dp), DIMENSION(nuv) :: det
409 
410 ! Start executable code.
411  ALLOCATE(gssf(nuv,ns_i), guuf(nuv,ns_i), gvvf(nuv,ns_i), &
412  gsuf(nuv,ns_i), gsvf(nuv,ns_i), guvf(nuv,ns_i), &
413  hss(nuv,ns_i), huu(nuv,ns_i), hvv(nuv,ns_i), &
414  hsu(nuv,ns_i), hsv(nuv,ns_i), huv(nuv,ns_i), stat=istat)
415 
416  CALL to_full_mesh(gss, gssf)
417  CALL to_full_mesh(guu, guuf)
418  CALL to_full_mesh(gvv, gvvf)
419  CALL to_full_mesh(gsu, gsuf)
420  CALL to_full_mesh(gsv, gsvf)
421  CALL to_full_mesh(guv, guvf)
422 
423 ! Compute upper metric elements (inverse of lower matrix) on full mesh.
424 ! NOTE: DET == 1/det|Gij| == 1/sqrtg(full)**2
425  DO js = 1, ns_i
426  det(:) = gssf(:,js)*(guuf(:,js)*gvvf(:,js) - guvf(:,js)*guvf(:,js)) &
427  + gsuf(:,js)*(guvf(:,js)*gsvf(:,js) - gsuf(:,js)*gvvf(:,js)) &
428  + gsvf(:,js)*(gsuf(:,js)*guvf(:,js) - gsvf(:,js)*guuf(:,js))
429  IF (any(det .LE. zero)) THEN
430  IF (iam .EQ. 0) THEN
431  WRITE (*,1000) js
432  END IF
433  IF (js .EQ. ns_i) THEN
434  det(:) = sqrtg(:,ns_i)**2
435  gssf(:,ns_i) = gss(:,ns_i)
436  guuf(:,ns_i) = guu(:,ns_i)
437  gvvf(:,ns_i) = gvv(:,ns_i)
438  gsuf(:,ns_i) = gsu(:,ns_i)
439  guvf(:,ns_i) = guv(:,ns_i)
440  gsvf(:,ns_i) = gsv(:,ns_i)
441  ELSE
442  det(:) = ((sqrtg(:,js) + sqrtg(:,js + 1))/2)**2
443  END IF
444  END IF
445  det = one/det
446  hss(:,js) = det(:)* &
447  (guuf(:,js)*gvvf(:,js) - guvf(:,js)*guvf(:,js))
448  hsu(:,js) = det(:)* &
449  (guvf(:,js)*gsvf(:,js) - gsuf(:,js)*gvvf(:,js))
450  hsv(:,js) = det(:)* &
451  (gsuf(:,js)*guvf(:,js) - gsvf(:,js)*guuf(:,js))
452  huu(:,js) = det(:)* &
453  (gssf(:,js)*gvvf(:,js) - gsvf(:,js)*gsvf(:,js))
454  huv(:,js) = det(:)* &
455  (gsvf(:,js)*gsuf(:,js) - gssf(:,js)*guvf(:,js))
456  hvv(:,js) = det(:)* &
457  (gssf(:,js)*guuf(:,js) - gsuf(:,js)*gsuf(:,js))
458 
459  END DO
460 
461 1000 FORMAT('Determinant |gijf| <= 0 at js = ',i4)
462 
463  END SUBROUTINE full_mesh_metrics
464 
465 !-------------------------------------------------------------------------------
483 !-------------------------------------------------------------------------------
484  SUBROUTINE toupper(xsubsij, xsubuij, xsubvij, &
485  xsupsij, xsupuij, xsupvij, nsmin, nsmax)
486  USE stel_kinds
487  USE island_params, nuv=>nuv_i
488 
489  IMPLICIT NONE
490 
491 ! Declare Arguments
492  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(in) :: xsubsij
493  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(in) :: xsubuij
494  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(in) :: xsubvij
495  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(out) :: xsupsij
496  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(out) :: xsupuij
497  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(out) :: xsupvij
498  INTEGER, INTENT(in) :: nsmin
499  INTEGER, INTENT(in) :: nsmax
500 
501 ! Local variables
502  INTEGER :: js
503 
504 ! Start executable code.
505  DO js = nsmin, nsmax
506  xsupsij(:,js) = hss(:,js)*xsubsij(:,js) &
507  + hsu(:,js)*xsubuij(:,js) &
508  + hsv(:,js)*xsubvij(:,js)
509  xsupuij(:,js) = hsu(:,js)*xsubsij(:,js) &
510  + huu(:,js)*xsubuij(:,js) &
511  + huv(:,js)*xsubvij(:,js)
512  xsupvij(:,js) = hsv(:,js)*xsubsij(:,js) &
513  + huv(:,js)*xsubuij(:,js) &
514  + hvv(:,js)*xsubvij(:,js)
515  END DO
516 
517  END SUBROUTINE
518 
519 !-------------------------------------------------------------------------------
537 !-------------------------------------------------------------------------------
538  SUBROUTINE tolowerh(xsupsij, xsupuij, xsupvij, &
539  xsubsij, xsubuij, xsubvij, nsmin, nsmax)
540  USE stel_kinds
541  USE island_params, ONLY: nuv=>nuv_i
542 
543  IMPLICIT NONE
544 
545 ! Declare Arguments
546  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(in) :: xsupsij
547  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(in) :: xsupuij
548  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(in) :: xsupvij
549  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(out) :: xsubsij
550  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(out) :: xsubuij
551  REAL(dp), DIMENSION(nuv,nsmin:nsmax), INTENT(out) :: xsubvij
552  INTEGER, INTENT(IN) :: nsmin
553  INTEGER, INTENT(IN) :: nsmax
554 
555 ! Start executable code.
556  xsubsij = gss(:,nsmin:nsmax)*xsupsij &
557  + gsu(:,nsmin:nsmax)*xsupuij &
558  + gsv(:,nsmin:nsmax)*xsupvij
559 
560  xsubuij = gsu(:,nsmin:nsmax)*xsupsij &
561  + guu(:,nsmin:nsmax)*xsupuij &
562  + guv(:,nsmin:nsmax)*xsupvij
563 
564  xsubvij = gsv(:,nsmin:nsmax)*xsupsij &
565  + guv(:,nsmin:nsmax)*xsupuij &
566  + gvv(:,nsmin:nsmax)*xsupvij
567 
568  END SUBROUTINE tolowerh
569 
570 !-------------------------------------------------------------------------------
588 !-------------------------------------------------------------------------------
589  SUBROUTINE tolowerf(xsupsij, xsupuij, xsupvij, &
590  xsubsij, xsubuij, xsubvij, nsmin, nsmax)
591  USE stel_kinds
592  USE island_params, ONLY: ns=>ns_i, nuv=>nuv_i
593 
594  IMPLICIT NONE
595 
596 ! Declare Arguments
597  REAL(dp), DIMENSION(nuv,ns), INTENT(in) :: xsupsij ! MRC: Why ns and not nsmin:nsmax here?
598  REAL(dp), DIMENSION(nuv,ns), INTENT(in) :: xsupuij
599  REAL(dp), DIMENSION(nuv,ns), INTENT(in) :: xsupvij
600  REAL(dp), DIMENSION(nuv,ns), INTENT(out) :: xsubsij
601  REAL(dp), DIMENSION(nuv,ns), INTENT(out) :: xsubuij
602  REAL(dp), DIMENSION(nuv,ns), INTENT(out) :: xsubvij
603  INTEGER, INTENT(in) :: nsmin
604  INTEGER, INTENT(in) :: nsmax
605 
606 ! Local variables
607  INTEGER :: js
608  INTEGER :: js1
609 
610 ! Start executable code.
611  DO js = nsmin, nsmax
612  js1 = js - nsmin + 1
613  xsubsij(:,js1) = gssf(:,js)*xsupsij(:,js1) &
614  + gsuf(:,js)*xsupuij(:,js1) &
615  + gsvf(:,js)*xsupvij(:,js1)
616 
617  xsubuij(:,js1) = gsuf(:,js)*xsupsij(:,js1) &
618  + guuf(:,js)*xsupuij(:,js1) &
619  + guvf(:,js)*xsupvij(:,js1)
620 
621  xsubvij(:,js1) = gsvf(:,js)*xsupsij(:,js1) &
622  + guvf(:,js)*xsupuij(:,js1) &
623  + gvvf(:,js)*xsupvij(:,js1)
624  END DO
625 
626 ! FIXME: CHECK THIS
627  IF (nsmin .eq. 1) THEN
628  xsubuij(:,1) = 0
629  END IF
630 
631  END SUBROUTINE tolowerf
632 
633 !-------------------------------------------------------------------------------
640 !-------------------------------------------------------------------------------
641  SUBROUTINE cleanup_metric_elements
643  USE vmec_info
644 
645  IMPLICIT NONE
646 
647 ! Local variables
648  INTEGER :: istat
649 
650 ! Start executable code.
651  DEALLOCATE(gss, gsu, gsv, guu, guv, gvv, &
652  hss, hsu, hsv, huu, huv, hvv, &
653  phipf_i, chipf_i, presf_i, stat=istat)
654 
655  DEALLOCATE (r1_i, z1_i, ru_i, zu_i, rv_i, zv_i)
656 
657  DEALLOCATE(fourier_context)
658  CALL vmec_info_destruct_island
659 
660  END SUBROUTINE cleanup_metric_elements
661 
662 !-------------------------------------------------------------------------------
664 !-------------------------------------------------------------------------------
665  SUBROUTINE dealloc_full_lower_metrics
666  USE stel_kinds
667 
668  IMPLICIT NONE
669 
670 ! Local variables
671  INTEGER :: istat
672 
673 ! Start executable code.
674  DEALLOCATE(gssf, guuf, gvvf, gsuf, gsvf, guvf, stat = istat)
675 ! IF (istat .NE. 0) CALL ErrorAbort('Problem dealloc. in LOWER_METRIC_FULL')
676 
677  END SUBROUTINE dealloc_full_lower_metrics
678 
679  END MODULE
shared_data::r1_i
real(dp), dimension(:,:,:), allocatable r1_i
R coordinates of the computational grid.
Definition: shared_data.f90:153
shared_data::lverbose
logical lverbose
Use verbose screen output.
Definition: shared_data.f90:234
shared_data::zv_i
real(dp), dimension(:,:,:), allocatable zv_i
dZ/dv coordinates of the computational grid.
Definition: shared_data.f90:163
metrics::hvv
real(dp), dimension(:,:), allocatable hvv
Upper metric tensor. e^v . e^v.
Definition: metrics.f90:62
metrics::cleanup_metric_elements
subroutine cleanup_metric_elements
Deallocate memory containing metric elements on the half mesh.
Definition: metrics.f90:642
siesta_namelist::l_vessel
logical l_vessel
If extended grid is to be used using an available vessel file.
Definition: siesta_namelist.f90:140
metrics::zmax
real(dp) zmax
Maximum of the grid Z inside the last closed flux surface.
Definition: metrics.f90:81
metrics::guuf
real(dp), dimension(:,:), allocatable guuf
Lower metric tensor full grid. e_u . e_u.
Definition: metrics.f90:70
metrics::sqrtg
real(dp), dimension(:,:), allocatable sqrtg
Jacobian on half grid.
Definition: metrics.f90:38
metrics::gsv
real(dp), dimension(:,:), allocatable gsv
Lower metric tensor half grid. e_s . e_v.
Definition: metrics.f90:44
metrics::dealloc_full_lower_metrics
subroutine dealloc_full_lower_metrics
Deallocate memory containing metric elements on the full mesh.
Definition: metrics.f90:666
metrics::gsu
real(dp), dimension(:,:), allocatable gsu
Lower metric tensor half grid. e_s . e_u.
Definition: metrics.f90:42
shared_data::torflux
real(dp), dimension(:), allocatable torflux
Toroidal flux profile.
Definition: shared_data.f90:169
metrics::gvvf
real(dp), dimension(:,:), allocatable gvvf
Lower metric tensor full grid. e_v . e_v.
Definition: metrics.f90:74
island_params::fourier_context
type(fourier_class), pointer fourier_context
Fourier transform object.
Definition: island_params.f90:76
fourier
Module contains subroutines for computing FFT on parallel radial intervals. Converts quantities betwe...
Definition: fourier.f90:13
metrics::skston
real(dp) skston
Start timer.
Definition: metrics.f90:86
shared_data::jsupvdota
real(dp) jsupvdota
FIXME: UNKNOWN.
Definition: shared_data.f90:167
fourier::f_sin
integer, parameter f_sin
Sine parity.
Definition: fourier.f90:27
metrics::skstoff
real(dp) skstoff
Stop timer.
Definition: metrics.f90:88
metrics::hsu
real(dp), dimension(:,:), allocatable hsu
Upper metric tensor. e^s . e^u.
Definition: metrics.f90:54
metrics::gsuf
real(dp), dimension(:,:), allocatable gsuf
Lower metric tensor full grid. e_s . e_u.
Definition: metrics.f90:66
metrics::guu
real(dp), dimension(:,:), allocatable guu
Lower metric tensor half grid. e_u . e_u.
Definition: metrics.f90:46
metrics::gvv
real(dp), dimension(:,:), allocatable gvv
Lower metric tensor half grid. e_v . e_v.
Definition: metrics.f90:50
island_params
This file contains fix parameters related to the computational grids.
Definition: island_params.f90:10
metrics::guvf
real(dp), dimension(:,:), allocatable guvf
Lower metric tensor full grid. e_u . e_v.
Definition: metrics.f90:72
shared_data::zu_i
real(dp), dimension(:,:,:), allocatable zu_i
dZ/du coordinates of the computational grid.
Definition: shared_data.f90:159
metrics::rmax
real(dp) rmax
Maximum of the grid R inside the last closed flux surface.
Definition: metrics.f90:77
shared_data::polflux
real(dp), dimension(:), allocatable polflux
Poloidal flux profile.
Definition: shared_data.f90:171
metrics::huv
real(dp), dimension(:,:), allocatable huv
Upper metric tensor. e^u . e^v.
Definition: metrics.f90:60
island_params::gamma
real(dp), parameter gamma
Adiabatic constant.
Definition: island_params.f90:66
shared_data::lasym
logical lasym
Use non-stellarator symmetry.
Definition: shared_data.f90:230
metrics::rmin
real(dp) rmin
Minimum of the grid R inside the last closed flux surface.
Definition: metrics.f90:79
island_params::hs_i
real(dp) hs_i
Radial grid spacing. hs = s_i+1 - s_i.
Definition: island_params.f90:47
metrics::toupper
subroutine toupper(xsubsij, xsubuij, xsubvij,
Converts to contravariant component full grid.
Definition: metrics.f90:485
metrics::hsv
real(dp), dimension(:,:), allocatable hsv
Upper metric tensor. e^s . e^v.
Definition: metrics.f90:56
metrics::loadgrid
subroutine loadgrid(istat)
Load the R, Z and lambda arrays from VMEC.
Definition: metrics.f90:218
metrics::gss
real(dp), dimension(:,:), allocatable gss
Lower metric tensor half grid. e_s . e_s.
Definition: metrics.f90:40
metrics::init_metric_elements
subroutine init_metric_elements()
Initialize the metric elements.
Definition: metrics.f90:103
shared_data::rv_i
real(dp), dimension(:,:,:), allocatable rv_i
dR/dv coordinates of the computational grid.
Definition: shared_data.f90:161
metrics::tolowerf
subroutine tolowerf(xsupsij, xsupuij, xsupvij,
Converts to covariant component full grid.
Definition: metrics.f90:590
metrics::gsvf
real(dp), dimension(:,:), allocatable gsvf
Lower metric tensor full grid. e_s . e_v.
Definition: metrics.f90:68
metrics::half_mesh_metrics
subroutine half_mesh_metrics(r1_i, ru_i, rv_i, z1_i, zu_i, zv_i)
Compute the metric elements on the half mesh.
Definition: metrics.f90:308
fourier::f_du
integer, parameter f_du
Poloidal derivative flag.
Definition: fourier.f90:36
shared_data::nsp
integer nsp
Total radial grid size in the VMEC region.
Definition: shared_data.f90:64
metrics::zmin
real(dp) zmin
Minimum of the grid Z inside the last closed flux surface.
Definition: metrics.f90:83
metrics::rz_to_ijsp
subroutine rz_to_ijsp(rmn, zmn, lasym)
Transform from fourier R Z to real space qauntities.
Definition: metrics.f90:257
metrics::tolowerh
subroutine tolowerh(xsupsij, xsupuij, xsupvij,
Converts to covariant component half grid.
Definition: metrics.f90:539
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
shared_data::ndims
integer ndims
Number of independent variables.
Definition: shared_data.f90:58
metrics::surfavg
subroutine surfavg(average, q3d, nsmin, nsmax)
Surface average a quantity.
Definition: metrics.f90:183
island_params::ohs_i
real(dp) ohs_i
Radial grid derivative factor. d X/ ds = (x_i+1 - x_i)/(s_i+1 - s_i) (1) Where ohs = 1/(s_i+1 - s_i) ...
Definition: island_params.f90:51
metrics::gssf
real(dp), dimension(:,:), allocatable gssf
Lower metric tensor full grid. e_s . e_s.
Definition: metrics.f90:64
shared_data::mblk_size
integer mblk_size
Block size. (mpol + 1)*(2*ntor + 1)*ndims.
Definition: shared_data.f90:62
metrics::huu
real(dp), dimension(:,:), allocatable huu
Upper metric tensor. e^u . e^u.
Definition: metrics.f90:58
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
fourier::f_cos
integer, parameter f_cos
Cosine parity.
Definition: fourier.f90:25
siesta_namelist
This file contains all the variables and maximum sizes of the inputs for a SIESTA namelist input file...
Definition: siesta_namelist.f90:103
shared_data::ru_i
real(dp), dimension(:,:,:), allocatable ru_i
dR/du coordinates of the computational grid.
Definition: shared_data.f90:157
metrics::hss
real(dp), dimension(:,:), allocatable hss
Upper metric tensor half grid. e^s . e^s.
Definition: metrics.f90:52
shared_data
This file contains variables and parameters used by many modules in SIESTA.
Definition: shared_data.f90:10
metrics::set_grid_sizes
subroutine set_grid_sizes(mpol_in, ntor_in)
Set grid sizes.
Definition: metrics.f90:127
fourier::f_dv
integer, parameter f_dv
Toroidal derivative flag.
Definition: fourier.f90:38
metrics::full_mesh_metrics
subroutine full_mesh_metrics
Gets the full grid metric elements.
Definition: metrics.f90:401
metrics::guv
real(dp), dimension(:,:), allocatable guv
Lower metric tensor half grid. e_u . e_v.
Definition: metrics.f90:48
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