MOM6
MOM_checksum_packages.F90
1 !> Provides routines that do checksums of groups of MOM variables
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 ! This module provides several routines that do check-sums of groups
7 ! of variables in the various dynamic solver routines.
8 
9 use mom_coms, only : min_across_pes, max_across_pes, reproducing_sum
10 use mom_debugging, only : hchksum, uvchksum
11 use mom_error_handler, only : mom_mesg, is_root_pe
12 use mom_grid, only : ocean_grid_type
16 
17 implicit none ; private
18 
19 public mom_state_chksum, mom_thermo_chksum, mom_accel_chksum
20 public mom_state_stats, mom_surface_chksum
21 
22 !> Write out checksums of the MOM6 state variables
24  module procedure mom_state_chksum_5arg
25  module procedure mom_state_chksum_3arg
26 end interface
27 
28 #include <MOM_memory.h>
29 
30 !> A type for storing statistica about a variable
31 type :: stats ; private
32  real :: minimum = 1.e34 !< The minimum value
33  real :: maximum = -1.e34 !< The maximum value
34  real :: average = 0. !< The average value
35 end type stats
36 
37 contains
38 
39 ! =============================================================================
40 
41 !> Write out chksums for the model's basic state variables, including transports.
42 subroutine mom_state_chksum_5arg(mesg, u, v, h, uh, vh, G, GV, US, haloshift, symmetric, vel_scale)
43  character(len=*), &
44  intent(in) :: mesg !< A message that appears on the chksum lines.
45  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
46  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
47  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
48  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1] or other units.
49  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
50  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1] or other units.
51  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
52  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
53  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
54  intent(in) :: uh !< Volume flux through zonal faces = u*h*dy
55  !! [H L2 T-1 ~> m3 s-1 or kg s-1].
56  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
57  intent(in) :: vh !< Volume flux through meridional faces = v*h*dx
58  !! [H L2 T-1 ~> m3 s-1 or kg s-1].
59  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
60  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
61  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully symmetric
62  !! computational domain.
63  real, optional, intent(in) :: vel_scale !< The scaling factor to convert velocities to [m s-1]
64 
65  real :: scale_vel ! The scaling factor to convert velocities to [m s-1]
66  logical :: sym
67  integer :: is, ie, js, je, nz, hs
68  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
69 
70  ! Note that for the chksum calls to be useful for reproducing across PE
71  ! counts, there must be no redundant points, so all variables use is..ie
72  ! and js...je as their extent.
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
76 
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
82 
83 ! =============================================================================
84 
85 !> Write out chksums for the model's basic state variables.
86 subroutine mom_state_chksum_3arg(mesg, u, v, h, G, GV, US, haloshift, symmetric)
87  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
88  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
89  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
90  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
91  intent(in) :: u !< Zonal velocity [L T-1 ~> m s-1] or [m s-1].
92  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
93  intent(in) :: v !< Meridional velocity [L T-1 ~> m s-1] or [m s-1]..
94  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
95  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
96  type(unit_scale_type), optional, intent(in) :: US !< A dimensional unit scaling type, which is
97  !! used to rescale u and v if present.
98  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
99  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully
100  !! symmetric computational domain.
101  real :: L_T_to_m_s ! A rescaling factor for velocities [m T s-1 L-1 ~> nondim] or [nondim]
102  integer :: is, ie, js, je, nz, hs
103  logical :: sym
104 
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
107 
108  ! Note that for the chksum calls to be useful for reproducing across PE
109  ! counts, there must be no redundant points, so all variables use is..ie
110  ! and js...je as their extent.
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
116 
117 ! =============================================================================
118 
119 !> Write out chksums for the model's thermodynamic state variables.
120 subroutine mom_thermo_chksum(mesg, tv, G, US, haloshift)
121  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
122  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
123  !! thermodynamic variables.
124  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
125  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
126  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
127 
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
131 
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)
138 
139 end subroutine mom_thermo_chksum
140 
141 ! =============================================================================
142 
143 !> Write out chksums for the ocean surface variables.
144 subroutine mom_surface_chksum(mesg, sfc_state, G, US, haloshift, symmetric)
145  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
146  type(surface), intent(inout) :: sfc_state !< transparent ocean surface state structure
147  !! shared with the calling routine data in this
148  !! structure is intent out.
149  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
150  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
151  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
152  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully symmetric
153  !! computational domain.
154 
155  integer :: hs
156  logical :: sym
157 
158  sym = .false. ; if (present(symmetric)) sym = symmetric
159  hs = 1 ; if (present(haloshift)) hs = haloshift
160 
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, &
166  scale=us%Z_to_m)
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, &
169  scale=us%L_T_to_m_s)
170 ! if (allocated(sfc_state%salt_deficit)) &
171 ! call hchksum(sfc_state%salt_deficit, mesg//" salt deficit", G%HI, haloshift=hs, scale=US%RZ_to_kg_m2)
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)
174 
175 end subroutine mom_surface_chksum
176 
177 ! =============================================================================
178 
179 !> Write out chksums for the model's accelerations
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 !< A message that appears on the chksum lines.
183  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
184  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
185  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
186  intent(in) :: cau !< Zonal acceleration due to Coriolis
187  !! and momentum advection terms [L T-2 ~> m s-2].
188  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
189  intent(in) :: cav !< Meridional acceleration due to Coriolis
190  !! and momentum advection terms [L T-2 ~> m s-2].
191  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
192  intent(in) :: pfu !< Zonal acceleration due to pressure gradients
193  !! (equal to -dM/dx) [L T-2 ~> m s-2].
194  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
195  intent(in) :: pfv !< Meridional acceleration due to pressure gradients
196  !! (equal to -dM/dy) [L T-2 ~> m s-2].
197  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
198  intent(in) :: diffu !< Zonal acceleration due to convergence of the
199  !! along-isopycnal stress tensor [L T-2 ~> m s-2].
200  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
201  intent(in) :: diffv !< Meridional acceleration due to convergence of
202  !! the along-isopycnal stress tensor [L T-2 ~> m s-2].
203  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
204  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
205  optional, intent(in) :: pbce !< The baroclinic pressure anomaly in each layer
206  !! due to free surface height anomalies
207  !! [L2 T-2 H-1 ~> m s-2 or m4 s-2 kg-1].
208  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
209  optional, intent(in) :: u_accel_bt !< The zonal acceleration from terms in the
210  !! barotropic solver [L T-2 ~> m s-2].
211  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
212  optional, intent(in) :: v_accel_bt !< The meridional acceleration from terms in
213  !! the barotropic solver [L T-2 ~> m s-2].
214  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully symmetric
215  !! computational domain.
216 
217  integer :: is, ie, js, je, nz
218  logical :: sym
219 
220  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
221  sym=.false.; if (present(symmetric)) sym=symmetric
222 
223  ! Note that for the chksum calls to be useful for reproducing across PE
224  ! counts, there must be no redundant points, so all variables use is..ie
225  ! and js...je as their extent.
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)
229  if (present(pbce)) &
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
235 
236 ! =============================================================================
237 
238 !> Monitor and write out statistics for the model's state variables.
239 subroutine mom_state_stats(mesg, u, v, h, Temp, Salt, G, GV, US, allowChange, permitDiminishing)
240  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
241  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
242  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
243  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
244  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
245  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
246  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
247  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
248  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
249  real, pointer, dimension(:,:,:), &
250  intent(in) :: temp !< Temperature [degC].
251  real, pointer, dimension(:,:,:), &
252  intent(in) :: salt !< Salinity [ppt].
253  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
254  logical, optional, intent(in) :: allowchange !< do not flag an error
255  !! if the statistics change.
256  logical, optional, intent(in) :: permitdiminishing !< do not flag error if the
257  !! extrema are diminishing.
258 
259  ! Local variables
260  real, dimension(G%isc:G%iec, G%jsc:G%jec) :: &
261  tmp_a, & ! The area per cell [m2] (unscaled to permit reproducing sum).
262  tmp_v, & ! The column-integrated volume [m3] (unscaled to permit reproducing sum)
263  tmp_t, & ! The column-integrated temperature [degC m3]
264  tmp_s ! The column-integrated salinity [ppt m3]
265  real :: vol, dv ! The total ocean volume and its change [m3] (unscaled to permit reproducing sum).
266  real :: area ! The total ocean surface area [m2] (unscaled to permit reproducing sum).
267  real :: h_minimum ! The minimum layer thicknesses [H ~> m or kg m-2]
268  logical :: do_ts ! If true, evaluate statistics for temperature and salinity
269  type(stats) :: t, s, delt, dels
270 
271  ! NOTE: save data is not normally allowed but we use it for debugging purposes here on the
272  ! assumption we will not turn this on with threads
273  type(stats), save :: oldt, olds
274  logical, save :: firstcall = .true.
275  real, save :: oldvol ! The previous total ocean volume [m3]
276 
277  character(len=80) :: lmsg
278  integer :: is, ie, js, je, nz, i, j, k
279 
280  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
281  do_ts = associated(temp) .and. associated(salt)
282 
283  tmp_a(:,:) = 0.0
284  tmp_v(:,:) = 0.0
285  tmp_t(:,:) = 0.0
286  tmp_s(:,:) = 0.0
287 
288  ! First collect local stats
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)
291  enddo ; enddo
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)
304  endif
305  if (h_minimum > h(i,j,k)) h_minimum = h(i,j,k)
306  endif
307  enddo ; enddo ; enddo
308  area = reproducing_sum( tmp_a ) ; vol = reproducing_sum( tmp_v )
309  if (do_ts) then
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 )
312  t%average = reproducing_sum( tmp_t ) ; s%average = reproducing_sum( tmp_s )
313  t%average = t%average / vol ; s%average = s%average / vol
314  endif
315  if (is_root_pe()) then
316  if (.not.firstcall) then
317  dv = vol - oldvol
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))
324  if (do_ts) then
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))
333  endif
334  else
335  write(lmsg(1:80),'(a,es12.4)') 'Mean thickness =', vol/area
336  call mom_mesg(lmsg//trim(mesg))
337  if (do_ts) then
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))
342  endif
343  endif
344  endif
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
348 
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'
354  do k = 1, nz
355  write(0,'(i3,3es12.4)') k, h(i,j,k), temp(i,j,k), salt(i,j,k)
356  enddo
357  stop 'Extremum detected'
358  endif
359  enddo ; enddo
360  endif
361 
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'
367  do k = 1, nz
368  write(0,'(i3,3es12.4)') k, h(i,j,k), temp(i,j,k), salt(i,j,k)
369  enddo
370  stop 'Negative thickness detected'
371  endif
372  enddo ; enddo
373  endif
374 
375 end subroutine mom_state_stats
376 
377 end module mom_checksum_packages
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.
Definition: MOM_grid.F90:26
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
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.
Definition: MOM_coms.F90:53
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.