Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_obs / da_transform_xtoy_adj.inc
blobdbbe9ddd151368d6e0d9a70d9a2442e22c4c3abe
1 subroutine da_transform_xtoy_adj(cv_size, cv, grid, iv, jo_grad_y, jo_grad_x &
2 #if (WRF_CHEM == 1)
3            , jo_grad_xchem &
4 #endif
7    !--------------------------------------------------------------------------
8    ! Purpose: TBD
9    !    Updated for Analysis on Arakawa-C grid
10    !    Author: Syed RH Rizvi,  MMM/ESSL/NCAR,  Date: 10/22/2008
11    !--------------------------------------------------------------------------
12    
13    implicit none
14    
15    integer, intent(in)           :: cv_size         ! Size of cv array.
16    real, intent(inout)           :: cv(1:cv_size)   ! control variables.
17    type (domain),  intent(inout) :: grid
18    type (iv_type), intent(inout) :: iv          ! obs. inc vector (o-b).
19    type (y_type),  intent(inout) :: jo_grad_y   ! grad_y(jo)
20    type (x_type),  intent(inout) :: jo_grad_x   ! grad_x(jo)
22 #if (WRF_CHEM == 1)
23    type (xchem_type),  optional, intent(inout) :: jo_grad_xchem   ! grad_x(jo)
24 #endif
26    integer :: i,j,k
27    real, dimension(:), allocatable :: adj_ref
28    real, dimension(:), allocatable :: adj_ref_tot
30    if (trace_use) call da_trace_entry("da_transform_xtoy_adj")
31   
32    !--------------------------------------------------------------------------
33    ! [1.0] observation operator y = H(x):
34    !--------------------------------------------------------------------------
35   
36    if (iv%info(sound)%nlocal    > 0) call da_transform_xtoy_sound_adj    (grid, iv, jo_grad_y, jo_grad_x)
37    if (iv%info(sonde_sfc)%nlocal  > 0) call da_transform_xtoy_sonde_sfc_adj(grid, iv, jo_grad_y, jo_grad_x)
38    if (iv%info(mtgirs)%nlocal   > 0) call da_transform_xtoy_mtgirs_adj   (grid, iv, jo_grad_y, jo_grad_x)
39    if (iv%info(tamdar)%nlocal   > 0) call da_transform_xtoy_tamdar_adj   (grid, iv, jo_grad_y, jo_grad_x)
40    if (iv%info(tamdar_sfc)%nlocal > 0) call da_transform_xtoy_tamdar_sfc_adj(grid, iv, jo_grad_y, jo_grad_x)
41    if (iv%info(synop)%nlocal    > 0) call da_transform_xtoy_synop_adj    (grid, iv, jo_grad_y, jo_grad_x)
42    if (iv%info(geoamv)%nlocal   > 0) call da_transform_xtoy_geoamv_adj   (grid, iv, jo_grad_y, jo_grad_x)
43    if (iv%info(polaramv)%nlocal > 0) call da_transform_xtoy_polaramv_adj (grid, iv, jo_grad_y, jo_grad_x)   
44    if (iv%info(airep)%nlocal    > 0) call da_transform_xtoy_airep_adj    (grid, iv, jo_grad_y, jo_grad_x)
45    if (iv%info(metar)%nlocal    > 0) call da_transform_xtoy_metar_adj    (grid, iv, jo_grad_y, jo_grad_x)
46    if (iv%info(ships)%nlocal    > 0) call da_transform_xtoy_ships_adj    (grid, iv, jo_grad_y, jo_grad_x)
47    if (iv%info(gpspw)%nlocal    > 0) then
48       if (use_gpspwobs) then
49          call da_transform_xtoy_gpspw_adj    (grid, iv, jo_grad_y, jo_grad_x)
50       else if (use_gpsztdobs) then
51          call da_transform_xtoy_gpsztd_adj   (grid, iv, jo_grad_y, jo_grad_x)
52       endif
53    end if
54    if (iv%info(ssmi_tb)%nlocal  > 0) call da_transform_xtoy_ssmi_tb_adj  (grid, iv, jo_grad_y, jo_grad_x)
55    if (iv%info(ssmi_rv)%nlocal  > 0) call da_transform_xtoy_ssmi_rv_adj  (grid, iv, jo_grad_y, jo_grad_x)
56    if (iv%info(pilot)%nlocal    > 0) call da_transform_xtoy_pilot_adj    (grid, iv, jo_grad_y, jo_grad_x)
57    if (iv%info(satem)%nlocal    > 0) call da_transform_xtoy_satem_adj    (grid, iv, jo_grad_y, jo_grad_x)
58    if (iv%info(ssmt1)%nlocal    > 0) call da_transform_xtoy_ssmt1_adj    (iv, jo_grad_y, jo_grad_x)
59    if (iv%info(ssmt2)%nlocal    > 0) call da_transform_xtoy_ssmt2_adj    (iv, jo_grad_y, jo_grad_x)
60    if (iv%info(qscat)%nlocal    > 0) call da_transform_xtoy_qscat_adj    (grid, iv, jo_grad_y, jo_grad_x)
61    if (iv%info(profiler)%nlocal > 0) call da_transform_xtoy_profiler_adj (grid, iv, jo_grad_y, jo_grad_x)
62    if (iv%info(buoy)%nlocal     > 0) call da_transform_xtoy_buoy_adj     (grid, iv, jo_grad_y, jo_grad_x)
63    if (iv%info(gpsref)%nlocal   > 0) call da_transform_xtoy_gpsref_adj   (iv, jo_grad_y, jo_grad_x)
64    if (iv%info(gpseph)%nlocal   > 0) call da_transform_xtoy_gpseph_adj   (iv, jo_grad_y, jo_grad_x)
66 #if (WRF_CHEM == 1)
67    if( present(jo_grad_xchem) ) then
68      if (iv%info(chemic_surf)%nlocal      > 0) then
69          call da_transform_xtoy_chem_sfc_adj (grid, iv, jo_grad_y, jo_grad_xchem)
70      end if
71    end if
72 #endif
74    if ( use_gpsephobs ) then
75 #ifdef DM_PARALLEL
76       allocate ( adj_ref ((ide-ids+1)*(jde-jds+1)*(kde-kds+1)) )
77       allocate ( adj_ref_tot ((ide-ids+1)*(jde-jds+1)*(kde-kds+1)) )
79       adj_ref = 0.0
80       adj_ref_tot = 0.0
82       adj_ref = reshape ( global_adj_ref, (/(ide-ids+1)*(jde-jds+1)*(kde-kds+1)/) )
83       global_adj_ref = 0.0
85       call wrf_dm_sum_reals ( adj_ref, adj_ref_tot)
86       global_adj_ref = reshape ( adj_ref_tot, (/(ide-ids+1),(jde-jds+1),(kde-kds+1)/) )
88       do k = kts, kte
89          do j = jts, jte
90             do i = its, ite
91                jo_grad_x%ref(i,j,k) = jo_grad_x%ref(i,j,k) + global_adj_ref(i,j,k)
92             enddo
93          enddo
94       enddo
95       global_adj_ref = 0.0
97       deallocate (adj_ref)
98       deallocate (adj_ref_tot)
99 #else
100       do k = kts, kte
101          do j = jts, jte
102             do i = its, ite
103                jo_grad_x%ref(i,j,k) = jo_grad_x%ref(i,j,k) + global_adj_ref(i,j,k)
104             enddo
105          enddo
106       enddo
107       global_adj_ref = 0.0
108 #endif
109    endif
110    if (iv%info(radar)%nlocal    > 0) call da_transform_xtoy_radar_adj    (grid, iv, jo_grad_y, jo_grad_x)
111    if (iv%info(bogus)%nlocal    > 0) call da_transform_xtoy_bogus_adj    (grid, iv, jo_grad_y, jo_grad_x)
112    if (iv%info(airsr)%nlocal    > 0) call da_transform_xtoy_airsr_adj    (iv, jo_grad_y, jo_grad_x)
113    if (iv%info(pseudo)%nlocal   > 0) call da_transform_xtoy_pseudo_adj   (iv, jo_grad_y, jo_grad_x)
115 #if defined(CRTM) || defined(RTTOV)
116    if (use_rad) then
117       if (rtm_option == rtm_option_rttov) then
118 #ifdef RTTOV
119          call da_transform_xtoy_rttov_adj (iv, jo_grad_y, jo_grad_x)
120 #endif
121       elseif (rtm_option == rtm_option_crtm) then
122 #ifdef CRTM
123          call da_transform_xtoy_crtm_adj (cv_size, cv, iv, jo_grad_y, jo_grad_x)
124 #endif
125       else
126          call da_warning(__FILE__,__LINE__,(/"Unknown radiative transfer model"/))
127       end if
128    end if
129 #endif
131    if (trace_use) call da_trace_exit("da_transform_xtoy_adj")
133 end subroutine da_transform_xtoy_adj