ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
mesh_accumulation.f90
Go to the documentation of this file.
1 subroutine mesh_accumulation(nrtmp, nv, s_acc, sig, weights, &
2  equidistant, sg, dsg, ddsg, in_psi)
3 !-----------------------------------------------------------------------
4 ! subroutine to construct non-equidistant radial mesh in helena
5 ! for in_psi == .true. accumulation in psi, else in s
6 !-----------------------------------------------------------------------
7 
8  use itm_types
9  use mod_helena_io, only : iu6
10  use mod_profiles
11  use helena_spline
12 
13  implicit none
14 
15  integer(itm_i4), intent(in) :: nrtmp, nv
16  real(r8), dimension(nv), intent(in) :: s_acc, sig, weights
17  real(r8), intent(in) :: equidistant
18  logical, intent(in) :: in_psi
19 
20  real(r8), dimension(nrtmp), intent(out) :: sg, dsg, ddsg
21 
22  integer(itm_i4), parameter :: nmax = 1001
23 
24  type(spline_coefficients) :: f_spline
25  real(r8), allocatable :: s1(:), fsum(:)
26  real(r8), dimension(3) :: abltg
27  real(r8) :: alfa, beta
28  real(r8) :: fi
29  real(r8) :: dx
30  real(r8) :: weight, gauss_integral
31  integer(itm_i4) :: typ
32  integer(itm_i4) :: i, n_int
33 
34  real(r8) :: fgauss
35 
36  if (any(s_acc(1 : nv) < 0._r8) .or. any(s_acc(1 : nv) > 1._r8)) then
37  write(iu6, 2)
38  stop
39  else if (any(sig(1 : nv) <= 0._r8)) then
40  write(iu6, 3)
41  stop
42  else if (any(weights(1 : nv) < 0._r8) .or. any(weights(1 : nv) > 1)) then
43  write(iu6, 4)
44  stop
45  end if
46 
47  weight = 0._r8
48  do i = 1, nv
49  weight = weight + weights(i)
50  end do
51  if (weight <= 0._r8) then
52  write(iu6, 5)
53  stop
54  end if
55 
56 ! --- first integrate the complete Gauss function for normalization
57  n_int = max(nint(1._r8 / minval(sig(1 : nv))), nmax)
58 
59  allocate(s1(n_int))
60  call allocate_spline_coefficients(f_spline, n_int)
61  allocate(fsum(n_int))
62 
63  dx = 1._r8 / dble(n_int - 1)
64  gauss_integral = 0._r8
65  fsum(1) = 0._r8
66  do i = 1, n_int - 1
67  gauss_integral = gauss_integral + fgauss(dble(i - 0.5_r8) * dx, &
68  nv, equidistant, s_acc(1 : nv), weights(1 : nv), sig(1 : nv)) * dx
69  fsum(i + 1) = gauss_integral
70  end do
71 
72  do i = 1, n_int
73  s1(i) = dble(i - 1) / dble(n_int - 1)
74  fsum(i) = fsum(i) / fsum(n_int)
75  end do
76 
77  if (in_psi) fsum = fsum**2
78 
79  alfa = 0._r8
80  beta = 0._r8
81  typ = 3
82  call spline(n_int, fsum, s1, alfa, beta, typ, f_spline)
83 
84  do i = 1, nrtmp
85  fi = dble(i - 1) / dble(nrtmp - 1)
86  sg(i) = spwert(n_int, fi, f_spline, fsum, abltg, 2)
87  dsg(i) = abltg(1)
88  ddsg(i) = abltg(2)
89  end do
90 
91  deallocate(s1)
92  call deallocate_spline_coefficients(f_spline)
93  deallocate(fsum)
94 
95  2 format(//, 1x, 'ERROR:', /, 1x, &
96  'accumulation points s_acc must be >= 0 and <= 1', /, 1x, &
97  'stop in subroutine mesh_accumulation')
98  3 format(//, 1x, 'ERROR:', /, 1x, &
99  'Gaussian widths sig must be > 0', /, 1x, &
100  'stop in subroutine mesh_accumulation')
101  4 format(//, 1x, 'ERROR:', /, 1x, &
102  'weighting factors weights must be >= 0 and <= 1', /, 1x, &
103  'stop in subroutine mesh_accumulation')
104  5 format(//, 1x, 'ERROR:', /, 1x, &
105  'at least one weighting factor must be > 0', /, 1x, &
106  'stop in subroutine mesh_accumulation')
107  return
108 
109 end subroutine mesh_accumulation
real(r8) function fgauss(x, nv, equidistant, s_acc, weights, sig)
Definition: fgauss.f90:1
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 mesh_accumulation(nrtmp, nv, s_acc, sig, weights, equidistant, sg, dsg, ddsg, in_psi)