da_allocate_background_errors.inc

References to this file elsewhere.
1 subroutine da_allocate_background_errors  (jy, kz, l, e, be_eval_loc, &
2                                            be_evec_loc, be_sub)
3 
4    !---------------------------------------------------------------------------
5    ! Purpose: Allocate components of wrfvar background errors.
6    !---------------------------------------------------------------------------
7 
8    implicit none
9 
10    integer, intent(in)              :: jy                 ! i/y dimension.
11    integer, intent(in)              :: kz                 ! k/z dimension.
12    real, intent(in)                 :: l(:)               ! Global eigenvalue.
13    real, intent(in)                 :: e(:,:)             ! Global eigenvectors.
14    real, intent(in)                 :: be_eval_loc(:,:)   ! Std dev/local evalue.
15    real, intent(in)                 :: be_evec_loc(:,:,:) ! Local eigenvectors.
16    type (be_subtype), intent(inout) :: be_sub             ! Backgrd error struct.
17     
18    integer                          :: mz                 ! Vertical truncation.
19    integer                          :: j, m, k            ! Loop counter.
20 
21    if (trace_use) call da_trace_entry("da_allocate_background_errors")
22 
23    !--------------------------------------------------------------------------
24    ! [1.0] Initialise:
25    !--------------------------------------------------------------------------
26 
27    mz = be_sub % mz
28    
29    !--------------------------------------------------------------------------
30    ! [2.0] Allocate components of be_sub structure:
31    !--------------------------------------------------------------------------
32 
33    if (mz > 0) then
34       allocate  (be_sub % val(1:jy,1:mz))
35       
36       if (vert_corr == 2) then
37 
38          !--------------------------------------------------------------------
39          ! [3.0] Allocate eigenvalues of vertical error covariance matrix:
40          !--------------------------------------------------------------------
41 
42          if (vert_evalue == 1) then
43             ! use global eigenvalues:
44             do m = 1, mz
45                be_sub % val(1:jy,m) = sqrt (l(m))
46             end do
47          else
48             ! use eigenvalues varying with j-direction.
49             do j = 1, jy
50                do k = 1, mz
51                   if (be_eval_loc(j,k) <=0) then
52                      write (unit=message(1),fmt='(A,I5,A,I5,A,F10.2)') &
53                         "At lat= ",j," For mode = ",k," local eigen value= ",be_eval_loc(j,k)
54                      call da_error(__FILE__,__LINE__,message(1:1))
55                   end if
56                end do
57             end do
58             be_sub % val(1:jy,1:mz) = sqrt (be_eval_loc(1:jy,1:mz))            
59          end if
60  
61          if (print_detail_be) then
62             write(unit=stdout,fmt='(/A,A)')' j*k Eigenvalues for ', be_sub % name
63             call da_array_print(2, be_sub % val(1:jy,1:mz), 'j*k Eigenvalues')
64          end if
65 
66          !----------------------------------------------------------------------- 
67          ! [4.0] Allocate global eigenvectors of vertical error cov.:
68          !-----------------------------------------------------------------------
69 
70          allocate  (be_sub % evec(1:jy,1:kz,1:mz))
71          
72          if (vert_evalue == 1) then
73             ! use global eigenvectors:
74             do j = 1, jy
75                be_sub % evec(j,1:kz,1:mz) = e(1:kz,1:mz)
76             end do
77          else
78             ! use eigenvectors varying with i-direction.
79             be_sub % evec(1:jy,1:kz,1:mz) =  be_evec_loc(1:jy,1:kz,1:mz)
80          end if
81          
82          if (print_detail_be) then      
83             write(unit=stdout,fmt='(/A,A)')' k*k Eigenvectors for j = 1', be_sub % name
84             call da_array_print (2, be_sub % evec(1,1:kz,1:mz), &
85                                  'Eigenvectors for j=1')
86          
87             write(unit=stdout,fmt='(/A,A)')' k*k Eigenvectors for j = jy', be_sub % name
88             call da_array_print (2, be_sub % evec(jy,1:kz,1:mz), &
89                                  'Eigenvectors for j=jy')
90          end if
91 
92          allocate (be_sub%val_g(1:mz))
93          allocate (be_sub%evec_g(1:kz,1:mz))
94   
95          be_sub % val_g(1:mz) = l(1:mz)
96          be_sub % evec_g(1:kz,1:mz) = e(1:kz,1:mz)
97       else if (vert_corr == 1) then
98          if (print_detail_be) then
99            write(unit=stdout,fmt='(A)') &
100               'Change be std dev to variance in NMC code'
101          end if
102          be_sub % val(1:jy,1:mz) = be_eval_loc(1:jy,1:mz)
103       end if
104 
105       !-----------------------------------------------------------------------
106       ! [2.2] Allocate recursive filter lengthscales and variance scaling factors:
107       !-----------------------------------------------------------------------
108 
109       allocate (be_sub % rf_alpha(1:mz))
110 
111       be_sub % rf_alpha(1:mz) = 1.0    
112    end if
113 
114    if (trace_use) call da_trace_exit("da_allocate_background_errors")
115 
116 end subroutine da_allocate_background_errors
117 
118