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
!=== 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 ===
!==================================================================
!
! Discretisation module of the model problem
!
!
! -omega2 * (d^2/dx^2 + d^2/dy^2 + lambda2 * d^2/dz^2 ) u
! + delta u = RHS
! [Cartesian]
!
! or
!
! -omega2 * (laplace_{2d} + lambda2/r^2 d/dr (r^2 d/dr)) u
! + delta u = RHS
! [Spherical]
!
! We use a cell centered finite volume discretisation with
! The equation is discretised either in a unit cube or on 1/6th
! of a cubed sphere grid.
!
! The vertical grid spacing is not necessarily uniform and can
! be chosen by specifying the vertical grid in a vector.
!
! The following boundary conditions are used:
!
! * Dirichlet in the horizontal
! * Neumann in the vertical
!
! For delta = 0 the operator reduces to the Poisson operator.
!
! Eike Mueller, University of Bath, Feb 2012
!
!==================================================================
module discretisation
use parameters
use messages
use datatypes
use communication

Juan Escobar
committed
#ifndef MNH
use mpi
#else
use modd_mpif

Juan Escobar
committed
#endif
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
implicit none
private
type model_parameters
real(kind=rl) :: omega2 ! omega^2
real(kind=rl) :: lambda2 ! lambda^2
real(kind=rl) :: delta ! delta
end type model_parameters
! --- Stencil ---
!
! Grid traversal direction in SOR
integer, parameter :: DIRECTION_FORWARD = 1
integer, parameter :: DIRECTION_BACKWARD = 2
! Ordering in SOR
! Lexicographic ordering
integer, parameter :: ORDERING_LEX = 1
! Red-black ordering
integer, parameter :: ORDERING_RB = 2
type smoother_parameters
! smoother
integer :: smoother
! relaxation parameter
real(kind=rl) :: rho
! ordering of degrees of freedom
integer :: ordering
end type smoother_parameters
! Allowed smoothers
integer, parameter :: SMOOTHER_LINE_SOR = 3
integer, parameter :: SMOOTHER_LINE_SSOR = 4
integer, parameter :: SMOOTHER_LINE_JAC = 6
! Number of levels
integer :: nlev
! Grid parameters
type(grid_parameters) :: grid_param
! Model parameters
type(model_parameters) :: model_param
! Smoother parameters
type(smoother_parameters) :: smoother_param
! Arrays for measuring the residual reduction
real(kind=rl), allocatable :: log_resreduction(:)
integer, allocatable :: nsmooth_total(:)
! Data structure for storing the vertical discretisation
type vertical_coefficients

ESCOBAR MUNOZ Juan
committed
real(kind=rl), pointer , contiguous :: a(:)
real(kind=rl), pointer , contiguous :: b(:)
real(kind=rl), pointer , contiguous :: c(:)
real(kind=rl), pointer , contiguous :: d(:)
end type vertical_coefficients
! Stoarge for vertical coefficients
type(vertical_coefficients) :: vert_coeff

Juan Escobar
committed
type Temp_jacobi
real(kind=rl), pointer , contiguous :: r(:)
real(kind=rl), pointer , contiguous :: c(:), utmp(:)
real(kind=rl), pointer , contiguous :: u0(:,:,:)
real(kind=rl), pointer , contiguous :: ut0(:,:,:)
type(scalar3d) , pointer :: Sr,Sc,Sut0,Sutmp
end type Temp_jacobi
integer , parameter :: max_lev = 128
type (Temp_jacobi) , save , dimension(max_lev) :: Tjacobi
real(kind=rl), pointer , contiguous , dimension(:) :: zt1d_discretisation
INTEGER :: zt1d_discretisation_size = 10000000

ESCOBAR MUNOZ Juan
committed
INTEGER , parameter :: zt1d_discretisation_size_factor = 4 * 8

Juan Escobar
committed
INTEGER :: nt1d_discretisation_top = 0 , nt1d_discretisation_max = 0
INTEGER , ALLOCATABLE, DIMENSION (:) :: nt1d_discretisation_pointer , nt1d_discretisation_pointer_ksize
INTEGER , parameter :: zt1d_discretisation_pointer_size = 128
INTEGER :: nt1d_discretisation_pointer_top = 0 , nt1d_discretisation_pointer_max = 0
logical :: gfirst_call_zt1d_discretisation = .true.

Juan Escobar
committed
public::discretisation_initialise_mnh
public::discretisation_initialise
public::discretisation_finalise
public::smooth_mnh
public::smooth
public::line_SOR
public::line_SSOR
public::line_jacobi_mnh
public::line_jacobi

Juan Escobar
committed
public::calculate_residual_mnh
public::calculate_residual

Juan Escobar
committed
public::apply_mnh
public::apply
public::model_parameters
public::smoother_parameters
public::volume_of_element
public::SMOOTHER_LINE_SOR
public::SMOOTHER_LINE_SSOR
public::SMOOTHER_LINE_JAC
public::DIRECTION_FORWARD
public::DIRECTION_BACKWARD
public::ORDERING_LEX
public::ORDERING_RB

Juan Escobar
committed
public::zt1d_discretisation_init
public::zt1d_discretisation_allocate3d
contains
!==================================================================
! Initialise module
!==================================================================

Juan Escobar
committed
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
subroutine zt1d_discretisation_init(kiu,kju,kku)
implicit none
integer , optional , intent(in) :: kiu,kju,kku ! size of 1 3D MNH array
if (gfirst_call_zt1d_discretisation) then
!
gfirst_call_zt1d_discretisation = .false.
!
if (present(kiu)) then
zt1d_discretisation_size = kiu*kju*kku * zt1d_discretisation_size_factor
end if
allocate(zt1d_discretisation(zt1d_discretisation_size))
zt1d_discretisation = 0.0
!$acc enter data copyin(zt1d_discretisation)
!
allocate(nt1d_discretisation_pointer(zt1d_discretisation_pointer_size))
nt1d_discretisation_pointer = 0.0
allocate(nt1d_discretisation_pointer_ksize(zt1d_discretisation_pointer_size))
nt1d_discretisation_pointer_ksize = 0.0
!
end if
end subroutine zt1d_discretisation_init
function zt1d_discretisation_allocate3d(ptab3d,kib,kie,kjb,kje,kkb,kke) result(kindex)
implicit none
real(kind=rl), pointer , contiguous , dimension(:,:,:) :: ptab3d
!real(kind=rl), pointer , dimension(:,:,:) :: ptab3d
integer , intent(in) :: kib,kie,kjb,kje,kkb,kke
integer :: kindex
! local var
integer :: ksize, kbeg, kend
if ( nt1d_discretisation_pointer_top == zt1d_discretisation_pointer_size ) then
print*,"ERROR zt1d_discretisation_allocate3d:: zt1d_discretisation_pointer_size to small=",zt1d_discretisation_pointer_size
print*,"ERROR zt1d_discretisation_allocate3d:: Augment it "

Juan Escobar
committed
stop
end if
ksize = (kie-kib+1)*(kje-kjb+1)*(kke-kkb+1)
kbeg = nt1d_discretisation_top + 1
kend = nt1d_discretisation_top + ksize
if ( kend > zt1d_discretisation_size ) then
print*,"ERROR zt1d_discretisation_allocate3d:: zt1d_discretisation_size to small=",zt1d_discretisation_size
print*,"ERROR zt1d_discretisation_allocate3d:: Augment it "
stop
end if
ptab3d(kib:kie,kjb:kje,kkb:kke) => zt1d_discretisation(kbeg:)
nt1d_discretisation_pointer_top = nt1d_discretisation_pointer_top + 1
nt1d_discretisation_top = kend
nt1d_discretisation_pointer(nt1d_discretisation_pointer_top) = kend
nt1d_discretisation_pointer_ksize(nt1d_discretisation_pointer_top) = ksize
if ( nt1d_discretisation_pointer_top > nt1d_discretisation_pointer_max ) then
nt1d_discretisation_pointer_max = nt1d_discretisation_pointer_top
!!$ print*,"zt1d_discretisation_allocate3d:: nt1d_discretisation_pointer_max=",nt1d_discretisation_pointer_max
!!$ call flush(6)

Juan Escobar
committed
endif
if ( nt1d_discretisation_top > nt1d_discretisation_max ) then
nt1d_discretisation_max = nt1d_discretisation_top
!!$ print*,"zt1d_discretisation_allocate3d:: nt1d_discretisation_max=",nt1d_discretisation_max
!!$ call flush(6)

Juan Escobar
committed
endif

ESCOBAR MUNOZ Juan
committed
kindex = nt1d_discretisation_pointer_top

Juan Escobar
committed
end function zt1d_discretisation_allocate3d

Juan Escobar
committed
subroutine discretisation_initialise_mnh(grid_param_in, &
model_param_in, &
smoother_param_in, &
nlev_in, &
PA_K,PB_K,PC_K,PD_K)

Juan Escobar
committed
implicit none
type(grid_parameters), intent(in) :: grid_param_in
type(model_parameters), intent(in) :: model_param_in
type(smoother_parameters), intent(in) :: smoother_param_in
integer, intent(in) :: nlev_in
real(kind=rl) , optional , intent (in) :: PA_K(:),PB_K(:),PC_K(:),PD_K(:)
! local var

Juan Escobar
committed
integer :: k

Juan Escobar
committed
grid_param = grid_param_in
model_param = model_param_in
smoother_param = smoother_param_in
nlev = nlev_in
allocate(log_resreduction(nlev))
allocate(nsmooth_total(nlev))
log_resreduction(:) = 0.0_rl
nsmooth_total(:) = 0
allocate(r_grid(grid_param%nz+1))
if (grid_param%graded) then
do k=1,grid_param%nz+1
r_grid(k) = grid_param%H*(1.0_rl*(k-1.0_rl)/grid_param%nz)**2
end do
else
do k=1,grid_param%nz+1
r_grid(k) = grid_param%H*(1.0_rl*(k-1.0_rl)/grid_param%nz)
end do
end if
! Allocate arrays for vertical discretisation matrices
! and calculate matrix entries
allocate(vert_coeff%a(grid_param%nz))
allocate(vert_coeff%b(grid_param%nz))
allocate(vert_coeff%c(grid_param%nz))
allocate(vert_coeff%d(grid_param%nz))
call construct_vertical_coeff_mnh(PA_K,PB_K,PC_K,PD_K)

ESCOBAR MUNOZ Juan
committed
!$acc enter data copyin(vert_coeff%a,vert_coeff%b,vert_coeff%c,vert_coeff%d)

Juan Escobar
committed
end subroutine discretisation_initialise_mnh
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
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
subroutine discretisation_initialise(grid_param_in, &
model_param_in, &
smoother_param_in, &
nlev_in)
implicit none
type(grid_parameters), intent(in) :: grid_param_in
type(model_parameters), intent(in) :: model_param_in
type(smoother_parameters), intent(in) :: smoother_param_in
integer, intent(in) :: nlev_in
integer :: k
grid_param = grid_param_in
model_param = model_param_in
smoother_param = smoother_param_in
nlev = nlev_in
allocate(log_resreduction(nlev))
allocate(nsmooth_total(nlev))
log_resreduction(:) = 0.0_rl
nsmooth_total(:) = 0
allocate(r_grid(grid_param%nz+1))
if (grid_param%graded) then
do k=1,grid_param%nz+1
r_grid(k) = grid_param%H*(1.0_rl*(k-1.0_rl)/grid_param%nz)**2
end do
else
do k=1,grid_param%nz+1
r_grid(k) = grid_param%H*(1.0_rl*(k-1.0_rl)/grid_param%nz)
end do
end if
#ifdef CARTESIANGEOMETRY
#else
r_grid(:) = 1.0_rl + r_grid(:)
#endif
! Allocate arrays for vertical discretisation matrices
! and calculate matrix entries
allocate(vert_coeff%a(grid_param%nz))
allocate(vert_coeff%b(grid_param%nz))
allocate(vert_coeff%c(grid_param%nz))
allocate(vert_coeff%d(grid_param%nz))
call construct_vertical_coeff()
end subroutine discretisation_initialise
!==================================================================
! Finalise module
!==================================================================
subroutine discretisation_finalise()
implicit none
integer :: level
real(kind=rl) :: rho_avg
#ifdef MEASURESMOOTHINGRATE
if (i_am_master_mpi) then
write(STDOUT,'("Average smoothing rates:")')
do level=nlev,1,-1
if (nsmooth_total(level) > 0) then
rho_avg = exp(log_resreduction(level)/nsmooth_total(level))
else
rho_avg = 1.0_rl
end if
write(STDOUT,'("rho_{avg}(",I3,") = ",E10.4," ( ",I5," x )")') &
level, rho_avg, nsmooth_total(level)
end do
end if
#endif
deallocate(log_resreduction)
deallocate(nsmooth_total)
deallocate(r_grid)
! Deallocate storage for vertical discretisation matrices
deallocate(vert_coeff%a)
deallocate(vert_coeff%b)
deallocate(vert_coeff%c)
deallocate(vert_coeff%d)
end subroutine discretisation_finalise
!==================================================================
! Construct alpha_{i',j'} and |T_{ij}| needed for the
! horizontal stencil
! ( alpha_{i+1,j},
! alpha_{i-1,j},
! alpha_{i,j+1},
! alpha_{i,j-1},
! alpha_{ij})
! (ix,iy) are LOCAL indices of the grid boxes, which are
! converted to global indices

Juan Escobar
committed
!==================================================================
subroutine construct_alpha_T_mnh(grid_param,ix,iy,alpha_T,Tij)
implicit none
type(grid_parameters), intent(in) :: grid_param
integer, intent(in) :: ix
integer, intent(in) :: iy
real(kind=rl), intent(inout), dimension(5) :: alpha_T
real(kind=rl), intent(out) :: Tij

Juan Escobar
committed
!local var
real(kind=rl) :: h, rho_i, sigma_j

Juan Escobar
committed
real(kind=rl) :: xcoef

Juan Escobar
committed
logical :: l2nd

Juan Escobar
committed
h = grid_param%L/grid_param%n
! Cartesian coefficients
Tij = h**2

ESCOBAR Juan
committed
! optimisation for newman MNH case = all coef constant
alpha_T(1:4) = 1.0_rl
alpha_T(5) = 4.0_rl
return

Juan Escobar
committed
xcoef = 0.5_rl ! 0.0

Juan Escobar
committed
l2nd = .false. ! .true. ! .false.

Juan Escobar
committed
alpha_T(1) = 1.0
alpha_T(2) = 1.0

Juan Escobar
committed
if (ix == grid_param%n) then

Juan Escobar
committed
alpha_T(1) = xcoef * 2.0_rl
if (l2nd) alpha_T(2) = 2.0_rl

Juan Escobar
committed
end if
if (ix == 1) then

Juan Escobar
committed
alpha_T(2) = xcoef * 2.0_rl
if (l2nd) alpha_T(1) = 2.0_rl

Juan Escobar
committed
end if

Juan Escobar
committed
alpha_T(3) = 1.0
alpha_T(4) = 1.0

Juan Escobar
committed
if (iy == grid_param%n) then

Juan Escobar
committed
alpha_T(3) = xcoef * 2.0_rl
if (l2nd) alpha_T(4) = 2.0

Juan Escobar
committed
end if
if (iy == 1) then

Juan Escobar
committed
alpha_T(4) = xcoef * 2.0_rl
if (l2nd) alpha_T(3) = 2.0

Juan Escobar
committed
end if

Juan Escobar
committed
alpha_T(5) = alpha_T(1) + alpha_T(2) + alpha_T(3) + alpha_T(4)
end subroutine construct_alpha_T_mnh

Juan Escobar
committed
! constant coef for MNH
subroutine construct_alpha_T_cst_mnh(grid_param,alpha_T,Tij)
implicit none
type(grid_parameters), intent(in) :: grid_param
real(kind=rl), intent(inout), dimension(5) :: alpha_T
real(kind=rl), intent(out) :: Tij
!local var
real(kind=rl) :: h
h = grid_param%L/grid_param%n
! Cartesian coefficients
Tij = h**2
! optimisation for newman MNH case = all coef constant
alpha_T(1:4) = 1.0_rl
alpha_T(5) = 4.0_rl
end subroutine construct_alpha_T_cst_mnh
!==================================================================
subroutine construct_alpha_T(grid_param,ix,iy,alpha_T,Tij)
implicit none
type(grid_parameters), intent(in) :: grid_param
integer, intent(in) :: ix
integer, intent(in) :: iy
real(kind=rl), intent(inout), dimension(5) :: alpha_T
real(kind=rl), intent(out) :: Tij
real(kind=rl) :: h, rho_i, sigma_j
#ifdef CARTESIANGEOMETRY
h = grid_param%L/grid_param%n
! Cartesian coefficients
Tij = h**2
if (ix == grid_param%n) then

Juan Escobar
committed
alpha_T(1) = 2.0_rl
else
alpha_T(1) = 1.0_rl
end if
if (ix == 1) then

Juan Escobar
committed
alpha_T(2) = 2.0_rl
else
alpha_T(2) = 1.0_rl
end if
if (iy == grid_param%n) then

Juan Escobar
committed
alpha_T(3) = 2.0_rl
else
alpha_T(3) = 1.0_rl
end if
if (iy == 1) then

Juan Escobar
committed
alpha_T(4) = 2.0_rl
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
else
alpha_T(4) = 1.0_rl
end if
#else
! Coefficients in cubed sphere geometry
! (rho_i,sigma_j) \in [-1,1] x [-1,1] are the coordinates of the
! cubed sphere segment
h = 2.0_rl/grid_param%n
Tij = volume_of_element(ix,iy,grid_param)
rho_i = 2.0_rl*(1.0_rl*ix-0.5_rl)/grid_param%n-1.0_rl
sigma_j = 2.0_rl*(1.0_rl*iy-0.5_rl)/grid_param%n-1.0_rl
! alpha_{i+1,j}
if (ix == grid_param%n) then
alpha_T(1) = 2.0_rl*DSQRT((1.0_rl+(rho_i+0.25_rl*h)**2)/(1.0_rl+sigma_j**2))
else
alpha_T(1) = DSQRT((1.0_rl+(rho_i+0.5_rl*h)**2)/(1.0_rl+sigma_j**2))
end if
! alpha_{i-1,j}
if (ix == 1) then
alpha_T(2) = 2.0_rl*DSQRT((1.0_rl+(rho_i-0.25_rl*h)**2)/(1.0_rl+sigma_j**2))
else
alpha_T(2) = DSQRT((1.0_rl+(rho_i-0.5_rl*h)**2)/(1.0_rl+sigma_j**2))
end if
! alpha_{i,j+1}
if (iy == grid_param%n) then
alpha_T(3) = 2.0_rl*DSQRT((1.0_rl+(sigma_j+0.25_rl*h)**2)/(1.0_rl+rho_i**2))
else
alpha_T(3) = DSQRT((1.0_rl+(sigma_j+0.5_rl*h)**2)/(1.0_rl+rho_i**2))
end if
! alpha_{i,j-1}
if (iy == 1) then
alpha_T(4) = 2.0_rl*DSQRT((1.0_rl+(sigma_j-0.25_rl*h)**2)/(1.0_rl+rho_i**2))
else
alpha_T(4) = DSQRT((1.0_rl+(sigma_j-0.5_rl*h)**2)/(1.0_rl+rho_i**2))
end if
#endif
alpha_T(5) = alpha_T(1) + alpha_T(2) + alpha_T(3) + alpha_T(4)
end subroutine construct_alpha_T
!==================================================================
! Construct coefficients of tridiagonal matrix A_T
! describing the coupling in the vertical direction and the
! diagonal matrix diag(d)
!==================================================================
subroutine construct_vertical_coeff_mnh(PA_K,PB_K,PC_K,PD_K)

Juan Escobar
committed
implicit none
real(kind=rl) , optional , intent (in) :: PA_K(:),PB_K(:),PC_K(:),PD_K(:)
!local var

Juan Escobar
committed
real(kind=rl) :: a_k_tmp, b_k_tmp, c_k_tmp, d_k_tmp
real(kind=rl) :: omega2, lambda2, delta, vol_r, surface_k, surface_kp1
integer :: k
IF (.NOT. PRESENT(PA_K)) THEN

Juan Escobar
committed
omega2 = model_param%omega2
lambda2 = model_param%lambda2
delta = model_param%delta
do k = 1, grid_param%nz

Juan Escobar
committed
vol_r = r_grid(k+1)-r_grid(k)
surface_k = 1.0_rl
surface_kp1 = 1.0_rl

Juan Escobar
committed
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
! Diagonal element
a_k_tmp = delta*vol_r
! off diagonal elements
! Boundary conditions
! Top
if (k == grid_param%nz) then
if (grid_param%vertbc == VERTBC_DIRICHLET) then
b_k_tmp = - 2.0_rl * omega2*lambda2 &
* surface_kp1/(r_grid(k+1)-r_grid(k))
else
b_k_tmp = 0.0_rl
end if
else
b_k_tmp = - 2.0_rl*omega2*lambda2 &
* surface_kp1/(r_grid(k+2)-r_grid(k))
end if
! Bottom
if (k == 1) then
if (grid_param%vertbc == VERTBC_DIRICHLET) then
c_k_tmp = - 2.0_rl * omega2*lambda2 &
* surface_k/(r_grid(k+1)-r_grid(k))
else
c_k_tmp = 0.0_rl
end if
else
c_k_tmp = - 2.0_rl * omega2 * lambda2 &
* surface_k/(r_grid(k+1)-r_grid(k-1))
end if
! Diagonal matrix d_k
d_k_tmp = - omega2 * (r_grid(k+1)-r_grid(k))
vert_coeff%a(k) = a_k_tmp/d_k_tmp
vert_coeff%b(k) = b_k_tmp/d_k_tmp
vert_coeff%c(k) = c_k_tmp/d_k_tmp
vert_coeff%d(k) = d_k_tmp
end do
ELSE
do k = 1, grid_param%nz
vert_coeff%a(k) = PA_K(k)
vert_coeff%b(k) = PB_K(k)
vert_coeff%c(k) = PC_K(k)
vert_coeff%d(k) = PD_K(k)
end do
ENDIF

Juan Escobar
committed
end subroutine construct_vertical_coeff_mnh
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
subroutine construct_vertical_coeff()
implicit none
real(kind=rl) :: a_k_tmp, b_k_tmp, c_k_tmp, d_k_tmp
real(kind=rl) :: omega2, lambda2, delta, vol_r, surface_k, surface_kp1
integer :: k
omega2 = model_param%omega2
lambda2 = model_param%lambda2
delta = model_param%delta
do k = 1, grid_param%nz
#ifdef CARTESIANGEOMETRY
vol_r = r_grid(k+1)-r_grid(k)
surface_k = 1.0_rl
surface_kp1 = 1.0_rl
#else
vol_r = (r_grid(k+1)**3 - r_grid(k)**3)/3.0_rl
surface_k = r_grid(k)**2
surface_kp1 = r_grid(k+1)**2
#endif
! Diagonal element
a_k_tmp = delta*vol_r
! off diagonal elements
! Boundary conditions
! Top
if (k == grid_param%nz) then
if (grid_param%vertbc == VERTBC_DIRICHLET) then
b_k_tmp = - 2.0_rl * omega2*lambda2 &
* surface_kp1/(r_grid(k+1)-r_grid(k))
else
b_k_tmp = 0.0_rl
end if
else
b_k_tmp = - 2.0_rl*omega2*lambda2 &
* surface_kp1/(r_grid(k+2)-r_grid(k))
end if
! Bottom
if (k == 1) then
if (grid_param%vertbc == VERTBC_DIRICHLET) then
c_k_tmp = - 2.0_rl * omega2*lambda2 &
* surface_k/(r_grid(k+1)-r_grid(k))
else
c_k_tmp = 0.0_rl
end if
else
c_k_tmp = - 2.0_rl * omega2 * lambda2 &
* surface_k/(r_grid(k+1)-r_grid(k-1))
end if
! Diagonal matrix d_k
d_k_tmp = - omega2 * (r_grid(k+1)-r_grid(k))
vert_coeff%a(k) = a_k_tmp/d_k_tmp
vert_coeff%b(k) = b_k_tmp/d_k_tmp
vert_coeff%c(k) = c_k_tmp/d_k_tmp
vert_coeff%d(k) = d_k_tmp
end do
end subroutine construct_vertical_coeff
!==================================================================
! Calculate local residual r = b - A.u
!==================================================================

Juan Escobar
committed
subroutine calculate_residual_mnh(level,m,b,u,r)
implicit none
integer, intent(in) :: level
integer, intent(in) :: m
type(scalar3d), intent(in) :: b
type(scalar3d), intent(inout) :: u
type(scalar3d), intent(inout) :: r
integer :: ix,iy,iz

Juan Escobar
committed
integer :: iib,iie,ijb,ije,ikb,ike

Juan Escobar
committed

ESCOBAR MUNOZ Juan
committed
real , dimension(:,:,:) , pointer , contiguous :: zr_st , zb_st

ESCOBAR MUNOZ Juan
committed
integer :: nib,nie,njb,nje,nzb,nze

Juan Escobar
committed
! r <- A.u

Juan Escobar
committed
!call boundary_mnh(u)

Juan Escobar
committed
call apply_mnh(u,r)
! r <- b - r = b - A.u
if (LUseO) then
do ix=u%icompx_min,u%icompx_max
do iy=u%icompy_min,u%icompy_max
do iz=1,u%grid_param%nz
r%s(iz,iy,ix) = b%s(iz,iy,ix) - r%s(iz,iy,ix)
end do
end do
end do
endif
if (LUseT) then

Juan Escobar
committed
!!$ do iz=1,u%grid_param%nz
!!$ do iy=u%icompy_min,u%icompy_max
!!$ do ix=u%icompx_min,u%icompx_max
!!$ r%st(ix,iy,iz) = b%st(ix,iy,iz) - r%st(ix,iy,iz)
!!$ end do
!!$ end do
!!$ end do
!-----------------
iib=u%icompx_min
iie=u%icompx_max
ijb=u%icompy_min
ije=u%icompy_max
ikb=1
ike=u%grid_param%nz
zr_st => r%st
zb_st => b%st

ESCOBAR MUNOZ Juan
committed
nib = Lbound(zr_st,1) ; nie = Ubound(zr_st,1)
njb = Lbound(zr_st,2) ; nje = Ubound(zr_st,2)
nzb = Lbound(zr_st,3) ; nze = Ubound(zr_st,3)
call calculate_residual_mnh_dim(zr_st,zb_st)
endif

ESCOBAR MUNOZ Juan
committed
contains
subroutine calculate_residual_mnh_dim(pzr_st,pzb_st)
implicit none
real :: pzr_st(nib:nie,njb:nje,nzb:nze), &
pzb_st(nib:nie,njb:nje,nzb:nze)
!$acc kernels present(pzr_st,pzb_st)
pzr_st(iib:iie,ijb:ije,ikb:ike) = pzb_st(iib:iie,ijb:ije,ikb:ike) - pzr_st(iib:iie,ijb:ije,ikb:ike)
!$acc end kernels
end subroutine calculate_residual_mnh_dim

Juan Escobar
committed
end subroutine calculate_residual_mnh
subroutine calculate_residual(level,m,b,u,r)
implicit none
integer, intent(in) :: level
integer, intent(in) :: m
type(scalar3d), intent(in) :: b
type(scalar3d), intent(inout) :: u
type(scalar3d), intent(inout) :: r
integer :: ix,iy,iz

Juan Escobar
committed
! r <- A.u
call apply(u,r)
! r <- b - r = b - A.u
do ix=u%icompx_min,u%icompx_max
do iy=u%icompy_min,u%icompy_max
do iz=1,u%grid_param%nz
r%s(iz,iy,ix) = b%s(iz,iy,ix) - r%s(iz,iy,ix)
end do
end do
end do
end subroutine calculate_residual
!==================================================================
! Apply operator v = A.u
!==================================================================

Juan Escobar
committed
subroutine apply_mnh(u,v)
implicit none

Juan Escobar
committed
type(scalar3d), intent(inout) :: u

Juan Escobar
committed
type(scalar3d), intent(inout) :: v

Juan Escobar
committed
! local var

Juan Escobar
committed
real(kind=rl), dimension(5) :: alpha_T
real(kind=rl) :: Tij
real(kind=rl) :: a_k, b_k, c_k, d_k
integer :: ix,iy,iz
real(kind=rl) :: tmp

Juan Escobar
committed
integer :: iib,iie,ijb,ije

Juan Escobar
committed

ESCOBAR MUNOZ Juan
committed
real(kind=rl), dimension(:,:,:) , pointer , contiguous :: zv_st , zu_st
real(kind=rl), dimension(:) , pointer , contiguous :: za_k, zb_k, zc_k, zd_k
integer :: ii,ij

Juan Escobar
committed
integer :: ize

ESCOBAR MUNOZ Juan
committed
integer :: nib,nie,njb,nje,nzb,nze

Juan Escobar
committed
call boundary_mnh(u)
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
if (LUseO) then
do ix=u%icompx_min,u%icompx_max
do iy=u%icompy_min,u%icompy_max
! Construct horizontal part of stencil
call construct_alpha_T_mnh(u%grid_param, &
ix+u%ix_min-1, &
iy+u%iy_min-1, &
alpha_T,Tij)
do iz=1,u%grid_param%nz
a_k = vert_coeff%a(iz)
b_k = vert_coeff%b(iz)
c_k = vert_coeff%c(iz)
d_k = vert_coeff%d(iz)
tmp = ((a_k-b_k-c_k)*Tij ) * u%s(iz,iy,ix)
if (iz < grid_param%nz) then
tmp = tmp + b_k*Tij * u%s(iz+1,iy,ix)
end if
if (iz > 1) then
tmp = tmp + c_k*Tij * u%s(iz-1,iy,ix)
end if
if ((iz > 1) .and. (iz < grid_param%nz)) then
tmp = tmp - alpha_T(5) * u%s(iz, iy ,ix ) &
+ alpha_T(1) * u%s(iz, iy ,ix+1) &
+ alpha_T(2) * u%s(iz, iy ,ix-1) &
+ alpha_T(3) * u%s(iz, iy+1,ix ) &
+ alpha_T(4) * u%s(iz, iy-1,ix )
end if
v%s(iz,iy,ix) = d_k*tmp
end do
end do
end do
endif
if (LUseT) then

Juan Escobar
committed
call construct_alpha_T_cst_mnh(u%grid_param,alpha_T,Tij)
!-----------------------------------------------------------
iib=u%icompx_min
iie=u%icompx_max
ijb=u%icompy_min
ije=u%icompy_max

Juan Escobar
committed
ize=u%grid_param%nz

Juan Escobar
committed
!
zv_st => v%st
zu_st => u%st
zb_k => vert_coeff%b
zc_k => vert_coeff%c
zd_k => vert_coeff%d

ESCOBAR MUNOZ Juan
committed
nib = Lbound(zv_st,1) ; nie = Ubound(zv_st,1)
njb = Lbound(zv_st,2) ; nje = Ubound(zv_st,2)
nzb = Lbound(zv_st,3) ; nze = Ubound(zv_st,3)
call apply_mnh_dim(zv_st,zu_st,zb_k,zc_k,zd_k)

Juan Escobar
committed

ESCOBAR MUNOZ Juan
committed
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
endif
contains
subroutine apply_mnh_dim(pzv_st,pzu_st,pzb_k,pzc_k,pzd_k)
implicit none
real :: pzv_st(nib:nie,njb:nje,nzb:nze), &
pzu_st(nib:nie,njb:nje,nzb:nze), &
pzb_k(ize),pzc_k(ize),pzd_k(ize)
!$acc kernels present_cr(pzb_k,pzv_st,pzu_st)
iz=1
!$mnh_do_concurrent( ii=iib:iie , ij=ijb:ije )
pzv_st(ii,ij,iz) = pzd_k(iz)* ( (-pzb_k(iz)-pzc_k(iz))*Tij * pzu_st(ii,ij,iz ) &
+pzb_k(iz) *Tij * pzu_st(ii,ij,iz+1) )
!$mnh_end_do()
!
!$mnh_do_concurrent( ii=iib:iie , ij=ijb:ije , iz=2:ize-1 )
pzv_st(ii,ij,iz) = pzd_k(iz)* ( ((-pzb_k(iz)-pzc_k(iz))*Tij - 4.0_rl ) * pzu_st(ii,ij,iz) &
+pzb_k(iz) *Tij * pzu_st(ii,ij,iz+1) &
+pzc_k(iz) *Tij * pzu_st(ii,ij,iz-1) &
+ pzu_st(ii+1,ij,iz) &
+ pzu_st(ii-1,ij,iz) &
+ pzu_st(ii,ij+1,iz) &
+ pzu_st(ii,ij-1,iz) &

ESCOBAR MUNOZ Juan
committed
!$mnh_end_do()

Juan Escobar
committed
!

Juan Escobar
committed
iz=ize

ESCOBAR MUNOZ Juan
committed
!$mnh_do_concurrent( ii=iib:iie , ij=ijb:ije )

ESCOBAR MUNOZ Juan
committed
pzv_st(ii,ij,iz) = pzd_k(iz)* ( (-pzb_k(iz)-pzc_k(iz))*Tij * pzu_st(ii,ij,iz) &
+pzc_k(iz) *Tij * pzu_st(ii,ij,iz-1) )

ESCOBAR MUNOZ Juan
committed
!$mnh_end_do()

ESCOBAR MUNOZ Juan
committed
!$acc end kernels

ESCOBAR MUNOZ Juan
committed
end subroutine apply_mnh_dim

Juan Escobar
committed
end subroutine apply_mnh
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
subroutine apply(u,v)
implicit none
type(scalar3d), intent(in) :: u
type(scalar3d), intent(inout) :: v
real(kind=rl), dimension(5) :: alpha_T
real(kind=rl) :: Tij
real(kind=rl) :: a_k, b_k, c_k, d_k
integer :: ix,iy,iz
real(kind=rl) :: tmp
do ix=u%icompx_min,u%icompx_max
do iy=u%icompy_min,u%icompy_max
! Construct horizontal part of stencil
call construct_alpha_T(u%grid_param, &
ix+u%ix_min-1, &
iy+u%iy_min-1, &
alpha_T,Tij)
do iz=1,u%grid_param%nz
a_k = vert_coeff%a(iz)
b_k = vert_coeff%b(iz)
c_k = vert_coeff%c(iz)
d_k = vert_coeff%d(iz)
tmp = ((a_k-b_k-c_k)*Tij - alpha_T(5)) * u%s(iz,iy,ix)
if (iz < grid_param%nz) then
tmp = tmp + b_k*Tij * u%s(iz+1,iy,ix)
end if
if (iz > 1) then
tmp = tmp + c_k*Tij * u%s(iz-1,iy,ix)
end if
tmp = tmp + alpha_T(1) * u%s(iz, iy ,ix+1) &
+ alpha_T(2) * u%s(iz, iy ,ix-1) &
+ alpha_T(3) * u%s(iz, iy+1,ix ) &
+ alpha_T(4) * u%s(iz, iy-1,ix )
v%s(iz,iy,ix) = d_k*tmp
end do
end do
end do
end subroutine apply
!==================================================================
!==================================================================
!
! S M O O T H E R S
!
!==================================================================
!==================================================================
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
!==================================================================
! Perform nsmooth smoother iterations
!==================================================================
subroutine smooth_mnh(level,m,nsmooth,direction,b,u)
implicit none
integer, intent(in) :: level
integer, intent(in) :: m
integer, intent(in) :: nsmooth ! Number of smoothing steps
integer, intent(in) :: direction ! Direction
type(scalar3d), intent(inout) :: b ! RHS
type(scalar3d), intent(inout) :: u ! solution vector
integer :: i
real(kind=rl) :: log_res_initial, log_res_final
type(scalar3d) :: r
integer :: halo_size
integer :: nlocal, nz
#ifdef MEASURESMOOTHINGRATE
r%ix_min = u%ix_min
r%ix_max = u%ix_max
r%iy_min = u%iy_min
r%iy_max = u%iy_max
r%icompx_min = u%icompx_min
r%icompx_max = u%icompx_max
r%icompy_min = u%icompy_min
r%icompy_max = u%icompy_max
r%halo_size = u%halo_size
r%isactive = u%isactive
r%grid_param = u%grid_param
nlocal = r%ix_max-r%ix_min+1
halo_size = r%halo_size
nz = r%grid_param%nz
if (LUseO) then
allocate(r%s(0:nz+1, &
1-halo_size:nlocal+halo_size, &
1-halo_size:nlocal+halo_size))
end if
if (LUseT) then
allocate(r%st(1-halo_size:nlocal+halo_size, &
1-halo_size:nlocal+halo_size, &
0:nz+1))

ESCOBAR MUNOZ Juan
committed
!$acc enter data create (r%st)
end if
call calculate_residual(level,m,b,u,r)
log_res_initial = log(l2norm(r))
#endif
! Carry out nsmooth iterations of the smoother
if (smoother_param%smoother == SMOOTHER_LINE_SOR) then
do i=1,nsmooth

Juan Escobar
committed
call line_SOR_mnh(level,m,direction,b,u)
end do
else if (smoother_param%smoother == SMOOTHER_LINE_SSOR) then
do i=1,nsmooth

Juan Escobar
committed
call line_SSOR_mnh(level,m,direction,b,u)
end do
else if (smoother_param%smoother == SMOOTHER_LINE_JAC) then
do i=1,nsmooth
call line_jacobi_mnh(level,m,b,u)
end do
end if
#ifdef MEASURESMOOTHINGRATE

Juan Escobar
committed
call calculate_residual_mnh(level,m,b,u,r)
log_res_final = log(l2norm(r))
log_resreduction(level) = log_resreduction(level) &
+ (log_res_final - log_res_initial)
nsmooth_total(level) = nsmooth_total(level) + nsmooth
if (LUseO) deallocate(r%s)

ESCOBAR MUNOZ Juan
committed
if (LUseT) then
!$acc exit data delete(r%st)
deallocate(r%st)
end if
#endif
end subroutine smooth_mnh
!==================================================================
! Perform nsmooth smoother iterations
!==================================================================
subroutine smooth(level,m,nsmooth,direction,b,u)
implicit none
integer, intent(in) :: level
integer, intent(in) :: m
integer, intent(in) :: nsmooth ! Number of smoothing steps
integer, intent(in) :: direction ! Direction
type(scalar3d), intent(inout) :: b ! RHS
type(scalar3d), intent(inout) :: u ! solution vector
integer :: i
real(kind=rl) :: log_res_initial, log_res_final