ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
check_solution.f90
Go to the documentation of this file.
1 subroutine check_solution(equilibrium_out)
2 !-----------------------------------------------------------------------
3 ! - checks the error in the grad-shafranov equation using the
4 ! straight field line system (not called by Helena)
5 !-----------------------------------------------------------------------
6 
7  use itm_types
8  use mod_map
9  use mod_mesh, only : nr
10  use mod_parameter
11  use mod_profiles
12  use mod_spline
13  use mod_helena_io, only : out_he
14  use mod_output
15  use helena_spline
16 
17  use euitm_schemas
18 
19  implicit none
20 
21  type(type_equilibrium), intent(in) :: equilibrium_out
22 
23  type(spline_coefficients) :: gp_spline
24  integer(itm_i4) :: i, j, im, jm
25  real(r8) :: maxerr, err
26  real(r8) :: dgp0, dgp1, dgp
27  real(r8) :: sps2
28  real(r8) :: dp_axis, dp_boundary
29  real(r8) :: drbphi_axis, drbphi_boundary
30  real(r8), dimension(npsi - 1) :: c1
31  real(r8), dimension(3) :: dummy
32  real(r8), dimension(nchi) :: dgt
33  real(r8), dimension(npsi) :: gp
34  real(r8), dimension(npsi - 1, nchi) :: error, error2
35 ! real :: ff(nchi), fp(nchi), zc(10)
36 
37 !-- derivatives at plasma boundaries (all zero except for dp_boundary)
38  dp_axis = 0._r8
39  drbphi_axis = 0._r8
40  drbphi_boundary = 0._r8
41  dp_boundary = equilibrium_out%profiles_1d%pprime(nr) &
42  * 2 * equilibrium_out%profiles_1d%psi(nr) * mu0 &
43  / equilibrium_out%global_param%mag_axis%bphi**2
44 
45 !TODO: is this really supposed to change npsi?
46  npsi = nr
47 
48  do i = 1, nchi
49  c1(:) = g12_hel(2 : npsi - 1, i)
50  call spline(npsi - 1, cs(2), c1, 0._r8, 0._r8, 2, q_spline)
51  g12_hel(1, i) = spwert(npsi - 1, 0._r8, q_spline, cs(2), dummy, 0)
52  end do
53 
54  call spline(npsi, cs, p0, dp_axis, dp_boundary, 1, p_spline)
55  call spline(npsi, cs, rbphi, drbphi_axis, drbphi_boundary, 1, rbphi_spline)
56 
57  maxerr = -1.e20_r8
58 
59  call allocate_spline_coefficients(gp_spline, npsi)
60 
61  do j = 1, nchi
62  do i = 1, npsi
63  gp(i) = g11_hel(i, j) * qs(i) / rbphi(i)
64  end do
65  dgp0 = (gp(2) - gp(1)) / (cs(2) - cs(1))
66  dgp1 = (gp(npsi) - gp(npsi - 1)) / (cs(npsi) - cs(npsi - 1))
67  call spline(npsi, cs, gp, dgp0, dgp1, 2, gp_spline)
68 ! call spline(npsi, cs, gp, 0._r8, 0._r8, 2, gp_spline)
69 
70  do i = 2, npsi - 1
71  sps2 = 2._r8 * cs(i) * cpsurf
72  call theta_derivative(g12_hel(i, 1), dgt, nchi)
73  dgp = (g11_hel(i + 1, j) * qs(i + 1) / rbphi(i + 1) &
74  - g11_hel(i - 1, j) * qs(i - 1) / rbphi(i - 1)) &
75  / (cs(i + 1) - cs(i - 1))
76  err = abs(rbphi(i) / qs(i) * dgp + dgt(j) * sps2 &
77  + (rbphi(i) * rbphi_spline%sp2(i) + g33_hel(i, j) * p_spline%sp2(i)))
78  if (err > maxerr) then
79  maxerr = err
80  im = i
81  jm = j
82  end if
83 ! if (standard_output) &
84 ! write(out_he, 21) i, j, cs(i), chi(j), err, log(err)
85  error(i - 1, j) = log(err)
86  error2(i - 1, j) = err
87  end do
88  end do
89  if (standard_output) then
90  write(out_he, *)
91  write(out_he, 23) maxerr, im, jm
92  write(out_he, *)
93  end if
94 
95 ! if (standard_output) then
96 ! write(out_he, *)
97 ! write(out_he, *) ' error in gs : harmonics'
98 ! write(out_he, *)
99 ! do j = 1, nchi
100 ! fp(j) = float(j - 1)
101 ! end do
102 ! do i = 2, npsi
103 ! do j = 1, nchi
104 ! ff(j) = error2(i - 1, j)
105 ! end do
106 ! call rft2(ff, nchi, 1)
107 ! do m = 1, nchi / 2
108 ! ff(m) = 2 * ff(m) / real(nchi)
109 ! write(out_he, 22) i, m, cs(i), ff(m)
110 ! end do
111 ! end do
112 ! end if
113 
114  call deallocate_spline_coefficients(gp_spline)
115 
116  return
117 
118 ! 21 format(2i3, 2f8.3, 3e12.4)
119 ! 22 format(2i3, f8.3, 3e12.4)
120  23 format(' max error in gs after mapping : ', e12.4, 2i5)
121 
122 end subroutine check_solution
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
subroutine deallocate_spline_coefficients(spline)
subroutine theta_derivative(arrin, darr, nchi)
subroutine check_solution(equilibrium_out)