Update version info for release v4.6.1 (#2122)
[WRF.git] / var / da / da_radiance / da_biascorr.inc
blob05c5761836ece97565e46dd49331a0b986416147
1 subroutine da_biascorr ( i, ob, iv)
3    !---------------------------------------------------------------------------
4    ! Purpose: perform bias correction for radiance data.
5    !
6    ! METHOD:  omb(corrected)=omb-scanbias-airbias
7    !---------------------------------------------------------------------------
9    implicit none
11    integer,        intent(in)    :: i       ! sensor index.
12    type (y_type),  intent(in)    :: ob      ! Observation structure.
13    type (iv_type), intent(inout) :: iv      ! O-B structure.
15    ! Local variables
16    integer   :: k,iband,iscan, n,j,npred,nlevels, num_rad
17    real      :: pred(4),airbias
18    real,allocatable      :: q(:), temp(:), hum(:), pf(:)
20    num_rad = iv%instid(i)%info%n2-iv%instid(i)%info%n1+1
21    if (num_rad < 1) return
23    if (trace_use) call da_trace_entry("da_biascorr")
25    npred=4
26    nlevels=iv%instid(i)%nlevels-1
27    allocate(temp(nlevels))
28    allocate(hum(nlevels))
29    allocate(pf(0:nlevels))
31    do n=iv%instid(i)%info%n1,iv%instid(i)%info%n2
32       ! get airmass predictors
33       !-------------------------
34       if (rtm_option==rtm_option_rttov) then
35 #ifdef RTTOV
36          q(1:nlevels) = iv%instid(i)%mr(1:nlevels,n)/q2ppmv
37          call da_predictor_rttov(i, pred(1:npred), npred, nlevels, &
38               iv%instid(i)%t(1:nlevels,n), q(1:nlevels), iv%instid(i)%ts(n))
39 #endif
40       else if (rtm_option==rtm_option_crtm) then
41 #ifdef CRTM
42 ! FIX? problems with IBM AIX COMPILER
43          temp(1:nlevels) = iv%instid(i)%tm(1:nlevels,n)
44          hum(1:nlevels) = iv%instid(i)%qm(1:nlevels,n)
45          pf(0:nlevels) = iv%instid(i)%pf(0:nlevels,n)
46          call da_predictor_crtm(pred(1:npred), npred, nlevels,temp, &
47             hum, iv%instid(i)%ts(n), pf)
48 #endif
49       end if
50       iscan = iv%instid(i)%scanpos(n)
51       iband = floor(iv%instid(i)%info%lat(1,n)/10.0001) + 10
52       do k=1,iv%instid(i)%nchan
53          ! only do bias correction for used channels
54          if ( (satinfo(i)%iuse(k) == 1) .and. (iv%instid(i)%tb_inv(k,n) > missing_r) ) then
55             ! scan bias correction
56             !-----------------------
57             if (global) then
58                iv%instid(i)%tb_inv(k,n) = iv%instid(i)%tb_inv(k,n) - satinfo(i)%scanbias_b(k,iscan,iband)
59             else
60                iv%instid(i)%tb_inv(k,n) = iv%instid(i)%tb_inv(k,n) - satinfo(i)%scanbias(k,iscan) 
61             end if
62             ! airmass bias correction
63             !----------------------------
64             airbias = satinfo(i)%bcoef0(k)
65             do j=1,npred
66                airbias= airbias + satinfo(i)%bcoef(k,j)*pred(j)
67             end do
68             iv%instid(i)%tb_inv(k,n) = iv%instid(i)%tb_inv(k,n)-airbias
69          end if
70       end do
71    end do
73    deallocate(q)
74    deallocate(temp)
75    deallocate(hum)
76    deallocate(pf)
78    if (trace_use) call da_trace_exit("da_biascorr")
80 end subroutine da_biascorr