ROMS
Loading...
Searching...
No Matches
shapiro_mod Module Reference

Functions/Subroutines

subroutine shapiro2d_tile (ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, amask, a)
 
subroutine shapiro3d_tile (ng, tile, model, lbi, ubi, lbj, ubj, lbk, ubk, imins, imaxs, jmins, jmaxs, amask, a)
 

Function/Subroutine Documentation

◆ shapiro2d_tile()

subroutine shapiro_mod::shapiro2d_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
real(r8), dimension(lbi:,lbj:), intent(in) amask,
real(r8), dimension(lbi:,lbj:), intent(inout) a )

Definition at line 26 of file shapiro.F.

33!***********************************************************************
34!
35 USE mod_param
36!
37! Imported variable declarations.
38!
39 integer, intent(in) :: ng, tile, model
40 integer, intent(in) :: LBi, UBi, LBj, UBj
41 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
42
43#ifdef ASSUMED_SHAPE
44# ifdef MASKING
45 real(r8), intent(in) :: Amask(LBi:,LBj:)
46# endif
47 real(r8), intent(inout) :: A(LBi:,LBj:)
48#else
49# ifdef MASKING
50 real(r8), intent(in) :: Amask(LBi:UBi,LBj:UBj)
51# endif
52 real(r8), intent(inout) :: A(LBi:UBi,LBj:UBj)
53#endif
54!
55! Local variable declarations.
56!
57 integer :: i, j
58!
59 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Awrk1
60 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Awrk2
61
62#include "set_bounds.h"
63!
64!-----------------------------------------------------------------------
65! Shapiro filter requested 2D field.
66!-----------------------------------------------------------------------
67!
68! This subroutine will apply a Shapiro filter of order 2 (defined
69! as twice the order in Shapiro (1970), with N even) to an array, A.
70! The order of the filter is reduced at the boundaries and at the
71! mask edges, if any.
72!
73! Initialize filter in the Y-direction.
74!
75 DO j=jstr,jend
76 DO i=istr-1,iend+1
77#ifdef MASKING
78 awrk1(i,j)=0.25_r8* &
79 & (a(i,j-1)*amask(i,j-1)+ &
80 & a(i,j+1)*amask(i,j+1)- &
81 & 2.0_r8*a(i,j)*amask(i,j))* &
82 & amask(i,j-1)*amask(i,j+1)*amask(i,j)
83#else
84 awrk1(i,j)=0.25_r8* &
85 & (a(i,j-1)+a(i,j+1)-2.0_r8*a(i,j))
86#endif
87 END DO
88 END DO
89!
90! Add the changes to the field.
91!
92 DO j=jstr,jend
93 DO i=istr-1,iend+1
94 awrk2(i,j)=a(i,j)+awrk1(i,j)
95 END DO
96 END DO
97!
98! Initialize filter in the X-direction.
99!
100 DO j=jstr,jend
101 DO i=istr,iend
102#ifdef MASKING
103 awrk1(i,j)=0.25_r8* &
104 & (awrk2(i-1,j)*amask(i-1,j)+ &
105 & awrk2(i+1,j)*amask(i+1,j)- &
106 & 2.0_r8*awrk2(i,j)*amask(i,j))* &
107 & amask(i-1,j)*amask(i+1,j)*amask(i,j)
108#else
109 awrk1(i,j)=0.25_r8* &
110 & (awrk2(i-1,j)+awrk2(i+1,j)-2.0_r8*awrk2(i,j))
111#endif
112 END DO
113 END DO
114!
115! Add changes to field.
116!
117 DO j=jstr,jend
118 DO i=istr,iend
119 a(i,j)=awrk2(i,j)+awrk1(i,j)
120 END DO
121 END DO
122!
123 RETURN

Referenced by lmd_bkpp_tile(), lmd_skpp_tile(), regrid_mod::regrid_nf90(), and regrid_mod::regrid_pio().

Here is the caller graph for this function:

◆ shapiro3d_tile()

subroutine shapiro_mod::shapiro3d_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbk,
integer, intent(in) ubk,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
real(r8), dimension(lbi:,lbj:), intent(in) amask,
real(r8), dimension(lbi:,lbj:,lbk:), intent(inout) a )

Definition at line 129 of file shapiro.F.

136!***********************************************************************
137!
138 USE mod_param
139!
140! Imported variable declarations.
141!
142 integer, intent(in) :: ng, tile, model
143 integer, intent(in) :: LBi, UBi, LBj, UBj, LBk, UBk
144 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
145!
146# ifdef ASSUMED_SHAPE
147# ifdef MASKING
148 real(r8), intent(in) :: Amask(LBi:,LBj:)
149# endif
150 real(r8), intent(inout) :: A(LBi:,LBj:,LBk:)
151# else
152# ifdef MASKING
153 real(r8), intent(in) :: Amask(LBi:UBi,LBj:UBj)
154# endif
155 real(r8), intent(inout) :: A(LBi:UBi,LBj:UBj,LBk:UBk)
156# endif
157!
158! Local variable declarations.
159!
160 integer :: i, j, k
161!
162 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Awrk1
163 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Awrk2
164
165# include "set_bounds.h"
166!
167!-----------------------------------------------------------------------
168! Shapiro filter requested 3D field.
169!-----------------------------------------------------------------------
170!
171! This subroutine will apply a Shapiro filter of order 2 (defined
172! as twice the order in Shapiro (1970), with N even) to an array, A.
173! The order of the filter is reduced at the boundaries and at the
174! mask edges, if any.
175!
176! Initialize filter in the Y-direction.
177!
178 DO k=lbk,ubk
179 DO j=jstr,jend
180 DO i=istr-1,iend+1
181# ifdef MASKING
182 awrk1(i,j)=0.25_r8* &
183 & (a(i,j-1,k)*amask(i,j-1)+ &
184 & a(i,j+1,k)*amask(i,j+1)- &
185 & 2.0_r8*a(i,j,k)*amask(i,j))* &
186 & amask(i,j-1)*amask(i,j+1)*amask(i,j)
187# else
188 awrk1(i,j)=0.25_r8* &
189 & (a(i,j-1,k)+a(i,j+1,k)-2.0_r8*a(i,j,k))
190# endif
191 END DO
192 END DO
193!
194! Add the changes to the field.
195!
196 DO j=jstr,jend
197 DO i=istr-1,iend+1
198 awrk2(i,j)=a(i,j,k)+awrk1(i,j)
199 END DO
200 END DO
201!
202! Initialize filter in the X-direction.
203!
204 DO j=jstr,jend
205 DO i=istr,iend
206# ifdef MASKING
207 awrk1(i,j)=0.25_r8* &
208 & (awrk2(i-1,j)*amask(i-1,j)+ &
209 & awrk2(i+1,j)*amask(i+1,j)- &
210 & 2.0_r8*awrk2(i,j)*amask(i,j))* &
211 & amask(i-1,j)*amask(i+1,j)*amask(i,j)
212# else
213 awrk1(i,j)=0.25_r8* &
214 & (awrk2(i-1,j)+awrk2(i+1,j)-2.0_r8*awrk2(i,j))
215# endif
216 END DO
217 END DO
218!
219! Add changes to field.
220!
221 DO j=jstr,jend
222 DO i=istr,iend
223 a(i,j,k)=awrk2(i,j)+awrk1(i,j)
224 END DO
225 END DO
226 END DO
227!
228 RETURN