MOM6
mom_checksum_packages Module Reference

Data Types

interface  mom_state_chksum
 
type  stats
 

Functions/Subroutines

subroutine mom_state_chksum_5arg (mesg, u, v, h, uh, vh, G, GV, haloshift, symmetric)
 
subroutine mom_state_chksum_3arg (mesg, u, v, h, G, GV, haloshift, symmetric)
 
subroutine, public mom_thermo_chksum (mesg, tv, G, haloshift)
 
subroutine, public mom_accel_chksum (mesg, CAu, CAv, PFu, PFv, diffu, diffv, G, GV, pbce, u_accel_bt, v_accel_bt, symmetric)
 
subroutine, public mom_state_stats (mesg, u, v, h, Temp, Salt, G, allowChange, permitDiminishing)
 

Function/Subroutine Documentation

◆ mom_accel_chksum()

subroutine, public mom_checksum_packages::mom_accel_chksum ( character(len=*), intent(in)  mesg,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in)  CAu,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in)  CAv,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in)  PFu,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in)  PFv,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in)  diffu,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in)  diffv,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in), optional  pbce,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in), optional  u_accel_bt,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in), optional  v_accel_bt,
logical, intent(in), optional  symmetric 
)
Parameters
[in]mesgA message that appears on the chksum lines.
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]cauZonal acceleration due to Coriolis
[in]cavMeridional acceleration due to Coriolis
[in]pfuZonal acceleration due to pressure gradients
[in]pfvMeridional acceleration due to pressure gradients
[in]diffuZonal acceleration due to convergence of the
[in]diffvMeridional acceleration due to convergence of
[in]pbceThe baroclinic pressure anomaly in each layer
[in]u_accel_btThe zonal acceleration from terms in the
[in]v_accel_btThe meridional acceleration from terms in

Definition at line 161 of file MOM_checksum_packages.F90.

Referenced by mom_dynamics_legacy_split::step_mom_dyn_legacy_split(), mom_dynamics_split_rk2::step_mom_dyn_split_rk2(), mom_dynamics_unsplit::step_mom_dyn_unsplit(), and mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2().

161  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
162  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
163  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
164  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
165  intent(in) :: cau !< Zonal acceleration due to Coriolis
166  !! and momentum advection terms, in m s-2.
167  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
168  intent(in) :: cav !< Meridional acceleration due to Coriolis
169  !! and momentum advection terms, in m s-2.
170  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
171  intent(in) :: pfu !< Zonal acceleration due to pressure gradients
172  !! (equal to -dM/dx) in m s-2.
173  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
174  intent(in) :: pfv !< Meridional acceleration due to pressure gradients
175  !! (equal to -dM/dy) in m s-2.
176  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
177  intent(in) :: diffu !< Zonal acceleration due to convergence of the
178  !! along-isopycnal stress tensor, in m s-2.
179  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
180  intent(in) :: diffv !< Meridional acceleration due to convergence of
181  !! the along-isopycnal stress tensor, in m s-2.
182  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
183  optional, intent(in) :: pbce !< The baroclinic pressure anomaly in each layer
184  !! due to free surface height anomalies, in
185  !! m2 s-2 H-1. !! NULL.
186  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
187  optional, intent(in) :: u_accel_bt !< The zonal acceleration from terms in the
188  !! barotropic solver,in m s-2.
189  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
190  optional, intent(in) :: v_accel_bt !< The meridional acceleration from terms in
191  !! the barotropic solver,in m s-2.
192  logical, optional, intent(in) :: symmetric
193 
194 ! This subroutine writes out chksums for the model's accelerations.
195 ! Arguments: mesg - A message that appears on the chksum lines.
196 ! (in) CAu - Zonal acceleration due to Coriolis and momentum
197 ! advection terms, in m s-2.
198 ! (in) CAv - Meridional acceleration due to Coriolis and
199 ! momentum advection terms, in m s-2.
200 ! (in) PFu - Zonal acceleration due to pressure gradients
201 ! (equal to -dM/dx) in m s-2.
202 ! (in) PFv - Meridional acceleration due to pressure
203 ! gradients (equal to -dM/dy) in m s-2.
204 ! (in) diffu - Zonal acceleration due to convergence of the
205 ! along-isopycnal stress tensor, in m s-2.
206 ! (in) diffv - Meridional acceleration due to convergence of
207 ! the along-isopycnal stress tensor, in m s-2.
208 ! (in) G - The ocean's grid structure.
209 ! (in) GV - The ocean's vertical grid structure.
210 ! (in) pbce - the baroclinic pressure anomaly in each layer
211 ! due to free surface height anomalies, in m2 s-2 H-1.
212 ! pbce points to a space with nz layers or NULL.
213 ! (in) u_accel_bt - The zonal acceleration from terms in the barotropic
214 ! solver, in m s-2.
215 ! (in) v_accel_bt - The meridional acceleration from terms in the
216 ! barotropic solver, in m s-2.
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)
227  call uvchksum(mesg//" PF[uv]", pfu, pfv, g%HI, haloshift=0, symmetric=sym)
228  call uvchksum(mesg//" diffu", diffu, diffv, g%HI,haloshift=0, symmetric=sym)
229  if (present(pbce)) &
230  call hchksum(pbce, mesg//" pbce",g%HI,haloshift=0, scale=gv%m_to_H)
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)
Here is the caller graph for this function:

◆ mom_state_chksum_3arg()

subroutine mom_checksum_packages::mom_state_chksum_3arg ( character(len=*), intent(in)  mesg,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in)  u,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in)  v,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  symmetric 
)
private
Parameters
[in]mesgA message that appears on the chksum lines.
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]uZonal velocity, in m s-1.
[in]vMeridional velocity, in m s-1.
[in]hLayer thicknesses, in H (usually m or kg m-2).

Definition at line 99 of file MOM_checksum_packages.F90.

99  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
100  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
101  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
102  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
103  intent(in) :: u !< Zonal velocity, in m s-1.
104  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
105  intent(in) :: v !< Meridional velocity, in m s-1.
106  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
107  intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2).
108  integer, optional, intent(in) :: haloshift
109  logical, optional, intent(in) :: symmetric
110 ! This subroutine writes out chksums for the model's basic state variables.
111 ! Arguments: mesg - A message that appears on the chksum lines.
112 ! (in) u - Zonal velocity, in m s-1.
113 ! (in) v - Meridional velocity, in m s-1.
114 ! (in) h - Layer thickness, in m.
115 ! (in) uh - Volume flux through zonal faces = u*h*dy, m3 s-1.
116 ! (in) vh - Volume flux through meridional faces = v*h*dx, in m3 s-1.
117 ! (in) G - The ocean's grid structure.
118 ! (in) GV - The ocean's vertical grid structure.
119  integer :: is, ie, js, je, nz, hs
120  logical :: sym
121  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
122 
123  ! Note that for the chksum calls to be useful for reproducing across PE
124  ! counts, there must be no redundant points, so all variables use is..ie
125  ! and js...je as their extent.
126  hs=1; if (present(haloshift)) hs=haloshift
127  sym=.false.; if (present(symmetric)) sym=symmetric
128  call uvchksum(mesg//" u", u, v, g%HI,haloshift=hs, symmetric=sym)
129  call hchksum(h, mesg//" h",g%HI, haloshift=hs, scale=gv%H_to_m)

◆ mom_state_chksum_5arg()

subroutine mom_checksum_packages::mom_state_chksum_5arg ( character(len=*), intent(in)  mesg,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in)  u,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in)  v,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in)  uh,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in)  vh,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  symmetric 
)
private
Parameters
[in]mesgA message that appears on the chksum lines.
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]uThe zonal velocity, in m s-1.
[in]vThe meridional velocity, in m s-1.
[in]hLayer thicknesses, in H (usually m or kg m-2).
[in]uhVolume flux through zonal faces = u*h*dy, m3 s-1.
[in]vhVolume flux through meridional

Definition at line 55 of file MOM_checksum_packages.F90.

55  character(len=*), &
56  intent(in) :: mesg !< A message that appears on the chksum lines.
57  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
58  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
59  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
60  intent(in) :: u !< The zonal velocity, in m s-1.
61  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
62  intent(in) :: v !< The meridional velocity, in m s-1.
63  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
64  intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2).
65  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
66  intent(in) :: uh !< Volume flux through zonal faces = u*h*dy, m3 s-1.
67  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
68  intent(in) :: vh !< Volume flux through meridional
69  !! faces = v*h*dx, in m3 s-1.
70  integer, optional, intent(in) :: haloshift
71  logical, optional, intent(in) :: symmetric
72 ! This subroutine writes out chksums for the model's basic state variables.
73 ! Arguments: mesg - A message that appears on the chksum lines.
74 ! (in) u - Zonal velocity, in m s-1.
75 ! (in) v - Meridional velocity, in m s-1.
76 ! (in) h - Layer thickness, in m.
77 ! (in) uh - Volume flux through zonal faces = u*h*dy, m3 s-1.
78 ! (in) vh - Volume flux through meridional faces = v*h*dx, in m3 s-1.
79 ! (in) G - The ocean's grid structure.
80 ! (in) GV - The ocean's vertical grid structure.
81  integer :: is, ie, js, je, nz, hs
82  logical :: sym
83  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
84 
85  ! Note that for the chksum calls to be useful for reproducing across PE
86  ! counts, there must be no redundant points, so all variables use is..ie
87  ! and js...je as their extent.
88  hs=1; if (present(haloshift)) hs=haloshift
89  sym=.false.; if (present(symmetric)) sym=symmetric
90  call uvchksum(mesg//" [uv]", u, v, g%HI, haloshift=hs, symmetric=sym)
91  call hchksum(h, mesg//" h", g%HI, haloshift=hs, scale=gv%H_to_m)
92  call uvchksum(mesg//" [uv]h", uh, vh, g%HI, haloshift=hs, &
93  symmetric=sym, scale=gv%H_to_m)

◆ mom_state_stats()

subroutine, public mom_checksum_packages::mom_state_stats ( character(len=*), intent(in)  mesg,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in)  u,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in)  v,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
real, dimension(:,:,:), intent(in), pointer  Temp,
real, dimension(:,:,:), intent(in), pointer  Salt,
type(ocean_grid_type), intent(in)  G,
logical, intent(in), optional  allowChange,
logical, intent(in), optional  permitDiminishing 
)
Parameters
[in]gThe ocean's grid structure.
[in]mesgA message that appears on the chksum lines.
[in]uThe zonal velocity, in m s-1.
[in]vThe meridional velocity, in m s-1.
[in]hLayer thicknesses, in H (usually m or kg m-2).
[in]tempTemperature in degree C.
[in]saltSalinity, in ppt.
[in]allowchangedo not flag an error if the statistics change.
[in]permitdiminishingdo not flag error

Definition at line 238 of file MOM_checksum_packages.F90.

References mom_error_handler::is_root_pe(), and mom_error_handler::mom_mesg().

Referenced by mom_diabatic_driver::diabatic().

238  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
239  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
240  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
241  intent(in) :: u !< The zonal velocity, in m s-1.
242  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
243  intent(in) :: v !< The meridional velocity, in m s-1.
244  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
245  intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2).
246  real, pointer, dimension(:,:,:), &
247  intent(in) :: temp !< Temperature in degree C.
248  real, pointer, dimension(:,:,:), &
249  intent(in) :: salt !< Salinity, in ppt.
250 
251  logical, optional, intent(in) :: allowchange !< do not flag an error
252  !! if the statistics change.
253  logical, optional, &
254  intent(in) :: permitdiminishing !< do not flag error
255  !!if the extrema are diminishing.
256 ! This subroutine monitors statistics for the model's state variables.
257 ! Arguments: mesg - A message that appears on the chksum lines.
258 ! (in) u - Zonal velocity, in m s-1.
259 ! (in) v - Meridional velocity, in m s-1.
260 ! (in) h - Layer thickness, in m.
261 ! (in) T - Temperature, in degree C.
262 ! (in) S - Salinity, in ppt.
263 ! (in) G - The ocean's grid structure.
264 ! (in) allowChange - do not flag an error if the statistics change
265 ! (in) permitDiminishing - do not flag an error if the extrema are diminishing
266  integer :: is, ie, js, je, nz, i, j, k
267  real :: vol, dv, area, h_minimum
268  type(stats) :: t, s, delt, dels
269  type(stats), save :: oldt, olds ! NOTE: save data is not normally allowed but
270  logical, save :: firstcall = .true. ! we use it for debugging purposes here on the
271  logical :: do_ts
272  real, save :: oldvol ! assumption we will not turn this on with threads
273  character(len=80) :: lmsg
274  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
275 
276  do_ts = associated(temp) .and. associated(salt)
277 
278  ! First collect local stats
279  area = 0. ; vol = 0.
280  do j = js, je ; do i = is, ie
281  area = area + g%areaT(i,j)
282  enddo ; enddo
283  t%minimum = 1.e34 ; t%maximum = -1.e34 ; t%average = 0.
284  s%minimum = 1.e34 ; s%maximum = -1.e34 ; s%average = 0.
285  h_minimum = 1.e34
286  do k = 1, nz ; do j = js, je ; do i = is, ie
287  if (g%mask2dT(i,j)>0.) then
288  dv = g%areaT(i,j)*h(i,j,k) ; vol = vol + dv
289  if (do_ts .and. h(i,j,k)>0.) then
290  t%minimum = min( t%minimum, temp(i,j,k) ) ; t%maximum = max( t%maximum, temp(i,j,k) )
291  t%average = t%average + dv*temp(i,j,k)
292  s%minimum = min( s%minimum, salt(i,j,k) ) ; s%maximum = max( s%maximum, salt(i,j,k) )
293  s%average = s%average + dv*salt(i,j,k)
294  endif
295  if (h_minimum > h(i,j,k)) h_minimum = h(i,j,k)
296  endif
297  enddo ; enddo ; enddo
298  call sum_across_pes( area ) ; call sum_across_pes( vol )
299  if (do_ts) then
300  call min_across_pes( t%minimum ) ; call max_across_pes( t%maximum ) ; call sum_across_pes( t%average )
301  call min_across_pes( s%minimum ) ; call max_across_pes( s%maximum ) ; call sum_across_pes( s%average )
302  t%average = t%average / vol ; s%average = s%average / vol
303  endif
304  if (is_root_pe()) then
305  if (.not.firstcall) then
306  dv = vol - oldvol
307  delt%minimum = t%minimum - oldt%minimum ; delt%maximum = t%maximum - oldt%maximum
308  delt%average = t%average - oldt%average
309  dels%minimum = s%minimum - olds%minimum ; dels%maximum = s%maximum - olds%maximum
310  dels%average = s%average - olds%average
311  write(lmsg(1:80),'(2(a,es12.4))') 'Mean thickness =',vol/area,' frac. delta=',dv/vol
312  call mom_mesg(lmsg//trim(mesg))
313  if (do_ts) then
314  write(lmsg(1:80),'(a,3es12.4)') 'Temp min/mean/max =',t%minimum,t%average,t%maximum
315  call mom_mesg(lmsg//trim(mesg))
316  write(lmsg(1:80),'(a,3es12.4)') 'delT min/mean/max =',delt%minimum,delt%average,delt%maximum
317  call mom_mesg(lmsg//trim(mesg))
318  write(lmsg(1:80),'(a,3es12.4)') 'Salt min/mean/max =',s%minimum,s%average,s%maximum
319  call mom_mesg(lmsg//trim(mesg))
320  write(lmsg(1:80),'(a,3es12.4)') 'delS min/mean/max =',dels%minimum,dels%average,dels%maximum
321  call mom_mesg(lmsg//trim(mesg))
322  endif
323  else
324  write(lmsg(1:80),'(a,es12.4)') 'Mean thickness =',vol/area
325  call mom_mesg(lmsg//trim(mesg))
326  if (do_ts) then
327  write(lmsg(1:80),'(a,3es12.4)') 'Temp min/mean/max =',t%minimum,t%average,t%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  endif
332  endif
333  endif
334  firstcall = .false. ; oldvol = vol
335  oldt%minimum = t%minimum ; oldt%maximum = t%maximum ; oldt%average = t%average
336  olds%minimum = s%minimum ; olds%maximum = s%maximum ; olds%average = s%average
337 
338  if (do_ts .and. t%minimum<-5.0) then
339  do j = js, je ; do i = is, ie
340  if (minval(temp(i,j,:)) == t%minimum) then
341  write(0,'(a,2f12.5)') 'x,y=',g%geoLonT(i,j),g%geoLatT(i,j)
342  write(0,'(a3,3a12)') 'k','h','Temp','Salt'
343  do k = 1, nz
344  write(0,'(i3,3es12.4)') k,h(i,j,k),temp(i,j,k),salt(i,j,k)
345  enddo
346  stop 'Extremum detected'
347  endif
348  enddo ; enddo
349  endif
350 
351  if (h_minimum<0.0) then
352  do j = js, je ; do i = is, ie
353  if (minval(h(i,j,:)) == h_minimum) then
354  write(0,'(a,2f12.5)') 'x,y=',g%geoLonT(i,j),g%geoLatT(i,j)
355  write(0,'(a3,3a12)') 'k','h','Temp','Salt'
356  do k = 1, nz
357  write(0,'(i3,3es12.4)') k,h(i,j,k),temp(i,j,k),salt(i,j,k)
358  enddo
359  stop 'Negative thickness detected'
360  endif
361  enddo ; enddo
362  endif
363 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mom_thermo_chksum()

subroutine, public mom_checksum_packages::mom_thermo_chksum ( character(len=*), intent(in)  mesg,
type(thermo_var_ptrs), intent(in)  tv,
type(ocean_grid_type), intent(in)  G,
integer, intent(in), optional  haloshift 
)
Parameters
[in]mesgA message that appears on the chksum lines.
[in]tvA structure pointing to various thermodynamic variables.
[in]gThe ocean's grid structure.

Definition at line 135 of file MOM_checksum_packages.F90.

Referenced by mom::step_mom_thermo().

135  character(len=*), intent(in) :: mesg !< A message that appears on the chksum lines.
136  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
137  !! thermodynamic variables.
138  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
139  integer, optional, intent(in) :: haloshift
140 ! This subroutine writes out chksums for the model's thermodynamic state
141 ! variables.
142 ! Arguments: mesg - A message that appears on the chksum lines.
143 ! (in) tv - A structure containing pointers to any thermodynamic
144 ! fields that are in use.
145 ! (in) G - The ocean's grid structure.
146  integer :: is, ie, js, je, nz, hs
147  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
148  hs=1; if (present(haloshift)) hs=haloshift
149 
150  if (associated(tv%T)) call hchksum(tv%T, mesg//" T",g%HI,haloshift=hs)
151  if (associated(tv%S)) call hchksum(tv%S, mesg//" S",g%HI,haloshift=hs)
152  if (associated(tv%frazil)) call hchksum(tv%frazil, mesg//" frazil",g%HI,haloshift=hs)
153  if (associated(tv%salt_deficit)) call hchksum(tv%salt_deficit, mesg//" salt deficit",g%HI,haloshift=hs)
154 
Here is the caller graph for this function: