81
82
87
88# ifdef DISTRIBUTE
89
91# endif
92
93 implicit none
94
95
96
97 integer, intent(in) :: ng, tile
98 integer, intent(in) :: LBi, UBi, LBj, UBj
99 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
100 integer, intent(in) :: nstp, kstp
101
102# ifdef ASSUMED_SHAPE
103 real(r8), intent(in) :: h(LBi:,LBj:)
104 real(r8), intent(in) :: omn(LBi:,LBj:)
105# ifdef SOLVE3D
106 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
107 real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:)
108 real(r8), intent(in) :: tl_z_w(LBi:,LBj:,0:)
109 real(r8), intent(in) :: tl_rho(LBi:,LBj:,:)
110 real(r8), intent(in) :: tl_u(LBi:,LBj:,:,:)
111 real(r8), intent(in) :: tl_v(LBi:,LBj:,:,:)
112# endif
113 real(r8), intent(in) :: tl_ubar(LBi:,LBj:,:)
114 real(r8), intent(in) :: tl_vbar(LBi:,LBj:,:)
115 real(r8), intent(in) :: tl_zeta(LBi:,LBj:,:)
116# else
117 real(r8), intent(in) :: h(LBi:UBi,LBj:UBj)
118 real(r8), intent(in) :: omn(LBi:UBi,LBj:UBj)
119# ifdef SOLVE3D
120 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
121 real(r8), intent(in) :: tl_z_r(LBi:UBi,LBj:UBj,N(ng))
122 real(r8), intent(in) :: tl_z_w(LBi:UBi,LBj:UBj,0:N(ng))
123 real(r8), intent(in) :: tl_rho(LBi:UBi,LBj:UBj,N(ng))
124 real(r8), intent(in) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
125 real(r8), intent(in) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
126# endif
127 real(r8), intent(in) :: tl_ubar(LBi:UBi,LBj:UBj,:)
128 real(r8), intent(in) :: tl_vbar(LBi:UBi,LBj:UBj,:)
129 real(r8), intent(in) :: tl_zeta(LBi:UBi,LBj:UBj,:)
130# endif
131
132
133
134 integer :: NSUB, i, j, k, trd
135 integer :: idia, istep
136
137# ifdef DISTRIBUTE
138 real(r8), dimension(3) :: rbuffer
139 character (len=3), dimension(3) :: op_handle
140# else
141 integer :: my_threadnum
142# endif
143
144 real(r8) :: cff, my_avgke, my_avgpe, my_volume
145
146 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ke2d
147 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: pe2d
148
149 character (len=8 ) :: kechar, pechar
150 character (len=22) :: DateTime
151
152# include "set_bounds.h"
153
154
155
156
157
158
159
160
161
162
163 istep=int(mod(real(
iic(ng)-1,r8),1.0e+10_r8))
164# ifdef SOLVE3D
165 idia=nstp
166# else
167 idia=kstp
168# endif
170
171
172
173 IF (mod(istep,
ninfo(ng)).eq.0)
THEN
174 DO j=jstr,jend
175# ifdef SOLVE3D
176 DO i=istr,iend
177 ke2d(i,j)=0.0_r8
178 pe2d(i,j)=0.5_r8*
g*tl_z_w(i,j,
n(ng))*tl_z_w(i,j,
n(ng))
179 END DO
182 DO i=istr,iend
183 ke2d(i,j)=ke2d(i,j)+ &
184 & tl_hz(i,j,k)* &
185 & 0.25_r8*(tl_u(i ,j,k,idia)*tl_u(i ,j,k,idia)+ &
186 & tl_u(i+1,j,k,idia)*tl_u(i+1,j,k,idia)+ &
187 & tl_v(i,j ,k,idia)*tl_v(i,j ,k,idia)+ &
188 & tl_v(i,j+1,k,idia)*tl_v(i,j+1,k,idia))
189 pe2d(i,j)=pe2d(i,j)+ &
190 & cff*tl_hz(i,j,k)*(tl_rho(i,j,k)+1000.0_r8)* &
191 & (tl_z_r(i,j,k)-tl_z_w(i,j,0))
192 END DO
193 END DO
194# else
196 DO i=istr,iend
197 ke2d(i,j)=(tl_zeta(i,j,idia)+h(i,j))* &
198 & 0.25_r8*(tl_ubar(i ,j,idia)*tl_ubar(i ,j,idia)+ &
199 & tl_ubar(i+1,j,idia)*tl_ubar(i+1,j,idia)+ &
200 & tl_vbar(i,j ,idia)*tl_vbar(i,j ,idia)+ &
201 & tl_vbar(i,j+1,idia)*tl_vbar(i,j+1,idia))
202 pe2d(i,j)=cff*tl_zeta(i,j,idia)*tl_zeta(i,j,idia)
203 END DO
204# endif
205 END DO
206
207
208
209
210
211
212
213
214
215 DO i=istr,iend
216 pe2d(i,jend+1)=0.0_r8
217 pe2d(i,jstr-1)=0.0_r8
218 ke2d(i,jstr-1)=0.0_r8
219 END DO
220 DO j=jstr,jend
221 DO i=istr,iend
222# ifdef SOLVE3D
223
224
225# else
226
227
228# endif
229 pe2d(i,jend+1)=pe2d(i,jend+1)+omn(i,j)*h(i,j)
230 pe2d(i,jstr-1)=pe2d(i,jstr-1)+omn(i,j)*pe2d(i,j)
231 ke2d(i,jstr-1)=ke2d(i,jstr-1)+omn(i,j)*ke2d(i,j)
232 END DO
233 END DO
234 my_volume=0.0_r8
235 my_avgpe=0.0_r8
236 my_avgke=0.0_r8
237 DO i=istr,iend
238 my_volume=my_volume+pe2d(i,jend+1)
239 my_avgpe =my_avgpe +pe2d(i,jstr-1)
240 my_avgke =my_avgke +ke2d(i,jstr-1)
241 END DO
242
243
244
245
246
247
248
249# ifdef DISTRIBUTE
250 nsub=1
251# else
252 IF (
domain(ng)%SouthWest_Corner(tile).and. &
253 &
domain(ng)%NorthEast_Corner(tile))
THEN
254 nsub=1
255 ELSE
257 END IF
258# endif
259
264 END IF
271# ifdef DISTRIBUTE
275 op_handle(1)='SUM'
276 op_handle(2)='SUM'
277 op_handle(3)='SUM'
283# else
285# endif
289
290
291
294 IF (
master.and.(ng.eq.1))
THEN
295 WRITE (
stdout,10)
'TIME-STEP',
'YYYY-MM-DD hh:mm:ss.ss', &
296 & 'KINETIC_ENRG', 'POTEN_ENRG', &
297# ifdef NESTING
298 & 'TOTAL_ENRG', 'NET_VOLUME', 'Grid'
299# else
300 & 'TOTAL_ENRG', 'NET_VOLUME'
301# endif
302# ifdef NESTING
303 10 FORMAT (/,1x,a,1x,a,2x,a,3x,a,4x,a,4x,a,2x,a)
304# else
305 10 FORMAT (/,1x,a,1x,a,2x,a,3x,a,4x,a,4x,a)
306# endif
307 END IF
308 END IF
309
311 WRITE (
stdout,20) istep, datetime, &
312# ifdef NESTING
314# else
316# endif
317# ifdef NESTING
318 20 FORMAT (i10,1x,a,4(1pe14.6),2x,i2.2)
319# else
320 20 FORMAT (i10,1x,a,4(1pe14.6))
321# endif
323 END IF
324
325
326
327 WRITE (kechar,
'(1pe8.1)')
avgke
328 WRITE (pechar,
'(1pe8.1)')
avgpe
329 DO i=1,8
330 IF ((kechar(i:i).eq.'N').or.(pechar(i:i).eq.'N').or. &
331 & (kechar(i:i).eq.'n').or.(pechar(i:i).eq.'n').or. &
332 & (kechar(i:i).eq.'*').or.(pechar(i:i).eq.'*')) THEN
334 END IF
335 END DO
336 END IF
337
338 END IF
339
340 RETURN
integer function my_threadnum()
integer, dimension(:), allocatable n
integer, dimension(:), allocatable ntilex
type(t_domain), dimension(:), allocatable domain
integer, dimension(:), allocatable ntilee
integer, dimension(:), allocatable ninfo
integer, dimension(:), allocatable iic
integer, dimension(:), allocatable first_time
character(len=22), dimension(:), allocatable time_code