-
Notifications
You must be signed in to change notification settings - Fork 23
/
psmoo.F90
297 lines (294 loc) · 8.69 KB
/
psmoo.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
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
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
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
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
subroutine psmooth(a,margin_a,margin_smooth, imask, w)
use mod_xc ! HYCOM communication interface
use mod_cb_arrays ! HYCOM saved arrays
implicit none
!
integer margin_a,margin_smooth
integer imask(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
real a(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy), &
w(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
!
! --- ragged boundary version of basic 9-point smoothing routine.
! --- this routine is set up to smooth data carried at -p- points.
!
! --- margin_a is the margin on entry, and
! --- margin_smooth is the margin on exit.
! --- imask is ip or ishlf.
!
! --- see also psmooth_ice and psmooth_dif.
!
! --- w used as workspace, overwritten on exit
!
integer i,ismth,j,jsmth,msmth
real qc,sh
!
real c(-1:1,-1:1)
save c
data c / 1.0, 2.0, 1.0, &
2.0, 4.0, 2.0, &
1.0, 2.0, 1.0 /
!
qc = 1.0/sum(c(:,:))
!
msmth = min(margin_smooth,nbdy-1)
!
if (margin_a.lt.msmth+1) then
! --- update the halo
call xctilr(a,1,1, msmth+1,msmth+1, halo_ps)
endif
!
!$OMP PARALLEL DO PRIVATE(j,i) &
!$OMP SCHEDULE(STATIC)
do j=0-msmth,jj+msmth+1
do i=0-msmth,ii+msmth+1
w(i,j) = a(i,j)
enddo
enddo
!$OMP END PARALLEL DO
!
!$OMP PARALLEL DO PRIVATE(j,i,sh,jsmth,ismth) &
!$OMP SCHEDULE(STATIC)
do j=1-msmth,jj+msmth
do i=1-msmth,ii+msmth
if (imask(i,j).eq.1) then
sh = 0.0
do jsmth= -1,1
do ismth= -1,1
if (imask(i+ismth,j+jsmth).eq.1) then
sh = sh + c(ismth,jsmth)*w(i+ismth,j+jsmth)
else
sh = sh + c(ismth,jsmth)*w(i, j)
endif
enddo
enddo
a(i,j) = sh*qc
endif !imask.eq.1
enddo
enddo
!$OMP END PARALLEL DO
return
end subroutine psmooth
subroutine psmooth_new(a,b,margin_a,margin_smooth, imask, w)
use mod_xc ! HYCOM communication interface
use mod_cb_arrays ! HYCOM saved arrays
implicit none
!
integer margin_a,margin_smooth
integer imask(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
real a(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy), &
b(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy), &
w(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
!
! --- ragged boundary version of basic 9-point smoothing routine.
! --- this routine is set up to smooth data carried at -p- points.
! --- input in a, output in b.
!
! --- margin_a is the margin on entry, and
! --- margin_smooth is the margin on exit.
! --- imask is ip or ishlf.
!
! --- see also psmooth.
!
! --- a and b must not be the same array.
!
integer i,ismth,j,jsmth,msmth
real qc,sh
!
real c(-1:1,-1:1)
save c
data c / 1.0, 2.0, 1.0, &
2.0, 4.0, 2.0, &
1.0, 2.0, 1.0 /
!
qc = 1.0/sum(c(:,:))
!
msmth = min(margin_smooth,nbdy-1)
!
if (margin_a.lt.msmth+1) then
! --- update the halo
call xctilr(a,1,1, msmth+1,msmth+1, halo_ps)
endif
!
!$OMP PARALLEL DO PRIVATE(j,i,sh,jsmth,ismth) &
!$OMP SCHEDULE(STATIC)
do j=1-msmth,jj+msmth
do i=1-msmth,ii+msmth
if (imask(i,j).eq.1) then
sh = 0.0
do jsmth= -1,1
do ismth= -1,1
if (imask(i+ismth,j+jsmth).eq.1) then
sh = sh + c(ismth,jsmth)*a(i+ismth,j+jsmth)
else
sh = sh + c(ismth,jsmth)*a(i, j)
endif
enddo
enddo
b(i,j) = sh*qc
endif !imask.eq.1
enddo
enddo
!$OMP END PARALLEL DO
return
end subroutine psmooth_new
subroutine psmooth_ice(a,margin_a,margin_smooth, imask, w)
use mod_xc ! HYCOM communication interface
use mod_cb_arrays ! HYCOM saved arrays
implicit none
!
integer margin_a,margin_smooth
integer imask(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
real a(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy), &
w(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
!
! --- ragged boundary version of basic 9-point smoothing routine.
! --- this routine is set up to smooth data carried at -p- points.
! --- it also smooths covice=1.0 and covice=0.0 regions separately,
! --- leaving areas with fractional covice untouched. note that
! --- covice must be valid in the halo out to margin_smooth.
!
! --- margin_a is the margin on entry, and
! --- margin_smooth is the margin on exit.
! --- imask is ip or ishlf.
!
! --- see also psmooth and psmooth_dif
!
! --- array a can't be covice
! --- w used as workspace, overwritten on exit
!
integer i,ismth,j,jsmth,msmth
real qc,sh,ci
!
real c(-1:1,-1:1)
save c
data c / 1.0, 2.0, 1.0, &
2.0, 4.0, 2.0, &
1.0, 2.0, 1.0 /
!
qc = 1.0/sum(c(:,:))
!
msmth = min(margin_smooth,nbdy-1)
!
if (margin_a.lt.msmth+1) then
! --- update the halo
call xctilr(a,1,1, msmth+1,msmth+1, halo_ps)
endif
!
!$OMP PARALLEL DO PRIVATE(j,i) &
!$OMP SCHEDULE(STATIC)
do j=0-msmth,jj+msmth+1
do i=0-msmth,ii+msmth+1
w(i,j) = a(i,j)
enddo
enddo
!$OMP END PARALLEL DO
!
!$OMP PARALLEL DO PRIVATE(j,i,sh,ci,jsmth,ismth) &
!$OMP SCHEDULE(STATIC)
do j=1-msmth,jj+msmth
do i=1-msmth,ii+msmth
if (imask(i,j).eq.1) then
ci = covice(i,j)
if (ci.eq.0.0 .or. &
ci.eq.1.0 ) then !full sea or full ice
sh = 0.0
do jsmth= -1,1
do ismth= -1,1
if ( imask(i+ismth,j+jsmth).eq.1 .and. &
covice(i+ismth,j+jsmth).eq.ci ) then
sh = sh + c(ismth,jsmth)*w(i+ismth,j+jsmth)
else
sh = sh + c(ismth,jsmth)*w(i, j)
endif
enddo
enddo
a(i,j) = sh*qc
endif !full sea or full ice
endif !imask.eq.1
enddo
enddo
!$OMP END PARALLEL DO
return
end subroutine psmooth_ice
subroutine psmooth_dif(a,aklist,k,margin_a,margin_smooth, &
imask, w)
use mod_xc ! HYCOM communication interface
use mod_cb_arrays ! HYCOM saved arrays
implicit none
!
integer k,margin_a,margin_smooth
integer imask( 1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
real a( 1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy), &
aklist(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy), &
w( 1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
!
! --- ragged boundary version of basic 9-point smoothing routine.
! --- this routine is set up to smooth vcty carried at -p- points.
! --- it return the maximum of the original and smoothed value and
! --- ignores locations where k > aklist(i,j).
!
! --- margin_a is the margin on entry, and
! --- margin_smooth is the margin on exit.
! --- imask is ip or ishlf.
!
! --- see also psmooth and psmooth_ice.
!
! --- w used as workspace, overwritten on exit
! --- assumes that aklist's halo is valid out to msmth+1.
!
integer i,ismth,j,jsmth,msmth
real qc,sh
!
real c(-1:1,-1:1)
save c
data c / 1.0, 2.0, 1.0, &
2.0, 4.0, 2.0, &
1.0, 2.0, 1.0 /
!
qc = 1.0/sum(c(:,:))
!
msmth = min(margin_smooth,nbdy-1)
!
if (margin_a.lt.msmth+1) then
! --- update the halo
call xctilr(a,1,1, msmth+1,msmth+1, halo_ps)
endif
!
!$OMP PARALLEL DO PRIVATE(j,i) &
!$OMP SCHEDULE(STATIC)
do j=0-msmth,jj+msmth+1
do i=0-msmth,ii+msmth+1
w(i,j) = a(i,j)
enddo
enddo
!$OMP END PARALLEL DO
!
!$OMP PARALLEL DO PRIVATE(j,i,sh,jsmth,ismth) &
!$OMP SCHEDULE(STATIC)
do j=1-msmth,jj+msmth
do i=1-msmth,ii+msmth
if ( imask(i,j).eq.1 .and. &
aklist(i,j).ge.k ) then
sh = 0.0
do jsmth= -1,1
do ismth= -1,1
if ( imask(i+ismth,j+jsmth).eq.1 .and. &
aklist(i+ismth,j+jsmth).ge.k ) then
sh = sh + c(ismth,jsmth)*w(i+ismth,j+jsmth)
else
sh = sh + c(ismth,jsmth)*w(i, j)
endif
enddo
enddo
a(i,j) = max( a(i,j), sh*qc )
endif !imask.eq.1
enddo
enddo
!$OMP END PARALLEL DO
return
end subroutine psmooth_dif
!
!
!> Revision history:
!>
!> Apr. 2014 - added imask