82
83
88
89# ifdef DISTRIBUTE
90
92# endif
93
94 implicit none
95
96
97
98 integer, intent(in) :: ng, tile
99 integer, intent(in) :: LBi, UBi, LBj, UBj
100 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
101 integer, intent(in) :: nstp, kstp
102
103# ifdef ASSUMED_SHAPE
104 real(r8), intent(in) :: h(LBi:,LBj:)
105 real(r8), intent(in) :: omn(LBi:,LBj:)
106# ifdef SOLVE3D
107 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
108 real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:)
109 real(r8), intent(in) :: tl_z_w(LBi:,LBj:,0:)
110 real(r8), intent(in) :: tl_rho(LBi:,LBj:,:)
111 real(r8), intent(in) :: tl_u(LBi:,LBj:,:,:)
112 real(r8), intent(in) :: tl_v(LBi:,LBj:,:,:)
113# endif
114 real(r8), intent(in) :: tl_ubar(LBi:,LBj:,:)
115 real(r8), intent(in) :: tl_vbar(LBi:,LBj:,:)
116 real(r8), intent(in) :: tl_zeta(LBi:,LBj:,:)
117# else
118 real(r8), intent(in) :: h(LBi:UBi,LBj:UBj)
119 real(r8), intent(in) :: omn(LBi:UBi,LBj:UBj)
120# ifdef SOLVE3D
121 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
122 real(r8), intent(in) :: tl_z_r(LBi:UBi,LBj:UBj,N(ng))
123 real(r8), intent(in) :: tl_z_w(LBi:UBi,LBj:UBj,0:N(ng))
124 real(r8), intent(in) :: tl_rho(LBi:UBi,LBj:UBj,N(ng))
125 real(r8), intent(in) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
126 real(r8), intent(in) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
127# endif
128 real(r8), intent(in) :: tl_ubar(LBi:UBi,LBj:UBj,:)
129 real(r8), intent(in) :: tl_vbar(LBi:UBi,LBj:UBj,:)
130 real(r8), intent(in) :: tl_zeta(LBi:UBi,LBj:UBj,:)
131# endif
132
133
134
135 integer :: NSUB, i, j, k, trd
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
151# include "set_bounds.h"
152
153
154
155
156
157
158
159 IF (mod(
iic(ng)-1,
ninfo(ng)).eq.0)
THEN
160 DO j=jstr,jend
161# ifdef SOLVE3D
162 DO i=istr,iend
163 ke2d(i,j)=0.0_r8
164 pe2d(i,j)=0.5_r8*
g*tl_z_w(i,j,
n(ng))*tl_z_w(i,j,
n(ng))
165 END DO
168 DO i=istr,iend
169 ke2d(i,j)=ke2d(i,j)+ &
170 & tl_hz(i,j,k)* &
171 & 0.25_r8*(tl_u(i ,j,k,nstp)*tl_u(i ,j,k,nstp)+ &
172 & tl_u(i+1,j,k,nstp)*tl_u(i+1,j,k,nstp)+ &
173 & tl_v(i,j ,k,nstp)*tl_v(i,j ,k,nstp)+ &
174 & tl_v(i,j+1,k,nstp)*tl_v(i,j+1,k,nstp))
175 pe2d(i,j)=pe2d(i,j)+ &
176 & cff*tl_hz(i,j,k)*(tl_rho(i,j,k)+1000.0_r8)* &
177 & (tl_z_r(i,j,k)-tl_z_w(i,j,0))
178 END DO
179 END DO
180# else
182 DO i=istr,iend
183 ke2d(i,j)=(tl_zeta(i,j,kstp)+h(i,j))* &
184 & 0.25_r8*(tl_ubar(i ,j,kstp)*tl_ubar(i ,j,kstp)+ &
185 & tl_ubar(i+1,j,kstp)*tl_ubar(i+1,j,kstp)+ &
186 & tl_vbar(i,j ,kstp)*tl_vbar(i,j ,kstp)+ &
187 & tl_vbar(i,j+1,kstp)*tl_vbar(i,j+1,kstp))
188 pe2d(i,j)=cff*tl_zeta(i,j,kstp)*tl_zeta(i,j,kstp)
189 END DO
190# endif
191 END DO
192
193
194
195
196
197
198
199
200
201 DO i=istr,iend
202 pe2d(i,jend+1)=0.0_r8
203 pe2d(i,jstr-1)=0.0_r8
204 ke2d(i,jstr-1)=0.0_r8
205 END DO
206 DO j=jstr,jend
207 DO i=istr,iend
208# ifdef SOLVE3D
209
210
211# else
212
213
214# endif
215 pe2d(i,jend+1)=pe2d(i,jend+1)+omn(i,j)*h(i,j)
216 pe2d(i,jstr-1)=pe2d(i,jstr-1)+omn(i,j)*pe2d(i,j)
217 ke2d(i,jstr-1)=ke2d(i,jstr-1)+omn(i,j)*ke2d(i,j)
218 END DO
219 END DO
220 my_volume=0.0_r8
221 my_avgpe=0.0_r8
222 my_avgke=0.0_r8
223 DO i=istr,iend
224 my_volume=my_volume+pe2d(i,jend+1)
225 my_avgpe =my_avgpe +pe2d(i,jstr-1)
226 my_avgke =my_avgke +ke2d(i,jstr-1)
227 END DO
228
229
230
231
232
233
234
235# ifdef DISTRIBUTE
236 nsub=1
237# else
238 IF (
domain(ng)%SouthWest_Corner(tile).and. &
239 &
domain(ng)%NorthEast_Corner(tile))
THEN
240 nsub=1
241 ELSE
243 END IF
244# endif
245
250 END IF
257# ifdef DISTRIBUTE
261 op_handle(1)='SUM'
262 op_handle(2)='SUM'
263 op_handle(3)='SUM'
269# else
271# endif
277 IF (
master.and.(ng.eq.1))
THEN
278 WRITE (
stdout,10)
'TIME-STEP',
'YYYY-MM-DD hh:mm:ss.ss', &
279 & 'KINETIC_ENRG', 'POTEN_ENRG', &
280# ifdef NESTING
281 & 'TOTAL_ENRG', 'NET_VOLUME', 'Grid'
282# else
283 & 'TOTAL_ENRG', 'NET_VOLUME'
284# endif
285# ifdef NESTING
286 10 FORMAT (/,1x,a,1x,a,2x,a,3x,a,4x,a,4x,a,2x,a)
287# else
288 10 FORMAT (/,1x,a,1x,a,2x,a,3x,a,4x,a,4x,a)
289# endif
290 END IF
291 END IF
293 WRITE (
stdout,20) int(mod(real(
iic(ng)-1,r8),1.0e+10_r8)), &
295# ifdef NESTING
297# else
299# endif
300# ifdef NESTING
301 20 FORMAT (i10,1x,a,4(1pe14.6),2x,i2.2)
302# else
303 20 FORMAT (i10,1x,a,4(1pe14.6))
304# endif
306 END IF
307
308
309
310 WRITE (kechar,
'(1pe8.1)')
avgke
311 WRITE (pechar,
'(1pe8.1)')
avgpe
312 DO i=1,8
313 IF ((kechar(i:i).eq.'N').or.(pechar(i:i).eq.'N').or. &
314 & (kechar(i:i).eq.'n').or.(pechar(i:i).eq.'n').or. &
315 & (kechar(i:i).eq.'*').or.(pechar(i:i).eq.'*')) THEN
317 END IF
318 END DO
319 END IF
320
321 END IF
322
323 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