1 subroutine da_allocate_background_errors (jy, kz, l, e, be_eval_loc, &
4 !---------------------------------------------------------------------------
5 ! Purpose: Allocate components of wrfvar background errors.
6 !---------------------------------------------------------------------------
10 integer, intent(in) :: jy ! i/y dimension.
11 integer, intent(in) :: kz ! k/z dimension.
12 real*8, intent(in) :: l(:) ! Global eigenvalue.
13 real*8, intent(in) :: e(:,:) ! Global eigenvectors.
14 real*8, intent(in) :: be_eval_loc(:,:) ! Std dev/local evalue.
15 real*8, intent(in) :: be_evec_loc(:,:,:) ! Local eigenvectors.
16 type (be_subtype), intent(inout) :: be_sub ! Backgrd error struct.
18 integer :: mz ! Vertical truncation.
19 integer :: j, m, k ! Loop counter.
21 if (trace_use_dull) call da_trace_entry("da_allocate_background_errors")
23 !--------------------------------------------------------------------------
25 !--------------------------------------------------------------------------
29 !--------------------------------------------------------------------------
30 ! [2.0] Allocate components of be_sub structure:
31 !--------------------------------------------------------------------------
34 allocate (be_sub % val(1:jy,1:mz))
36 if (vert_corr == vert_corr_2) then
38 !--------------------------------------------------------------------
39 ! [3.0] Allocate eigenvalues of vertical error covariance matrix:
40 !--------------------------------------------------------------------
42 if (vert_evalue == vert_evalue_global) then
43 ! use global eigenvalues:
45 be_sub % val(1:jy,m) = sqrt (l(m))
48 ! use eigenvalues varying with j-direction.
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))
58 be_sub % val(1:jy,1:mz) = sqrt (be_eval_loc(1:jy,1:mz))
61 if (print_detail_be) then
62 write(unit=message(1),fmt='(A,A)') 'j*k Eigenvalues for ', be_sub % name
63 call da_array_print(2, be_sub % val(1:jy,1:mz), message(1))
66 !-----------------------------------------------------------------------
67 ! [4.0] Allocate global eigenvectors of vertical error cov.:
68 !-----------------------------------------------------------------------
70 allocate (be_sub % evec(1:jy,1:kz,1:mz))
72 if (vert_evalue == vert_evalue_global) then
73 ! use global eigenvectors:
75 be_sub % evec(j,1:kz,1:mz) = e(1:kz,1:mz)
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)
82 if (print_detail_be) then
83 write(unit=message(1),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), message(1))
86 write(unit=message(1),fmt='(A,A)') 'k*k Eigenvectors for j = jy ', be_sub % name
87 call da_array_print (2, be_sub % evec(jy,1:kz,1:mz), message(1))
90 allocate (be_sub%val_g(1:mz))
91 allocate (be_sub%evec_g(1:kz,1:mz))
93 be_sub % val_g(1:mz) = l(1:mz)
94 be_sub % evec_g(1:kz,1:mz) = e(1:kz,1:mz)
95 else if (vert_corr == vert_corr_1) then
96 if (print_detail_be) then
97 write(unit=message(1),fmt='(A)') 'Change BE std dev to variance in NMC code'
98 call da_message(message(1:1))
100 if (vert_evalue == vert_evalue_global) then
101 ! use global eigenvalues:
103 be_sub % val(1:jy,m) = l(m)
106 be_sub % val(1:jy,1:mz) = be_eval_loc(1:jy,1:mz)
110 !-----------------------------------------------------------------------
111 ! [2.2] Allocate recursive filter lengthscales and variance scaling factors:
112 !-----------------------------------------------------------------------
114 allocate (be_sub % rf_alpha(1:mz))
116 be_sub % rf_alpha(1:mz) = 1.0
119 if (trace_use_dull) call da_trace_exit("da_allocate_background_errors")
121 end subroutine da_allocate_background_errors