MOM6
mom_debugging::check_redundant_t Interface Reference

Detailed Description

Definition at line 62 of file MOM_debugging.F90.

Private functions

subroutine check_redundant_st3d (mesg, array, G, is, ie, js, je)
 
subroutine check_redundant_st2d (mesg, array, G, is, ie, js, je)
 
subroutine check_redundant_vt3d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 
subroutine check_redundant_vt2d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 

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

Definition at line 431 of file MOM_debugging.F90.

431  character(len=*), intent(in) :: mesg
432  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
433  real, dimension(G%isd:,G%jsd:), intent(in) :: array
434  integer, optional, intent(in) :: is, ie, js, je
435 ! Arguments: array - The array being checked.
436 ! (in) mesg - A message indicating what is being checked.
437 ! (in) G - The ocean's grid structure.
438 ! (in/opt) is, ie, js, je - the i- and j- range of indices to check.
439 
440  real :: a_nonsym(g%isd:g%ied,g%jsd:g%jed)
441  character(len=128) :: mesg2
442 
443  integer :: i, j, is_ch, ie_ch, js_ch, je_ch
444  integer :: isq, ieq, jsq, jeq, isd, ied, jsd, jed
445  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
446 
447  is_ch = g%isc ; ie_ch = g%iec ; js_ch = g%jsc ; je_ch = g%jec
448  if (present(is)) is_ch = is ; if (present(ie)) ie_ch = ie
449  if (present(js)) js_ch = js ; if (present(js)) je_ch = je
450 
451  ! This only works on points outside of the standard computational domain.
452  if ((is_ch == g%isc) .and. (ie_ch == g%iec) .and. &
453  (js_ch == g%jsc) .and. (je_ch == g%jec)) return
454 
455  do i=isd,ied ; do j=jsd,jed
456  a_nonsym(i,j) = array(i,j)
457  enddo ; enddo
458 
459  call pass_var(a_nonsym, g%Domain)
460 
461  do i=is_ch,ie_ch ; do j=js_ch,je_ch
462  if (a_nonsym(i,j) /= array(i,j) .and. &
463  redundant_prints(1) < max_redundant_prints) then
464  write(mesg2,'(" Redundant points",2(1pe12.4)," differ by ", &
465  & 1pe12.4," at i,j = ",2i4," on pe ",i4)') &
466  array(i,j), a_nonsym(i,j),array(i,j)-a_nonsym(i,j),i,j,pe_here()
467  write(0,'(A130)') trim(mesg)//trim(mesg2)
468  redundant_prints(1) = redundant_prints(1) + 1
469  endif
470  enddo ; enddo
471 

◆ 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
Parameters
[in,out]gThe ocean's grid structure

Definition at line 406 of file MOM_debugging.F90.

406  character(len=*), intent(in) :: mesg
407  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
408  real, dimension(G%isd:,G%jsd:,:), intent(in) :: array
409  integer, optional, intent(in) :: is, ie, js, je
410 ! Arguments: array - The array being checked.
411 ! (in) mesg - A message indicating what is being checked.
412 ! (in) G - The ocean's grid structure.
413 ! (in/opt) is, ie, js, je - the i- and j- range of indices to check.
414 
415  character(len=24) :: mesg_k
416  integer :: k
417 
418  do k=1,size(array,3)
419  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
420  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
421  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
422  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
423 
424  call check_redundant_st2d(trim(mesg)//trim(mesg_k), array(:,:,k), &
425  g, is, ie, js, je)
426  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
Parameters
[in,out]gThe ocean's grid structure

Definition at line 506 of file MOM_debugging.F90.

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

◆ 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
Parameters
[in,out]gThe ocean's grid structure

Definition at line 477 of file MOM_debugging.F90.

477  character(len=*), intent(in) :: mesg
478  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
479  real, dimension(G%isd:,G%jsd:,:), intent(in) :: u_comp
480  real, dimension(G%isd:,G%jsd:,:), intent(in) :: v_comp
481  integer, optional, intent(in) :: is, ie, js, je
482  integer, optional, intent(in) :: direction
483 ! Arguments: u_comp - The u-component of the vector being checked.
484 ! (in) v_comp - The v-component of the vector being checked.
485 ! (in) mesg - A message indicating what is being checked.
486 ! (in) G - The ocean's grid structure.
487 ! (in/opt) is, ie, js, je - the i- and j- range of indices to check.
488 ! (in/opt) direction - the direction flag to be passed to pass_vector.
489 
490  character(len=24) :: mesg_k
491  integer :: k
492 
493  do k=1,size(u_comp,3)
494  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
495  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
496  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
497  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
498 
499  call check_redundant_vt2d(trim(mesg)//trim(mesg_k), u_comp(:,:,k), &
500  v_comp(:,:,k), g, is, ie, js, je, direction)
501  enddo

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