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, haloshift, symmetric)
 Write out chksums for the model's basic state variables, including transports. More...
 
subroutine mom_state_chksum_3arg (mesg, u, v, h, G, GV, 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,
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 [m s-1].
[in]vMeridional velocity [m s-1].
[in]hLayer thicknesses [H ~> m or kg m-2].
[in]haloshiftThe width of halos to check (default 0).
[in]symmetricIf true, do checksums on the fully symmetric computationoal domain.

Definition at line 82 of file MOM_checksum_packages.F90.

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)

◆ 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,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  symmetric 
)
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 [m s-1].
[in]vThe meridional velocity [m s-1].
[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]haloshiftThe width of halos to check (default 0).
[in]symmetricIf true, do checksums on the fully symmetric computationoal domain.

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 [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)

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