da_transform_xtoy_crtm_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_crtm_adj ( iv, jo_grad_y, jo_grad_x )
2 
3    !---------------------------------------------------------------------------
4    ! PURPOSE: transform gradient from obs space to model grid space.
5    !
6    ! METHOD:  jo_grad_x = H^T jo_grad_y =  - H^T R^-1 ( d - H delta_x )
7    !           1. input gradient in obs space and reference state of RTTOV
8    !           2. call adjoint of RTM
9    !           3. adjoint of interpolation from model grid to obs loc
10    !
11    !  HISTORY: 11/19/2006 - Creation            Zhiquan Liu
12    !
13    !---------------------------------------------------------------------------
14 
15    implicit none
16 
17    type (x_type), intent(inout)   :: jo_grad_x ! 
18    type (y_type),  intent(in)     :: jo_grad_y ! H' delta_x
19    type (iv_type), intent(in)     :: iv        ! O-B structure.
20 
21 #ifdef CRTM
22 
23    integer                        :: k  ! Index dimension.
24    integer                        :: num_rad  ! Number of radiance obs
25    integer                        :: inst, nchanl, n
26    real, allocatable              :: q_ad(:,:)
27    real, allocatable              :: t(:,:)
28    real, allocatable              :: p(:)
29 !   type(infa_type), pointer       :: info                
30 !   integer                        :: i,j
31 !   real                           :: dx,dy,dxm,dym
32 
33 
34    ! CRTM local varaibles and types
35    integer :: wmo_sensor_id, Error_Status, Allocate_Status
36    type (CRTM_RTSolution_type ), allocatable :: RTSolution(:,:),RTSolution_AD(:,:)
37    type (CRTM_Atmosphere_type )   :: Atmosphere(1), Atmosphere_AD(1)
38    type (CRTM_Surface_type )      :: Surface(1), Surface_AD(1)
39    type (CRTM_GeometryInfo_type ) :: GeometryInfo(1)
40    !type(CRTM_Options_type )      :: Options
41 
42 !---------------------------------------------------------
43 
44    if ( iv%num_inst < 1 ) return
45 
46    if (trace_use) call da_trace_entry("da_transform_xtoy_crtm_adj")
47 
48 !----------------------------------------------------------------------------
49 ! CRTM allocation
50 !
51 ! Atmosphere structure
52     Atmosphere(1)%n_Layers=(kte-kts)+1   ! number of vertical levels
53     Atmosphere(1)%n_Absorbers=2
54     Atmosphere(1)%n_Clouds=0
55     Atmosphere(1)%n_Aerosols=0
56     if (crtm_cloud) Atmosphere(1)%n_Clouds=6
57  
58     Error_Status = CRTM_Allocate_Atmosphere( Atmosphere(1)%n_Layers, &
59                                              Atmosphere(1)%n_Absorbers, &
60                                              Atmosphere(1)%n_Clouds, &
61                                              Atmosphere(1)%n_Aerosols, &
62                                              Atmosphere)
63     if ( Error_Status /= 0 ) then
64        call da_error(__FILE__,__LINE__, &
65          (/"Error in allocatting CRTM Atmosphere Structure"/))
66     end if
67     Atmosphere(1)%Absorber_ID(1)=H2O_ID
68     Atmosphere(1)%Absorber_ID(2)=O3_ID
69  
70     if (crtm_cloud) then
71        Atmosphere(1)%Cloud(1)%Type=WATER_CLOUD
72        Atmosphere(1)%Cloud(2)%Type=ICE_CLOUD
73        Atmosphere(1)%Cloud(3)%Type=RAIN_CLOUD
74        Atmosphere(1)%Cloud(4)%Type=SNOW_CLOUD
75        Atmosphere(1)%Cloud(5)%Type=GRAUPEL_CLOUD
76        Atmosphere(1)%Cloud(6)%Type=HAIL_CLOUD
77     end if
78 
79 !-------------------------------------------------------------------------------
80 
81    do inst = 1, iv%num_inst                 ! loop for sensor
82       if ( iv%instid(inst)%num_rad < 1 ) cycle
83       num_rad = iv%instid(inst)%num_rad
84       allocate (q_ad(kts:kte,num_rad))
85       allocate (t(kts:kte,num_rad))
86       allocate (p(num_rad))
87 
88 !      info => iv%instid(inst)%info
89 
90   ! CRTM channel information structure
91       ! Error_Status = CRTM_Set_ChannelInfo(Sensor_Descriptor(inst),ChannelInfo)
92       ! if ( Error_Status /= 0) then
93       !    call da_error(__FILE__,__LINE__, &
94       !     (/"Error in calling CRTM_Set_ChannelInfo"/))
95       ! end if
96       nchanl    = ChannelInfo(inst)%n_channels
97  
98   ! Allocate forward model solution RTSolution array to number of channels
99       allocate( RTSolution( ChannelInfo(inst)%n_Channels, 1 )   , &
100                 RTSolution_AD( ChannelInfo(inst)%n_Channels, 1 ), &
101                STAT = Allocate_Status )
102       if ( Allocate_Status /= 0 ) then
103          call da_error(__FILE__,__LINE__, &
104           (/"Error in allocatting RTSolution"/))
105       END IF
106  
107   ! CRTM Surface Structure
108       if (trim(crtm_sensor_name(rtminit_sensor(inst))) =='amsua') then
109          wmo_sensor_id=WMO_AMSUA
110       elseif (trim(crtm_sensor_name(rtminit_sensor(inst))) =='amsub') then
111          wmo_sensor_id=WMO_AMSUB
112       elseif (trim(crtm_sensor_name(rtminit_sensor(inst))) =='amsre') then
113          wmo_sensor_id=WMO_AMSRE
114       elseif (trim(crtm_sensor_name(rtminit_sensor(inst))) =='ssmi') then
115          wmo_sensor_id=WMO_SSMI
116       else
117          wmo_sensor_id=INVALID_WMO_SENSOR_ID
118       end if
119  
120       Error_Status = CRTM_Allocate_Surface( nchanl,     &  ! Input
121                                    Surface)  ! Output
122       if ( Error_Status /= 0 ) then
123         call da_error(__FILE__,__LINE__, &
124           (/"Error in allocatting CRTM Surface Structure"/))
125       end if
126  
127   ! CRTM Options structure
128       !Options%n_channels = nchanl
129       !Error_Status = CRTM_Allocate_Options( nchanl,     &  ! Input
130       !                             Options)  ! Output
131       !if ( Error_Status /= 0 ) then
132       !  call da_error(__FILE__,__LINE__, &
133       !    (/"Error in allocatting CRTM Options Structure"/))
134       !endif
135 
136       do n=1,num_rad
137 
138       ! [1.0] Extract base state Atmosphere variables
139             Atmosphere(1)%level_pressure(0) = iv%instid(inst)%pf(0,n)
140            do k=1,Atmosphere(1)%n_layers
141             Atmosphere(1)%pressure(k) = iv%instid(inst)%pm(k,n)
142             Atmosphere(1)%level_pressure(k) = iv%instid(inst)%pf(k,n)
143             Atmosphere(1)%temperature(k) = iv%instid(inst)%tm(k,n)
144             Atmosphere(1)%absorber(k,1) = iv%instid(inst)%qm(k,n)
145            if (crtm_cloud) then
146             Atmosphere(1)%cloud(1)%water_content(k)=iv%instid(inst)%qcw(k,n)
147             Atmosphere(1)%cloud(2)%water_content(k)=iv%instid(inst)%qci(k,n)
148             Atmosphere(1)%cloud(3)%water_content(k)=iv%instid(inst)%qrn(k,n)
149             Atmosphere(1)%cloud(4)%water_content(k)=iv%instid(inst)%qsn(k,n)
150             Atmosphere(1)%cloud(5)%water_content(k)=iv%instid(inst)%qgr(k,n)
151             Atmosphere(1)%cloud(6)%water_content(k)=0.0
152             Atmosphere(1)%cloud(1)%effective_radius(k)=10.0
153             Atmosphere(1)%cloud(2)%effective_radius(k)=200.0
154             Atmosphere(1)%cloud(3)%effective_radius(k)=200.0
155             Atmosphere(1)%cloud(4)%effective_radius(k)=200.0
156             Atmosphere(1)%cloud(5)%effective_radius(k)=200.0
157             Atmosphere(1)%cloud(6)%effective_radius(k)=200.0
158            end if
159            end do
160 
161       ! [1.1] User-supplied emissivity
162           ! Options%emissivity_switch = 1
163           ! Options%emissivity(1:Options%n_channels) = &
164           !     iv%instid(inst)%emiss(1:Options%n_channels,n)
165 
166       ! [1.1] CRTM Surface parameter data
167      Surface(1)%Land_Coverage=iv%instid(inst)%land_coverage(n)
168      Surface(1)%Water_Coverage=iv%instid(inst)%water_coverage(n)
169      Surface(1)%Snow_Coverage=iv%instid(inst)%snow_coverage(n)
170      Surface(1)%Ice_Coverage=iv%instid(inst)%ice_coverage(n)
171 
172      if (Surface(1)%Land_Coverage > 0.0) then
173        Surface(1)%Land_Type=GRASS_SOIL           ! land type (User guide appendix 3)
174        Surface(1)%Land_Temperature=iv%instid(inst)%ts(n)      ! K
175        Surface(1)%Soil_Moisture_Content=iv%instid(inst)%smois(n)  !0.05    ! volumetric water content (g/cm**3)
176        !Surface(1)%Canopy_Water_Content=0.05      ! gravimetric water content
177        Surface(1)%Vegetation_Fraction=iv%instid(inst)%vegtyp(n)
178        Surface(1)%Soil_Temperature=iv%instid(inst)%tslb(n)
179      end if
180      if (Surface(1)%Water_Coverage > 0.0) then
181        !Surface(1)%Water_Type=SEA_WATER          ! (Currently NOT used)
182        Surface(1)%Water_Temperature=iv%instid(inst)%ts(n)     ! K
183        Surface(1)%Wind_Speed=sqrt((iv%instid(inst)%u10(n))**2+ &
184                                (iv%instid(inst)%v10(n))**2)  ! m/sec
185        !surface(1)%Wind_Direction=0.0            ! NOT used
186        Surface(1)%Salinity=33.0                  ! ppmv
187      end if
188      if (Surface(1)%Snow_Coverage > 0.0) then
189        Surface(1)%Snow_Type=NEW_SNOW             ! User guide appendix 3
190        Surface(1)%Snow_Temperature=iv%instid(inst)%ts(n)      ! K
191        Surface(1)%Snow_Depth=iv%instid(inst)%snowh(n)         ! mm
192        !Surface(1)%Snow_Density=0.2               ! g/cm**3
193        !Surface(1)%Snow_Grain_Size=2.0            ! mm
194      end if
195      if (Surface(1)%Ice_Coverage > 0.0) then
196        !Surface(1)%Ice_Type=FRESH_ICE             ! NO Table offered, single example is FRESH_ICE
197        Surface(1)%Ice_Temperature=iv%instid(inst)%ts(n)       ! K
198        Surface(1)%Ice_Thickness=10.0              ! mm
199        !Surface(1)%Ice_Density=0.9                ! g/cm**3
200        !Surface(1)%Ice_Roughness=0.0               ! NO Table offered, single example is ZERO
201      end if
202      Surface(1)%SensorData%n_channels = nchanl
203      Surface(1)%SensorData%Sensor_ID  = wmo_sensor_id
204      Surface(1)%SensorData%Tb(1:nchanl) = iv%instid(inst)%tb_inv(1:nchanl,n) + &
205                                           iv%instid(inst)%tb_xb(1:nchanl,n)
206      
207       ! -- Copy the adjoint atmosphere structure
208       Error_Status = CRTM_Assign_Atmosphere( Atmosphere, Atmosphere_AD )
209 
210       if ( Error_Status /= 0 ) then
211         call da_error(__FILE__,__LINE__, &
212           (/"Error copying Atmosphere_AD structure"/))
213       END IF
214 
215       ! -- Copy the adjoint surface structure
216       Error_Status = CRTM_Assign_Surface( Surface, Surface_AD )
217 
218       if ( Error_Status /= 0 ) then
219         call da_error(__FILE__,__LINE__, &
220           (/"Error copying Surface_AD structure"/))
221       END IF
222 
223     ! -- Zero the Adjoint outputs 
224     ! Important: adjoint variables must be initialized
225        call CRTM_Zero_Atmosphere( Atmosphere_AD )
226        call CRTM_Zero_Surface( Surface_AD )
227 
228       ! [1.2] CRTM GeometryInfo Structure
229         GeometryInfo(1)%Sensor_Zenith_Angle=iv%instid(inst)%satzen(n)
230         GeometryInfo(1)%Source_Zenith_Angle=iv%instid(inst)%solzen(n)
231    !     GeometryInfo(1)%Satellite_Height=830.0
232    !     GeometryInfo(1)%Sensor_Scan_Angle=
233    !     GeometryInfo(1)%Sensor_Zenith_Angle=
234    !     GeometryInfo(1)%Sensor_Scan_Angle=
235    !     GeometryInfo(1)%Source_Zenith_Angle=
236 
237 
238       ! [1.3] assign tb = R^-1 Re :
239 
240          RTSolution_AD(:,1)%brightness_temperature = jo_grad_y%instid(inst)%tb(:,n)
241          RTSolution_AD(:,1)%radiance = 0.0 ! must assign zero, since each call of AD model will return non-zero value
242        !if (n <= 10) write(6,'(15f8.3)') RTSolution_AD(:,1)%brightness_temperature
243 
244       ! [1.4] Call CRTM_AD model
245 
246          call da_crtm_ad (1, nchanl, 1, Atmosphere,   &
247                             Surface,      &
248                             RTSolution_AD,&
249                             GeometryInfo, &
250                             ChannelInfo(inst:inst),  &
251                             Atmosphere_AD,&
252                             Surface_AD,   &
253                             RTSolution) !,   &
254                             !Options = Options)
255 
256          !Error_Status = CRTM_Adjoint(Atmosphere,   &
257          !                   Surface,      &
258          !                   RTSolution_AD,&
259          !                   GeometryInfo, &
260          !                   ChannelInfo,  &
261          !                   Atmosphere_AD,&
262          !                   Surface_AD,   &
263          !                   RTSolution) !,   &
264          !                   !Options = Options)
265          !if ( Error_Status /= 0 ) then
266          !     call da_error(__FILE__,__LINE__, &
267          !        (/"Error in calling CRTM_Adjoint"/))
268          !endif
269 
270       ! [1.5] Scale transformation and fill zero for no-control variable
271 
272          Atmosphere_AD(1)%Level_Pressure(Atmosphere(1)%n_layers) = &
273                  Atmosphere_AD(1)%Level_Pressure(Atmosphere(1)%n_layers)*0.01 
274          Atmosphere_AD(1)%Absorber(:,1) = 1000.0*Atmosphere_AD(1)%Absorber(:,1) ! in g/kg
275 
276             ! [1.7] Adjoint of Interpolate horizontally from ob to grid:
277 
278             do k=kts,kte ! from bottom to top
279                if (atmosphere(1)%pressure(kte-k+1) < 75.0) then
280                   q_ad(k,n) = 0.0
281                else
282                   q_ad(k,n) = Atmosphere_AD(1)%Absorber(kte-k+1,1)
283                end if
284                t(k,n)=Atmosphere_AD(1)%Temperature(kte-k+1)
285             end do
286 
287             p(n) = Atmosphere_AD(1)%Level_Pressure(atmosphere(1)%n_layers)
288 
289                !if (n <=10 ) then
290                !   write(6,'(15f8.2)') RTSolution(:,1)%brightness_temperature
291                !   write(6,'(e15.5)') jo_grad_x% psfc(i,j)
292                !  do k=kts,kte
293                !   write(6,'(2e15.5)') jo_grad_x%t(i,j,k), jo_grad_x%q(i,j,k)
294                !  end do
295                !endif
296 
297             Error_Status = CRTM_Destroy_Atmosphere( Atmosphere_AD )
298             if ( Error_Status /= 0 ) then
299               call da_error(__FILE__,__LINE__, &
300                  (/"Error in deallocatting CRTM Atmosphere_AD Structure"/))
301             end if
302 
303             Error_Status = CRTM_Destroy_Surface(Surface_AD)
304             if ( Error_Status /= 0 ) then
305               call da_error(__FILE__,__LINE__, &
306                  (/"Error in deallocatting CRTM Surface_AD Structure"/))
307             end if
308 
309          end do       !  end loop for pixels
310 
311          call da_interp_lin_2d_adj_partial(jo_grad_x%t(:,:,kts:kte),    iv%instid(inst)%info, kts,kte, t)
312          call da_interp_lin_2d_adj_partial(jo_grad_x%q(:,:,kts:kte),    iv%instid(inst)%info, kts,kte, q_ad)
313          call da_interp_lin_2d_adj_partial(jo_grad_x%psfc, iv%instid(inst)%info, 1,1, p)
314 
315          deallocate (q_ad)
316          deallocate (t)
317          deallocate (p)
318 
319       !-------------------------------------------------------------------
320       ! [2.0] Deallocating CRTM structures
321       !-------------------------------------------------------------------
322          deallocate( RTSolution, RTSolution_AD, STAT = Allocate_Status )
323          if ( Allocate_Status /= 0 ) then
324             call da_error(__FILE__,__LINE__, &
325               (/"Error in deallocatting RTSolution"/))
326          END IF
327                                                                                                                           
328          Error_Status = CRTM_Destroy_Surface(Surface)
329          if ( Error_Status /= 0 ) then
330             call da_error(__FILE__,__LINE__, &
331                (/"Error in deallocatting CRTM Surface Structure"/))
332          end if
333    end do        ! end loop for sensor
334 
335       !-------------------------------------------------------------------
336       ! [3.0] Deallocating CRTM Atmosphere structures
337       !-------------------------------------------------------------------
338     Error_Status = CRTM_Destroy_Atmosphere( Atmosphere )
339     if ( Error_Status /= 0 ) then
340        call da_error(__FILE__,__LINE__, &
341          (/"Error in deallocatting CRTM Atmosphere Structure"/))
342     end if
343        
344    if (trace_use) call da_trace_exit("da_transform_xtoy_crtm_adj")
345 #else
346     call da_error(__FILE__,__LINE__, &
347        (/"Must compile with $CRTM option for radiances"/))
348 #endif
349  
350 end subroutine da_transform_xtoy_crtm_adj
351