Merge branch 'release-v4.6.0'
[WPS.git] / geogrid / src / bitarray_module.F
blob31a4d23be2c31e493f8f541d8f1cec090e42f003
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 ! Module: bitarray_module
4 ! Purpose: This module provides a two-dimensional bit array and a set of 
5 !   routines to manipulate and examine the bits of the array.
6 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
7 module bitarray_module
9    use module_debug
11    type bitarray
12       integer, pointer, dimension(:,:) :: iarray ! Storage array
13       integer :: nx, ny                 ! Number of bits in the x and y directions
14       integer :: x_int_dim, y_int_dim   ! Number of integers in the x and y directions
15       integer :: integer_size           ! Number of bits in an integer
16    end type bitarray
18    contains
20    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
21    ! Name: bitarray_create
22    !
23    ! Purpose: Allocate and initialize a bit array so that all bits are FALSE 
24    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
25    subroutine bitarray_create(b, i, j)
26    
27       implicit none
28   
29       ! Arguments
30       integer, intent(in) :: i, j   
31       type (bitarray), intent(out) :: b
32   
33       b%integer_size = bit_size(b%integer_size) 
34   
35       b%nx = i
36       b%ny = j
37   
38       b%x_int_dim = ceiling(real(b%nx)/real(b%integer_size))
39       b%y_int_dim = b%ny 
41       nullify(b%iarray)
42       allocate(b%iarray(b%x_int_dim, b%y_int_dim))
43       b%iarray = 0
45    end subroutine bitarray_create
48    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
49    ! Name: bitarray_copy
50    !
51    ! Purpose: Duplicate a bitarray.
52    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
53    subroutine bitarray_copy(src, dst)
54    
55       implicit none
56   
57       ! Arguments
58       type (bitarray), intent(in) :: src
59       type (bitarray), intent(out) :: dst
60   
61       dst%integer_size = src%integer_size
62   
63       dst%nx = src%nx
64       dst%ny = src%ny
65   
66       dst%x_int_dim = src%x_int_dim
67       dst%y_int_dim = src%y_int_dim
69       allocate(dst%iarray(dst%x_int_dim, dst%y_int_dim))
70       dst%iarray = src%iarray
72    end subroutine bitarray_copy
75    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
76    ! Name: bitarray_set
77    !
78    ! Purpose: Set the bit located at (i,j) to TRUE
79    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
80    subroutine bitarray_set(b, i, j)
82       implicit none
83   
84       ! Arguments
85       integer, intent(in) :: i, j
86       type (bitarray), intent(inout) :: b
87   
88       ! Local variables
89       integer :: n_integer, n_bit
90   
91       n_integer = ((i-1) / b%integer_size) + 1
92       n_bit = mod((i-1), b%integer_size)
93   
94       b%iarray(n_integer, j) = ibset(b%iarray(n_integer, j), n_bit) 
96    end subroutine bitarray_set
99    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
100    ! Name: bitarray_clear
101    !
102    ! Purpose: Set the bit located at (i,j) to FALSE 
103    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
104    subroutine bitarray_clear(b, i, j)
106       implicit none
107   
108       ! Arguments
109       integer, intent(in) :: i, j
110       type (bitarray), intent(inout) :: b
111   
112       ! Local variables
113       integer :: n_integer, n_bit
114   
115       n_integer = ((i-1) / b%integer_size) + 1
116       n_bit = mod((i-1), b%integer_size)
117   
118       b%iarray(n_integer, j) = ibclr(b%iarray(n_integer, j), n_bit) 
120    end subroutine bitarray_clear
123    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
124    ! Name: bitarray_test
125    !
126    ! Purpose: To return the value of the bit located at (i,j)
127    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
128    function bitarray_test(b, i, j)
130       implicit none
131   
132       ! Arguments
133       integer, intent(in) :: i, j
134       type (bitarray), intent(in) :: b
135   
136       ! Local variables
137       logical :: bitarray_test
138       integer :: n_integer, n_bit
139   
140       n_integer = ((i-1) / b%integer_size) + 1
141       n_bit = mod((i-1), b%integer_size)
142   
143       bitarray_test = btest(b%iarray(n_integer,j), n_bit) 
145    end function bitarray_test
148    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
149    ! Name: bitarray_merge
150    !
151    ! Purpose: The first bitarray argument, b1, is set to the union of the .TRUE. 
152    !   bits in b1 and b2. That is, after returning, a bit x in b1 is set if
153    !   either x was set in b1 or x was set in b2. Thus, b1 AND b2 MUST BE BIT 
154    !   ARRAYS OF THE SAME DIMENSIONS.
155    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
156    subroutine bitarray_merge(b1, b2)
158       implicit none
159   
160       ! Arguments
161       type (bitarray), intent(inout) :: b1, b2
162   
163       ! Local variables
164       integer :: i, j
165   
166       if (b1%x_int_dim /= b2%x_int_dim .or. b1%y_int_dim /= b2%y_int_dim) then
167          call mprintf(.true.,ERROR,'In bitarray_merge(), b1 and b2 have different dimensions.')
168       end if
169   
170       do i=1,b1%x_int_dim
171          do j=1,b1%y_int_dim
172             b1%iarray(i,j) = ior(b1%iarray(i,j), b2%iarray(i,j))  
173          end do
174       end do
176    end subroutine bitarray_merge
177    
179    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
180    ! Name: bitarray_destroy
181    !
182    ! Purpose: To deallocate all allocated memory associated with the bit array
183    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
184    subroutine bitarray_destroy(b)
186       implicit none
187   
188       ! Arguments
189       type (bitarray), intent(inout) :: b
190   
191       if (associated(b%iarray)) then
192          deallocate(b%iarray)
193       else
194          call mprintf(.true.,WARN,'In bitarray_destroy(), b is not allocated.')
195       end if
197    end subroutine bitarray_destroy
199 end module bitarray_module