MOM6
mom_debugging::check_redundant_t Interface Reference

Detailed Description

Check for consistency between the duplicated points of an A-grid vector or scalar.

Definition at line 46 of file MOM_debugging.F90.

Private functions

subroutine check_redundant_st3d (mesg, array, G, is, ie, js, je)
 Check for consistency between the duplicated points of a 3-D scalar at tracer points. More...
 
subroutine check_redundant_st2d (mesg, array, G, is, ie, js, je)
 Check for consistency between the duplicated points of a 2-D scalar at tracer points. More...
 
subroutine check_redundant_vt3d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 Check for consistency between the duplicated points of a 3-D A-grid vector. More...
 
subroutine check_redundant_vt2d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 Check for consistency between the duplicated points of a 2-D A-grid vector. More...
 

Functions and subroutines

◆ check_redundant_st2d()

subroutine mom_debugging::check_redundant_t::check_redundant_st2d ( character(len=*), intent(in)  mesg,
real, dimension(g%isd:,g%jsd:), 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 tracer 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 426 of file MOM_debugging.F90.

426  character(len=*), intent(in) :: mesg !< An identifying message
427  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
428  real, dimension(G%isd:,G%jsd:), intent(in) :: array !< The array to be checked for consistency
429  integer, optional, intent(in) :: is !< The starting i-index to check
430  integer, optional, intent(in) :: ie !< The ending i-index to check
431  integer, optional, intent(in) :: js !< The starting j-index to check
432  integer, optional, intent(in) :: je !< The ending j-index to check
433  ! Local variables
434  real :: a_nonsym(G%isd:G%ied,G%jsd:G%jed)
435  character(len=128) :: mesg2
436 
437  integer :: i, j, is_ch, ie_ch, js_ch, je_ch
438  integer :: Isq, Ieq, Jsq, Jeq, isd, ied, jsd, jed
439  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
440 
441  is_ch = g%isc ; ie_ch = g%iec ; js_ch = g%jsc ; je_ch = g%jec
442  if (present(is)) is_ch = is ; if (present(ie)) ie_ch = ie
443  if (present(js)) js_ch = js ; if (present(js)) je_ch = je
444 
445  ! This only works on points outside of the standard computational domain.
446  if ((is_ch == g%isc) .and. (ie_ch == g%iec) .and. &
447  (js_ch == g%jsc) .and. (je_ch == g%jec)) return
448 
449  do i=isd,ied ; do j=jsd,jed
450  a_nonsym(i,j) = array(i,j)
451  enddo ; enddo
452 
453  call pass_var(a_nonsym, g%Domain)
454 
455  do i=is_ch,ie_ch ; do j=js_ch,je_ch
456  if (a_nonsym(i,j) /= array(i,j) .and. &
457  redundant_prints(1) < max_redundant_prints) then
458  write(mesg2,'(" Redundant points",2(1pe12.4)," differ by ", &
459  & 1pe12.4," at i,j = ",2i4," on pe ",i4)') &
460  array(i,j), a_nonsym(i,j),array(i,j)-a_nonsym(i,j),i,j,pe_here()
461  write(0,'(A130)') trim(mesg)//trim(mesg2)
462  redundant_prints(1) = redundant_prints(1) + 1
463  endif
464  enddo ; enddo
465 

◆ check_redundant_st3d()

subroutine mom_debugging::check_redundant_t::check_redundant_st3d ( character(len=*), intent(in)  mesg,
real, dimension(g%isd:,g%jsd:,:), 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 tracer 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 401 of file MOM_debugging.F90.

401  character(len=*), intent(in) :: mesg !< An identifying message
402  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
403  real, dimension(G%isd:,G%jsd:,:), intent(in) :: array !< The array to be checked for consistency
404  integer, optional, intent(in) :: is !< The starting i-index to check
405  integer, optional, intent(in) :: ie !< The ending i-index to check
406  integer, optional, intent(in) :: js !< The starting j-index to check
407  integer, optional, intent(in) :: je !< The ending j-index to check
408  ! Local variables
409  character(len=24) :: mesg_k
410  integer :: k
411 
412  do k=1,size(array,3)
413  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
414  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
415  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
416  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
417 
418  call check_redundant_st2d(trim(mesg)//trim(mesg_k), array(:,:,k), &
419  g, is, ie, js, je)
420  enddo

◆ check_redundant_vt2d()

subroutine mom_debugging::check_redundant_t::check_redundant_vt2d ( character(len=*), intent(in)  mesg,
real, dimension(g%isd:,g%jsd:), intent(in)  u_comp,
real, dimension(g%isd:,g%jsd:), 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 A-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 501 of file MOM_debugging.F90.

501  character(len=*), intent(in) :: mesg !< An identifying message
502  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
503  real, dimension(G%isd:,G%jsd:), intent(in) :: u_comp !< The u-component of the vector
504  !! to be checked for consistency
505  real, dimension(G%isd:,G%jsd:), intent(in) :: v_comp !< The v-component of the vector
506  !! to be checked for consistency
507  integer, optional, intent(in) :: is !< The starting i-index to check
508  integer, optional, intent(in) :: ie !< The ending i-index to check
509  integer, optional, intent(in) :: js !< The starting j-index to check
510  integer, optional, intent(in) :: je !< The ending j-index to check
511  integer, optional, intent(in) :: direction !< the direction flag to be
512  !! passed to pass_vector
513  ! Local variables
514  real :: u_nonsym(G%isd:G%ied,G%jsd:G%jed)
515  real :: v_nonsym(G%isd:G%ied,G%jsd:G%jed)
516  character(len=128) :: mesg2
517 
518  integer :: i, j, is_ch, ie_ch, js_ch, je_ch
519  integer :: Isq, Ieq, Jsq, Jeq, isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
520  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
521  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
522  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
523 
524  is_ch = g%isc ; ie_ch = g%iec ; js_ch = g%jsc ; je_ch = g%jec
525  if (present(is)) is_ch = is ; if (present(ie)) ie_ch = ie
526  if (present(js)) js_ch = js ; if (present(js)) je_ch = je
527 
528  ! This only works on points outside of the standard computational domain.
529  if ((is_ch == g%isc) .and. (ie_ch == g%iec) .and. &
530  (js_ch == g%jsc) .and. (je_ch == g%jec)) return
531 
532  do i=isd,ied ; do j=jsd,jed
533  u_nonsym(i,j) = u_comp(i,j) ; v_nonsym(i,j) = v_comp(i,j)
534  enddo ; enddo
535 
536  call pass_vector(u_nonsym, v_nonsym, g%Domain, direction, stagger=agrid)
537 
538  do i=is_ch,ie_ch ; do j=js_ch+1,je_ch
539  if (u_nonsym(i,j) /= u_comp(i,j) .and. &
540  redundant_prints(1) < max_redundant_prints) then
541  write(mesg2,'(" redundant u-components",2(1pe12.4)," differ by ", &
542  & 1pe12.4," at i,j = ",2i4," on pe ",i4)') &
543  u_comp(i,j), u_nonsym(i,j),u_comp(i,j)-u_nonsym(i,j),i,j,pe_here()
544  write(0,'(A130)') trim(mesg)//trim(mesg2)
545  redundant_prints(1) = redundant_prints(1) + 1
546  endif
547  enddo ; enddo
548  do i=is_ch+1,ie_ch ; do j=js_ch,je_ch
549  if (v_nonsym(i,j) /= v_comp(i,j) .and. &
550  redundant_prints(1) < max_redundant_prints) then
551  write(mesg2,'(" redundant v-comps",2(1pe12.4)," differ by ", &
552  & 1pe12.4," at i,j = ",2i4," x,y = ",2(1pe12.4)" on pe ",i4)') &
553  v_comp(i,j), v_nonsym(i,j),v_comp(i,j)-v_nonsym(i,j),i,j, &
554  g%geoLonBu(i,j), g%geoLatBu(i,j), pe_here()
555  write(0,'(A155)') trim(mesg)//trim(mesg2)
556  redundant_prints(1) = redundant_prints(1) + 1
557  endif
558  enddo ; enddo
559 

◆ check_redundant_vt3d()

subroutine mom_debugging::check_redundant_t::check_redundant_vt3d ( character(len=*), intent(in)  mesg,
real, dimension(g%isd:,g%jsd:,:), intent(in)  u_comp,
real, dimension(g%isd:,g%jsd:,:), 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 A-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 471 of file MOM_debugging.F90.

471  character(len=*), intent(in) :: mesg !< An identifying message
472  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
473  real, dimension(G%isd:,G%jsd:,:), intent(in) :: u_comp !< The u-component of the vector
474  !! to be checked for consistency
475  real, dimension(G%isd:,G%jsd:,:), intent(in) :: v_comp !< The v-component of the vector
476  !! to be checked for consistency
477  integer, optional, intent(in) :: is !< The starting i-index to check
478  integer, optional, intent(in) :: ie !< The ending i-index to check
479  integer, optional, intent(in) :: js !< The starting j-index to check
480  integer, optional, intent(in) :: je !< The ending j-index to check
481  integer, optional, intent(in) :: direction !< the direction flag to be
482  !! passed to pass_vector
483  ! Local variables
484  character(len=24) :: mesg_k
485  integer :: k
486 
487  do k=1,size(u_comp,3)
488  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
489  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
490  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
491  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
492 
493  call check_redundant_vt2d(trim(mesg)//trim(mesg_k), u_comp(:,:,k), &
494  v_comp(:,:,k), g, is, ie, js, je, direction)
495  enddo

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