ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
current_calculations.f90
Go to the documentation of this file.
2 
3  use itm_types
4  use mod_helena_io
5 
6  implicit none
7 
8  real(r8), public :: j0 ! constant if alfa constant
9  real(r8), public :: cscale ! used for calculation of cur0
10 
11  private
12 
13  real(r8), allocatable :: j_tor_hat(:)
14  real(r8), allocatable :: j_tor_hat_av(:), I_tor_hat(:)
15 
16  real(r8), allocatable, public :: j_tor_loc(:), j_tor_av(:), I_tor(:)
17 
18  public :: current_densities
21  public :: total_current, si_currents
22  public :: q_calculation
23 
24  contains
25 
26 !--------------------------------------------------------------------------
28 !--------------------------------------------------------------------------
29 
30  implicit none
31 
32  integer(itm_i4), intent(in) :: nr, np
33 
34  allocate(j_tor_hat((nr - 1) * (np - 1)))
35  allocate(j_tor_hat_av(nr))
36  allocate(i_tor_hat(nr))
37  allocate(j_tor_loc((nr - 1) * (np - 1)))
38  allocate(j_tor_av(nr))
39  allocate(i_tor(nr))
40 
41 !-- initialize I_tor_hat
42  i_tor_hat = 0._r8
43 
45 
46 !--------------------------------------------------------------------------
48 !--------------------------------------------------------------------------
49 
50  implicit none
51 
52  deallocate(j_tor_hat, j_tor_hat_av)
53  deallocate(j_tor_loc, j_tor_av)
54  deallocate(i_tor_hat, i_tor)
55 
56  end subroutine current_calculations_destructor
57 
58 !--------------------------------------------------------------------------
59  subroutine current_densities(a, xaxis, type)
60 !--------------------------------------------------------------------------
61 ! Calculates averaged normalized current densities
62 ! The averaging process depends on the type:
63 ! type = 'line ' : line average in poloidal cross section
64 ! = 'area ' : surface average over poloidal cross section between
65 ! two flux surfaces
66 ! = 'volume' : volume average between two flux surfaces
67 ! = 'theta ' : average over surface area of a flux surface
68 !--------------------------------------------------------------------------
69 
70  use mod_dat, only : eps
71  use mod_mesh, only : nr, np
73  use phys_functions, only : rh_side, identity
74 
75  implicit none
76 
77  real(r8), intent(in) :: a, xaxis
78  character(len = *), intent(in) :: type
79 
80  real(r8) :: dl(np - 1), df(np - 1)
81  integer(itm_i4) :: i, j, n
82 
83 !-- initialize dl and df
84  dl = 0._r8
85  df = 0._r8
86 
87 !-- current density calculation excluding magnetic axis (i = 1)
88  n = 0
89  do i = 2, nr
90  call element_average(i, type, xaxis, a, 0._r8, rh_side, df, dl)
91  do j = 1, np - 1
92  n = n + 1
93 !-- normalized current density
94  j_tor_hat(n) = -eps * df(j)
95  end do
96  end do
97 !-- flux surface averaged current density (including magnetic axis, i = 1)
98  do i = 1, nr
99  j_tor_hat_av(i) = -eps * flux_surface_average(i, xaxis, a, 0._r8, &
100  type, rh_side)
101  end do
102 
103  return
104  end subroutine current_densities
105 
106 !--------------------------------------------------------------------------
107  subroutine total_current(a, xaxis, toroidal_current)
108 !--------------------------------------------------------------------------
109 ! calculates the normalized current profile I_tor_hat(i) but returns the
110 ! total plasma toroidal current in SI units in toroidal_current
111 ! Attention: recalculates the current densities for averaging type 'area '
112 !--------------------------------------------------------------------------
113 
114  use mod_dat, only : rvac, bvac, eps
115  use mod_mesh, only : nr, np, ias
116  use mod_helena_io, only : out_he
117  use mod_output
118  use mod_parameter
119  use math_functions, only : element_average
120  use phys_functions, only : identity
121 
122  implicit none
123 
124  character(len = 6), parameter :: type = 'area '
125 
126  real(r8), intent(in) :: a, xaxis
127  real(r8), intent(out) :: toroidal_current
128 
129  real(r8) :: df(np - 1), da(np - 1)
130  real(r8) :: da_sum, i_sum, factas
131  integer(itm_i4) :: i, j
132 
133 !-- initialize dA and df
134  da = 0._r8
135  df = 0._r8
136 
137  call current_densities(a, xaxis, type)
138 
139  factas = 2._r8
140  if (ias == 1) factas = 1._r8
141 
142  i_sum = 0._r8
143  do i = 1, nr
144  da_sum = 0._r8
145  call element_average(i, type, xaxis, 1._r8, 0._r8, identity, df, da)
146  do j = 1, np -1
147  da_sum = da_sum + da(j)
148  end do
149  i_sum = i_sum + factas * j_tor_hat_av(i) * da_sum
150  i_tor_hat(i) = i_sum
151  end do
152 
153  toroidal_current = i_tor_hat(nr) * rvac * bvac * eps / (4.e-7_r8 * pi)
154 
155  if (standard_output) &
156  write(out_he, 4) toroidal_current
157  4 format(' total current (alfa = 1) : ', e12.4)
158 
159  return
160  end subroutine total_current
161 
162 !--------------------------------------------------------------------------
163  subroutine si_currents
164 !--------------------------------------------------------------------------
165 ! denormalizes the currents and current densities
166 !--------------------------------------------------------------------------
167 
168  use mod_parameter
169  use mod_dat
170 
171  implicit none
172 
173  real(r8) :: j0, i0
174 
175 !-- normalization constants
176  j0 = bvac / (mu0 * eps * rvac * alfa)
177  i0 = j0 * (eps * rvac)**2
178 
179 !-- denormalization
180  j_tor_loc = j0 * j_tor_hat
181  j_tor_av = j0 * j_tor_hat_av
182 
183  i_tor = i_tor_hat * i0
184 
185  end subroutine si_currents
186 
187 !--------------------------------------------------------------------------
188  subroutine q_calculation(xaxis, cx, cy, a, q)
189 !--------------------------------------------------------------------------
190 ! calculated the q profile
191 ! requires flux surface aligned coordinates, i.e. prior call to remesh
192 !--------------------------------------------------------------------------
193 
194  use mod_parameter
195  use mod_dat
196  use mod_gauss
198  use mod_map
199  use mod_mesh
200 
201  implicit none
202 
203  real(r8), intent(in) :: xaxis, cx, cy
204  real(r8), intent(in) :: a
205  real(r8), intent(inout) :: q(nr)
206  real(r8) :: factas
207  real(r8) :: r, s
208  real(r8) :: sumq
209  real(r8) :: x, xr, xs, xrs, xrr, xss
210  real(r8) :: y, yr, ys, yrs, yrr, yss
211  real(r8) :: ps, psr, pss, psrs, psrr, psss
212  real(r8) :: jacobian, bigr, dl, grad_psi2
213  integer(itm_i4) :: i, j, ngs
214  integer(itm_i4) :: n1, n2, n3, n4
215 
216 ! profiles yields the pressure normalized to eps * Bvac^2 / (mu_0 * alfa^2)
217 ! and R Bphi normalized to Rvac * Bvac
218  call profiles(p0, rbphi, dp0, drbphi, a)
219 
220  factas = 2._r8
221  if (ias == 1) factas = 1._r8
222 
223 !------------------------------------- nelm elements ------------------
224  r = -1._r8
225  do i = 1, nr - 1
226  sumq = 0._r8
227  do j = 1, np - 1
228  n1 = (i - 1) * np + j
229  n2 = n1 + 1
230  n3 = n2 + np
231  n4 = n1 + np
232 !------------------------------------- 4 point gaussian int. in s -----
233  do ngs = 1, 4
234  s = xgauss(ngs)
235  call interpolation(2, xx(1, n1), xx(1, n2), xx(1, n3), &
236  xx(1, n4), r, s, x, xr, xs, xrs, xrr, xss)
237  call interpolation(2, yy(1, n1), yy(1, n2), yy(1, n3), &
238  yy(1, n4), r, s, y, yr, ys, yrs, yrr, yss)
239  call interpolation(2, psi(4 * (n1 - 1) + 1), &
240  psi(4 * (n2 - 1) + 1), psi(4 * (n3 - 1) + 1), &
241  psi(4 * (n4 - 1) + 1), r, s, ps, psr, pss, psrs, psrr, psss)
242  jacobian = xr * ys - xs * yr
243  bigr = 1._r8 + eps * x
244  dl = sqrt(xs**2 + ys**2)
245  grad_psi2= psr**2 * (xs**2 + ys**2) / jacobian**2
246  sumq = sumq + wgauss(ngs) / (bigr * sqrt(grad_psi2)) * dl
247  end do
248  end do
249  q(nr - i + 1) = 0.5_r8 * factas * sumq * rbphi(nr - i + 1)
250  end do
251  q(1) = rbphi(1) * pi / (2._r8 * sqrt(cx * cy) * (1. + eps * xaxis))
252  do i = 1, nr
253  q(i) = q(i) * alfa / pi
254  end do
255 
256  end subroutine q_calculation
257 
258 end module current_calculations
subroutine profiles(p0, rbphi, dp0, drbphi, a)
Definition: profiles.f90:1
subroutine, public q_calculation(xaxis, cx, cy, a, q)
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)
real(r8) function rh_side(a, x, xr, xs, yr, ys, psi, psir, F_dia)
real(r8) function, public flux_surface_average(i, xaxis, a, F_dia, type, func)
subroutine, public si_currents
subroutine, public current_densities(a, xaxis, type)
subroutine, public element_average(i, type, xaxis, a, F_dia, func, df, dl)
subroutine, public total_current(a, xaxis, toroidal_current)
real(r8) function identity(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine, public current_calculations_destructor
subroutine, public current_calculations_constructor(nr, np)