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

Functions/Subroutines

subroutine, public tl_diag (ng, tile)
 
subroutine tl_diag_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nstp, kstp, h, omn, tl_hz, tl_z_r, tl_z_w, tl_rho, tl_u, tl_v, tl_ubar, tl_vbar, tl_zeta)
 

Function/Subroutine Documentation

◆ tl_diag()

subroutine, public tl_diag_mod::tl_diag ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 24 of file tl_diag.F.

25!***********************************************************************
26!
27 USE mod_param
28 USE mod_grid
29 USE mod_ocean
30 USE mod_stepping
31!
32! Imported variable declarations.
33!
34 integer, intent(in) :: ng, tile
35!
36! Local variable declarations.
37!
38 character (len=*), parameter :: MyFile = &
39 & __FILE__
40!
41# include "tile.h"
42!
43# ifdef PROFILE
44 CALL wclock_on (ng, itlm, 7, __line__, myfile)
45# endif
46 CALL tl_diag_tile (ng, tile, &
47 & lbi, ubi, lbj, ubj, &
48 & imins, imaxs, jmins, jmaxs, &
49 & nstp(ng), kstp(ng), &
50 & grid(ng) % h, &
51 & grid(ng) % omn, &
52# ifdef SOLVE3D
53 & grid(ng) % tl_Hz, &
54 & grid(ng) % tl_z_r, &
55 & grid(ng) % tl_z_w, &
56 & ocean(ng) % tl_rho, &
57 & ocean(ng) % tl_u, &
58 & ocean(ng) % tl_v, &
59# endif
60 & ocean(ng) % tl_ubar, &
61 & ocean(ng) % tl_vbar, &
62 & ocean(ng) % tl_zeta)
63# ifdef PROFILE
64 CALL wclock_off (ng, itlm, 7, __line__, myfile)
65# endif
66!
67 RETURN
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable nstp
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_grid::grid, mod_param::itlm, mod_stepping::kstp, mod_stepping::nstp, mod_ocean::ocean, tl_diag_tile(), wclock_off(), and wclock_on().

Referenced by tl_main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ tl_diag_tile()

subroutine tl_diag_mod::tl_diag_tile ( integer, intent(in) ng,
integer, intent(in) tile,
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,
integer, intent(in) nstp,
integer, intent(in) kstp,
real(r8), dimension(lbi:,lbj:), intent(in) h,
real(r8), dimension(lbi:,lbj:), intent(in) omn,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_hz,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_z_r,
real(r8), dimension(lbi:,lbj:,0:), intent(in) tl_z_w,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_rho,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_v,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_zeta )
private

Definition at line 71 of file tl_diag.F.

81!***********************************************************************
82!
83 USE mod_param
84 USE mod_parallel
85 USE mod_iounits
86 USE mod_scalars
87
88# ifdef DISTRIBUTE
89!
90 USE distribute_mod, ONLY : mp_reduce
91# endif
92!
93 implicit none
94!
95! Imported variable declarations.
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! Local variable declarations.
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! Compute and report out volume averaged kinetic, potential total
156! energy, and volume of tangent linear fields. Notice that the
157! proper adjoint of these quantities are not coded.
158!-----------------------------------------------------------------------
159!
160! Set time timestep counter and time level indices to process. Restart
161! counter after 10 billion steps.
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
169 datetime=time_code(ng)
170!
171! Compute kinetic and potential energy.
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
180 cff=g/rho0
181 DO k=n(ng),1,-1
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
195 cff=0.5_r8*g
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! Integrate horizontally within one tile. In order to reduce the
208! round-off errors, the summation is performed in two stages. First,
209! the index j is collapsed and then the accumulation is carried out
210! along index i. In this order, the partial sums consist on much
211! fewer number of terms than in a straightforward two-dimensional
212! summation. Thus, adding numbers which are orders of magnitude
213! apart is avoided.
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!! pe2d(i,Jend+1)=pe2d(i,Jend+1)+ &
224!! & omn(i,j)*(z_w(i,j,N(ng))-z_w(i,j,0))
225# else
226!! pe2d(i,Jend+1)=pe2d(i,Jend+1)+ &
227!! & omn(i,j)*(zeta(i,j,idia)+h(i,j))
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! Perform global summation: whoever gets first to the critical region
244! resets global sums before global summation starts; after the global
245! summation is completed, thread, which is the last one to enter the
246! critical region, finalizes the computation of diagnostics and prints
247! them out.
248!
249# ifdef DISTRIBUTE
250 nsub=1 ! distributed-memory
251# else
252 IF (domain(ng)%SouthWest_Corner(tile).and. &
253 & domain(ng)%NorthEast_Corner(tile)) THEN
254 nsub=1 ! non-tiled application
255 ELSE
256 nsub=ntilex(ng)*ntilee(ng) ! tiled application
257 END IF
258# endif
259!$OMP CRITICAL (TL_DIAGNOSTICS)
260 IF (tile_count.eq.0) THEN
261 volume=0.0_r8
262 avgke=0.0_r8
263 avgpe=0.0_r8
264 END IF
265 volume=volume+my_volume
266 avgke=avgke+my_avgke
267 avgpe=avgpe+my_avgpe
269 IF (tile_count.eq.nsub) THEN
270 tile_count=0
271# ifdef DISTRIBUTE
272 rbuffer(1)=volume
273 rbuffer(2)=avgke
274 rbuffer(3)=avgpe
275 op_handle(1)='SUM'
276 op_handle(2)='SUM'
277 op_handle(3)='SUM'
278 CALL mp_reduce (ng, itlm, 3, rbuffer, op_handle)
279 volume=rbuffer(1)
280 avgke=rbuffer(2)
281 avgpe=rbuffer(3)
282 trd=mymaster
283# else
284 trd=my_threadnum()
285# endif
289!
290! Report global run diagnotics values for the tangent linear kernel.
291!
292 IF (first_time(ng).eq.0) THEN
293 first_time(ng)=1
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!
310 IF (master) THEN
311 WRITE (stdout,20) istep, datetime, &
312# ifdef NESTING
313 & avgke, avgpe, avgkp, volume, ng
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
322 FLUSH (stdout)
323 END IF
324!
325! If blowing-up, set exit_flag to stop computations.
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
333 exit_flag=1
334 END IF
335 END DO
336 END IF
337!$OMP END CRITICAL (TL_DIAGNOSTICS)
338 END IF
339!
340 RETURN
integer function my_threadnum()
integer stdout
logical master
integer mymaster
integer tile_count
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer, dimension(:), allocatable ntilex
Definition mod_param.F:685
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable ntilee
Definition mod_param.F:686
integer, dimension(:), allocatable ninfo
integer, dimension(:), allocatable iic
real(dp) avgke
real(dp) avgpe
real(dp) volume
integer, dimension(:), allocatable first_time
character(len=22), dimension(:), allocatable time_code
integer exit_flag
real(dp) avgkp
real(dp) g
real(dp) rho0

References mod_scalars::avgke, mod_scalars::avgkp, mod_scalars::avgpe, mod_param::domain, mod_scalars::exit_flag, mod_scalars::first_time, mod_scalars::g, mod_scalars::iic, mod_param::itlm, mod_parallel::master, mod_parallel::mymaster, mod_scalars::ninfo, mod_param::ntilee, mod_param::ntilex, mod_scalars::rho0, mod_iounits::stdout, mod_parallel::tile_count, mod_scalars::time_code, and mod_scalars::volume.

Referenced by tl_diag().

Here is the caller graph for this function: