da_biasprep.inc

References to this file elsewhere.
1 subroutine da_biasprep(inst,ob,iv)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: Output information files for bias correction progs
5    !-----------------------------------------------------------------------
6 
7    implicit none
8 
9    integer       ,  intent(in)      :: inst
10    type (y_type) ,  intent(in)      :: ob         ! O structure.
11    type (ob_type),  intent(in)      :: iv         ! O-B structure.
12 
13    integer  :: n,jx,npred,nchan,num_rad,nlevels
14    character(len=80)  :: filename
15    character(len=1)   :: s1
16    real               :: pred(6), q(43)
17    type (bias_type)   :: radbias
18 
19    num_rad         = iv%ob_numb(iv%current_ob_time)%radiance(inst)- &
20                       iv%ob_numb(iv%current_ob_time-1)%radiance(inst)
21 
22    if (num_rad < 1) return
23 
24    if (trace_use) call da_trace_entry("da_biasprep")
25 
26 #ifdef DM_PARALLEL
27    write(filename, '(a,i2.2)') 'biasprep_'//trim(iv%instid(inst)%rttovid_string)//'.', myproc
28 #else
29    write(filename, '(a)') 'biasprep_'//trim(iv%instid(inst)%rttovid_string)
30 #endif
31 
32    call da_get_unit(biasprep_unit)
33    open(unit=biasprep_unit,FILE=filename,FORM='unformatted')
34 
35    !---------------------------------------------------------------------------
36    npred = 4
37    nchan = iv%instid(inst)%nchan 
38    nlevels = iv%instid(inst)%nlevels-1
39 
40    allocate(radbias%tb(nchan))
41    allocate(radbias%omb(nchan))
42    allocate(radbias%bias(nchan))
43    allocate(radbias%qc_flag(nchan))
44    allocate(radbias%cloud_flag(nchan))
45    allocate(radbias%pred(npred))
46 
47    do n= iv%ob_numb(iv%current_ob_time-1)%radiance(inst)+1, &
48       iv%ob_numb(iv%current_ob_time)%radiance(inst)
49       if (iv%instid(inst)%proc_domain(n)) then 
50 
51         if (rtm_option==rtm_option_rttov) then
52          q(1:43) = iv%instid(inst)%mr(1:43,n)/q2ppmv
53          call da_predictor_rttov(pred(1:npred), npred, iv%instid(inst)%t(1:43,n), &
54             q(1:43), iv%instid(inst)%ts(n))
55         else if (rtm_option==rtm_option_crtm) then
56 #ifdef CRTM
57          call da_predictor_crtm(pred(1:npred), npred, nlevels,iv%instid(inst)%tm(1:nlevels,n), &
58             iv%instid(inst)%qm(1:nlevels,n), iv%instid(inst)%ts(n), &
59             iv%instid(inst)%pm(1:nlevels,n),iv%instid(inst)%pf(0:nlevels,n))
60 #endif
61         end if
62 
63          ! transfer information to bias structure
64          radbias%platform_id  = iv%instid(inst)%platform_id
65          radbias%satellite_id = iv%instid(inst)%satellite_id
66          radbias%sensor_id    = iv%instid(inst)%sensor_id
67 
68          read(iv%instid(inst)%info(n)%date_char,'(i4,5(a1,i2))') &
69                                    radbias%year,s1, radbias%month,s1, radbias%day, &
70                                    s1,radbias%hour, s1,radbias%min, s1,radbias%sec
71 
72          radbias%scanline     = iv%instid(inst)%scanline(n)    ! not available
73          radbias%scanpos      = iv%instid(inst)%scanpos(n)
74          radbias%landmask     = iv%instid(inst)%landsea_mask(n)
75          radbias%elevation    = iv%instid(inst)%info(n)%elv
76          radbias%lat          = iv%instid(inst)%info(n)%lat
77          radbias%lon          = iv%instid(inst)%info(n)%lon
78          radbias%surf_flag    = iv%instid(inst)%isflg(n)
79          radbias%ps           = iv%instid(inst)%ps(n)
80          radbias%t2m          = iv%instid(inst)%t2m(n)
81          radbias%q2m          = iv%instid(inst)%mr2m(n)/q2ppmv
82          radbias%tsk          = iv%instid(inst)%ts(n)
83          radbias%clwp         = iv%instid(inst)%clwp(n)  ! in mm
84 
85          radbias%nchan        = nchan 
86          radbias%tb(1:nchan)  = ob%instid(inst)%tb(1:nchan,n)
87          radbias%omb(1:nchan) = iv%instid(inst)%tb_inv(1:nchan,n)
88          radbias%bias(1:nchan) = 0.
89 
90          radbias%npred         = npred
91          radbias%pred(1:npred) = pred(1:npred)
92 
93          radbias%qc_flag(1:nchan)= iv%instid(inst)%tb_qc(1:nchan,n)
94          radbias%cloud_flag(1:nchan)= iv%instid(inst)%cloud_flag(1:nchan,n)
95 
96          ! set missing data and bad data to missing
97          do jx=1,nchan   
98             if (radbias%tb(jx) < 150. .or. radbias%tb(jx) > 400. ) then
99                radbias%tb(jx)   = missing_r
100                radbias%omb(jx)  = missing_r 
101             end if
102          end do
103 
104          !write(unit=biasprep_unit) radbias ! can not compiled with pointer
105 
106          call da_write_biasprep(radbias)
107 
108       end if
109    end do
110   
111    close(unit=biasprep_unit)
112    call da_free_unit(biasprep_unit)
113 
114    deallocate(radbias%tb)
115    deallocate(radbias%omb)
116    deallocate(radbias%bias)
117    deallocate(radbias%qc_flag)
118    deallocate(radbias%cloud_flag)
119    deallocate(radbias%pred)
120 
121    if (trace_use) call da_trace_exit("da_biasprep")
122 
123 end subroutine da_biasprep
124 
125