Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
!=== COPYRIGHT AND LICENSE STATEMENT ===
!
! This file is part of the TensorProductMultigrid code.
!
! (c) The copyright relating to this work is owned jointly by the
! Crown, Met Office and NERC [2014]. However, it has been created
! with the help of the GungHo Consortium, whose members are identified
! at https://puma.nerc.ac.uk/trac/GungHo/wiki .
!
! Main Developer: Eike Mueller
!
! TensorProductMultigrid is free software: you can redistribute it and/or
! modify it under the terms of the GNU Lesser General Public License as
! published by the Free Software Foundation, either version 3 of the
! License, or (at your option) any later version.
!
! TensorProductMultigrid is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU Lesser General Public License for more details.
!
! You should have received a copy of the GNU Lesser General Public License
! along with TensorProductMultigrid (see files COPYING and COPYING.LESSER).
! If not, see <http://www.gnu.org/licenses/>.
!
!=== COPYRIGHT AND LICENSE STATEMENT ===
!==================================================================
!
! Analytical forms of RHS vectors
!
! Eike Mueller, University of Bath, Feb 2012
!
!==================================================================
module profiles
use communication
use parameters
use datatypes
use discretisation
implicit none
public::initialise_rhs
public::analytical_solution
private
contains
!==================================================================
! Initialise RHS vector
!==================================================================
subroutine initialise_rhs(grid_param,model_param,b)
implicit none
type(grid_parameters), intent(in) :: grid_param
type(model_parameters), intent(in) :: model_param
type(scalar3d), intent(inout) :: b
integer :: ix, iy, iz, ix_min, ix_max, iy_min, iy_max
real(kind=rl) :: x, y, z
real(kind=rl) :: rho, sigma, theta, phi, r, b_low, b_up, pi
#ifdef TESTCONVERGENCE
real(kind=rl) :: px,py,pz
#endif
ix_min = b%ix_min
ix_max = b%ix_max
iy_min = b%iy_min
iy_max = b%iy_max
b_low = 1.0_rl+0.25*b%grid_param%H
b_up = 1.0_rl+0.75*b%grid_param%H
pi = 4.0_rl*atan2(1.0_rl,1.0_rl)
! Initialise RHS
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,b%grid_param%nz
x = 1.0_rl*((ix-0.5_rl)/(1.0_rl*b%grid_param%n))
y = 1.0_rl*((iy-0.5_rl)/(1.0_rl*b%grid_param%n))
z = 1.0_rl*((iz-0.5_rl)/(1.0_rl*b%grid_param%nz))
#ifdef TESTCONVERGENCE
! RHS for analytical solution x*(1-x)*y*(1-y)*z*(1-z)
if (grid_param%vertbc == VERTBC_DIRICHLET) then
px = x*(1.0_rl-x)
py = y*(1.0_rl-y)
pz = z*(1.0_rl-z)
b%s(iz,iy-iy_min+1,ix-ix_min+1) = &
( 2.0_rl*model_param%omega2*((px+py)*pz &
+ model_param%lambda2*px*py)+model_param%delta*px*py*pz)
else
px = x*(1.0_rl-x)
py = y*(1.0_rl-y)
pz = 1.0_rl
b%s(iz,iy-iy_min+1,ix-ix_min+1) = &
( 2.0_rl*model_param%omega2*((px+py)*pz)+model_param%delta*px*py*pz)
end if
#else
b%s(iz,iy-iy_min+1,ix-ix_min+1) = 0.0_rl
#ifdef CARTESIANGEOMETRY
if ( ( x .ge. 0.1_rl ) .and. ( x .le. 0.4_rl ) &
.and. (y .ge. 0.3_rl ) .and. ( y .le. 0.6_rl ) &
.and. (z .ge. 0.2_rl ) .and. ( z .le. 0.7_rl ) ) &
then
b%s(iz,iy-iy_min+1,ix-ix_min+1) = 1.0_rl
end if
#else
rho = 2.0_rl*(1.0_rl*ix-0.5_rl)/grid_param%n-1.0_rl
sigma = 2.0_rl*(1.0_rl*iy-0.5_rl)/grid_param%n-1.0_rl
phi = atan(sigma)
theta = atan(rho/sqrt(1.0_rl+sigma**2))
x = sin(theta)
y = cos(theta)*sin(phi)
z = cos(theta)*cos(phi)
phi = atan2(x,y)
theta = atan2(sqrt(x**2+y**2),z)
r = 0.5_rl*(r_grid(iz)+r_grid(iz+1))
if (( (r > b_low) .and. (r < b_up) ) .and. &
(((theta>pi/10.0_rl) .and. (theta<pi/5.0_rl )) .or. &
((theta>3.0_rl*pi/8.0_rl) .and. (theta<5.0_rl*pi/8.0_rl )) .or. &
((theta>4.0_rl*pi/5.0_rl) .and. (theta<9.0_rl*pi/10.0_rl)))) then
b%s(iz,iy-iy_min+1,ix-ix_min+1) = 1.0_rl
end if
! RHS used in GPU code:
! if ( (r > b_low) .and. (r < b_up) .and. &
! (rho > -0.5) .and. (rho < 0.5) .and. &
! (sigma > -0.5).and. (sigma < 0.5) ) then
! b%s(iz,iy-iy_min+1,ix-ix_min+1) = 1.0_rl
! end if
#endif
#endif
end do
end do
end do
end subroutine initialise_rhs
!==================================================================
! Exact solution for test problem
! u(x,y,z) = x*(1-x)*y*(1-y)*z*(1-z)
!==================================================================
subroutine analytical_solution(grid_param,u)
implicit none
type(grid_parameters), intent(in) :: grid_param
type(scalar3d), intent(inout) :: u
integer :: ix, iy, iz, ix_min, ix_max, iy_min, iy_max
real(kind=rl) :: x, y, z
ix_min = u%ix_min
ix_max = u%ix_max
iy_min = u%iy_min
iy_max = u%iy_max
! Initialise RHS
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,u%grid_param%nz
x = u%grid_param%L*((ix-0.5_rl)/(1.0_rl*u%grid_param%n))
y = u%grid_param%L*((iy-0.5_rl)/(1.0_rl*u%grid_param%n))
z = u%grid_param%H*((iz-0.5_rl)/(1.0_rl*u%grid_param%nz))
if (grid_param%vertbc == VERTBC_DIRICHLET) then
u%s(iz,iy-iy_min+1,ix-ix_min+1) &
= x*(1.0_rl-x) &
* y*(1.0_rl-y) &
* z*(1.0_rl-z)
else
u%s(iz,iy-iy_min+1,ix-ix_min+1) &
= x*(1.0_rl-x) &
* y*(1.0_rl-y)
end if
end do
end do
end do
end subroutine analytical_solution
end module profiles