|
V3FIT
|
Go to the documentation of this file.
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,
22 rmns_vmec=>rmns, zmnc_vmec=>zmnc, lmnc_vmec=>lmnc,
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,
27 USE descriptor_mod,
ONLY: iam
29 USE utilities,
ONLY: to_full_mesh
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
103 USE timer_mod,
ONLY: init_timers
133 INTEGER,
INTENT(IN) :: mpol_in
134 INTEGER,
INTENT(IN) :: ntor_in
156 ELSE IF (mod(nu_i, 2) .ne. 1)
THEN
161 nv_i = nv_i + mod(nv_i, 2)
163 IF (ntor_i .EQ. 0)
THEN
167 mnmax_i = (mpol_i + 1)*(2*ntor_i + 1)
170 CALL assert(
mblk_size .NE. 0,
'mblk_size = 0 in set_grid_sizes')
182 SUBROUTINE surfavg(average, q3d, nsmin, nsmax)
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
196 REAL(dp),
ALLOCATABLE,
DIMENSION(:,:,:) :: wint
199 ALLOCATE(wint(nu_i,nv_i,nsmin:nsmax))
204 average(js) = sum(q3d(:,:,js)*wint(:,:,js))
223 INTEGER,
INTENT(out) :: istat
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),
263 REAL(dp),
DIMENSION(:,:,:),
INTENT(in) :: rmn
264 REAL(dp),
DIMENSION(:,:,:),
INTENT(in) :: zmn
265 LOGICAL,
INTENT(in) :: lasym
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
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
336 REAL(dp),
PARAMETER :: p5 = 1.d0/2.d0
337 REAL(dp),
PARAMETER :: zero = 0
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)
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)
376 mintest = minval(
sqrtg(:,2:))
377 maxtest = maxval(
sqrtg(:,2:))
379 IF (mintest*maxtest .LE. zero)
THEN
380 imin = minloc(
sqrtg(:,2:))
381 imax = maxloc(
sqrtg(:,2:))
383 WRITE(*,1000) mintest, imin(2), maxtest, imax(2)
385 CALL assert(mintest*maxtest.GT.zero,
386 ' Jacobian changed sign in half_mesh_metrics')
389 1000
FORMAT(
' MIN-G: ',1p,e12.4,
' AT JS: ',i5,
390 ' MAX-G: ',1p,e12.4,
' AT JS: ',i5)
408 REAL (dp),
DIMENSION(nuv) :: det
411 ALLOCATE(
gssf(nuv,ns_i),
guuf(nuv,ns_i),
gvvf(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
429 IF (any(det .LE. zero))
THEN
433 IF (js .EQ. ns_i)
THEN
434 det(:) =
sqrtg(:,ns_i)**2
442 det(:) = ((
sqrtg(:,js) +
sqrtg(:,js + 1))/2)**2
461 1000
FORMAT(
'Determinant |gijf| <= 0 at js = ',i4)
484 SUBROUTINE toupper(xsubsij, xsubuij, xsubvij, &
485 xsupsij, xsupuij, xsupvij, nsmin, nsmax)
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
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)
538 SUBROUTINE tolowerh(xsupsij, xsupuij, xsupvij, &
539 xsubsij, xsubuij, xsubvij, nsmin, nsmax)
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
556 xsubsij =
gss(:,nsmin:nsmax)*xsupsij
557 +
gsu(:,nsmin:nsmax)*xsupuij
558 +
gsv(:,nsmin:nsmax)*xsupvij
560 xsubuij =
gsu(:,nsmin:nsmax)*xsupsij
561 +
guu(:,nsmin:nsmax)*xsupuij
562 +
guv(:,nsmin:nsmax)*xsupvij
564 xsubvij =
gsv(:,nsmin:nsmax)*xsupsij
565 +
guv(:,nsmin:nsmax)*xsupuij
566 +
gvv(:,nsmin:nsmax)*xsupvij
589 SUBROUTINE tolowerf(xsupsij, xsupuij, xsupvij, &
590 xsubsij, xsubuij, xsubvij, nsmin, nsmax)
597 REAL(dp),
DIMENSION(nuv,ns),
INTENT(in) :: xsupsij
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
613 xsubsij(:,js1) =
gssf(:,js)*xsupsij(:,js1)
614 +
gsuf(:,js)*xsupuij(:,js1)
615 +
gsvf(:,js)*xsupvij(:,js1)
617 xsubuij(:,js1) =
gsuf(:,js)*xsupsij(:,js1)
618 +
guuf(:,js)*xsupuij(:,js1)
619 +
guvf(:,js)*xsupvij(:,js1)
621 xsubvij(:,js1) =
gsvf(:,js)*xsupsij(:,js1)
622 +
guvf(:,js)*xsupuij(:,js1)
623 +
gvvf(:,js)*xsupvij(:,js1)
627 IF (nsmin .eq. 1)
THEN
653 phipf_i, chipf_i, presf_i, stat=istat)
658 CALL vmec_info_destruct_island
real(dp), dimension(:,:,:), allocatable r1_i
R coordinates of the computational grid.
logical lverbose
Use verbose screen output.
real(dp), dimension(:,:,:), allocatable zv_i
dZ/dv coordinates of the computational grid.
real(dp), dimension(:,:), allocatable hvv
Upper metric tensor. e^v . e^v.
subroutine cleanup_metric_elements
Deallocate memory containing metric elements on the half mesh.
logical l_vessel
If extended grid is to be used using an available vessel file.
real(dp) zmax
Maximum of the grid Z inside the last closed flux surface.
real(dp), dimension(:,:), allocatable guuf
Lower metric tensor full grid. e_u . e_u.
real(dp), dimension(:,:), allocatable sqrtg
Jacobian on half grid.
real(dp), dimension(:,:), allocatable gsv
Lower metric tensor half grid. e_s . e_v.
subroutine dealloc_full_lower_metrics
Deallocate memory containing metric elements on the full mesh.
real(dp), dimension(:,:), allocatable gsu
Lower metric tensor half grid. e_s . e_u.
real(dp), dimension(:), allocatable torflux
Toroidal flux profile.
real(dp), dimension(:,:), allocatable gvvf
Lower metric tensor full grid. e_v . e_v.
type(fourier_class), pointer fourier_context
Fourier transform object.
Module contains subroutines for computing FFT on parallel radial intervals. Converts quantities betwe...
real(dp) skston
Start timer.
real(dp) jsupvdota
FIXME: UNKNOWN.
integer, parameter f_sin
Sine parity.
real(dp) skstoff
Stop timer.
real(dp), dimension(:,:), allocatable hsu
Upper metric tensor. e^s . e^u.
real(dp), dimension(:,:), allocatable gsuf
Lower metric tensor full grid. e_s . e_u.
real(dp), dimension(:,:), allocatable guu
Lower metric tensor half grid. e_u . e_u.
real(dp), dimension(:,:), allocatable gvv
Lower metric tensor half grid. e_v . e_v.
This file contains fix parameters related to the computational grids.
real(dp), dimension(:,:), allocatable guvf
Lower metric tensor full grid. e_u . e_v.
real(dp), dimension(:,:,:), allocatable zu_i
dZ/du coordinates of the computational grid.
real(dp) rmax
Maximum of the grid R inside the last closed flux surface.
real(dp), dimension(:), allocatable polflux
Poloidal flux profile.
real(dp), dimension(:,:), allocatable huv
Upper metric tensor. e^u . e^v.
real(dp), parameter gamma
Adiabatic constant.
logical lasym
Use non-stellarator symmetry.
real(dp) rmin
Minimum of the grid R inside the last closed flux surface.
real(dp) hs_i
Radial grid spacing. hs = s_i+1 - s_i.
subroutine toupper(xsubsij, xsubuij, xsubvij,
Converts to contravariant component full grid.
real(dp), dimension(:,:), allocatable hsv
Upper metric tensor. e^s . e^v.
subroutine loadgrid(istat)
Load the R, Z and lambda arrays from VMEC.
real(dp), dimension(:,:), allocatable gss
Lower metric tensor half grid. e_s . e_s.
subroutine init_metric_elements()
Initialize the metric elements.
real(dp), dimension(:,:,:), allocatable rv_i
dR/dv coordinates of the computational grid.
subroutine tolowerf(xsupsij, xsupuij, xsupvij,
Converts to covariant component full grid.
real(dp), dimension(:,:), allocatable gsvf
Lower metric tensor full grid. e_s . e_v.
subroutine half_mesh_metrics(r1_i, ru_i, rv_i, z1_i, zu_i, zv_i)
Compute the metric elements on the half mesh.
integer, parameter f_du
Poloidal derivative flag.
integer nsp
Total radial grid size in the VMEC region.
real(dp) zmin
Minimum of the grid Z inside the last closed flux surface.
subroutine rz_to_ijsp(rmn, zmn, lasym)
Transform from fourier R Z to real space qauntities.
subroutine tolowerh(xsupsij, xsupuij, xsupvij,
Converts to covariant component half grid.
integer, parameter f_sum
Sum fouier real space transformed quantity. This is used when a non-stellarator symmetric parity is b...
integer ndims
Number of independent variables.
subroutine surfavg(average, q3d, nsmin, nsmax)
Surface average a quantity.
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) ...
real(dp), dimension(:,:), allocatable gssf
Lower metric tensor full grid. e_s . e_s.
integer mblk_size
Block size. (mpol + 1)*(2*ntor + 1)*ndims.
real(dp), dimension(:,:), allocatable huu
Upper metric tensor. e^u . e^u.
real(dp), dimension(:,:,:), allocatable z1_i
Z coordinates of the computational grid.
Module for reading in VMEC input and computing metric elements on the SIESTA (sqrt-flux) radial grid....
integer, parameter f_cos
Cosine parity.
This file contains all the variables and maximum sizes of the inputs for a SIESTA namelist input file...
real(dp), dimension(:,:,:), allocatable ru_i
dR/du coordinates of the computational grid.
real(dp), dimension(:,:), allocatable hss
Upper metric tensor half grid. e^s . e^s.
This file contains variables and parameters used by many modules in SIESTA.
subroutine set_grid_sizes(mpol_in, ntor_in)
Set grid sizes.
integer, parameter f_dv
Toroidal derivative flag.
subroutine full_mesh_metrics
Gets the full grid metric elements.
real(dp), dimension(:,:), allocatable guv
Lower metric tensor half grid. e_u . e_v.
integer, parameter f_none
Do not sum fouier real space transformed quantity. This is used when a stellarator symmetric parity i...