updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_varbc_tamdar / da_varbc_tamdar_init.inc
blob733c468440b30a1572fa63c8f008a4da3bcd389e
1    subroutine da_varbc_tamdar_init (iv)
3    !-----------------------------------------------------!
4    !  Initialize Varbc for TAMDAR from VARBC_TAMDAR.tbl  !
5    !-----------------------------------------------------!
7    implicit none
9    type (iv_type), intent(inout)   :: iv
11    integer, parameter              :: nmaxpred = 5  ! 1,w,(optional: Mach,dT/dt,T)
12    integer, parameter              :: nphase = 3    ! descent/ascent/cruise
14    integer                         :: npred,nair
15    integer                         :: iunit,i,id,iflt,iobs,iphase
17    logical                         :: ifexist_table
18    character(len=120)              :: filename
19    character(len=5), allocatable   :: tail_id(:)
22    if (trace_use) call da_trace_entry("da_varbc_tamdar_init")
23    
24    if (rootproc) then
25        write(unit=varbc_tamdar_unit,fmt='(//5X,A/)')      'VARBC_TAMDAR namelist options'
26        write(unit=varbc_tamdar_unit,fmt='(10X,A,7X,L)')   'use_varbc_tamdar :         ', use_varbc_tamdar
27        write(unit=varbc_tamdar_unit,fmt='(10X,A,7X,I1)')  'varbc_tamdar_bm :          ', varbc_tamdar_bm       ! bias models
28        write(unit=varbc_tamdar_unit,fmt='(10X,A,F8.3)')   'varbc_tamdar_pred0 :       ', varbc_tamdar_pred0    ! predictor(1) 
29        write(unit=varbc_tamdar_unit,fmt='(10X,A,I8)')     'varbc_tamdar_nbgerr :      ', varbc_tamdar_nbgerr
30        write(unit=varbc_tamdar_unit,fmt='(10X,A,I8)')     'varbc_tamdar_nobsmin :     ', varbc_tamdar_nobsmin
31    end if
33    filename = 'VARBC_TAMDAR.tbl'
34    inquire(file=trim(adjustl(filename)), exist=ifexist_table)
36    if (ifexist_table) then
38        if (rootproc) write(unit=varbc_tamdar_unit,fmt='(//5X,A/)') &
39           'Reading in VARBC_TAMDAR.tbl and initializating'
41        iv%varbc_tamdar%nmaxpred = nmaxpred
42        iv%varbc_tamdar%nphase   = nphase
44        call da_get_unit(iunit)
45        open(unit=iunit,file=filename, status='old')
46        do i = 1, 11; read(iunit,*); end do
47        read(iunit,*) nair
48        read(iunit,*)
50        iv%varbc_tamdar%nair = nair
52        npred = 2 
53        if (varbc_tamdar_bm == 2) npred = 5
55        iv%varbc_tamdar%npred = npred
57        allocate ( iv%varbc_tamdar%tail_id  (nair) )
58        allocate ( iv%varbc_tamdar%ifuse    (nphase,nair) )
59        allocate ( iv%varbc_tamdar%nobs     (nphase+1,nair) )
60        allocate ( iv%varbc_tamdar%nobs_sum (nphase+1,nair) )
61        allocate ( iv%varbc_tamdar%param    (nmaxpred,nphase,nair) )
62        allocate ( iv%varbc_tamdar%pred     (nmaxpred,nphase,nair) )
63        allocate ( iv%varbc_tamdar%bgerr    (npred,nphase,nair) )
64        allocate ( iv%varbc_tamdar%index    (npred,nphase,nair) )
65        allocate ( iv%varbc_tamdar%vtox     (npred,npred,nphase,nair) )
67        iv%varbc_tamdar%tail_id(:)    = 0
68        iv%varbc_tamdar%ifuse(:,:)    = 0
69        iv%varbc_tamdar%nobs(:,:)     = 0
70        iv%varbc_tamdar%nobs_sum(:,:) = 0
71        iv%varbc_tamdar%param(:,:,:)  = 0.
72        iv%varbc_tamdar%pred(:,:,:)   = 0.
73        iv%varbc_tamdar%pred(1,:,:)   = varbc_tamdar_pred0
74        iv%varbc_tamdar%bgerr(:,:,:)  = 0.
75        iv%varbc_tamdar%index(:,:,:)  = 0
76        iv%varbc_tamdar%vtox(:,:,:,:) = 0.
79        write(iv%varbc_tamdar%fmt_param,*) '(1X,A5,1X,3(1X,I1),15F7.3))'
81        allocate( tail_id (nair) )
82        tail_id = ''
83        do iflt = 1, nair
84           read(iunit, fmt=trim(adjustl(iv%varbc_tamdar%fmt_param))) &
85                tail_id(iflt), iv%varbc_tamdar%ifuse(1:nphase,iflt),     &
86               (iv%varbc_tamdar%param(1:nmaxpred,iphase,iflt), iphase=1,nphase) 
88           read(tail_id(iflt),'(I5)') iv%varbc_tamdar%tail_id(iflt)
90           do iobs = 1, iv%info(tamdar)%nlocal
91              read(iv%info(tamdar)%id(iobs),'(I5)') id
92              if (iv%varbc_tamdar%tail_id(iflt) .eq. id) &
93                  iv%varbc_tamdar%nobs(nphase+1,iflt) = iv%varbc_tamdar%nobs(nphase+1,iflt) + 1
94           end do
96           iv%varbc_tamdar%nobs_sum(nphase+1,iflt) = wrf_dm_sum_integer(iv%varbc_tamdar%nobs(nphase+1,iflt))
97        end do
99        iv%varbc_tamdar%nmaxobs = MAXVAL( iv%varbc_tamdar%nobs_sum(nphase+1,:) )
101        allocate ( iv%varbc_tamdar%obs_sn(iv%varbc_tamdar%nmaxobs,nphase+1,nair) )
102        iv%varbc_tamdar%obs_sn(:,:,:) = 0
104        do iflt = 1, nair
105           i = 0
106           do iobs = 1, iv%info(tamdar)%nlocal
107              read(iv%info(tamdar)%id(iobs),'(I5)') id
108              if (iv%varbc_tamdar%tail_id(iflt) .eq. id) then
109                  i = i + 1
110                  iv%varbc_tamdar%obs_sn(i,nphase+1,iflt) = iobs
111              end if
112           end do
113        end do
115        if (rootproc) then
116            write(unit=varbc_tamdar_unit,fmt='(10X,A,I4)') &
117                     'Number of aircrafts in table : ', nair
118            write(unit=varbc_tamdar_unit,fmt='(10X,A,I4)') &
119                     'Max. Obs. number in aircraft : ', iv%varbc_tamdar%nmaxobs
120        end if
122        close(iunit)
123        call da_free_unit(iunit)
125        deallocate (tail_id)
127    else
129        use_varbc_tamdar = .false.
130        if (rootproc) write(unit=varbc_tamdar_unit,fmt='(/5X,A/)') &
131                     'Could not find VARBC_TAMDAR.tbl file. VARBC_TAMDAR is switched off.'
133    end if
135    if (trace_use) call da_trace_exit("da_varbc_tamdar_init")
137    end subroutine da_varbc_tamdar_init