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
!=== 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_u_mnh
public::get_u_mnh
public::initialise_rhs_mnh
public::initialise_rhs
public::analytical_solution
private
contains
!==================================================================
! Initialise U vector
!==================================================================
subroutine initialise_u_mnh(grid_param,model_param,u,KIB,KIE,KIU,KJB,KJE,KJU,KKU,PU)
implicit none
type(grid_parameters), intent(in) :: grid_param
type(model_parameters), intent(in) :: model_param
type(scalar3d), intent(inout) :: u
integer :: ix, iy, iz, ix_min, ix_max, iy_min, iy_max, nz
real(kind=rl) :: x, y, z
real(kind=rl) :: rho, sigma, theta, phi, r, b_low, b_up, pi
integer , optional, intent(in) :: KIB,KIE,KIU,KJB,KJE,KJU,KKU
real(kind=rl) , optional, intent(in) :: PU(:,:,:)
real , dimension(:,:,:) , pointer , contiguous :: zu_st
ix_min = u%ix_min
ix_max = u%ix_max
iy_min = u%iy_min
iy_max = u%iy_max
nz = u%grid_param%nz
IF (.NOT. PRESENT(KIB) ) THEN
! Initialise RHS
if (LUseO) then
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,nz
u%s(iz,iy-iy_min+1,ix-ix_min+1) = 0.0_rl
end do
end do
end do
end if
if (LUseT) then
zu_st => u%st

WAUTELET Philippe
committed
!$acc kernels
!$acc loop independent collapse(3)
do iz=1,nz
do iy=iy_min, iy_max
do ix=ix_min, ix_max
zu_st(ix-ix_min+1,iy-iy_min+1,iz) = 0.0_rl
end do
end do
end do
!$acc end kernels
end if
ELSE
! Initialise RHS

Juan Escobar
committed
if (LUseO) then
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,nz

Juan Escobar
committed
u%s(iz,iy-iy_min+1,ix-ix_min+1) = PU(IX-ix_min+KIB,IY-iy_min+KJB,IZ)
end do
end do
end do
end if
if (LUseT) then

WAUTELET Philippe
committed
!$acc kernels
!$acc loop independent collapse(3)
do iz=1,nz

Juan Escobar
committed
do iy=iy_min, iy_max
do ix=ix_min, ix_max
zu_st(ix-ix_min+1,iy-iy_min+1,iz) = PU(IX-ix_min+KIB,IY-iy_min+KJB,IZ)

Juan Escobar
committed
end do
end do
end do
!$acc end kernels

Juan Escobar
committed
end if
END IF
end subroutine initialise_u_mnh
!==================================================================
! Get U vector
!==================================================================
subroutine get_u_mnh(grid_param,model_param,u,KIB,KIE,KIU,KJB,KJE,KJU,KKU,PU)
implicit none
type(grid_parameters), intent(in) :: grid_param
type(model_parameters), intent(in) :: model_param
type(scalar3d), intent(inout) :: u
integer :: ix, iy, iz, ix_min, ix_max, iy_min, iy_max, nz
real(kind=rl) :: x, y, z
real(kind=rl) :: rho, sigma, theta, phi, r, b_low, b_up, pi
integer , optional, intent(in) :: KIB,KIE,KIU,KJB,KJE,KJU,KKU
real(kind=rl) , optional, intent(inout) :: PU(:,:,:)

ESCOBAR MUNOZ Juan
committed
real , dimension(:,:,:) , pointer , contiguous :: zu_st
ix_min = u%ix_min
ix_max = u%ix_max
iy_min = u%iy_min
iy_max = u%iy_max
nz = u%grid_param%nz
IF (.NOT. PRESENT(KIB) ) THEN

Juan Escobar
committed
!
ELSE

Juan Escobar
committed
! Get PU
if (LUseO) then
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,nz

Juan Escobar
committed
PU(IX-ix_min+KIB,IY-iy_min+KJB,IZ) = u%s(iz,iy-iy_min+1,ix-ix_min+1)
end do
end do
end do
else

ESCOBAR MUNOZ Juan
committed
zu_st => u%st

WAUTELET Philippe
committed
!$acc kernels
!$acc loop independent collapse(3)
do iz=1,nz

Juan Escobar
committed
do iy=iy_min, iy_max
do ix=ix_min, ix_max

ESCOBAR MUNOZ Juan
committed
PU(IX-ix_min+KIB,IY-iy_min+KJB,IZ) = zu_st(ix-ix_min+1,iy-iy_min+1,iz)

Juan Escobar
committed
end do
end do
end do

ESCOBAR MUNOZ Juan
committed
!$acc end kernels

Juan Escobar
committed
end if
END IF
end subroutine get_u_mnh
!==================================================================
! Initialise RHS vector
!==================================================================
subroutine initialise_rhs_mnh(grid_param,model_param,b,KIB,KIE,KIU,KJB,KJE,KJU,KKU,PY)
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, nz, n
real(kind=rl) :: x, y, z
real(kind=rl) :: rho, sigma, theta, phi, r, b_low, b_up, pi
integer , optional, intent(in) :: KIB,KIE,KIU,KJB,KJE,KJU,KKU
real(kind=rl) , optional, intent(in) :: PY(:,:,:)
real , dimension(:,:,:) , pointer , contiguous :: zb_st
ix_min = b%ix_min
ix_max = b%ix_max
iy_min = b%iy_min
iy_max = b%iy_max
nz = b%grid_param%nz
n = b%grid_param%n
!
IF (.NOT. PRESENT(KIB) ) THEN
! Initialise RHS
if (LUseO) then
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,nz
x = 1.0_rl*((ix-0.5_rl)/(1.0_rl*n))
y = 1.0_rl*((iy-0.5_rl)/(1.0_rl*n))
z = 1.0_rl*((iz-0.5_rl)/(1.0_rl*nz))
b%s(iz,iy-iy_min+1,ix-ix_min+1) = 0.0_rl
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
end do
end do
end do
endif
if (LUseT) then
zb_st => b%st

WAUTELET Philippe
committed
!$acc kernels
!$acc loop independent collapse(3)
do iz=1,nz
do iy=iy_min, iy_max
do ix=ix_min, ix_max
x = 1.0_rl*((ix-0.5_rl)/(1.0_rl*n))
y = 1.0_rl*((iy-0.5_rl)/(1.0_rl*n))
z = 1.0_rl*((iz-0.5_rl)/(1.0_rl*nz))
zb_st(ix-ix_min+1,iy-iy_min+1,iz) = 0.0_rl
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
zb_st(ix-ix_min+1,iy-iy_min+1,iz) = 1.0_rl
end if
end do
end do
end do
!$acc end kernels
end if
ELSE
! Initialise RHS
if (LUseO) then
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,nz
b%s(iz,iy-iy_min+1,ix-ix_min+1) = PY(IX-ix_min+KIB,IY-iy_min+KJB,IZ)
end do
end do
end do
end if
if (LUseT) then
zb_st => b%st

WAUTELET Philippe
committed
!$acc kernels
!$acc loop independent collapse(3)
do iz=1,nz
do iy=iy_min, iy_max
do ix=ix_min, ix_max
zb_st(ix-ix_min+1,iy-iy_min+1,iz) = PY(IX-ix_min+KIB,IY-iy_min+KJB,IZ)
end do
end do
end do
!$acc end kernels
end if
END IF
end subroutine initialise_rhs_mnh
!==================================================================
! 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, nz, n
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
nz = b%grid_param%nz
n = b%grid_param%n
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,nz
x = 1.0_rl*((ix-0.5_rl)/(1.0_rl*n))
y = 1.0_rl*((iy-0.5_rl)/(1.0_rl*n))
z = 1.0_rl*((iz-0.5_rl)/(1.0_rl*nz))
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
#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, nz, n
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
nz = u%grid_param%nz
n = u%grid_param%n
! Initialise RHS
do ix=ix_min, ix_max
do iy=iy_min, iy_max
do iz=1,nz
x = u%grid_param%L*((ix-0.5_rl)/(1.0_rl*n))
y = u%grid_param%L*((iy-0.5_rl)/(1.0_rl*n))
z = u%grid_param%H*((iz-0.5_rl)/(1.0_rl*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