MOM6
mom_debugging::check_redundant_c Interface Reference

Detailed Description

Check for consistency between the duplicated points of a C-grid vector.

Definition at line 37 of file MOM_debugging.F90.

Private functions

subroutine check_redundant_vc3d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 Check for consistency between the duplicated points of a 3-D C-grid vector. More...
 
subroutine check_redundant_vc2d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 Check for consistency between the duplicated points of a 2-D C-grid vector. More...
 

Functions and subroutines

◆ check_redundant_vc2d()

subroutine mom_debugging::check_redundant_c::check_redundant_vc2d ( character(len=*), intent(in)  mesg,
real, dimension(g%isdb:,g%jsd:), intent(in)  u_comp,
real, dimension(g%isd:,g%jsdb:), intent(in)  v_comp,
type(ocean_grid_type), intent(inout)  G,
integer, intent(in), optional  is,
integer, intent(in), optional  ie,
integer, intent(in), optional  js,
integer, intent(in), optional  je,
integer, intent(in), optional  direction 
)
private

Check for consistency between the duplicated points of a 2-D C-grid vector.

Parameters
[in]mesgAn identifying message
[in,out]gThe ocean's grid structure
[in]u_compThe u-component of the vector to be checked for consistency
[in]v_compThe u-component of the vector to be checked for consistency
[in]isThe starting i-index to check
[in]ieThe ending i-index to check
[in]jsThe starting j-index to check
[in]jeThe ending j-index to check
[in]directionthe direction flag to be passed to pass_vector

Definition at line 136 of file MOM_debugging.F90.

136  character(len=*), intent(in) :: mesg !< An identifying message
137  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
138  real, dimension(G%IsdB:,G%jsd:), intent(in) :: u_comp !< The u-component of the vector
139  !! to be checked for consistency
140  real, dimension(G%isd:,G%JsdB:), intent(in) :: v_comp !< The u-component of the vector
141  !! to be checked for consistency
142  integer, optional, intent(in) :: is !< The starting i-index to check
143  integer, optional, intent(in) :: ie !< The ending i-index to check
144  integer, optional, intent(in) :: js !< The starting j-index to check
145  integer, optional, intent(in) :: je !< The ending j-index to check
146  integer, optional, intent(in) :: direction !< the direction flag to be
147  !! passed to pass_vector
148  ! Local variables
149  real :: u_nonsym(G%isd:G%ied,G%jsd:G%jed)
150  real :: v_nonsym(G%isd:G%ied,G%jsd:G%jed)
151  real :: u_resym(G%IsdB:G%IedB,G%jsd:G%jed)
152  real :: v_resym(G%isd:G%ied,G%JsdB:G%JedB)
153  character(len=128) :: mesg2
154  integer :: i, j, is_ch, ie_ch, js_ch, je_ch
155  integer :: Isq, Ieq, Jsq, Jeq, isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
156 
157  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
158  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
159  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
160 
161  if (.not.(present(is) .or. present(ie) .or. present(js) .or. present(je))) then
162  ! This only works with symmetric memory, so otherwise return.
163  if ((isd == isdb) .and. (jsd == jsdb)) return
164  endif
165 
166  do i=isd,ied ; do j=jsd,jed
167  u_nonsym(i,j) = u_comp(i,j) ; v_nonsym(i,j) = v_comp(i,j)
168  enddo ; enddo
169 
170  if (.not.associated(g%Domain_aux)) call mom_error(fatal," check_redundant"//&
171  " called with a non-associated auxiliary domain the grid type.")
172  call pass_vector(u_nonsym, v_nonsym, g%Domain_aux, direction)
173 
174  do i=isdb,iedb ; do j=jsd,jed ; u_resym(i,j) = u_comp(i,j) ; enddo ; enddo
175  do i=isd,ied ; do j=jsdb,jedb ; v_resym(i,j) = v_comp(i,j) ; enddo ; enddo
176  do i=isd,ied ; do j=jsd,jed
177  u_resym(i,j) = u_nonsym(i,j) ; v_resym(i,j) = v_nonsym(i,j)
178  enddo ; enddo
179  call pass_vector(u_resym, v_resym, g%Domain, direction)
180 
181  is_ch = isq ; ie_ch = ieq ; js_ch = jsq ; je_ch = jeq
182  if (present(is)) is_ch = is ; if (present(ie)) ie_ch = ie
183  if (present(js)) js_ch = js ; if (present(js)) je_ch = je
184 
185  do i=is_ch,ie_ch ; do j=js_ch+1,je_ch
186  if (u_resym(i,j) /= u_comp(i,j) .and. &
187  redundant_prints(3) < max_redundant_prints) then
188  write(mesg2,'(" redundant u-components",2(1pe12.4)," differ by ", &
189  & 1pe12.4," at i,j = ",2i4," on pe ",i4)') &
190  u_comp(i,j), u_resym(i,j),u_comp(i,j)-u_resym(i,j),i,j,pe_here()
191  write(0,'(A130)') trim(mesg)//trim(mesg2)
192  redundant_prints(3) = redundant_prints(3) + 1
193  endif
194  enddo ; enddo
195  do i=is_ch+1,ie_ch ; do j=js_ch,je_ch
196  if (v_resym(i,j) /= v_comp(i,j) .and. &
197  redundant_prints(3) < max_redundant_prints) then
198  write(mesg2,'(" redundant v-comps",2(1pe12.4)," differ by ", &
199  & 1pe12.4," at i,j = ",2i4," x,y = ",2(1pe12.4)" on pe ",i4)') &
200  v_comp(i,j), v_resym(i,j),v_comp(i,j)-v_resym(i,j),i,j, &
201  g%geoLonBu(i,j), g%geoLatBu(i,j), pe_here()
202  write(0,'(A155)') trim(mesg)//trim(mesg2)
203  redundant_prints(3) = redundant_prints(3) + 1
204  endif
205  enddo ; enddo
206 

◆ check_redundant_vc3d()

subroutine mom_debugging::check_redundant_c::check_redundant_vc3d ( character(len=*), intent(in)  mesg,
real, dimension(g%isdb:,g%jsd:,:), intent(in)  u_comp,
real, dimension(g%isd:,g%jsdb:,:), intent(in)  v_comp,
type(ocean_grid_type), intent(inout)  G,
integer, intent(in), optional  is,
integer, intent(in), optional  ie,
integer, intent(in), optional  js,
integer, intent(in), optional  je,
integer, intent(in), optional  direction 
)
private

Check for consistency between the duplicated points of a 3-D C-grid vector.

Parameters
[in]mesgAn identifying message
[in,out]gThe ocean's grid structure
[in]u_compThe u-component of the vector to be checked for consistency
[in]v_compThe u-component of the vector to be checked for consistency
[in]isThe starting i-index to check
[in]ieThe ending i-index to check
[in]jsThe starting j-index to check
[in]jeThe ending j-index to check
[in]directionthe direction flag to be passed to pass_vector

Definition at line 106 of file MOM_debugging.F90.

106  character(len=*), intent(in) :: mesg !< An identifying message
107  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
108  real, dimension(G%IsdB:,G%jsd:,:), intent(in) :: u_comp !< The u-component of the vector
109  !! to be checked for consistency
110  real, dimension(G%isd:,G%JsdB:,:), intent(in) :: v_comp !< The u-component of the vector
111  !! to be checked for consistency
112  integer, optional, intent(in) :: is !< The starting i-index to check
113  integer, optional, intent(in) :: ie !< The ending i-index to check
114  integer, optional, intent(in) :: js !< The starting j-index to check
115  integer, optional, intent(in) :: je !< The ending j-index to check
116  integer, optional, intent(in) :: direction !< the direction flag to be
117  !! passed to pass_vector
118  ! Local variables
119  character(len=24) :: mesg_k
120  integer :: k
121 
122  do k=1,size(u_comp,3)
123  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
124  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
125  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
126  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
127 
128  call check_redundant_vc2d(trim(mesg)//trim(mesg_k), u_comp(:,:,k), &
129  v_comp(:,:,k), g, is, ie, js, je, direction)
130  enddo

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