MOM6
mom_debugging::check_redundant_b Interface Reference

Detailed Description

Check for consistency between the duplicated points of a B-grid vector or scalar.

Definition at line 41 of file MOM_debugging.F90.

Private functions

subroutine check_redundant_vb3d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 Check for consistency between the duplicated points of a 3-D B-grid vector. More...
 
subroutine check_redundant_vb2d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 Check for consistency between the duplicated points of a 2-D B-grid vector. More...
 
subroutine check_redundant_sb3d (mesg, array, G, is, ie, js, je)
 Check for consistency between the duplicated points of a 3-D scalar at corner points. More...
 
subroutine check_redundant_sb2d (mesg, array, G, is, ie, js, je)
 Check for consistency between the duplicated points of a 2-D scalar at corner points. More...
 

Functions and subroutines

◆ check_redundant_sb2d()

subroutine mom_debugging::check_redundant_b::check_redundant_sb2d ( character(len=*), intent(in)  mesg,
real, dimension(g%isdb:,g%jsdb:), intent(in)  array,
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 
)
private

Check for consistency between the duplicated points of a 2-D scalar at corner points.

Parameters
[in]mesgAn identifying message
[in,out]gThe ocean's grid structure
[in]arrayThe array 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

Definition at line 236 of file MOM_debugging.F90.

236  character(len=*), intent(in) :: mesg !< An identifying message
237  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
238  real, dimension(G%IsdB:,G%JsdB:), intent(in) :: array !< The array to be checked for consistency
239  integer, optional, intent(in) :: is !< The starting i-index to check
240  integer, optional, intent(in) :: ie !< The ending i-index to check
241  integer, optional, intent(in) :: js !< The starting j-index to check
242  integer, optional, intent(in) :: je !< The ending j-index to check
243  ! Local variables
244  real :: a_nonsym(G%isd:G%ied,G%jsd:G%jed)
245  real :: a_resym(G%IsdB:G%IedB,G%JsdB:G%JedB)
246  character(len=128) :: mesg2
247  integer :: i, j, is_ch, ie_ch, js_ch, je_ch
248  integer :: Isq, Ieq, Jsq, Jeq, isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
249 
250  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
251  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
252  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
253 
254  if (.not.(present(is) .or. present(ie) .or. present(js) .or. present(je))) then
255  ! This only works with symmetric memory, so otherwise return.
256  if ((isd == isdb) .and. (jsd == jsdb)) return
257  endif
258 
259  do i=isd,ied ; do j=jsd,jed
260  a_nonsym(i,j) = array(i,j)
261  enddo ; enddo
262 
263  if (.not.associated(g%Domain_aux)) call mom_error(fatal," check_redundant"//&
264  " called with a non-associated auxiliary domain the grid type.")
265  call pass_vector(a_nonsym, a_nonsym, g%Domain_aux, &
266  direction=to_all+scalar_pair, stagger=bgrid_ne)
267 
268  do i=isdb,iedb ; do j=jsdb,jedb ; a_resym(i,j) = array(i,j) ; enddo ; enddo
269  do i=isd,ied ; do j=jsd,jed
270  a_resym(i,j) = a_nonsym(i,j)
271  enddo ; enddo
272  call pass_vector(a_resym, a_resym, g%Domain, direction=to_all+scalar_pair, &
273  stagger=bgrid_ne)
274 
275  is_ch = isq ; ie_ch = ieq ; js_ch = jsq ; je_ch = jeq
276  if (present(is)) is_ch = is ; if (present(ie)) ie_ch = ie
277  if (present(js)) js_ch = js ; if (present(js)) je_ch = je
278 
279  do i=is_ch,ie_ch ; do j=js_ch,je_ch
280  if (a_resym(i,j) /= array(i,j) .and. &
281  redundant_prints(2) < max_redundant_prints) then
282  write(mesg2,'(" Redundant points",2(1pe12.4)," differ by ", &
283  & 1pe12.4," at i,j = ",2i4," on pe ",i4)') &
284  array(i,j), a_resym(i,j),array(i,j)-a_resym(i,j),i,j,pe_here()
285  write(0,'(A130)') trim(mesg)//trim(mesg2)
286  redundant_prints(2) = redundant_prints(2) + 1
287  endif
288  enddo ; enddo
289 

◆ check_redundant_sb3d()

subroutine mom_debugging::check_redundant_b::check_redundant_sb3d ( character(len=*), intent(in)  mesg,
real, dimension(g%isdb:,g%jsdb:,:), intent(in)  array,
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 
)
private

Check for consistency between the duplicated points of a 3-D scalar at corner points.

Parameters
[in]mesgAn identifying message
[in,out]gThe ocean's grid structure
[in]arrayThe array 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

Definition at line 211 of file MOM_debugging.F90.

211  character(len=*), intent(in) :: mesg !< An identifying message
212  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
213  real, dimension(G%IsdB:,G%JsdB:,:), intent(in) :: array !< The array to be checked for consistency
214  integer, optional, intent(in) :: is !< The starting i-index to check
215  integer, optional, intent(in) :: ie !< The ending i-index to check
216  integer, optional, intent(in) :: js !< The starting j-index to check
217  integer, optional, intent(in) :: je !< The ending j-index to check
218 
219  ! Local variables
220  character(len=24) :: mesg_k
221  integer :: k
222 
223  do k=1,size(array,3)
224  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
225  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
226  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
227  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
228 
229  call check_redundant_sb2d(trim(mesg)//trim(mesg_k), array(:,:,k), &
230  g, is, ie, js, je)
231  enddo

◆ check_redundant_vb2d()

subroutine mom_debugging::check_redundant_b::check_redundant_vb2d ( character(len=*), intent(in)  mesg,
real, dimension(g%isdb:,g%jsdb:), intent(in)  u_comp,
real, dimension(g%isdb:,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 B-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 v-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 325 of file MOM_debugging.F90.

325  character(len=*), intent(in) :: mesg !< An identifying message
326  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
327  real, dimension(G%IsdB:,G%JsdB:), intent(in) :: u_comp !< The u-component of the vector
328  !! to be checked for consistency
329  real, dimension(G%IsdB:,G%JsdB:), intent(in) :: v_comp !< The v-component of the vector
330  !! to be checked for consistency
331  integer, optional, intent(in) :: is !< The starting i-index to check
332  integer, optional, intent(in) :: ie !< The ending i-index to check
333  integer, optional, intent(in) :: js !< The starting j-index to check
334  integer, optional, intent(in) :: je !< The ending j-index to check
335  integer, optional, intent(in) :: direction !< the direction flag to be
336  !! passed to pass_vector
337  ! Local variables
338  real :: u_nonsym(G%isd:G%ied,G%jsd:G%jed)
339  real :: v_nonsym(G%isd:G%ied,G%jsd:G%jed)
340  real :: u_resym(G%IsdB:G%IedB,G%JsdB:G%JedB)
341  real :: v_resym(G%IsdB:G%IedB,G%JsdB:G%JedB)
342  character(len=128) :: mesg2
343  integer :: i, j, is_ch, ie_ch, js_ch, je_ch
344  integer :: Isq, Ieq, Jsq, Jeq, isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
345 
346  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
347  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
348  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
349 
350  if (.not.(present(is) .or. present(ie) .or. present(js) .or. present(je))) then
351  ! This only works with symmetric memory, so otherwise return.
352  if ((isd == isdb) .and. (jsd == jsdb)) return
353  endif
354 
355  do i=isd,ied ; do j=jsd,jed
356  u_nonsym(i,j) = u_comp(i,j) ; v_nonsym(i,j) = v_comp(i,j)
357  enddo ; enddo
358 
359  if (.not.associated(g%Domain_aux)) call mom_error(fatal," check_redundant"//&
360  " called with a non-associated auxiliary domain the grid type.")
361  call pass_vector(u_nonsym, v_nonsym, g%Domain_aux, direction, stagger=bgrid_ne)
362 
363  do i=isdb,iedb ; do j=jsdb,jedb
364  u_resym(i,j) = u_comp(i,j) ; v_resym(i,j) = v_comp(i,j)
365  enddo ; enddo
366  do i=isd,ied ; do j=jsd,jed
367  u_resym(i,j) = u_nonsym(i,j) ; v_resym(i,j) = v_nonsym(i,j)
368  enddo ; enddo
369  call pass_vector(u_resym, v_resym, g%Domain, direction, stagger=bgrid_ne)
370 
371  is_ch = isq ; ie_ch = ieq ; js_ch = jsq ; je_ch = jeq
372  if (present(is)) is_ch = is ; if (present(ie)) ie_ch = ie
373  if (present(js)) js_ch = js ; if (present(js)) je_ch = je
374 
375  do i=is_ch,ie_ch ; do j=js_ch,je_ch
376  if (u_resym(i,j) /= u_comp(i,j) .and. &
377  redundant_prints(2) < max_redundant_prints) then
378  write(mesg2,'(" redundant u-components",2(1pe12.4)," differ by ", &
379  & 1pe12.4," at i,j = ",2i4," on pe ",i4)') &
380  u_comp(i,j), u_resym(i,j),u_comp(i,j)-u_resym(i,j),i,j,pe_here()
381  write(0,'(A130)') trim(mesg)//trim(mesg2)
382  redundant_prints(2) = redundant_prints(2) + 1
383  endif
384  enddo ; enddo
385  do i=is_ch,ie_ch ; do j=js_ch,je_ch
386  if (v_resym(i,j) /= v_comp(i,j) .and. &
387  redundant_prints(2) < max_redundant_prints) then
388  write(mesg2,'(" redundant v-comps",2(1pe12.4)," differ by ", &
389  & 1pe12.4," at i,j = ",2i4," x,y = ",2(1pe12.4)" on pe ",i4)') &
390  v_comp(i,j), v_resym(i,j),v_comp(i,j)-v_resym(i,j),i,j, &
391  g%geoLonBu(i,j), g%geoLatBu(i,j), pe_here()
392  write(0,'(A155)') trim(mesg)//trim(mesg2)
393  redundant_prints(2) = redundant_prints(2) + 1
394  endif
395  enddo ; enddo
396 

◆ check_redundant_vb3d()

subroutine mom_debugging::check_redundant_b::check_redundant_vb3d ( character(len=*), intent(in)  mesg,
real, dimension(g%isdb:,g%jsdb:,:), intent(in)  u_comp,
real, dimension(g%isdb:,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 B-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 v-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 295 of file MOM_debugging.F90.

295  character(len=*), intent(in) :: mesg !< An identifying message
296  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
297  real, dimension(G%IsdB:,G%JsdB:,:), intent(in) :: u_comp !< The u-component of the vector
298  !! to be checked for consistency
299  real, dimension(G%IsdB:,G%JsdB:,:), intent(in) :: v_comp !< The v-component of the vector
300  !! to be checked for consistency
301  integer, optional, intent(in) :: is !< The starting i-index to check
302  integer, optional, intent(in) :: ie !< The ending i-index to check
303  integer, optional, intent(in) :: js !< The starting j-index to check
304  integer, optional, intent(in) :: je !< The ending j-index to check
305  integer, optional, intent(in) :: direction !< the direction flag to be
306  !! passed to pass_vector
307  ! Local variables
308  character(len=24) :: mesg_k
309  integer :: k
310 
311  do k=1,size(u_comp,3)
312  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
313  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
314  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
315  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
316 
317  call check_redundant_vb2d(trim(mesg)//trim(mesg_k), u_comp(:,:,k), &
318  v_comp(:,:,k), g, is, ie, js, je, direction)
319  enddo

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