ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
flxint.f90
Go to the documentation of this file.
1 subroutine flxint(xaxis, a, nmg, fscale)
2 !----------------------------------------------------------------------
3 ! subroutine to evaluate new fdf profile
4 !----------------------------------------------------------------------
5 
6  use itm_types
7  use mod_parameter
8  use mod_corners
9  use mod_dat
10  use mod_gauss
11  use mod_mesh
12  use mod_profiles
14  use mod_spline
15  use current_calculations, only : j0, cscale
16  use phys_profiles, only : j_phi
17  use phys_functions
20  use helena_spline
21 
22  implicit none
23 
24  real(r8), intent(in) :: a, xaxis
25  real(r8), intent(inout) :: fscale
26  integer(itm_i4), intent(in) :: nmg
27 
28  type(spline_coefficients) :: fdf_spline
29 
30  real(r8), dimension(nr) :: psi_new, fdf_new
31  real(r8), dimension(3) :: abltg
32  real(r8) :: factas, cur0
33  real(r8) :: r, s, ws
34  real(r8) :: ps, psr, pss, psrs, psrr, psss
35  integer(itm_i4) :: i, nelm
36  integer(itm_i4) :: n1, n2, n3, n4
37  integer(itm_i4) :: error
38 
39  real(r8) :: dp_dpsi
40 
41  real(r8) :: pnorm
42 
43  if (ias == 1) then
44  factas = 1._r8
45  else
46  factas = 2._r8
47  end if
48 
49 !TODO: adapt B such that cur0 = - j_tor(1) / (eps * A * j0) ???
50  if (nmg > 0) then
51  j0 = bvac / (mu0 * eps * rvac * alfa)
52  pnorm = eps * bvac**2 / (mu0 * alfa**2)
53  cur0 = -cscale / (eps * a * j0)
54 ! write(*, *) 'cur0 = ', cur0
55 ! write(*, *) 'old b = ', b
56 ! b = -pscale * sign(1._r8, cpsurfin) / (eps * pnorm * a)
57 ! if (hbt) then
58 ! b = (cur0 * (1._r8 + eps * xaxis) - 1._r8) / (xaxis * (1._r8 &
59 ! + eps * xaxis / 2._r8))
60 ! else
61 ! b = (cur0 * (1._r8 + eps * xaxis) - 1._r8) / (1._r8 + eps &
62 ! * xaxis)**2
63 ! end if
64  end if
65 ! write(*, *) 'new b = ', b
66 
67  if (hbt) then
68  cur0 = (1._r8 + b * xaxis * (1._r8 + eps * xaxis / 2._r8)) &
69  / (1._r8 + eps * xaxis)
70  else
71  cur0 = (1._r8 + b * (1._r8 + eps * xaxis)**2) / (1._r8 + eps &
72  * xaxis)
73  end if
74 
75  r = -1._r8
76  s = xgauss(1)
77  ws = wgauss(1)
78 
79  do i = 1, nr - 1
80 !-- first element along flux surface sufficient to determine ps
81  nelm = (i - 1) * (np - 1) + 1
82  call set_node_number(nelm, n1, n2, n3, n4)
83  call interpolation(2, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) + 1), &
84  psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, ps, psr, pss, &
85  psrs, psrr, psss)
86 !-- elements in HELENA are ordered in inverse radial direction, i.e.,
87 ! psikn(nr) = 0, psikn(1) = 1
88  psi_new(nr + 1 - i) = ps
89  if (hbt) then
90  fdf_new(nr + 1 - i) = (cur0 * j_phi(ps) - b * dp_dpsi(ps) &
91  * flux_surface_average(nr + 1 - i, xaxis, 1._r8, 0._r8, &
92  current_averaging, rquad_minus_1_over_2_r)) &
93  / flux_surface_average(nr + 1 - i, xaxis, 1._r8, 0._r8, &
94  current_averaging, one_over_r)
95  else
96  fdf_new(nr + 1 - i) = (cur0 * j_phi(ps) - b * dp_dpsi(ps) &
97  * flux_surface_average(nr + 1 - i, xaxis, 1._r8, 0._r8, &
98  current_averaging, r_major)) &
99  / flux_surface_average(nr + 1 - i, xaxis, 1._r8, 0._r8, &
100  current_averaging, one_over_r)
101  end if
102  end do
103 
104 !-- extrapolate fdf_new on axis
105  psi_new(1) = 0_r8
106  fdf_new(1) = (psi_new(3) - psi_new(1)) / (psi_new(3) - psi_new(2)) &
107  * fdf_new(2) - (psi_new(2) - psi_new(1)) / (psi_new(3) - psi_new(2)) &
108  * fdf_new(3)
109 !-- renomalize fdf_new to 1 on axis
110  fdf_new = fdf_new / fdf_new(1)
111 
112 !-- calculate new equidistant gamma profile fdf based on fdf_new
113  call allocate_spline_coefficients(fdf_spline, nr)
114  call spline(nr, psi_new, fdf_new, 0._r8, 0._r8, 2, fdf_spline)
115  write(*, *) 'before evaluate'
116  call evaluate_spline(fdf_spline, psi_new, psivec, fdf, error)
117  write(*, *) 'after evaluate'
118  call deallocate_spline_coefficients(fdf_spline)
119 
120  return
121 end subroutine flxint
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 flxint(xaxis, a, nmg, fscale)
Definition: flxint.f90:1
subroutine spline(N, X, Y, ALFA, BETA, TYP, A, B, C, D)
Definition: solution6.f90:655
real(r8) function, public flux_surface_average(i, xaxis, a, F_dia, type, func)
subroutine deallocate_spline_coefficients(spline)
real(r8) function rquad_minus_1_over_2_r(a, x, xr, xs, yr, ys, psi, psir, F_dia)
real(r8) function one_over_r(a, x, xr, xs, yr, ys, psi, psir, F_dia)
real(r8) function r_major(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine set_node_number(n_node, n1, n2, n3, n4)
real(r8) function dp_dpsi(flux)
Definition: dp_dpsi.f90:1
real(r8) function j_phi(flux)
subroutine evaluate_spline(p, x, x_new, f_new, error, fprime_new)
Definition: mod_spline.f90:12