ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
initial_grid.f90
Go to the documentation of this file.
1 subroutine initial_grid(fr)
2 !--------------------------------------------------------------------
3 ! on exit xx and yy are filled with the values of x, xr, xs, xrs and
4 ! y, yr, ys, yrs of every node :
5 ! xx(1,node) = x, xx(2,node) = xr, xx(3,node) = xs, xx(4,node) = xrs
6 ! the shape of the boundary is given by fr(m)
7 !--------------------------------------------------------------------
8 
9  use itm_types
10  use mod_parameter
11  use mod_dat
12  use mod_gauss
13  use mod_mesh
14  use mod_meshacc
15  use mod_output
16  use plot_data
17 
18  implicit none
19 
20  real(r8), dimension(mfm + 2), intent(in) :: fr
21 
22  real(r8), dimension(nr) :: sg_equi, density
23  real(r8), dimension(np) :: theta, dtc
24  real(r8) :: dt, dr
25  real(r8) :: thtj, radius
26  real(r8) :: rm, drm, drmt, drmtr
27  integer(itm_i4) :: i, j, m
28  integer(itm_i4) :: node
29 
30  dt = (1._r8 + ias) * pi / dble(np - 1)
31  dr = 1._r8 / dble(nr - 1)
32 
33 !--------------------------- change theta grid to constant arclength
34  if (iarc == 1) then
35  call arclength(fr, mfm, theta, dtc, xgauss, wgauss)
36  else
37  do j = 1, np
38  theta(j) = dt * (j - 1)
39  dtc(j) = 1._r8
40  end do
41  end if
42 
43 !TODO: mesh accumulation in psi or s for all imesh possible
44 !TODO: allow for mesh accumulation in s or in psi in general
45  if (imesh == 2 .or. imesh == 3) then
46  call mesh_accumulation(nr, n_acc_points, s_acc(1 : n_acc_points), &
47  sig(1 : n_acc_points), weights(1 : n_acc_points), equidistant, &
48  sg, dsg, ddsg, .false.)
49  else if (imesh == 3) then
50  call mesh_accumulation(nr, n_acc_points, s_acc(1 : n_acc_points), &
51  sig(1 : n_acc_points), weights(1 : n_acc_points), equidistant, &
52  sg, dsg, ddsg, .true.)
53  else
54  do i = 1, nr
55  sg(i) = dble(i - 1) / dble(nr - 1)
56  dsg(i) = 1._r8
57  ddsg(i) = 0._r8
58  end do
59  end if
60 
61  do i = 1, nr
62  sg_equi(i) = dble(i)
63  density(i) = dble(nr - 1) / dsg(i)
64  end do
65  if (xmgrace_output) then
66  call plot_data_1d('line', sg, sg_equi, nr, 'xmgr_radial_points')
67  call plot_data_1d('line', sg, density, nr, 'xmgr_s_density')
68  end if
69 
70  do i = 1, nr
71  do j = 1, np
72  node = np * (i - 1) + j
73  thtj = theta(j)
74  radius = sg(nr - i + 1)
75  xx(1,node) = radius * fr(1) * cos(thtj) / 2._r8
76  xx(2,node) = fr(1) * cos(thtj) / 2._r8
77  xx(3,node) = -radius * fr(1) * sin(thtj) / 2._r8
78  xx(4,node) = -fr(1) * sin(thtj) / 2._r8
79  yy(1,node) = radius * fr(1) * sin(thtj) / 2._r8
80  yy(2,node) = fr(1) * sin(thtj) / 2._r8
81  yy(3,node) = radius * fr(1) * cos(thtj) / 2._r8
82  yy(4,node) = fr(1) * cos(thtj) / 2._r8
83 !---------------------------- keep ellipticity on axis -----------
84  do m = 2, mharm
85  if (m == 2) then
86  rm = radius * (fr(2 * m - 1) * cos((m - 1) * thtj) &
87  + fr(2 * m) * sin((m - 1) * thtj))
88  drm = (fr(2 * m - 1) * cos((m - 1) * thtj) + fr(2 * m) &
89  * sin((m - 1) * thtj))
90  drmt = radius * (-fr(2 * m - 1) * (m - 1) * sin((m - 1) &
91  * thtj) + fr(2 * m) * (m - 1) * cos((m - 1) * thtj))
92  drmtr = (-fr(2 * m - 1) * (m - 1) * sin((m - 1) * thtj) &
93  + fr(2 * m) * (m - 1) * cos((m - 1) * thtj))
94  else
95  rm = radius**(m - 1) * (fr(2 * m - 1) * cos((m - 1) &
96  * thtj) + fr(2 * m) * sin((m - 1) * thtj))
97  drm = (m - 1) * radius**(m - 2) * (fr(2 * m - 1) &
98  * cos((m - 1) * thtj) + fr(2 * m) * sin((m - 1) * thtj))
99  drmt = radius**(m - 1) * (-fr(2 * m - 1) * (m - 1) &
100  * sin((m - 1) * thtj) + fr(2 * m) * (m - 1) &
101  * cos((m - 1) * thtj))
102  drmtr = (m - 1) * radius**(m - 2) * (-fr(2 * m - 1) &
103  * (m - 1) * sin((m - 1) * thtj) + fr(2 * m) * (m - 1) &
104  * cos((m - 1) * thtj))
105  end if
106  xx(1, node) = xx(1, node) + rm * cos(thtj)
107  yy(1, node) = yy(1, node) + rm * sin(thtj)
108  xx(2, node) = xx(2, node) + drm * cos(thtj)
109  yy(2, node) = yy(2, node) + drm * sin(thtj)
110  xx(3, node) = xx(3, node) - rm * sin(thtj) + drmt * cos(thtj)
111  yy(3, node) = yy(3, node) + rm * cos(thtj) + drmt * sin(thtj)
112  xx(4, node) = xx(4, node) - drm * sin(thtj) + drmtr * cos(thtj)
113  yy(4, node) = yy(4, node) + drm * cos(thtj) + drmtr * sin(thtj)
114  end do
115  xx(2, node) = -xx(2, node) * dr / 2._r8 * dsg(nr - i + 1)
116  xx(3, node) = xx(3, node) * dt / 2._r8 * dtc(j)
117  xx(4, node) = -xx(4, node) * dr / 2._r8 * dt / 2._r8 * dtc(j) &
118  * dsg(nr - i + 1)
119  yy(2, node) = -yy(2, node) * dr / 2._r8 * dsg(nr - i + 1)
120  yy(3, node) = yy(3, node) * dt / 2._r8 * dtc(j)
121  yy(4, node) = -yy(4, node) * dr / 2._r8 * dt / 2._r8 * dtc(j) &
122  * dsg(nr - i + 1)
123  psi(4 * (node - 1) + 1) = radius**2
124  psi(4 * (node - 1) + 2) = -2._r8* radius * dr / 2._r8 &
125  * dsg(nr - i + 1)
126  psi(4 * (node - 1) + 3) = 0._r8
127  psi(4 * (node - 1) + 4) = 0._r8
128  end do
129  end do
130 
131  return
132 end subroutine initial_grid
subroutine plot_data_1d(type_plot, x, y, np, printname, z)
Definition: plot_data.f90:11
subroutine initial_grid(fr)
Definition: initial_grid.f90:1
subroutine mesh_accumulation(nrtmp, nv, s_acc, sig, weights, equidistant, sg, dsg, ddsg, in_psi)
subroutine arclength(fr, mfm, theta, dtc, wr, ws)
Definition: arclength.f90:1