ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
triangularity.f90
Go to the documentation of this file.
1 subroutine triangularity(xaxis, xell, xtria)
2 !----------------------------------------------------------------------
3 ! subroutine to calculate the triangularity of flux surfaces as a function
4 ! of the minor radius (not called by Helena)
5 !----------------------------------------------------------------------
6 
7  use itm_types
8  use mod_mesh
9  use mod_helena_io, only : out_he
11  use mod_output
12  use math_functions
13 
14  implicit none
15 
16  real(r8), intent(in) :: xaxis
17 
18  real(r8), dimension(nr), intent(out) :: xell, xtria
19 
20  real(r8), dimension(4 * (np - 1)) :: rrt
21  real(r8), dimension(8 * (np - 1)) :: rr
22  real(r8) :: r, s
23  real(r8) :: x, y, ps
24  real(r8) :: xrad, xshift, xgeo
25  real(r8) :: xleft, xright
26  real(r8) :: xtop, ytop, dummy
27  real(r8) :: ell, tri
28  real(r8), parameter :: tol = 1._r8 + 1.e-8_r8
29  integer(itm_i4) :: i, j, ngs
30  integer(itm_i4) :: n1, n2, n3, n4
31  integer(itm_i4) :: npt
32  integer(itm_i4) :: ifail
33 
34  if (standard_output) then
35  write(out_he, *) '**************************************************'
36  write(out_he, *) '* ellipticity and triangularity (fourier co.) *'
37  write(out_he, *) '**************************************************'
38  write(out_he, *) '* index, s, radius, shift, ellip, tria *'
39  write(out_he, *) '**************************************************'
40  end if
41  do i = 1, nr - 1
42  do j = 1, np - 1
43  n1 = (i - 1) * np + j
44  n2 = n1 + 1
45  n3 = n2 + np
46  n4 = n1 + np
47  do ngs = 1, 4
48  r = -1._r8
49  s = -1._r8 + 0.5_r8 * dble(ngs - 1)
50  call interpolation(0, xx(1, n1), xx(1, n2), xx(1, n3), xx(1, n4), &
51  r, s, x)
52  call interpolation(0, yy(1, n1), yy(1, n2), yy(1, n3), yy(1, n4), &
53  r, s, y)
54  call interpolation(0, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) + 1), &
55  psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, ps)
56  rrt(4 * (j - 1) + ngs) = sqrt((x - xaxis)**2 + y * y)
57  end do
58  end do
59  rrt(4 * (np - 1) + 1) = sqrt((xx(1, i * np) - xaxis)**2)
60  do j = 1, 4 * (np - 1) + 1
61  rr(j) = rrt(4 * (np - 1) - j + 2)
62  end do
63  do j = 1, 4 * (np - 1) - 1
64  rr(4 * (np - 1) + j + 1) = rr(4 * (np - 1) - j + 1)
65  end do
66  npt = 8 * (np - 1)
67  call rft2(rr, npt, 1)
68  xrad = rr(1) / dble(npt)
69  xshift = -2._r8 * rr(3) / dble(npt) / sqrt(ps)
70  xell(i) = -4._r8 * rr(5) / dble(npt) / sqrt(ps)
71  xtria(i) = 8._r8 * rr(7) / dble(npt) / ps
72  if (standard_output) &
73  write(out_he, 11) i, sqrt(ps), xrad, xshift, xell(i), xtria(i)
74  end do
75 
76 !-----------------------------------------------------------------
77 ! alternative to calculate the ellipticity and triangularity by
78 ! finding the point dy/ds=0 on each flux surface.
79 !-----------------------------------------------------------------
80  if (standard_output) then
81  write(out_he, *) '**************************************************'
82  write(out_he, *) '* ellipticity and triangularity (geometric co.) *'
83  write(out_he, *) '**************************************************'
84  end if
85  do i = 1, nr - 1
86  do j = 1, np - 1
87  n1 = (i - 1) * np + j
88  n2 = (i - 1) * np + j + 1
89  if (yy(3, n1) * yy(3, n2) <= 0._r8) then
90 !------------------------------------- quad. eq for s value at minimum -
91  call quad_equation(yy(1, n1), yy(3, n1), yy(1, n2), yy(3, n2), &
92  tol, s, ifail)
93  call cubic_interp_1d(xx(1, n1), xx(3, n1), xx(1, n2), &
94  xx(3, n2), s, xtop, dummy)
95  call cubic_interp_1d(yy(1, n1), yy(3, n1), yy(1, n2), &
96  yy(3, n2), s, ytop, dummy)
97  xleft = xx(1, (i - 1) * np + 1)
98  xright = xx(1, (i - 1) * np + np)
99  xgeo = (xleft + xright) / 2._r8
100  xrad = (xright - xleft) / 2._r8
101  xshift = (xaxis - xgeo) / xrad
102  ell = ytop / xrad - 1._r8
103  tri = -(xtop - xgeo) / xrad**2
104  ps = psi(4 * (n1 - 1) + 1)
105  if (standard_output) &
106  write(out_he, 11) i, sqrt(ps), xrad, xshift, ell, tri
107  end if
108  end do
109  end do
110 
111  11 format(i3, 5e16.8)
112 
113  return
114 end subroutine triangularity
subroutine interpolation(type_interp, xn1, xn2, xn3, xn4, r, s, x, xr, xs, xrs, xrr, xss, yn1, yn2, yn3, yn4, pn1, pn2, pn3, pn4, yr, ys, ps)
subroutine rft2(data, nr, kr)
Definition: rft2.f90:1
subroutine cubic_interp_1d(x1, x1s, x2, x2s, s, x, xs)
subroutine, public quad_equation(p_m, p_mr, p_p, p_pr, tol, r, ifail)
subroutine triangularity(xaxis, xell, xtria)