-
Notifications
You must be signed in to change notification settings - Fork 3
/
rcfl_x.f90
119 lines (95 loc) · 4.42 KB
/
rcfl_x.f90
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
subroutine rcfl_x( im, jm, km, imax, al, bt, fld, inx, imx)
!---------------------------------------------------------------------------
! !
! Copyright 2018 Anna Teruzzi, OGS, Trieste !
! !
! This file is part of 3DVarBio.
! 3DVarBio is based on OceanVar (Dobricic, 2006) !
! !
! 3DVarBio is free software: you can redistribute it and/or modify. !
! it under the terms of the GNU General Public License as published by !
! the Free Software Foundation, either version 3 of the License, or !
! (at your option) any later version. !
! !
! 3DVarBio 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_rcx PARTICULAR PURPOSE. See the !
! GNU General Public License for more details. !
! !
! You should have received a_rcx copy of the GNU General Public License !
! along with OceanVar. If not, see <http://www.gnu.org/licenses/>. !
! !
!---------------------------------------------------------------------------
!-----------------------------------------------------------------------
! !
! Recursive filter in x direction !
! !
! Version 1: A. Teruzzi 2018 !
!-----------------------------------------------------------------------
use set_knd
use cns_str
use rcfl
use grd_str
use mpi_str
implicit none
INTEGER(i4) :: im, jm, km, imax
REAL(r8) :: fld(im,jm,km)
REAL(r8) :: al(jm,imax,km), bt(jm,imax,km)
INTEGER(i4) :: inx(im,jm,km), imx(km)
INTEGER(i4) :: i,j,k, ktr
INTEGER(i4) :: indSupWP
INTEGER nthreads, tid
integer :: OMP_GET_NUM_THREADS, OMP_GET_THREAD_NUM
tid=1
!$OMP PARALLEL &
!$OMP PRIVATE(k,j,i,ktr,indSupWP,tid)
!$ tid = OMP_GET_THREAD_NUM()+1
!$OMP DO
do k=1,km
a_rcx(:,:,tid) = 0.0
b_rcx(:,:,tid) = 0.0
c_rcx(:,:,tid) = 0.0
do j=1,jm
do i=1,im
a_rcx(j,inx(i,j,k),tid) = fld(i,j,k)
enddo
enddo
alp_rcx(:,:,tid) = al(:,:,k)
bta_rcx(:,:,tid) = bt(:,:,k)
do ktr = 1,rcf%ntr
! positive direction
if( ktr.eq.1 )then
b_rcx(:,1,tid) = (1.-alp_rcx(:,1,tid)) * a_rcx(:,1,tid)
elseif( ktr.eq.2 )then
b_rcx(:,1,tid) = a_rcx(:,1,tid) / (1.+alp_rcx(:,1,tid))
else
b_rcx(:,1,tid) = (1.-alp_rcx(:,1,tid)) * (a_rcx(:,1,tid)-alp_rcx(:,1,tid)**3 * a_rcx(:,2,tid)) / (1.-alp_rcx(:,1,tid)**2)**2
endif
do j=2,imx(k)
b_rcx(:,j,tid) = alp_rcx(:,j,tid)*b_rcx(:,j-1,tid) + (1.-alp_rcx(:,j,tid))*a_rcx(:,j,tid)
enddo
! negative direction
if( ktr.eq.1 )then
c_rcx(:,imx(k),tid) = b_rcx(:,imx(k),tid) / (1.+bta_rcx(:,imx(k),tid))
else
c_rcx(:,imx(k),tid) = (1.-bta_rcx(:,imx(k),tid)) * &
(b_rcx(:,imx(k),tid)-bta_rcx(:,imx(k),tid)**3 * b_rcx(:,imx(k)-1,tid)) / (1.-bta_rcx(:,imx(k),tid)**2)**2
endif
do j=imx(k)-1,1,-1
c_rcx(:,j,tid) = bta_rcx(:,j,tid)*c_rcx(:,j+1,tid) + (1.-bta_rcx(:,j,tid))*b_rcx(:,j,tid)
enddo
a_rcx(:,:,tid) = c_rcx(:,:,tid)
enddo
! This way fills land points with some values.
! We prefer not investigate at the mooment and use only the water points
do j=1,localCol
do i=1,GlobalRow
if(grd%global_msk(i,j + GlobalColOffset,1).eq.1) then
fld(i,j,k) = a_rcx(j,inx(i,j,k),tid)
end if
end do
end do
enddo
!$OMP END DO
!$OMP END PARALLEL
end subroutine rcfl_x