da_interpolate_regcoeff.inc

References to this file elsewhere.
1 subroutine da_interpolate_regcoeff(iy, iys, kz, kzs, meanl_stats, &
2                                     meanl_xb, meanp_stats, meanp_xb, &
3                                     pb_vert_reg_stats, pb_vert_reg)
4 
5    !---------------------------------------------------------------------------
6    ! Purpose: Interpolate statistical regression coefficient to new domain.
7    !
8    ! Method:  i,k,k Interpolation.
9    !---------------------------------------------------------------------------
10 
11    implicit none
12 
13    integer, intent(in) :: iy                       ! Number of rows in  xb.
14    integer, intent(in) :: iys                      ! Number of rows in stats.
15    integer, intent(in) :: kz                       ! Number of levels in xb.
16    integer, intent(in) :: kzs                      ! Number of levels in stats.
17    real, intent(in)    :: meanl_stats(:)           ! Mean latitude on stats rows.
18    real, intent(in)    :: meanl_xb(:)              ! Mean latitude on xb rows.
19    real, intent(in)    :: meanp_stats(:)           ! Mean pressure on stats levs.
20    real, intent(in)    :: meanp_xb(:)              ! Mean pressure on xb levs.
21    real, intent(in)    :: pb_vert_reg_stats(:,:,:) ! Coefficient on stats grid.
22    real, intent(out)   :: pb_vert_reg(:,:,:)       ! Coefficient on xb grid.
23      
24    integer             :: i, is, i_south           ! Loop counters.
25    integer             :: k1, k2, k, ks            ! Loop counters.
26    integer             :: k1s, k2s
27    real                :: lat_wgt
28 
29    integer             :: k_above(1:kz)
30    real                :: pb_vert_reg_temp(1:iys,1:kz,1:kz)
31    real                :: weight(1:kz)
32 
33    if (trace_use) call da_trace_entry("da_interpolate_regcoeff")
34 
35    pb_vert_reg = 0.0
36 
37    !------------------------------------------------------------------------
38    ! [1.0] Find xb levels/rows bounded by stats domain:
39    !------------------------------------------------------------------------
40 
41    do k = 1, kz
42       if (meanp_xb(k) <= meanp_stats(1)) then
43          weight(k) = 1.0e-6
44          k_above(k) = 1
45       else if (meanp_xb(k) >= meanp_stats(kzs)) then
46          weight(k) = 1.0-1.0e-6
47          k_above(k) = kzs-1
48       else
49          do ks = 1, kzs-1
50             if (meanp_xb(k) >= meanp_stats(ks) .AND. &
51                  meanp_xb(k) <= meanp_stats(ks+1)) then
52                weight(k) = (meanp_xb(k) - meanp_stats(ks)) / &
53                            (meanp_stats(ks+1) - meanp_stats(ks))
54                k_above(k) = ks
55                exit
56             end if
57          end do
58       end if
59    end do
60 
61    !------------------------------------------------------------------------   
62    ! [3.0] Interpolate regression coefficient from stats to xb levels:
63    !------------------------------------------------------------------------
64 
65    pb_vert_reg_temp(1:iys,1:kz,1:kz) = 0.0
66 
67    do is = 1, iys
68       do k1 = 1, kz
69          k1s = k_above(k1)
70          do k2 = 1, kz
71             k2s = k_above(k2)
72 
73             pb_vert_reg_temp(is,k1,k2) = (1.0-weight(k1)) * (1.0-weight(k2)) * &
74                                          pb_vert_reg_stats(is,k1s,k2s) + &
75                                          weight(k1) * (1.0-weight(k2)) * &
76                                          pb_vert_reg_stats(is,k1s+1,k2s) + &
77                                          weight(k2) * (1.0-weight(k1)) * &
78                                          pb_vert_reg_stats(is,k1s,k2s+1) + &
79                                          weight(k2) * weight(k1) * &
80                                          pb_vert_reg_stats(is,k1s+1,k2s+1)
81          end do
82       end do     
83    end do
84          
85    !------------------------------------------------------------------------
86    ! [4.0] Interpolate to from statistics latitudes to xb latitudes:
87    !------------------------------------------------------------------------
88 
89    i_south = 2
90 
91    do i = 1, iy
92    
93       ! Find position of xb latitude in statistics rows:
94 
95       if (meanl_xb(i) <= meanl_stats(   2)) then
96          i_south = 2
97          lat_wgt = 0.5
98       else if (meanl_xb(i) >= meanl_stats(iys-1)) then
99          i_south = iys-2
100          lat_wgt = 0.5
101       else
102          do is = 1, iys-1
103             if (meanl_xb(i) >= meanl_stats(is) .AND. &
104                  meanl_xb(i) <= meanl_stats(is+1)) then
105 
106                lat_wgt = (meanl_xb(i) - meanl_stats(is)) / &
107                          (meanl_stats(is+1) - meanl_stats(is))
108                i_south = is
109                exit
110             end if
111          end do
112       end if
113    
114       do k1 = 1, kz
115          do k2 = 1, kz
116             pb_vert_reg(i,k1,k2) = lat_wgt * pb_vert_reg_temp(i_south+1,k1,k2) + &
117                                    (1.0 - lat_wgt) * pb_vert_reg_temp(i_south,k1,k2)
118          end do
119       end do     
120    end do
121 
122    if (print_detail_regression) then
123       call da_array_print(1, pb_vert_reg_stats(1,:,:), &
124          'pb_vert_reg_stats(1,:,:)')
125       call da_array_print(1, pb_vert_reg(1,:,:), &
126          'pb_vert_reg(1,:,:)')
127       call da_array_print(1, pb_vert_reg_stats(:,1,:), &
128          'pb_vert_reg_stats(:,1,:)')
129       call da_array_print(1, pb_vert_reg(:,1,:), &
130          'pb_vert_reg(:,1,:)')
131    end if
132 
133    if (trace_use) call da_trace_exit("da_interpolate_regcoeff")
134    
135 end subroutine da_interpolate_regcoeff
136 
137