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