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_debugging, only : hchksum, uvchksum
10 use mom_domains, only : sum_across_pes, min_across_pes, max_across_pes
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, haloshift, symmetric)
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 [m s-1].
49  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
50  intent(in) :: v !< The meridional velocity [m s-1].
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 m2 s-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 m2 s-1 ~> m3 s-1 or kg s-1].
59  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
60  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully symmetric
61  !! computationoal domain.
62 
63  integer :: is, ie, js, je, nz, hs
64  logical :: sym
65  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
66 
67  ! Note that for the chksum calls to be useful for reproducing across PE
68  ! counts, there must be no redundant points, so all variables use is..ie
69  ! and js...je as their extent.
70  hs=1; if (present(haloshift)) hs=haloshift
71  sym=.false.; if (present(symmetric)) sym=symmetric
72  call uvchksum(mesg//" [uv]", u, v, g%HI, haloshift=hs, symmetric=sym)
73  call hchksum(h, mesg//" h", g%HI, haloshift=hs, scale=gv%H_to_m)
74  call uvchksum(mesg//" [uv]h", uh, vh, g%HI, haloshift=hs, &
75  symmetric=sym, scale=gv%H_to_m)
76 end subroutine mom_state_chksum_5arg
77 
78 ! =============================================================================
79 
80 !> Write out chksums for the model's basic state variables.
81 subroutine mom_state_chksum_3arg(mesg, u, v, h, G, GV, haloshift, symmetric)
82  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
83  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
84  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
85  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
86  intent(in) :: u !< Zonal velocity [m s-1].
87  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
88  intent(in) :: v !< Meridional velocity [m s-1].
89  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
90  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
91  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
92  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully symmetric
93  !! computationoal domain.
94 
95  integer :: is, ie, js, je, nz, hs
96  logical :: sym
97  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
98 
99  ! Note that for the chksum calls to be useful for reproducing across PE
100  ! counts, there must be no redundant points, so all variables use is..ie
101  ! and js...je as their extent.
102  hs=1; if (present(haloshift)) hs=haloshift
103  sym=.false.; if (present(symmetric)) sym=symmetric
104  call uvchksum(mesg//" u", u, v, g%HI,haloshift=hs, symmetric=sym)
105  call hchksum(h, mesg//" h",g%HI, haloshift=hs, scale=gv%H_to_m)
106 end subroutine mom_state_chksum_3arg
107 
108 ! =============================================================================
109 
110 !> Write out chksums for the model's thermodynamic state variables.
111 subroutine mom_thermo_chksum(mesg, tv, G, haloshift)
112  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
113  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
114  !! thermodynamic variables.
115  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
116  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
117 
118  integer :: is, ie, js, je, nz, hs
119  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
120  hs=1; if (present(haloshift)) hs=haloshift
121 
122  if (associated(tv%T)) call hchksum(tv%T, mesg//" T",g%HI,haloshift=hs)
123  if (associated(tv%S)) call hchksum(tv%S, mesg//" S",g%HI,haloshift=hs)
124  if (associated(tv%frazil)) call hchksum(tv%frazil, mesg//" frazil",g%HI,haloshift=hs)
125  if (associated(tv%salt_deficit)) call hchksum(tv%salt_deficit, mesg//" salt deficit",g%HI,haloshift=hs)
126 
127 end subroutine mom_thermo_chksum
128 
129 ! =============================================================================
130 
131 !> Write out chksums for the ocean surface variables.
132 subroutine mom_surface_chksum(mesg, sfc, G, haloshift, symmetric)
133  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
134  type(surface), intent(inout) :: sfc !< transparent ocean surface state
135  !! structure shared with the calling routine
136  !! data in this structure is intent out.
137  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
138  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0).
139  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully symmetric
140  !! computationoal domain.
141 
142  integer :: hs
143  logical :: sym
144 
145  sym = .false. ; if (present(symmetric)) sym = symmetric
146  hs = 1 ; if (present(haloshift)) hs = haloshift
147 
148  if (allocated(sfc%SST)) call hchksum(sfc%SST, mesg//" SST",g%HI,haloshift=hs)
149  if (allocated(sfc%SSS)) call hchksum(sfc%SSS, mesg//" SSS",g%HI,haloshift=hs)
150  if (allocated(sfc%sea_lev)) call hchksum(sfc%sea_lev, mesg//" sea_lev",g%HI,haloshift=hs)
151  if (allocated(sfc%Hml)) call hchksum(sfc%Hml, mesg//" Hml",g%HI,haloshift=hs)
152  if (allocated(sfc%u) .and. allocated(sfc%v)) &
153  call uvchksum(mesg//" SSU", sfc%u, sfc%v, g%HI, haloshift=hs, symmetric=sym)
154 ! if (allocated(sfc%salt_deficit)) call hchksum(sfc%salt_deficit, mesg//" salt deficit",G%HI,haloshift=hs)
155  if (associated(sfc%frazil)) call hchksum(sfc%frazil, mesg//" frazil",g%HI,haloshift=hs)
156 
157 end subroutine mom_surface_chksum
158 
159 ! =============================================================================
160 
161 !> Write out chksums for the model's accelerations
162 subroutine mom_accel_chksum(mesg, CAu, CAv, PFu, PFv, diffu, diffv, G, GV, US, pbce, &
163  u_accel_bt, v_accel_bt, symmetric)
164  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
165  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
166  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
167  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
168  intent(in) :: cau !< Zonal acceleration due to Coriolis
169  !! and momentum advection terms [m s-2].
170  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
171  intent(in) :: cav !< Meridional acceleration due to Coriolis
172  !! and momentum advection terms [m s-2].
173  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
174  intent(in) :: pfu !< Zonal acceleration due to pressure gradients
175  !! (equal to -dM/dx) [m s-2].
176  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
177  intent(in) :: pfv !< Meridional acceleration due to pressure gradients
178  !! (equal to -dM/dy) [m s-2].
179  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
180  intent(in) :: diffu !< Zonal acceleration due to convergence of the
181  !! along-isopycnal stress tensor [m s-1 T-1 ~> m s-2].
182  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
183  intent(in) :: diffv !< Meridional acceleration due to convergence of
184  !! the along-isopycnal stress tensor [m s-1 T-1 ~> m s-2].
185  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
186  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
187  optional, intent(in) :: pbce !< The baroclinic pressure anomaly in each layer
188  !! due to free surface height anomalies
189  !! [m2 s-2 H-1 ~> m s-2 or m4 s-2 kg-1].
190  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
191  optional, intent(in) :: u_accel_bt !< The zonal acceleration from terms in the
192  !! barotropic solver [m s-2].
193  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
194  optional, intent(in) :: v_accel_bt !< The meridional acceleration from terms in
195  !! the barotropic solver [m s-2].
196  logical, optional, intent(in) :: symmetric !< If true, do checksums on the fully symmetric
197  !! computationoal domain.
198 
199  integer :: is, ie, js, je, nz
200  logical :: sym
201 
202  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
203  sym=.false.; if (present(symmetric)) sym=symmetric
204 
205  ! Note that for the chksum calls to be useful for reproducing across PE
206  ! counts, there must be no redundant points, so all variables use is..ie
207  ! and js...je as their extent.
208  call uvchksum(mesg//" CA[uv]", cau, cav, g%HI, haloshift=0, symmetric=sym)
209  call uvchksum(mesg//" PF[uv]", pfu, pfv, g%HI, haloshift=0, symmetric=sym)
210  call uvchksum(mesg//" diffu", diffu, diffv, g%HI,haloshift=0, symmetric=sym, scale=us%s_to_T)
211  if (present(pbce)) &
212  call hchksum(pbce, mesg//" pbce",g%HI,haloshift=0, scale=gv%m_to_H*us%L_T_to_m_s**2)
213  if (present(u_accel_bt) .and. present(v_accel_bt)) &
214  call uvchksum(mesg//" [uv]_accel_bt", u_accel_bt, v_accel_bt, g%HI,haloshift=0, symmetric=sym)
215 end subroutine mom_accel_chksum
216 
217 ! =============================================================================
218 
219 !> Monitor and write out statistics for the model's state variables.
220 subroutine mom_state_stats(mesg, u, v, h, Temp, Salt, G, allowChange, permitDiminishing)
221  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
222  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
223  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
224  intent(in) :: u !< The zonal velocity [m s-1].
225  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
226  intent(in) :: v !< The meridional velocity [m s-1].
227  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
228  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
229  real, pointer, dimension(:,:,:), &
230  intent(in) :: temp !< Temperature [degC].
231  real, pointer, dimension(:,:,:), &
232  intent(in) :: salt !< Salinity [ppt].
233  logical, optional, intent(in) :: allowchange !< do not flag an error
234  !! if the statistics change.
235  logical, optional, intent(in) :: permitdiminishing !< do not flag error
236  !!if the extrema are diminishing.
237  ! Local variables
238  integer :: is, ie, js, je, nz, i, j, k
239  real :: vol, dv, area, h_minimum
240  type(stats) :: t, s, delt, dels
241  type(stats), save :: oldt, olds ! NOTE: save data is not normally allowed but
242  logical, save :: firstcall = .true. ! we use it for debugging purposes here on the
243  logical :: do_ts
244  real, save :: oldvol ! assumption we will not turn this on with threads
245  character(len=80) :: lmsg
246  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
247 
248  do_ts = associated(temp) .and. associated(salt)
249 
250  ! First collect local stats
251  area = 0. ; vol = 0.
252  do j = js, je ; do i = is, ie
253  area = area + g%areaT(i,j)
254  enddo ; enddo
255  t%minimum = 1.e34 ; t%maximum = -1.e34 ; t%average = 0.
256  s%minimum = 1.e34 ; s%maximum = -1.e34 ; s%average = 0.
257  h_minimum = 1.e34
258  do k = 1, nz ; do j = js, je ; do i = is, ie
259  if (g%mask2dT(i,j)>0.) then
260  dv = g%areaT(i,j)*h(i,j,k) ; vol = vol + dv
261  if (do_ts .and. h(i,j,k)>0.) then
262  t%minimum = min( t%minimum, temp(i,j,k) ) ; t%maximum = max( t%maximum, temp(i,j,k) )
263  t%average = t%average + dv*temp(i,j,k)
264  s%minimum = min( s%minimum, salt(i,j,k) ) ; s%maximum = max( s%maximum, salt(i,j,k) )
265  s%average = s%average + dv*salt(i,j,k)
266  endif
267  if (h_minimum > h(i,j,k)) h_minimum = h(i,j,k)
268  endif
269  enddo ; enddo ; enddo
270  call sum_across_pes( area ) ; call sum_across_pes( vol )
271  if (do_ts) then
272  call min_across_pes( t%minimum ) ; call max_across_pes( t%maximum ) ; call sum_across_pes( t%average )
273  call min_across_pes( s%minimum ) ; call max_across_pes( s%maximum ) ; call sum_across_pes( s%average )
274  t%average = t%average / vol ; s%average = s%average / vol
275  endif
276  if (is_root_pe()) then
277  if (.not.firstcall) then
278  dv = vol - oldvol
279  delt%minimum = t%minimum - oldt%minimum ; delt%maximum = t%maximum - oldt%maximum
280  delt%average = t%average - oldt%average
281  dels%minimum = s%minimum - olds%minimum ; dels%maximum = s%maximum - olds%maximum
282  dels%average = s%average - olds%average
283  write(lmsg(1:80),'(2(a,es12.4))') 'Mean thickness =',vol/area,' frac. delta=',dv/vol
284  call mom_mesg(lmsg//trim(mesg))
285  if (do_ts) then
286  write(lmsg(1:80),'(a,3es12.4)') 'Temp min/mean/max =',t%minimum,t%average,t%maximum
287  call mom_mesg(lmsg//trim(mesg))
288  write(lmsg(1:80),'(a,3es12.4)') 'delT min/mean/max =',delt%minimum,delt%average,delt%maximum
289  call mom_mesg(lmsg//trim(mesg))
290  write(lmsg(1:80),'(a,3es12.4)') 'Salt min/mean/max =',s%minimum,s%average,s%maximum
291  call mom_mesg(lmsg//trim(mesg))
292  write(lmsg(1:80),'(a,3es12.4)') 'delS min/mean/max =',dels%minimum,dels%average,dels%maximum
293  call mom_mesg(lmsg//trim(mesg))
294  endif
295  else
296  write(lmsg(1:80),'(a,es12.4)') 'Mean thickness =',vol/area
297  call mom_mesg(lmsg//trim(mesg))
298  if (do_ts) then
299  write(lmsg(1:80),'(a,3es12.4)') 'Temp min/mean/max =',t%minimum,t%average,t%maximum
300  call mom_mesg(lmsg//trim(mesg))
301  write(lmsg(1:80),'(a,3es12.4)') 'Salt min/mean/max =',s%minimum,s%average,s%maximum
302  call mom_mesg(lmsg//trim(mesg))
303  endif
304  endif
305  endif
306  firstcall = .false. ; oldvol = vol
307  oldt%minimum = t%minimum ; oldt%maximum = t%maximum ; oldt%average = t%average
308  olds%minimum = s%minimum ; olds%maximum = s%maximum ; olds%average = s%average
309 
310  if (do_ts .and. t%minimum<-5.0) then
311  do j = js, je ; do i = is, ie
312  if (minval(temp(i,j,:)) == t%minimum) then
313  write(0,'(a,2f12.5)') 'x,y=',g%geoLonT(i,j),g%geoLatT(i,j)
314  write(0,'(a3,3a12)') 'k','h','Temp','Salt'
315  do k = 1, nz
316  write(0,'(i3,3es12.4)') k,h(i,j,k),temp(i,j,k),salt(i,j,k)
317  enddo
318  stop 'Extremum detected'
319  endif
320  enddo ; enddo
321  endif
322 
323  if (h_minimum<0.0) then
324  do j = js, je ; do i = is, ie
325  if (minval(h(i,j,:)) == h_minimum) then
326  write(0,'(a,2f12.5)') 'x,y=',g%geoLonT(i,j),g%geoLatT(i,j)
327  write(0,'(a3,3a12)') 'k','h','Temp','Salt'
328  do k = 1, nz
329  write(0,'(i3,3es12.4)') k,h(i,j,k),temp(i,j,k),salt(i,j,k)
330  enddo
331  stop 'Negative thickness detected'
332  endif
333  enddo ; enddo
334  endif
335 
336 end subroutine mom_state_stats
337 
338 end module mom_checksum_packages
mom_checksum_packages::stats
A type for storing statistica about a variable.
Definition: MOM_checksum_packages.F90:31
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
mom_checksum_packages::mom_state_chksum
Write out checksums of the MOM6 state variables.
Definition: MOM_checksum_packages.F90:23
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:82
mom_checksum_packages
Provides routines that do checksums of groups of MOM variables.
Definition: MOM_checksum_packages.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25