ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
solver.f90
Go to the documentation of this file.
1 module solver
2 
3  use itm_types
4  use mod_corners
5  use mod_mesh
6  use mod_helena_io, only : out_he, iu6
7  use mod_output, only : standard_output
8 
9  implicit none
10 
11  real(r8), allocatable, private :: qq(:)
12  real(r8), allocatable, private :: gpx(:, :, :), gpjac(:, :, :)
13  integer(itm_i4), allocatable, private :: indices(:), inv_index(:)
14 
15 
16  contains
17 
18  subroutine allocate_solver(nr, np)
19 
20  use itm_types
21 
22  implicit none
23 
24  integer(itm_i4) :: nr, np
25 
26  allocate(qq(4 * nr * np))
27  allocate(gpx(4, 4, (nr - 1) * (np - 1)))
28  allocate(gpjac(4, 4, (nr - 1) * (np - 1)))
29  allocate(indices(nr * np))
30  allocate(inv_index(nr * np))
31 
32 !-- initialize
33  qq = 0._r8
34  gpx = 0._r8
35  gpjac = 0._r8
36  indices = 0
37  inv_index = 0
38 
39  return
40  end subroutine allocate_solver
41 
42  subroutine deallocate_solver
43 
44  implicit none
45 
46  deallocate(qq)
47  deallocate(gpx, gpjac)
48  deallocate(indices, inv_index)
49 
50  end subroutine deallocate_solver
51 
52  subroutine find_axis(psaxis, xaxis, yaxis, nax, rax, sax)
53 !-----------------------------------------------------------------------
54 ! subroutine to localize the position of the magnetic axis ; the minimum
55 ! of psi of all elements
56 !-----------------------------------------------------------------------
57 
60 
61  implicit none
62 
63  real(r8), intent(out) :: psaxis
64  real(r8), intent(inout) :: xaxis
65  real(r8), intent(inout) :: yaxis
66  real(r8), intent(inout) :: rax, sax
67  integer(itm_i4), intent(out) :: nax
68 
69  real(r8), dimension(2) :: x
70  real(r8) :: psmin
71  real(r8) :: x_err, f_err
72  real(r8), parameter :: tol_x=1.e-8_r8, tol_f=1.e-8_r8
73  real(r8), parameter :: tol_rs=1.0001_r8
74  real(r8) :: r, s, zpsi
75  real(r8) :: psi_00, psi_00_min
76  integer(itm_i4) :: i, n_start
77  integer(itm_i4) :: naxis
78  integer(itm_i4) :: n1, n2, n3, n4
79  integer(itm_i4) :: nelm
80  integer(itm_i4) :: neq2, ntrial, ifail
81  integer(itm_i4) :: naxis_00
82  logical :: minimum_found
83 
84  minimum_found = .false.
85 
86  psmin = 1.e20_r8
87  psi_00_min = psmin
88  naxis_00 = -1
89 
90  if (ias == 1) then
91 !----------------------------------- asymmetric part -----------
92  nelm = (nr - 1) * (np - 1)
93  neq2 = 2
94  ntrial = 150
95  do i = nelm, nelm / 2, -1
96  naxis = i
97  ifail = 1
98  x(1) = 0._r8
99  x(2) = 0._r8
100  call set_node_number(naxis, n1, n2, n3, n4)
101  call interpolation(0, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) + 1), &
102  psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), x(1), x(2), psi_00)
103  if (psi_00 .lt. psi_00_min) then
104  psi_00_min = psi_00
105  naxis_00 = naxis
106  end if
107  call mnewtax(ntrial, naxis, x, neq2, tol_x, tol_f, x_err, f_err)
108  if (x_err <= tol_x .or. f_err <= tol_f) then
109  ifail = 0
110  else
111 ! if (standard_output) &
112 ! write(out_he, *) ' accuracy not reached : ', x_err, f_err, naxis
113  end if
114  r = x(1)
115  s = x(2)
116  if (ifail == 0 .and. abs(r) <= tol_rs .and. abs(s) <= tol_rs) then
117  minimum_found = .true.
118  call set_node_number(naxis, n1, n2, n3, n4)
119  call interpolation(0, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) &
120  + 1), psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, zpsi)
121  if (zpsi < psmin) then
122  psmin = zpsi
123  nax = naxis
124  rax = r
125  sax = s
126  end if
127  end if
128  end do
129  if (naxis_00 .gt. 0) then
130 ! write(*,*) 'naxis_00, psi_00_min ', naxis_00, psi_00_min
131  end if
132  if (.not. minimum_found) then
133  write(iu6, *) 'fatal error: no axis found'
134  stop
135  end if
136  call set_node_number(nax, n1, n2, n3, n4)
137  call interpolation(0, xx(1, n1), xx(1, n2), xx(1, n3), xx(1, n4), &
138  rax, sax, xaxis)
139  call interpolation(0, yy(1, n1), yy(1, n2), yy(1, n3), yy(1, n4), &
140  rax, sax, yaxis)
141  psaxis = psmin
142  else
143 !**************************************** symmetric part ***************
144  n_start = 1
145  call axis_symm(psaxis, xaxis, yaxis, nax, rax, sax, n_start)
146  n_start = np - 1
147  call axis_symm(psaxis, xaxis, yaxis, nax, rax, sax, n_start)
148  yaxis = 0._r8
149  end if
150 
151  return
152  end subroutine find_axis
153 
154  subroutine axis_symm(psaxis, xaxis, yaxis, nax, rax, sax, n_start)
155 
156  use math_functions
157 
158  implicit none
159 
160  real(r8), intent(inout) :: psaxis
161  real(r8), intent(inout) :: xaxis
162  real(r8), intent(inout) :: yaxis
163  real(r8), intent(inout) :: rax, sax
164  integer(itm_i4), intent(inout) :: nax
165  integer(itm_i4), intent(in) :: n_start
166 
167  real(r8) :: r, dummy, psmin
168  real(r8), parameter :: tol = 1._r8 + 1.e-8_r8
169  integer(itm_i4) :: n, n1, n2, ifail
170 
171  psmin = 1.e20_r8
172 
173  do n = n_start, (nr - 1) * (np - 1), np - 1
174  if (n_start == 1) then
175  n1 = nodeno(n, 1)
176  n2 = nodeno(n, 4)
177  else
178  n1 = nodeno(n, 2)
179  n2 = nodeno(n, 3)
180  end if
181  if (psi(4 * (n1 - 1) + 2) * psi(4 * (n2 - 1) + 2) <= 0.) then
182 !------------------------------------- quad. eq for r value at minimum -
183  call quad_equation(psi(4 * (n1 - 1) + 1), psi(4 * (n1 - 1) + 2), &
184  psi(4 * (n2 - 1) + 1), psi(4 * (n2 - 1) + 2), tol, r, ifail)
185 !-------- the sign of r changes for elements on the left (see remesh) -
186  call cubic_interp_1d(xx(1, n1), xx(2, n1), xx(1, n2), &
187  xx(2, n2), r, xaxis, dummy)
188  call cubic_interp_1d(psi(4 * (n1 - 1) + 1), psi(4 * (n1 - 1) + 2), &
189  psi(4 * (n2 - 1) + 1), psi(4 * (n2 - 1) + 2),r, psaxis, dummy)
190  if (psaxis < psmin) then
191  psmin = psaxis
192  nax = n
193  rax = r
194  if (n_start == 1) then
195  sax = -1._r8
196  else
197  sax = 1._r8
198  end if
199  if (standard_output) &
200  write(out_he, *) 'local minimum at xaxis = ', xaxis
201  end if
202  end if
203  end do
204 
205  return
206  end subroutine axis_symm
207 
208  subroutine build_matrix(a, first)
209 !----------------------------------------------------------------------
210 ! subroutine to calculate the matrix kk and the right hand side array qq
211 ! no boundary conditions are yet used.
212 ! number of rows and columns : 4 * nr * np
213 ! nr : number of radial nodes
214 ! np : number of poloidal nodes
215 ! a : the total amplitude of the rhs (hbt definition)
216 ! b : measure of the total pressure (hbt definition)
217 ! eps : the inverse aspect ratio
218 ! hbt = .true. : hbt definition of gamma, .false. : ff' as gamma.
219 ! gpx, gpjac, indices used as static variables !
220 !----------------------------------------------------------------------
221 
222  use itm_types
223  use mod_corners
224  use mod_dat
225  use mod_gauss
226  use mod_mesh
227  use mod_solver
230  use phys_functions
231 
232  implicit none
233 
234  real(r8), intent(in) :: a
235  logical, intent(in) :: first
236 
237  real(r8) :: r, wr, s, ws, wrs
238  real(r8) :: x, xr, xs, yr, ys, ps
239  real(r8) :: xjac
240  real(r8) :: arhs
241  real(r8) :: sumq, sumk
242  real(r8) :: vx, vy
243  real(r8) :: psix, psiy
244  integer(itm_i4) :: i, j, k, l, n, ngr, ngs
245  integer(itm_i4) :: nelm
246  integer(itm_i4) :: n1, n2, n3, n4
247  integer(itm_i4) :: nrow, ncol, noff, nend
248  logical, parameter :: invert_index = .false.
249 
250  nelm = (np - 1) * (nr - 1)
251 
252 !-- initialize qq to zero
253  qq = 0._r8
254 
255  if (first) then
256  call set_index(indices, invert_index)
257 !------------------------------------- nelm elements ------------------
258  do i = 1, 4 * nr * np
259  do j = 1, kklda
260  kkbig(j, i) = 0._r8
261  end do
262  end do
263  do n = 1, nelm
264  call set_node_number(n, n1, n2, n3, n4)
265 !------------------------------------- 4 point gaussian int. in r -----
266  do ngr = 1, 4
267  r = xgauss(ngr)
268  wr = wgauss(ngr)
269 !------------------------------------- 4 point gaussian int. in s -----
270  do ngs = 1, 4
271  s = xgauss(ngs)
272  ws = wgauss(ngs)
273  wrs = wr * ws
274  call interpolation(3, xx(1, n1), xx(1, n2), xx(1, n3), &
275  xx(1, n4), r, s, x, xr, xs, yn1 = yy(1, n1), yn2 = yy(1, n2), &
276  yn3 = yy(1, n3), yn4 = yy(1, n4), pn1 = psi(4 * (n1- 1) + 1), &
277  pn2 = psi(4 * (n2 - 1) + 1), pn3 = psi(4 * (n3 - 1) + 1), &
278  pn4 = psi(4 * (n4 - 1) + 1), yr = yr, ys = ys, ps = ps)
279 
280  xjac = xr * ys - xs * yr
281  gpx(ngr, ngs, n) = x
282  gpjac(ngr, ngs, n) = xjac
283 !----------------------------------------- calculate right hand side ---
284  arhs = rh_side(a, x, 0._r8, 0._r8, 0._r8, 0._r8, ps, 0._r8, 0._r8)
285 !------------------------------------- 4 nodes per element ------------
286  do i = 1, 4
287 !------------------------------------- 4 functions v per node ---------
288  do j = 1, 4
289  nrow = 4 * (indices(nodeno(n, i)) - 1) + j
290  sumq = -arhs * h(j, i, ngr, ngs) * xjac
291  vx = ys * hr(j, i, ngr, ngs) - yr * hs(j, i, ngr, ngs)
292  vy = -xs * hr(j, i, ngr, ngs) + xr * hs(j, i, ngr, ngs)
293  qq(nrow) = qq(nrow) - wrs * sumq
294 !------------------------------------- 4 nodes of function psi --------
295  do k = 1, 4
296 !------------------------------------- 4 functions h in psi -----------
297  do l = 1, 4
298  ncol = 4 * (indices(nodeno(n, k)) - 1) + l
299  noff = nrow - ncol
300  if (noff >= 0) then
301  psix = ys * hr(l, k, ngr, ngs) - yr * hs(l, k, ngr, ngs)
302  psiy = -xs * hr(l, k, ngr, ngs) + xr * hs(l, k, ngr, ngs)
303  sumk = -1. / (1. + eps * x) * (psix * vx + psiy * vy) &
304  / xjac
305  kkbig(noff + 1, ncol) = kkbig(noff + 1, ncol) + wrs &
306  * sumk
307  end if
308  end do
309  end do
310  end do
311  end do
312  end do
313  end do
314  end do
315 !------------------------------------------- remove empty columns (bnd.)
316  nend = np - 1
317  if (ias == 0) nend = np
318  do j = 1, nend
319  kkbig(1, 4 * j - 3) = 1.e20_r8
320  kkbig(1, 4 * j - 1) = 1.e20_r8
321  end do
322 !------------------------------------------- if matrix exists then
323  else
324  do n = 1, nelm
325  call set_node_number(n, n1, n2, n3, n4)
326 !------------------------------------- 4 point gaussian int. in r -----
327  do ngr = 1, 4
328  r = xgauss(ngr)
329  wr = wgauss(ngr)
330 !------------------------------------- 4 point gaussian int. in s -----
331  do ngs = 1, 4
332  s = xgauss(ngs)
333  ws = wgauss(ngs)
334  x = gpx(ngr, ngs, n)
335  xjac = gpjac(ngr, ngs, n)
336  call interpolation(0, psi(4 * (n1 - 1) + 1), psi(4 * (n2 - 1) &
337  + 1), psi(4 * (n3 - 1) + 1), psi(4 * (n4 - 1) + 1), r, s, ps)
338 !----------------------------------------- calculate right hand side ---
339  arhs = rh_side(a, x, 0._r8, 0._r8, 0._r8, 0._r8, ps, 0._r8, 0._r8)
340 !------------------------------------- 4 nodes per element ------------
341  do i = 1, 4
342 !------------------------------------- 4 functions v per node ---------
343  do j = 1, 4
344  nrow = 4 * (indices(nodeno(n, i)) - 1) + j
345  sumq = -arhs * h(j, i, ngr, ngs) * xjac
346  qq(nrow) = qq(nrow) - wr * ws * sumq
347  end do
348  end do
349  end do
350  end do
351  end do
352  end if
353 
354  return
355  end subroutine build_matrix
356 
357  subroutine gaussian_elimination(first)
358 !-----------------------------------------------------------------------
359 ! subroutine to solve the system of equations using gaussian elimination
360 !-----------------------------------------------------------------------
361 
362  use itm_types
363  use mod_parameter
364  use mod_mesh
365  use mod_solver
366  use matrix_calculations, only : set_index
367 
368  implicit none
369 
370  logical, intent(in) :: first
371 
372  integer(itm_i4) :: nd, nrow, nvar
373  integer(itm_i4) :: i, j, k
374  integer(itm_i4) :: info
375  integer(itm_i4) :: ik, ik2, ikn
376  integer(itm_i4) :: nbase, nbase2
377  logical, parameter :: invert = .true.
378 
379  if (ias == 0) then
380  nd = 4 * nr * np
381  else
382  nd = 4 * nr * (np - 1)
383  end if
384  nrow = 4 * nr * np
385  nvar = 4 * np + 7
386 
387  if (first) then
388 !------------------------------inverse of indices in build_matrix to restore
389  if (ias == 1) then
390  inv_index = 0
391  call set_index(inv_index, invert)
392  end if
393 
394 !---------------------------------------------------- essl version
395 ! call dpbf(kkbig, kklda, nd, 4 * np + 8)
396 !---------------------------------------------------- lapack version
397  call dpbtrf('l', nd, nvar, kkbig, kklda, info)
398  end if
399 
400  psi(1 : 4 * nr * np) = qq
401 
402 !-------------------------------------------------- essl version
403 ! call dpbs(kkbig, kklda, nd, 4 * np + 8, psi)
404 !-------------------------------------------------- lapack version
405  call dpbtrs('l', nd, nvar, 1, kkbig, kklda, psi, nrow, info)
406 
407 !-------------- restore to simple clockwise numbering
408  if (ias == 1) then
409  qq = psi(1 : 4 * nr * np)
410  do i = 1, nr * (np - 1)
411  do k = 1, 4
412  ik = 4 * (i - 1) + k
413  ikn = 4 * (inv_index(i) - 1) + k
414  psi(ikn) = qq(ik)
415  end do
416  end do
417  do i = 1, nr
418  do k = 1, 4
419  ik = 4 * (i - 1) * np + k
420  ik2 = 4 * (i - 1) * np + 4 * (np - 1) + k
421  psi(ik2) = psi(ik)
422  end do
423  end do
424  end if
425 !-------------------------------- fill in boundary conditions
426  do i = 1, nr * np
427  psi(4 * (i - 1) + 1) = psi(4 * (i - 1) + 1) + 1._r8
428  end do
429  do j = 1, np
430  psi(4 * j - 3) = 1._r8
431  psi(4 * j - 1) = 0._r8
432  end do
433  if (ias == 1) then
434  do i = 1, nr
435  do k = 1, 4
436  nbase = 4 * (i - 1) * np + k
437  nbase2 = nbase + 4 * (np - 1)
438  psi(nbase2) = psi(nbase)
439  end do
440  end do
441  end if
442 
443  return
444  end subroutine gaussian_elimination
445 
446  subroutine normalize_psi(psaxis)
447 !-----------------------------------------------------------------------
448 ! subroutine to normalize psi to one on the boundary and zero on axis
449 !-----------------------------------------------------------------------
450 
451  use itm_types
452  use mod_mesh
453 
454  implicit none
455 
456  real(r8), intent(in) :: psaxis
457 
458  integer(itm_i4) :: i, l
459 
460  do i = 1, nr * np
461  psi(4 * (i - 1) + 1) = 1._r8 - (1._r8 - psi(4 * (i - 1) + 1)) &
462  / (1._r8 - psaxis)
463  do l = 2, 4
464  psi(4 * (i - 1) + l) = psi(4 * (i - 1) + l) / (1._r8 - psaxis)
465  end do
466  end do
467 
468  return
469  end subroutine normalize_psi
470 
471 end module solver
Definition: solver.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 deallocate_solver
Definition: solver.f90:42
real(r8) function rh_side(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine cubic_interp_1d(x1, x1s, x2, x2s, s, x, xs)
subroutine normalize_psi(psaxis)
Definition: solver.f90:446
subroutine, public quad_equation(p_m, p_mr, p_p, p_pr, tol, r, ifail)
subroutine set_index(indices, inv)
subroutine axis_symm(psaxis, xaxis, yaxis, nax, rax, sax, n_start)
Definition: solver.f90:154
subroutine set_node_number(n_node, n1, n2, n3, n4)
subroutine build_matrix(a, first)
Definition: solver.f90:208
subroutine allocate_solver(nr, np)
Definition: solver.f90:18
subroutine mnewtax(ntrial, naxis, x, n, tolx, tolf, errx, errf)
Definition: mnewtax.f90:1
subroutine find_axis(psaxis, xaxis, yaxis, nax, rax, sax)
Definition: solver.f90:52
subroutine gaussian_elimination(first)
Definition: solver.f90:357