MOM6
mom_checksum_packages::mom_state_chksum Interface Reference

Detailed Description

Write out checksums of the MOM6 state variables.

Definition at line 23 of file MOM_checksum_packages.F90.

Private functions

subroutine mom_state_chksum_5arg (mesg, u, v, h, uh, vh, G, GV, US, haloshift, symmetric, vel_scale)
 Write out chksums for the model's basic state variables, including transports. More...
 
subroutine mom_state_chksum_3arg (mesg, u, v, h, G, GV, US, haloshift, symmetric)
 Write out chksums for the model's basic state variables. More...
 

Functions and subroutines

◆ mom_state_chksum_3arg()

subroutine mom_checksum_packages::mom_state_chksum::mom_state_chksum_3arg ( 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,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in), optional  US,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  symmetric 
)
private

Write out chksums for the model's basic state variables.

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 [L T-1 ~> m s-1] or [m s-1].
[in]vMeridional velocity [L T-1 ~> m s-1] or [m s-1]..
[in]hLayer thicknesses [H ~> m or kg m-2].
[in]usA dimensional unit scaling type, which is used to rescale u and v if present.
[in]haloshiftThe width of halos to check (default 0).
[in]symmetricIf true, do checksums on the fully symmetric computational domain.

Definition at line 87 of file MOM_checksum_packages.F90.

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)

◆ mom_state_chksum_5arg()

subroutine mom_checksum_packages::mom_state_chksum::mom_state_chksum_5arg ( 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(szib_(g),szj_(g),szk_(g)), intent(in)  uh,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in)  vh,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  symmetric,
real, intent(in), optional  vel_scale 
)
private

Write out chksums for the model's basic state variables, including transports.

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 [L T-1 ~> m s-1] or other units.
[in]vThe meridional velocity [L T-1 ~> m s-1] or other units.
[in]hLayer thicknesses [H ~> m or kg m-2].
[in]uhVolume flux through zonal faces = u*h*dy
[in]vhVolume flux through meridional faces = v*h*dx
[in]usA dimensional unit scaling type
[in]haloshiftThe width of halos to check (default 0).
[in]symmetricIf true, do checksums on the fully symmetric computational domain.
[in]vel_scaleThe scaling factor to convert velocities to [m s-1]

Definition at line 43 of file MOM_checksum_packages.F90.

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)

The documentation for this interface was generated from the following file: