OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
remove_edge_scenes.f90
Go to the documentation of this file.
1  subroutine remove_edge_scenes_ext(&
2  CSR_results, &
3  xsize, ysize, &
4  status)
5 
10  use core_arrays
11 
12  implicit none
13 
14 
15 ! for each scene, if the pixel is not surrounded by cloud we eliminate it
16 ! as a retreival candidate
17 ! In effect we remove all edge pixels
18 
19  integer, intent(in) :: xsize, ysize
20  integer(integer_onebyte), dimension(xsize,ysize), intent(inout) :: CSR_results
21  integer, intent(out) :: status
22 
23  integer :: i, j, tempmask(3,3), lo_i, hi_i, lo_j, hi_j, sum_tempmask, main_tempmask, ii, jj
24  integer, dimension(:,:), allocatable :: mask_solution
25 
26 
27  status = 0
28 
29 #if MAS_INST | EMAS_INST
30 ! this code is modified to be using the CSR results from the hi-res CSR product MAS_PR06CSR
31 
32  do i=1, xsize
33  do j=1, ysize
34  if (restoral_pos(i,j) == 2) then
35  csr_results(i,j) = 1
36  endif
37  end do
38  end do
39 #endif
40 
41 
42  end subroutine remove_edge_scenes_ext
43 
44 
45 
46 subroutine remove_edge_scenes(cloudmask, &
47  CSR_results, &
48  xsize, ysize, &
49  status)
50 ! WDR
51 ! cloudmask I structure(?) with cloud state for each pixel in this chunk
52 ! CSR_results O indicator (= 1) that this pixel has a cloud-free
53 ! pixel adjacent to it
54 ! xsize, ysize I Size of the chunk of data in pix, lin
55 ! status O result = 0 always
56 
60  use generalauxtype
61 
62  implicit none
63 
64 
65 ! for each scene, if the pixel is not surrounded by cloud we eliminate it
66 ! as a retreival candidate
67 ! In effect we remove all edge pixels
68 
69  integer, intent(in) :: xsize, ysize
70  type(processflag), dimension(xsize,ysize), intent(in) :: cloudmask
71  integer(integer_onebyte), dimension(xsize,ysize), intent(inout) :: CSR_results
72  integer, intent(out) :: status
73 
74  integer :: i, j, tempmask(3,3), lo_i, hi_i, lo_j, hi_j, sum_tempmask, main_tempmask, ii, jj
75  integer, dimension(:,:), allocatable :: mask_solution
76 
77 
78  status = 0
79 
80 ! this code is reworked to incorporate the tests for edge lines. Basically it's quite
81 ! simple. If there is at least one point around the current point that is not cloudy,
82 ! we remove that point.
83 
84  i=1
85  j=1
86 
87  allocate(mask_solution(xsize, ysize))
88  mask_solution = 0
89 
90  do i=1, xsize
91  if (i==1) then
92  lo_i = i
93  hi_i = i+1
94  else if (i==xsize) then
95  lo_i = i-1
96  hi_i = i
97  else
98  lo_i = i-1
99  hi_i = i+1
100  endif
101 
102  if (i==1 .or. i==xsize) then
103  sum_tempmask = 6
104  else
105  sum_tempmask = 9
106  endif
107 
108  main_tempmask = sum_tempmask
109 
110  do j=1, ysize
111 
112  if (j==1) then
113  lo_j = j
114  hi_j = j+1
115  else if (j==ysize) then
116  lo_j = j-1
117  hi_j = j
118  else
119  lo_j = j-1
120  hi_j = j+1
121  endif
122 
123  if (j==1 .or. j==ysize) then
124  if (sum_tempmask == 6) then
125  sum_tempmask = 4
126  else
127  sum_tempmask = 6
128  endif
129  endif
130 
131 
132  tempmask = 0
133 
134  ! attempt edge removal only if a pixel is cloudy to begin with, otherwise, why would we
135  ! ever want to run the restoral in the first place. G.Wind 4.20.05
136  if (cloudmask(i,j)%cloudobserved) then
137 
138  do ii=lo_i, hi_i
139  do jj=lo_j, hi_j
140  if (cloudmask(ii,jj)%cloudobserved .and. csr_results(ii,jj) == 0) then
141  tempmask( ii-lo_i+1, jj-lo_j+1) = 1
142  end if
143  end do
144  end do
145 
146  if (sum(tempmask) < sum_tempmask) &
147  mask_solution(i,j) = 1
148 
149  endif
150 
151  sum_tempmask = main_tempmask
152 
153  enddo
154  enddo
155 
156  do i=1, xsize
157  do j=1, ysize
158 ! don't overwrite the pixels with other CSR results
159  if (mask_solution(i,j) == 1 .and. csr_results(i,j) == 0) then
160  csr_results(i,j) = 1
161  endif
162 
163  end do
164  end do
165 
166 
167  deallocate(mask_solution)
168 
169 end subroutine remove_edge_scenes
subroutine remove_edge_scenes_ext(CSR_results, xsize, ysize, status)
integer(integer_onebyte), dimension(:,:), allocatable restoral_pos
Definition: core_arrays.f90:84
subroutine remove_edge_scenes(cloudmask, CSR_results, xsize, ysize, status)