ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
profiles.f90
Go to the documentation of this file.
1 subroutine profiles(p0, rbphi, dp0, drbphi, a)
2 !-----------------------------------------------------------------------
3 ! subroutine to integrate the equilibrium profiles
4 ! using the hbt definitions
5 !-----------------------------------------------------------------------
6 
7  use itm_types
8  use mod_dat
9  use mod_mesh
10  use mod_nodes
11  use mod_helena_io, only : iu6
12  use mod_profiles
13  use mod_spline
14  use helena_spline
15 
16  implicit none
17 
18  integer(itm_i4) :: j
19 
20  real(r8), dimension(nr), intent(out) :: p0, rbphi, dp0, drbphi
21  real(r8), intent(in) :: a
22  real(r8) :: flux
23 
24  real(r8), dimension(3) :: ablt
25 
26  real(r8) :: dp_dpsi, dgamma_dpsi
27 
28 !-- spline p_int and gam_int with use of known derivatives
29  call allocate_spline_coefficients(p_spline, npts)
30  call allocate_spline_coefficients(rbphi_spline, npts)
31 
32 !-- sign of cpsurfin needed because psivec always positive
33  call spline(npts, psivec, p_int, dpres(1) * sign(1._r8, cpsurfin), &
34  dpres(npts) * sign(1._r8, cpsurfin), 1, p_spline)
35  call spline(npts, psivec, gam_int, dgam(1) * sign(1._r8, cpsurfin), &
36  dgam(npts) * sign(1._r8, cpsurfin), 1, rbphi_spline)
37 
38 !---------------------------- derivatives to r = sqrt(psi) !!!! --------
39 ! multiplication by sign(-alfa) * sign(bvac) compensates for loss of sign
40 ! due to normalization with pscale and fscale
41  do j = 1, nr
42  flux = psikn(nr - j + 1)
43  if (hbt) then
44  p0(j) = 0.5_r8 * a * b * spwert(npts, flux, p_spline, psivec, ablt, 0) &
45  * sign(1.0_r8, -alfa) * sign(1.0_r8, bvac)
46  dp0(j) = -a * b * dp_dpsi(flux) * sqrt(flux)
47  rbphi(j) = p0(j) + eps * a * spwert(npts, flux, rbphi_spline, psivec, &
48  ablt, 0) * sign(1.0_r8, -alfa) * sign(1.0_r8, bvac)
49  drbphi(j) = dp0(j) - eps * 2._r8 * a * dgamma_dpsi(flux) &
50  * sqrt(flux)
51  else
52  p0(j) = eps * a * b * spwert(npts, flux, p_spline, psivec, ablt, 0) &
53  * sign(1.0_r8, -alfa) * sign(1.0_r8, bvac)
54  dp0(j) = -2._r8 * eps * a * b * dp_dpsi(flux) * sqrt(flux)
55  rbphi(j) = eps * a * spwert(npts, flux, rbphi_spline, psivec, ablt, 0) &
56  * sign(1.0_r8, -alfa) * sign(1.0_r8, bvac)
57  drbphi(j) = -eps * 2._r8 * a * dgamma_dpsi(flux) * sqrt(flux)
58  end if
59 !-- rbphi still missing sign(bvac), added in fill_equilibrium
60  rbphi(j) = sqrt(1._r8 + 2._r8 * eps * rbphi(j) / alfa**2)
61  drbphi(j) = 1._r8 / (2._r8 * rbphi(j)) * 2._r8 * eps * drbphi(j) / alfa**2
62  end do
63 
64  if (verbosity > 2) write(iu6, *) 'profiles computed with A = ', a
65 
66  call deallocate_spline_coefficients(p_spline)
67  call deallocate_spline_coefficients(rbphi_spline)
68 
69  return
70 end subroutine profiles
subroutine profiles(p0, rbphi, dp0, drbphi, a)
Definition: profiles.f90:1
subroutine allocate_spline_coefficients(spline, n)
subroutine spline(N, X, Y, ALFA, BETA, TYP, A, B, C, D)
Definition: solution6.f90:655
subroutine flux(psitok, rk, zk, nk)
Definition: EQ1_m.f:786
REAL *8 function spwert(N, XWERT, A, B, C, D, X, ABLTG)
Definition: solution6.f90:155
real(r8) function dgamma_dpsi(flux)
Definition: dgamma_dpsi.f90:1
subroutine deallocate_spline_coefficients(spline)
real(r8) function dp_dpsi(flux)
Definition: dp_dpsi.f90:1