da_biascorr.inc

References to this file elsewhere.
1 subroutine da_biascorr ( i, ob, iv)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: perform bias correction for radiance data.
5    !
6    ! METHOD:  omb(corrected)=omb-scanbias-airbias
7    !---------------------------------------------------------------------------
8 
9    implicit none
10 
11    integer,        intent(in)    :: i       ! sensor index.
12    type (y_type),  intent(in)    :: ob      ! Observation structure.
13    type (ob_type), intent(inout) :: iv      ! O-B structure.
14 
15    ! Local variables
16    integer   :: k,iband,iscan, n1,n2,n,j,npred,nlevels
17    real      :: q(43),pred(4),airbias
18 
19    if (trace_use) call da_trace_entry("da_biascorr")
20 
21    npred=4
22    nlevels=iv%instid(i)%nlevels-1
23 
24    n1 = iv%ob_numb(iv%current_ob_time-1)%radiance(i)+1
25    n2 = iv%ob_numb(iv%current_ob_time)%radiance(i)
26 
27    do n=n1,n2
28       ! get airmass predictors
29       !-------------------------
30         if (rtm_option==rtm_option_rttov) then
31          q(1:43) = iv%instid(i)%mr(1:43,n)/q2ppmv
32          call da_predictor_rttov(pred(1:npred), npred, iv%instid(i)%t(1:43,n), &
33             q(1:43), iv%instid(i)%ts(n))
34         else if (rtm_option==rtm_option_crtm) then
35 #ifdef CRTM
36          call da_predictor_crtm(pred(1:npred), npred, nlevels,iv%instid(i)%tm(1:nlevels,n), &
37             iv%instid(i)%qm(1:nlevels,n), iv%instid(i)%ts(n), &
38             iv%instid(i)%pm(1:nlevels,n),iv%instid(i)%pf(0:nlevels,n))
39 #endif
40         end if
41         iscan = iv%instid(i)%scanpos(n)
42         iband = floor(iv%instid(i)%info(n)%lat/10.0001) + 10
43       do k=1,iv%instid(i)%nchan
44        ! scan bias correction
45        !-----------------------
46         if (global) then
47           iv%instid(i)%tb_inv(k,n) = iv%instid(i)%tb_inv(k,n) - satinfo(i)%scanbias_b(k,iscan,iband)
48         else
49           iv%instid(i)%tb_inv(k,n) = iv%instid(i)%tb_inv(k,n) - satinfo(i)%scanbias(k,iscan) 
50         end if
51        ! airmass bias correction
52        !----------------------------
53           airbias = satinfo(i)%bcoef0(k)
54          do j=1,npred
55           airbias= airbias + satinfo(i)%bcoef(k,j)*pred(j)
56          end do
57           iv%instid(i)%tb_inv(k,n) = iv%instid(i)%tb_inv(k,n)-airbias
58         !if (k==5) then
59         !  write(*,'(6f15.5)') satinfo(i)%scanbias(k,iscan), satinfo(i)%bcoef(k,1:npred),satinfo(i)%bcoef0(k)
60         !  write(*,'(5f15.5)') airbias,pred(1:npred)
61         !end if
62       end do
63    end do
64 
65    if (trace_use) call da_trace_exit("da_biascorr")
66 
67 end subroutine da_biascorr
68