NASA Logo
Ocean Color Science Software

ocssw V2022
specific_ancillary.f90
Go to the documentation of this file.
2 
3  implicit none
4 
5  INTEGER :: numberofshortwavelengths ! to change for OCI
6  INTEGER, PARAMETER :: numberoflongwavelengths = 2
7  INTEGER, DIMENSION(8) :: bandindexmapsw=(/ 1,2,3,4,5,6,7,8 /)
8  INTEGER, DIMENSION(NumberOfLongWavelengths) :: bandindexmaplw=(/ 7,8 /)
9 
10  integer, parameter :: trans_nband = 7
11  integer, parameter :: index_2way = trans_nband
12 
13 ! Ackerman 1971
14  real, parameter :: ozone_absorp_coeff = 5.56209e-5 ! 2.07e-21 * 2.687e16
15  real, parameter :: o3_coef_470 = 4.06e-22 * 2.687e16
16  real, parameter :: o3_coef_550 = 3.36e-21 * 2.687e16
17 
18 
19  logical :: do_37, have_fascode
20 
21  logical :: have_band(9)
22 
23 
24 contains
25 
26  subroutine setup_atm_corr
27 ! WDR it is of interest that have_band is only used in specific_ancillary.f90
28 ! These are not in the radiance measurement order, so...?
30  use mod06_run_settings, only : band_0370
31 
32  do_37 = .true.
33  if( c2_cmp_there(band_0370) == 0 ) do_37 = .false.
34  if( ( c2_sensor_id == oci_id ) .or. ( c2_sensor_id == ocis_id ) ) &
35  do_37 = .false.
36  have_fascode = .true.
37  ! make have_band include the 2.2 um band
38  ! 470, 550, 650, 0.86, 0.94, 1.2, 1.6, 2.1, 2.2
39  have_band = (/.true., .true., .true., .true., .true., .true., .true., &
40  .true., .true. /)
41  if( ( c2_sensor_id == oci_id ) .or. ( c2_sensor_id == ocis_id ) ) then
43  else
44  have_band(9) = .false.
46  endif
47 
48  end subroutine setup_atm_corr
49 !
50 ! *** remap_bands()
51 !
52  subroutine remap_bands(tau2way, temp_trans, sdev2way, temp_sdev)
54  use ch_xfr, only : c2_sensor_id, oci_id, ocis_id
55 
56  real, dimension(:), intent(in) :: temp_trans, temp_sdev
57  real, dimension(:), intent(inout) :: tau2way, sdev2way
58 
59  ! for the non OCI branch
60  if ( (c2_sensor_id /= oci_id) .and. (c2_sensor_id /= ocis_id) ) then
61  tau2way(band_0065) = temp_trans(1)
62  tau2way(band_0086) = temp_trans(2)
63  tau2way(band_0124) = temp_trans(4)
64  tau2way(band_0163) = temp_trans(5)
65  tau2way(band_0213) = temp_trans(6)
66  tau2way(band_0935) = temp_trans(3)
67  tau2way(band_0370) = temp_trans(7)
68 
69  if (temp_sdev(1) /= -10000.) then
70  sdev2way(band_0065) = temp_sdev(1)
71  sdev2way(band_0086) = temp_sdev(2)
72  sdev2way(band_0124) = temp_sdev(4)
73  sdev2way(band_0163) = temp_sdev(5)
74  sdev2way(band_0213) = temp_sdev(6)
75  sdev2way(band_0935) = temp_sdev(3)
76  sdev2way(band_0370) = temp_sdev(7)
77  endif
78  else ! for the OCI and OCIS
79  tau2way(band_0065) = temp_trans(1)
80  tau2way(band_0086) = temp_trans(2)
81  tau2way(band_0124) = temp_trans(4)
82  tau2way(band_0163) = temp_trans(6)
83  tau2way(band_0213) = temp_trans(7)
84  tau2way(band_0935) = temp_trans(3)
85  tau2way(band_0226) = temp_trans(8)
86 
87  if (temp_sdev(1) /= -10000.) then
88  sdev2way(band_0065) = temp_sdev(1)
89  sdev2way(band_0086) = temp_sdev(2)
90  sdev2way(band_0124) = temp_sdev(4)
91  sdev2way(band_0163) = temp_sdev(6)
92  sdev2way(band_0213) = temp_sdev(7)
93  sdev2way(band_0935) = temp_sdev(3)
94  sdev2way(band_0226) = temp_sdev(8)
95  endif
96  endif
97 
98  end subroutine remap_bands
99 
100 
101 ! This subroutine is basically a stub as far as MODIS is concerned
102 ! subroutine get_specific_ancillary( mod06input_filedata, &
103 ! pressure_name, &
104 ! temperature_name, &
105 ! ctm_name, &
106 ! sfc_temp_name, &
107 ! cpi_name, &
108 ! mod06_start, &
109 ! mod06_stride, &
110 ! mod06_edge, EXECUTE_STANDARD, REPLACE_ALBEDO)
111 !
112 ! integer, intent(inout),dimension(:) :: mod06_start, mod06_stride, mod06_edge, mod06input_filedata
113 ! logical, intent(inout) :: EXECUTE_STANDARD, REPLACE_ALBEDO
114 !
115 ! character(*), intent(inout) :: pressure_name, temperature_name, ctm_name, sfc_temp_name, cpi_name
116 !
117 !#ifdef CT_1KM
118 !
119 ! pressure_name = "cloud_top_pressure_1km"
120 ! temperature_name = "cloud_top_temperature_1km"
121 ! ctm_name = "cloud_top_method_1km"
122 ! sfc_temp_name = "surface_temperature_1km"
123 ! cpi_name = "Cloud_Phase_Infrared_1km"
124 !
125 ! mod06_start = mod06_start
126 ! mod06_edge = mod06_edge
127 !
128 !#else
129 !
130 ! pressure_name = "Cloud_Top_Pressure"
131 ! temperature_name = "Cloud_Top_Temperature"
132 ! ctm_name = "Cloud_Height_Method"
133 ! sfc_temp_name = "Surface_Temperature"
134 ! cpi_name = "Cloud_Phase_Infrared"
135 !
136 ! mod06_start = floor(real(mod06_start)/5.)
137 ! mod06_edge = floor(real(mod06_edge)/5.)
138 !
139 !
140 !#endif
141 !
142 ! mod06_stride = 1
143 !
144 ! EXECUTE_STANDARD = .true.
145 ! REPLACE_ALBEDO = .false.
146 !
147 !
148 ! end subroutine get_specific_ancillary
149 
150 
151 ! subroutine get_cloud_top_properties(mapi_filedata, &
152 ! pressure_sds, &
153 ! temperature_sds, &
154 ! ctm_sds, &
155 ! sfc_temp_sds, &
156 ! cpi_sds, &
157 ! start, &
158 ! stride, &
159 ! edge, &
160 ! start_1km, &
161 ! edge_1km, &
162 ! cloud_top_pressure, &
163 ! cloud_top_temperature, &
164 ! cloud_height_method, &
165 ! surface_temperature, &
166 ! cloud_phase_infrared, &
167 ! EXECUTE_STANDARD, &
168 ! status)
169 !
170 ! use GeneralAuxType
171 ! use nonscience_parameters
172 ! use mod06_run_settings
173 ! use general_array_io,only: read_int_array, read_byte_array
174 ! implicit none
175 !
176 ! character(*), intent(in) :: pressure_sds, temperature_sds, ctm_sds, sfc_temp_sds, cpi_sds
177 ! integer, intent(in) :: mapi_filedata(:)
178 ! integer, intent(inout), dimension(:) :: start, stride, edge, start_1km, edge_1km
179 ! real(single), dimension(:,:), intent(out) :: cloud_top_pressure, cloud_top_temperature, surface_temperature
180 ! integer*1, dimension(:,:), intent(out ) :: cloud_height_method, cloud_phase_infrared
181 ! integer, intent(out) :: status
182 ! logical, intent(in) :: EXECUTE_STANDARD
183 !
184 ! integer(integer_twobyte) :: fill_value
185 ! integer :: i, j, local_i, local_j, dimsizes(2)
186 ! real(single), dimension(:,:), allocatable :: ctp_temp_array, ctt_temp_array, sfc_temp_array
187 ! integer*1, dimension(:,:), allocatable :: ctm_temp_array, cpi_temp_array
188 ! logical :: offset
189 ! integer :: localstride(2)
190 !
191 !
192 ! cloud_top_pressure = 0.
193 ! cloud_top_temperature = 0.
194 ! cloud_height_method = 0
195 ! cloud_phase_infrared = 0.
196 !
197 ! offset = .true.
198 !
199 !#ifdef CT_1KM
200 !
201 ! call read_int_array(mapi_filedata, pressure_sds, start, stride, edge, offset, cloud_top_pressure, status)
202 ! call read_int_array(mapi_filedata, temperature_sds, start, stride, edge, offset, cloud_top_temperature, status)
203 ! call read_int_array(mapi_filedata, sfc_temp_sds, start, stride, edge, offset, surface_temperature, status)
204 ! call read_byte_array(mapi_filedata, ctm_sds, start, stride, edge, cloud_height_method, status)
205 ! call read_byte_array(mapi_filedata, cpi_sds, start, stride, edge, cloud_phase_infrared, status)
206 !
207 !
208 !#else
209 !
210 ! allocate(ctp_temp_array(edge(1), edge(2)))
211 ! allocate(ctt_temp_array(edge(1), edge(2)))
212 ! allocate(ctm_temp_array(edge(1), edge(2)))
213 ! allocate(sfc_temp_array(edge(1), edge(2)))
214 ! allocate(cpi_temp_array(edge(1), edge(2)))
215 !
216 ! call read_int_array(mapi_filedata, pressure_sds, start, stride, edge, offset, ctp_temp_array, status)
217 ! call read_int_array(mapi_filedata, temperature_sds, start, stride, edge, offset, ctt_temp_array, status)
218 ! call read_int_array(mapi_filedata, sfc_temp_sds, start, stride, edge, offset, sfc_temp_array, status)
219 ! call read_byte_array(mapi_filedata, ctm_sds, start, stride, edge, ctm_temp_array, status)
220 ! call read_byte_array(mapi_filedata, cpi_sds, start, stride, edge, cpi_temp_array, status)
221 !
222 ! dimsizes(1) = size(cloud_top_pressure, 1)
223 ! dimsizes(2) = size(cloud_top_pressure, 2)
224 !
225 ! if (status /= success) then
226 ! call local_message_handler('Problem reported, see earlier message/s',status,'get_cloud_top_properties')
227 ! endif
228 !
229 !
230 ! do i = 1, dimsizes(1)
231 ! local_i = (i - 1) / 5 + 1
232 ! if (local_i > edge(1) ) local_i = edge(1)
233 ! do j = 1, dimsizes(2)
234 ! local_j = (j - 1) / 5 + 1
235 ! if (local_j > edge(2)) local_j = edge(2)
236 ! cloud_top_pressure(i,j) = ctp_temp_array(local_i, local_j)
237 ! cloud_top_temperature(i,j) = ctt_temp_array(local_i, local_j)
238 ! cloud_height_method(i,j) = ctm_temp_array(local_i, local_j)
239 ! surface_temperature(i,j) = sfc_temp_array(local_i, local_j)
240 ! cloud_phase_infrared(i,j) = cpi_temp_array(local_i, local_j)
241 !
242 ! enddo
243 ! enddo
244 ! deallocate( ctt_temp_array, ctp_temp_array, ctm_temp_array, sfc_temp_array, cpi_temp_array )
245 !
246 !#endif
247 !
248 !
249 ! end subroutine get_cloud_top_properties
250 
251 
252 end module specific_ancillary
Definition: ch_xfr.f90:1
integer, parameter index_2way
integer, parameter band_0124
integer, dimension(8) bandindexmapsw
integer ocis_id
Definition: ch_xfr.f90:51
integer, dimension(numberoflongwavelengths) bandindexmaplw
real, parameter ozone_absorp_coeff
integer, parameter numberoflongwavelengths
integer, parameter trans_nband
integer, parameter band_0370
integer, parameter band_0226
integer, parameter band_0086
integer, parameter band_0213
real, parameter o3_coef_470
real, parameter o3_coef_550
integer, parameter band_0163
integer c2_sensor_id
Definition: ch_xfr.f90:50
integer, parameter band_0065
integer, parameter band_0935
integer oci_id
Definition: ch_xfr.f90:52
integer *1, dimension(:), allocatable c2_cmp_there
Definition: ch_xfr.f90:43
subroutine remap_bands(tau2way, temp_trans, sdev2way, temp_sdev)
logical, dimension(9) have_band