ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
safety_factor.f90
Go to the documentation of this file.
1 subroutine safety_factor(xaxis, a, cx, cy, qprf, q95out, q1out)
2 !----------------------------------------------------------------------
3 ! subroutine to evaluate the q-profile
4 !----------------------------------------------------------------------
5 
6  use itm_types
7  use mod_parameter
8  use mod_dat
9  use mod_gauss
10  use mod_mesh
11  use mod_profiles
13  use phys_functions
14  use mod_helena_io, only : iu6
15  use helena_spline
16 
17  implicit none
18 
19  type(spline_coefficients) :: dq_spline
20  real(r8), intent(in) :: xaxis, a, cx, cy
21  real(r8), dimension(npts), intent(out) :: qprf
22  real(r8), intent(out) :: q95out, q1out
23 
24  real(r8), dimension(nr) :: p0tmp, f2tmp, dptmp, dftmp
25  real(r8), dimension(nr) :: zps, qq, avc
26  real(r8), dimension(3) :: abltg
27 
28  integer(itm_i4) :: i, j, k
29  integer(itm_i4) :: n1, n2, n3, n4
30  real(r8) :: factas, sumq, avcur, xlength
31  real(r8) :: r, s, ws, x, xr, xs, y, yr, ys, ps, psr, pss
32  real(r8) :: xjac, bigr
33  real(r8) :: arhs
34  real(r8) :: ss
35 
36  if (ias == 1) then
37  factas = 1._r8
38  else
39  factas = 2._r8
40  endif
41  if (verbosity > 2) write(iu6, *) 'safety factor'
42  call profiles(p0tmp, f2tmp, dptmp, dftmp, a)
43 
44 !------------------------------------- nelm elements ------------------
45  r = -1._r8
46  do i = 1, nr - 1
47  sumq = 0._r8
48  avcur = 0._r8
49  xlength = 0._r8
50  do j = 1, np - 1
51  n1 = (i - 1) * np + j
52  n2 = n1 + 1
53  n3 = n2 + np
54  n4 = n1 + np
55 !------------------------------------- 4 point gaussian int. in s -----
56  do k = 1, 4
57  s = xgauss(k)
58  ws = wgauss(k)
59  call interpolation(1, xx(1, n1), xx(1, n2), xx(1, n3), xx(1, n4), &
60  r, s, x, xr, xs)
61  call interpolation(1, yy(1, n1), yy(1, n2), yy(1, n3), yy(1, n4), &
62  r, s, y, yr, ys)
63  call interpolation(1, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) + 1), &
64  psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, ps, &
65  psr, pss)
66  xjac = xr * ys - xs * yr
67  bigr = 1._r8 + eps * x
68  sumq = sumq - ws * xjac / (bigr * abs(psr))
69  xlength = xlength + sqrt(xs**2 + ys**2) * ws
70  arhs = rh_side(1._r8, x, 0._r8, 0._r8, 0._r8, 0._r8, ps, 0._r8, 0._r8)
71  avcur = avcur + ws * arhs * sqrt(xs**2 + ys**2)
72  end do
73  end do
74  zps(nr - i + 1) = ps
75  qq(nr - i + 1) = factas / 2._r8 * f2tmp(nr - i + 1) * sumq * alfa / pi
76  avc(nr - i + 1) = factas * a * eps * avcur / (alfa * xlength)
77  end do
78 !------------------------------------- calc equidistant q profile
79  qq(1) = f2tmp(1) * alfa / (2._r8 * sqrt(cx * cy) * (1._r8 + eps &
80  * xaxis))
81  zps(1) = 0._r8
82  call allocate_spline_coefficients(dq_spline, nr)
83  call spline(nr, zps, qq, 0._r8, 0._r8, 2, dq_spline)
84  do i = 1, npts
85  ss = dble(i - 1) / dble(npts - 1)
86  ps = ss * ss
87  qprf(i) = spwert(nr, ps, dq_spline, zps, abltg, 0)
88  end do
89  ps = 0.95_r8
90  q95out = spwert(nr, ps, dq_spline, zps, abltg, 0)
91  q1out = qq(nr)
92  if (verbosity > 2) write(iu6, *) ' q95 : ', q95out
93 !------------------------------------- calc equidistant current profile
94  avc(1) = avc(2) - (avc(3) - avc(2)) / (zps(3) - zps(2)) * zps(2)
95  call spline(nr, zps, avc, 0._r8, 0._r8, 2, dq_spline)
96  do i = 1, npts
97  ss = dble(i - 1) / dble(npts - 1)
98  ps = ss * ss
99  j_tor(i) = spwert(nr, ps, dq_spline, zps, abltg, 0)
100  end do
101 
102  call deallocate_spline_coefficients(dq_spline)
103 
104  return
105 end subroutine safety_factor
subroutine profiles(p0, rbphi, dp0, drbphi, a)
Definition: profiles.f90:1
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 allocate_spline_coefficients(spline, n)
subroutine spline(N, X, Y, ALFA, BETA, TYP, A, B, C, D)
Definition: solution6.f90:655
REAL *8 function spwert(N, XWERT, A, B, C, D, X, ABLTG)
Definition: solution6.f90:155
real(r8) function rh_side(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine deallocate_spline_coefficients(spline)
subroutine safety_factor(xaxis, a, cx, cy, qprf, q95out, q1out)