2D AMRCLAW
trimbd.f90
Go to the documentation of this file.
1 !
2 ! :::::::::::::::::::::::: TRIMBD :::::::::::::::::::::::::::;
3 ! if used array is completely set (=1.) then return set=true,
4 ! otherwise return false, along with the dimensions of the smallest
5 ! rectangle containing all unset points in il,ir,jb,jt.
6 ! ::::::::::::::::::::::::::::::::::::;::::::::::::::::::::::;
7 !
32 subroutine trimbd(used,nrow,ncol,set,unset_rect)
33 
34  implicit none
35 
36  ! Input
37  integer, intent(in) :: nrow, ncol
38  integer(kind=1), intent(in) :: used(nrow,ncol)
39 
40  ! Output
41  logical, intent(out) :: set
42 ! integer, intent(out) :: il, ir, jb, jt
43  integer, intent(out) :: unset_rect(4)
44 
45  ! Locals
46  integer :: i, j, utot
47  integer(kind=1) :: check
48 
49  utot = 0
50  do 100 j = 1,ncol
51  do 100 i = 1,nrow
52 100 utot = utot + used(i,j)
53 
54  if (utot .eq. nrow * ncol ) then
55  set = .true.
56  else
57  set = .false.
58 
59  check = 1
60  do i = 1,nrow
61  do j = 1,ncol
62  check = min(check,used(i,j))
63  enddo
64  unset_rect(1) = i
65  if (check == 0) exit
66  enddo
67 
68  check = 1
69  do i = 1,nrow
70  do j = 1,ncol
71  check = min(check,used(nrow - i + 1,j))
72  enddo
73  unset_rect(2) = nrow - i + 1
74  if (check == 0) exit
75  enddo
76 
77  check = 1
78  do j = 1,ncol
79  do i = 1,nrow
80  check = min(check,used(i,j))
81  enddo
82  unset_rect(3) = j
83  if (check == 0) exit
84  enddo
85 
86  check = 1
87  do j = 1,ncol
88  do i = 1,nrow
89  check = min(check,used(i,ncol - j + 1))
90  enddo
91  unset_rect(4) = ncol - j + 1
92  if (check == 0) exit
93  enddo
94 
95  endif
96 
97 end subroutine trimbd
subroutine trimbd(used, nrow, ncol, set, unset_rect)
Examine the setting status of a patch.
Definition: trimbd.f90:33