updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_varbc_tamdar / da_varbc_tamdar_update.inc
blob0b078f467cbf5548dee76c811712c65a11c0aa29
1    subroutine da_varbc_tamdar_update (cv_size, cv, iv)
3    !-----------------------------------!
4    !  Update and write out parameters  !
5    !-----------------------------------!
7    implicit none
9    integer,        intent(in)      :: cv_size
10    real,           intent(in)      :: cv(1:cv_size)
11    type (iv_type), intent(inout)   :: iv   
13    integer                         :: i,ipred,iflt,iphase,npred
14    integer                         :: iunit,iost
15    character(len=filename_len)     :: filename
16    character(len=99)               :: fmt_param
17    character(len=3)                :: cphz(3)
18    character(len=30)               :: stringn
20    character(len=5), allocatable   :: tail_id(:)
21    real, allocatable               :: sum_tl(:),varbc_param_tl(:,:,:),param(:,:,:)
22    integer, allocatable            :: ifuse(:,:)
25    if (trace_use) call da_trace_entry("da_varbc_tamdar_update")
27    if (rootproc) write(unit=varbc_tamdar_unit,fmt='(//5X,A)') 'Updating parameters'
28       
29    npred = iv%varbc_tamdar%npred 
30    cphz = (/'des','asc','cru'/)
32    write(stringn,'(I3)') iv%varbc_tamdar%npred 
33    stringn = '(9X,I5,2X,A,2X,A,'//trim(ADJUSTL(stringn))//'f9.3)'
34    stringn = trim(adjustl(stringn))
36    allocate( sum_tl(npred) )
37    allocate( varbc_param_tl(npred,iv%varbc_tamdar%nphase,iv%varbc_tamdar%nair) )
39    write(unit=varbc_tamdar_unit,fmt='(/10X,A)') &
40          " ID  Phase  Parameters updated (npred)"
42    varbc_param_tl(:,:,:) = 0.0
44    do iflt = 1, iv%varbc_tamdar%nair
45       do iphase = 1, iv%varbc_tamdar%nphase 
46          if (iv%varbc_tamdar%nobs_sum(iphase,iflt) >= varbc_tamdar_nobsmin) then
48              if (iv%varbc_tamdar%nobs(iphase,iflt) > 0 .and. iv%varbc_tamdar%ifuse(iphase,iflt) > 0) then
49                  do ipred = 1, npred
50                     varbc_param_tl(ipred,iphase,iflt) = &
51                               SUM( cv(iv%varbc_tamdar%index(1:npred,iphase,iflt)) * &
52                                    iv%varbc_tamdar%vtox(ipred,1:npred,iphase,iflt) )
53                  end do
54              end if
56              do ipred = 1, npred
57                 sum_tl(ipred) = wrf_dm_sum_real( varbc_param_tl(ipred,iphase,iflt) )
58                 iv%varbc_tamdar%param(ipred,iphase,iflt) = iv%varbc_tamdar%param(ipred,iphase,iflt) + sum_tl(ipred)
59              end do
61              if (rootproc .and. iv%varbc_tamdar%ifuse(iphase,iflt) > 0) then
62                  write(unit=varbc_tamdar_unit,fmt=stringn) &
63                        iv%varbc_tamdar%tail_id(iflt),cphz(iphase),' parameters updated : ',sum_tl(:)
64              end if
65          end if
66       end do
67    end do
69    if (.not. rootproc) then
70        deallocate(varbc_param_tl)
71        if (trace_use) call da_trace_exit("da_varbc_tamdar_update")
72        return
73    end if
75    !-----------------------------------!
76    !  Write VARBC_TAMDAR.tbl.out file  !
77    !-----------------------------------!
78    
79    allocate( tail_id (iv%varbc_tamdar%nair) )
80    allocate( ifuse   (iv%varbc_tamdar%nphase,iv%varbc_tamdar%nair) )
81    allocate( param   (iv%varbc_tamdar%nmaxpred,iv%varbc_tamdar%nphase,iv%varbc_tamdar%nair) )
83    call da_get_unit(iunit)
84    open(unit=iunit,file='VARBC_TAMDAR.tbl', form='formatted', status='old')
86    do i = 1, 13
87       read(iunit, *)
88    end do
89    do iflt = 1, iv%varbc_tamdar%nair
90       read(iunit, *) tail_id(iflt), & 
91                      ifuse(1:iv%varbc_tamdar%nphase,iflt), &
92                    ( param(1:iv%varbc_tamdar%nmaxpred,iphase,iflt), iphase=1,iv%varbc_tamdar%nphase )
93    end do
95    close(iunit)
96    call da_free_unit(iunit)
99  ! write updated parameters
101    call da_get_unit(iunit)
102    filename = 'VARBC_TAMDAR.tbl.out'
104    open(unit=iunit,file=trim(adjustl(filename)),iostat=iost,status='replace')
105    if (iost /= 0) then
106       message(1)="Cannot open TAMDAR bias correction file "//adjustl(filename)
107       call da_error(__FILE__,__LINE__,message(1:1))
108    end if
110    write(iunit, '(A)')'*'
111    write(iunit, '(A)')'  PARAMETER TABLE FOR TAMDAR VARBC'
112    write(iunit, '(A)')' '
113    write(iunit, '(A)')'  Table format:'
114    write(iunit, '(A)')'- ID (1X,A5)'
115    write(iunit, '(A)')'- Ifuse(ascent/desent/cruise) 1X,3(1X,I1)'
116    write(iunit, '(A)')'- Parameters 3(5F7.3)'
117    write(iunit, '(A)')' '
118    write(iunit, '(A)')'  Preditors(1:5): 1.0, w, (optional: Mach, dT/dt, T)'
119    write(iunit, '(A)')' '
120    write(iunit, '(A)')'  Number of aircrafts currently involved in this table'
121    write(iunit, '(I5)') iv%varbc_tamdar%nair
122    write(iunit, '(A)')'*'
124    do iflt = 1, iv%varbc_tamdar%nair
125       do iphase = 1, iv%varbc_tamdar%nphase
126          param(1:npred,iphase,iflt) = iv%varbc_tamdar%param(1:npred,iphase,iflt)
127       end do
129       write(iunit,fmt=trim(adjustl(iv%varbc_tamdar%fmt_param))) &
130             tail_id(iflt), ifuse(iflt,1:iv%varbc_tamdar%nphase), &
131           ( param(1:iv%varbc_tamdar%nmaxpred,iphase,iflt), iphase=1,iv%varbc_tamdar%nphase )
132    end do
134    close(iunit)
135    call da_free_unit(iunit)
137    deallocate(sum_tl, tail_id, ifuse, param, varbc_param_tl)
139    if (rootproc) write(unit=varbc_tamdar_unit,fmt='(//5X,A/)') 'VARBC_TAMDAR is done'
141    if (trace_use) call da_trace_exit("da_varbc_tamdar_update")
143    end subroutine da_varbc_tamdar_update