MOM6
mom_debugging::check_redundant_b Interface Reference

Detailed Description

Definition at line 58 of file MOM_debugging.F90.

Private functions

subroutine check_redundant_vb3d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 
subroutine check_redundant_vb2d (mesg, u_comp, v_comp, G, is, ie, js, je, direction)
 
subroutine check_redundant_sb3d (mesg, array, G, is, ie, js, je)
 
subroutine check_redundant_sb2d (mesg, array, G, is, ie, js, je)
 

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

Definition at line 242 of file MOM_debugging.F90.

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

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

Definition at line 217 of file MOM_debugging.F90.

217  character(len=*), intent(in) :: mesg
218  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
219  real, dimension(G%IsdB:,G%JsdB:,:), intent(in) :: array
220  integer, optional, intent(in) :: is, ie, js, je
221 ! Arguments: array - The array being checked.
222 ! (in) mesg - A message indicating what is being checked.
223 ! (in) G - The ocean's grid structure.
224 ! (in/opt) is, ie, js, je - the i- and j- range of indices to check.
225 
226  character(len=24) :: mesg_k
227  integer :: k
228 
229  do k=1,size(array,3)
230  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
231  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
232  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
233  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
234 
235  call check_redundant_sb2d(trim(mesg)//trim(mesg_k), array(:,:,k), &
236  g, is, ie, js, je)
237  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
Parameters
[in,out]gThe ocean's grid structure

Definition at line 331 of file MOM_debugging.F90.

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

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

Definition at line 302 of file MOM_debugging.F90.

302  character(len=*), intent(in) :: mesg
303  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
304  real, dimension(G%IsdB:,G%JsdB:,:), intent(in) :: u_comp
305  real, dimension(G%IsdB:,G%JsdB:,:), intent(in) :: v_comp
306  integer, optional, intent(in) :: is, ie, js, je
307  integer, optional, intent(in) :: direction
308 ! Arguments: u_comp - The u-component of the vector being checked.
309 ! (in) v_comp - The v-component of the vector being checked.
310 ! (in) mesg - A message indicating what is being checked.
311 ! (in) G - The ocean's grid structure.
312 ! (in/opt) is, ie, js, je - the i- and j- range of indices to check.
313 ! (in/opt) direction - the direction flag to be passed to pass_vector.
314 
315  character(len=24) :: mesg_k
316  integer :: k
317 
318  do k=1,size(u_comp,3)
319  if (k < 10) then ; write(mesg_k,'(" Layer",i2," ")') k
320  elseif (k < 100) then ; write(mesg_k,'(" Layer",i3," ")') k
321  elseif (k < 1000) then ; write(mesg_k,'(" Layer",i4," ")') k
322  else ; write(mesg_k,'(" Layer",i9," ")') k ; endif
323 
324  call check_redundant_vb2d(trim(mesg)//trim(mesg_k), u_comp(:,:,k), &
325  v_comp(:,:,k), g, is, ie, js, je, direction)
326  enddo

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