ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
shape_from_points.f90
Go to the documentation of this file.
1 subroutine shape_from_points(R_bnd, Z_bnd, n_bnd, mfm, Rgeo, amin, fr)
2 !---------------------------------------------------------------------------
3 ! subroutine interpolates the given n_bnd boundary points (R_bnd,Z_bnd)
4 ! onto a regular grid in polar theta with mfm points
5 !---------------------------------------------------------------------------
6 
7  use itm_types
8  use mod_dat, only : zgeo, verbosity
9  use mod_parameter, only : twopi
10  use mod_helena_io, only : iu6
11  use qsort_c_module
12 
13  implicit none
14 
15  real(r8), intent(in) :: r_bnd(n_bnd), z_bnd(n_bnd)
16  integer(itm_i4), intent(in) :: mfm
17  integer(itm_i4), intent(inout) :: n_bnd
18  real(r8), intent(out) :: fr(mfm)
19  real(r8), intent(out) :: rgeo, amin
20 
21  real(r8), parameter :: error = 1e-8_r8
22  real(r8), dimension(n_bnd) :: tht_tmp, fr_tmp
23  real(r8), dimension(3 * n_bnd + 6) :: work
24  real(r8), dimension(n_bnd + 2) :: tht_sort, fr_sort, dfr_sort
25  integer(itm_i4), dimension(n_bnd + 2) :: index_order
26  real(r8) :: reast, rwest
27  real(r8) :: tht, values(4)
28  integer(itm_i4) :: ieast(1), iwest(1), n_bnd_short
29  integer(itm_i4) :: i
30 
31  if (verbosity > 2) &
32  write(iu6, *) 'fshape : (R, Z) set given on ', n_bnd, ' points'
33 
34  reast = maxval(r_bnd(1 : n_bnd))
35  rwest = minval(r_bnd(1 : n_bnd))
36  ieast = maxloc(r_bnd(1 : n_bnd))
37  iwest = minloc(r_bnd(1 : n_bnd))
38  rgeo = (reast + rwest) / 2._r8
39  zgeo = (z_bnd(ieast(1)) + z_bnd(iwest(1))) / 2._r8
40  amin = (reast - rwest) / 2._r8
41 
42  if (verbosity > 3) &
43  write(iu6, '(A, 3f12.8)') 'Rgeo, Zgeo, a : ', rgeo, zgeo, amin
44 
45 !-- geometric axis defined by geometric center of boundary
46 
47  do i = 1, n_bnd
48  tht_tmp(i) = atan2(z_bnd(i) - zgeo, r_bnd(i) - rgeo)
49  fr_tmp(i) = sqrt((r_bnd(i) - rgeo)**2 + (z_bnd(i) - zgeo)**2)
50  end do
51 
52  if (abs(tht_tmp(n_bnd) - tht_tmp(1)) < error) n_bnd = n_bnd - 1
53 
54  do i = 1, n_bnd
55  index_order(i) = i
56  end do
57 
58  call qsortc(tht_tmp(1 : n_bnd), index_order(1 : n_bnd))
59 
60  tht_sort(1 : n_bnd) = tht_tmp(1 : n_bnd)
61  do i = 1, n_bnd
62  fr_sort(i) = fr_tmp(index_order(i))
63  end do
64 
65  n_bnd_short = n_bnd
66 
67  do i = 2, n_bnd
68 
69  if ((tht_sort(i) - tht_sort(1)) > twopi) then
70  n_bnd_short = i - 1
71  exit
72  end if
73 
74  end do
75 
76  if (abs(tht_sort(n_bnd_short) - tht_sort(1)) < error) then
77  tht_sort(n_bnd_short) = tht_sort(1) + twopi
78  fr_sort(n_bnd_short) = fr_sort(1)
79  else
80  tht_sort(n_bnd_short + 1) = tht_sort(1) + twopi
81  fr_sort(n_bnd_short + 1) = fr_sort(1)
82  n_bnd_short = n_bnd_short + 1
83  end if
84 
85  call tb15a(n_bnd_short, n_bnd, tht_sort, fr_sort, dfr_sort, work, iu6)
86 
87  do i = 1, mfm
88  tht = twopi * float(i - 1) / float(mfm)
89  if (tht < tht_sort(1)) tht = tht + twopi
90  if (tht > tht_sort(n_bnd_short)) tht = tht - twopi
91 
92  call tg02a(0, n_bnd_short, n_bnd, tht_sort, fr_sort, dfr_sort, tht, &
93  values)
94 
95  fr(i) = values(1)
96 
97  end do
98 
99  return
100 
101 end subroutine shape_from_points
subroutine tg02a(ix, n, n_bnd, u, s, d, x, v)
Definition: tg02a.f90:1
subroutine tb15a(n, n_bnd, x, f, d, w, lp)
Definition: tb15a.f90:1
recursive subroutine, public qsortc(A, indices)
subroutine shape_from_points(R_bnd, Z_bnd, n_bnd, mfm, Rgeo, amin, fr)