17 implicit none ;
private 20 public mom_state_stats, mom_surface_chksum
24 module procedure mom_state_chksum_5arg
25 module procedure mom_state_chksum_3arg
28 #include <MOM_memory.h> 32 real :: minimum = 1.e34
33 real :: maximum = -1.e34
42 subroutine mom_state_chksum_5arg(mesg, u, v, h, uh, vh, G, GV, US, haloshift, symmetric, vel_scale)
47 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
49 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
51 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
53 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
56 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
60 integer,
optional,
intent(in) :: haloshift
61 logical,
optional,
intent(in) :: symmetric
63 real,
optional,
intent(in) :: vel_scale
67 integer :: is, ie, js, je, nz, hs
68 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
73 hs = 1 ;
if (
present(haloshift)) hs=haloshift
74 sym = .false. ;
if (
present(symmetric)) sym=symmetric
75 scale_vel = us%L_T_to_m_s ;
if (
present(vel_scale)) scale_vel = vel_scale
77 call uvchksum(mesg//
" [uv]", u, v, g%HI, haloshift=hs, symmetric=sym, scale=scale_vel)
78 call hchksum(h, mesg//
" h", g%HI, haloshift=hs, scale=gv%H_to_m)
79 call uvchksum(mesg//
" [uv]h", uh, vh, g%HI, haloshift=hs, &
80 symmetric=sym, scale=gv%H_to_m*us%L_to_m**2*us%s_to_T)
81 end subroutine mom_state_chksum_5arg
86 subroutine mom_state_chksum_3arg(mesg, u, v, h, G, GV, US, haloshift, symmetric)
87 character(len=*),
intent(in) :: mesg
90 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
92 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
94 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
98 integer,
optional,
intent(in) :: haloshift
99 logical,
optional,
intent(in) :: symmetric
102 integer :: is, ie, js, je, nz, hs
105 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
106 l_t_to_m_s = 1.0 ;
if (
present(us)) l_t_to_m_s = us%L_T_to_m_s
111 hs=1;
if (
present(haloshift)) hs=haloshift
112 sym=.false.;
if (
present(symmetric)) sym=symmetric
113 call uvchksum(mesg//
" u", u, v, g%HI, haloshift=hs, symmetric=sym, scale=l_t_to_m_s)
114 call hchksum(h, mesg//
" h",g%HI, haloshift=hs, scale=gv%H_to_m)
115 end subroutine mom_state_chksum_3arg
120 subroutine mom_thermo_chksum(mesg, tv, G, US, haloshift)
121 character(len=*),
intent(in) :: mesg
126 integer,
optional,
intent(in) :: haloshift
128 integer :: is, ie, js, je, nz, hs
129 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
130 hs=1;
if (
present(haloshift)) hs=haloshift
132 if (
associated(tv%T))
call hchksum(tv%T, mesg//
" T", g%HI, haloshift=hs)
133 if (
associated(tv%S))
call hchksum(tv%S, mesg//
" S", g%HI, haloshift=hs)
134 if (
associated(tv%frazil))
call hchksum(tv%frazil, mesg//
" frazil", g%HI, haloshift=hs, &
135 scale=us%Q_to_J_kg*us%R_to_kg_m3*us%Z_to_m)
136 if (
associated(tv%salt_deficit)) &
137 call hchksum(tv%salt_deficit, mesg//
" salt deficit", g%HI, haloshift=hs, scale=us%RZ_to_kg_m2)
139 end subroutine mom_thermo_chksum
144 subroutine mom_surface_chksum(mesg, sfc_state, G, US, haloshift, symmetric)
145 character(len=*),
intent(in) :: mesg
146 type(
surface),
intent(inout) :: sfc_state
151 integer,
optional,
intent(in) :: haloshift
152 logical,
optional,
intent(in) :: symmetric
158 sym = .false. ;
if (
present(symmetric)) sym = symmetric
159 hs = 1 ;
if (
present(haloshift)) hs = haloshift
161 if (
allocated(sfc_state%SST))
call hchksum(sfc_state%SST, mesg//
" SST", g%HI, haloshift=hs)
162 if (
allocated(sfc_state%SSS))
call hchksum(sfc_state%SSS, mesg//
" SSS", g%HI, haloshift=hs)
163 if (
allocated(sfc_state%sea_lev))
call hchksum(sfc_state%sea_lev, mesg//
" sea_lev", g%HI, &
164 haloshift=hs, scale=us%Z_to_m)
165 if (
allocated(sfc_state%Hml))
call hchksum(sfc_state%Hml, mesg//
" Hml", g%HI, haloshift=hs, &
167 if (
allocated(sfc_state%u) .and.
allocated(sfc_state%v)) &
168 call uvchksum(mesg//
" SSU", sfc_state%u, sfc_state%v, g%HI, haloshift=hs, symmetric=sym, &
172 if (
allocated(sfc_state%frazil))
call hchksum(sfc_state%frazil, mesg//
" frazil", g%HI, &
173 haloshift=hs, scale=us%Q_to_J_kg*us%RZ_to_kg_m2)
175 end subroutine mom_surface_chksum
180 subroutine mom_accel_chksum(mesg, CAu, CAv, PFu, PFv, diffu, diffv, G, GV, US, pbce, &
181 u_accel_bt, v_accel_bt, symmetric)
182 character(len=*),
intent(in) :: mesg
185 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
188 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
191 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
194 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
197 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
200 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
204 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
205 optional,
intent(in) :: pbce
208 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
209 optional,
intent(in) :: u_accel_bt
211 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
212 optional,
intent(in) :: v_accel_bt
214 logical,
optional,
intent(in) :: symmetric
217 integer :: is, ie, js, je, nz
220 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
221 sym=.false.;
if (
present(symmetric)) sym=symmetric
226 call uvchksum(mesg//
" CA[uv]", cau, cav, g%HI, haloshift=0, symmetric=sym, scale=us%L_T2_to_m_s2)
227 call uvchksum(mesg//
" PF[uv]", pfu, pfv, g%HI, haloshift=0, symmetric=sym, scale=us%L_T2_to_m_s2)
228 call uvchksum(mesg//
" diffu", diffu, diffv, g%HI,haloshift=0, symmetric=sym, scale=us%L_T2_to_m_s2)
230 call hchksum(pbce, mesg//
" pbce",g%HI,haloshift=0, scale=gv%m_to_H*us%L_T_to_m_s**2)
231 if (
present(u_accel_bt) .and.
present(v_accel_bt)) &
232 call uvchksum(mesg//
" [uv]_accel_bt", u_accel_bt, v_accel_bt, g%HI,haloshift=0, symmetric=sym, &
233 scale=us%L_T2_to_m_s2)
234 end subroutine mom_accel_chksum
239 subroutine mom_state_stats(mesg, u, v, h, Temp, Salt, G, GV, US, allowChange, permitDiminishing)
242 character(len=*),
intent(in) :: mesg
243 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
245 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
247 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
249 real,
pointer,
dimension(:,:,:), &
251 real,
pointer,
dimension(:,:,:), &
254 logical,
optional,
intent(in) :: allowChange
256 logical,
optional,
intent(in) :: permitDiminishing
260 real,
dimension(G%isc:G%iec, G%jsc:G%jec) :: &
269 type(
stats) :: T, S, delT, delS
273 type(
stats),
save :: oldT, oldS
274 logical,
save :: firstCall = .true.
277 character(len=80) :: lMsg
278 integer :: is, ie, js, je, nz, i, j, k
280 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
281 do_ts =
associated(temp) .and.
associated(salt)
289 do j=js,je ;
do i=is,ie
290 tmp_a(i,j) = tmp_a(i,j) + us%L_to_m**2*g%areaT(i,j)
292 t%minimum = 1.e34 ; t%maximum = -1.e34 ; t%average = 0.
293 s%minimum = 1.e34 ; s%maximum = -1.e34 ; s%average = 0.
294 h_minimum = 1.e34*gv%m_to_H
295 do k=1,nz ;
do j=js,je ;
do i=is,ie
296 if (g%mask2dT(i,j)>0.)
then 297 dv = us%L_to_m**2*g%areaT(i,j)*gv%H_to_m*h(i,j,k)
298 tmp_v(i,j) = tmp_v(i,j) + dv
299 if (do_ts .and. h(i,j,k)>0.)
then 300 t%minimum = min( t%minimum, temp(i,j,k) ) ; t%maximum = max( t%maximum, temp(i,j,k) )
301 t%average = t%average + dv*temp(i,j,k)
302 s%minimum = min( s%minimum, salt(i,j,k) ) ; s%maximum = max( s%maximum, salt(i,j,k) )
303 s%average = s%average + dv*salt(i,j,k)
305 if (h_minimum > h(i,j,k)) h_minimum = h(i,j,k)
307 enddo ;
enddo ;
enddo 310 call min_across_pes( t%minimum ) ;
call max_across_pes( t%maximum )
311 call min_across_pes( s%minimum ) ;
call max_across_pes( s%maximum )
313 t%average = t%average / vol ; s%average = s%average / vol
315 if (is_root_pe())
then 316 if (.not.firstcall)
then 318 delt%minimum = t%minimum - oldt%minimum ; delt%maximum = t%maximum - oldt%maximum
319 delt%average = t%average - oldt%average
320 dels%minimum = s%minimum - olds%minimum ; dels%maximum = s%maximum - olds%maximum
321 dels%average = s%average - olds%average
322 write(lmsg(1:80),
'(2(a,es12.4))')
'Mean thickness =', vol/area,
' frac. delta=',dv/vol
323 call mom_mesg(lmsg//trim(mesg))
325 write(lmsg(1:80),
'(a,3es12.4)')
'Temp min/mean/max =',t%minimum,t%average,t%maximum
326 call mom_mesg(lmsg//trim(mesg))
327 write(lmsg(1:80),
'(a,3es12.4)')
'delT min/mean/max =',delt%minimum,delt%average,delt%maximum
328 call mom_mesg(lmsg//trim(mesg))
329 write(lmsg(1:80),
'(a,3es12.4)')
'Salt min/mean/max =',s%minimum,s%average,s%maximum
330 call mom_mesg(lmsg//trim(mesg))
331 write(lmsg(1:80),
'(a,3es12.4)')
'delS min/mean/max =',dels%minimum,dels%average,dels%maximum
332 call mom_mesg(lmsg//trim(mesg))
335 write(lmsg(1:80),
'(a,es12.4)')
'Mean thickness =', vol/area
336 call mom_mesg(lmsg//trim(mesg))
338 write(lmsg(1:80),
'(a,3es12.4)')
'Temp min/mean/max =', t%minimum, t%average, t%maximum
339 call mom_mesg(lmsg//trim(mesg))
340 write(lmsg(1:80),
'(a,3es12.4)')
'Salt min/mean/max =', s%minimum, s%average, s%maximum
341 call mom_mesg(lmsg//trim(mesg))
345 firstcall = .false. ; oldvol = vol
346 oldt%minimum = t%minimum ; oldt%maximum = t%maximum ; oldt%average = t%average
347 olds%minimum = s%minimum ; olds%maximum = s%maximum ; olds%average = s%average
349 if (do_ts .and. t%minimum<-5.0)
then 350 do j=js,je ;
do i=is,ie
351 if (minval(temp(i,j,:)) == t%minimum)
then 352 write(0,
'(a,2f12.5)')
'x,y=', g%geoLonT(i,j), g%geoLatT(i,j)
353 write(0,
'(a3,3a12)')
'k',
'h',
'Temp',
'Salt' 355 write(0,
'(i3,3es12.4)') k, h(i,j,k), temp(i,j,k), salt(i,j,k)
357 stop
'Extremum detected' 362 if (h_minimum<0.0)
then 363 do j=js,je ;
do i=is,ie
364 if (minval(h(i,j,:)) == h_minimum)
then 365 write(0,
'(a,2f12.5)')
'x,y=',g%geoLonT(i,j),g%geoLatT(i,j)
366 write(0,
'(a3,3a12)')
'k',
'h',
'Temp',
'Salt' 368 write(0,
'(i3,3es12.4)') k, h(i,j,k), temp(i,j,k), salt(i,j,k)
370 stop
'Negative thickness detected' 375 end subroutine mom_state_stats
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Ocean grid type. See mom_grid for details.
Provides the ocean grid type.
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Describes various unit conversion factors.
Provides routines that do checksums of groups of MOM variables.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
Describes the vertical ocean grid, including unit conversion factors.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Find an accurate and order-invariant sum of a distributed 2d or 3d field.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
A type for storing statistica about a variable.
Write out checksums of the MOM6 state variables.
Provides checksumming functions for debugging.