ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
betapol.f90
Go to the documentation of this file.
1 subroutine betapol(a, bpl)
2 !----------------------------------------------------------------------
3 ! subroutine to evaluate output quantities
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
13  use mod_profiles
14  use mod_spline
16  use phys_functions
17  use helena_spline
18 
19  implicit none
20 
21  real(r8), intent(in) :: a
22 
23  real(r8), intent(out) :: bpl
24 
25  real(r8), dimension(3) :: ablt
26  real(r8) :: factas
27  real(r8) :: parea, carea
28  real(r8) :: r, wr, s, ws
29  real(r8) :: x, xr, xs, y, yr, ys, ps, psr, pss
30  real(r8) :: xjac
31  real(r8) :: arhs
32  integer(itm_i4) :: i, j, ngr, ngs
33  integer(itm_i4) :: nelm
34  integer(itm_i4) :: n1, n2, n3, n4
35 
36  factas = 2._r8
37  if (ias == 1) factas = 1._r8
38 
39  parea = 0._r8
40  carea = 0._r8
41 
42 !-- spline p_int and gam_int with use of known derivatives
43  call allocate_spline_coefficients(p_spline, npts)
44  call spline(npts, psivec, p_int, dpres(1) * sign(1._r8, cpsurfin), &
45  dpres(npts) * sign(1._r8, cpsurfin), 1, p_spline)
46 
47 !------------------------------------- nelm elements ------------------
48  do i = nr - 1, 1, -1
49  do j = 1, np - 1
50  nelm = (i - 1) * (np - 1) + j
51  call set_node_number(nelm, n1, n2, n3, n4)
52 !------------------------------------- 4 point gaussian int. in r -----
53  do ngr = 1, 4
54  r = xgauss(ngr)
55  wr = wgauss(ngr)
56 !------------------------------------- 4 point gaussian int. in s -----
57  do ngs = 1, 4
58  s = xgauss(ngs)
59  ws = wgauss(ngs)
60  call interpolation(1, xx(1, n1), xx(1, n2), xx(1, n3), xx(1, n4), &
61  r, s, x, xr, xs)
62  call interpolation(1, yy(1, n1), yy(1, n2), yy(1, n3), yy(1, n4), &
63  r, s, y, yr, ys)
64  call interpolation(1, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) &
65  + 1), psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, ps, &
66  psr, pss)
67  xjac = xr * ys - xs * yr
68  parea = parea + wr * ws * xjac * spwert(npts, ps, p_spline, &
69  psivec, ablt, 0)
70  arhs = rh_side(a, x, 0._r8, 0._r8, 0._r8, 0._r8, ps, 0._r8, 0._r8)
71  carea = carea + wr * ws * xjac * arhs
72  end do
73  end do
74  end do
75  end do
76  if (hbt) then
77  parea = 0.5_r8 * factas * a * b * parea
78  carea = factas * a * carea
79  else
80  parea = factas * eps * a * b * parea
81  carea = factas * a * carea
82  end if
83  bpl = -(8._r8 * pi / eps) * parea / carea**2
84 
85  call deallocate_spline_coefficients(p_spline)
86 
87  return
88 end subroutine betapol
subroutine betapol(a, bpl)
Definition: betapol.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 set_node_number(n_node, n1, n2, n3, n4)