da_find_layer_tl.inc

References to this file elsewhere.
1 subroutine da_find_layer_tl(layer,tv,pre,pre_ma,tv_ma,ks,ke,     &
2                           TGL_tv,TGL_pre_ma,TGL_tv_ma)
3 
4    !-----------------------------------------------------------------------
5    ! Purpose: tangent-linear routine for da_find_layer
6    !-----------------------------------------------------------------------
7 
8    implicit none
9 
10    integer,intent(in)                       :: ks, ke
11    integer,intent(out)                      :: layer
12    real, dimension(ks-1:ke+1),intent(in)    :: tv_ma
13    real, dimension(ks-1:ke+1),intent(inout) :: pre_ma
14    real, dimension(ks-1:ke+1),intent(in)    :: TGL_tv_ma
15    real, dimension(ks-1:ke+1),intent(inout) :: TGL_pre_ma
16    real,intent(in)                          :: pre
17    real,intent(out)                         :: tv
18    real,intent(out)                         :: TGL_tv
19 
20    real                      :: TGL_alpha
21    integer                   :: k
22    real                      :: alpha, coef1, coef2
23 
24    ! coef1, coef2 are temporarily used in this routine
25 
26    if (pre >= pre_ma(ks)) then
27 
28       ! Below model bottom
29       layer = ks
30       coef1=log(pre/pre_ma(ks+1))/(pre_ma(ks)*     &
31             (log(pre_ma(ks)/pre_ma(ks+1)))**2)
32       coef2=log(pre_ma(ks)/pre)/(pre_ma(ks+1)*     &
33             (log(pre_ma(ks)/pre_ma(ks+1)))**2)
34       TGL_alpha = coef1 * TGL_pre_ma(ks) + coef2 * TGL_pre_ma(ks+1)
35       alpha = log(pre_ma(ks)/pre)/log(pre_ma(ks)/pre_ma(ks+1))
36 
37       TGL_tv = (1.-alpha)*TGL_tv_ma(ks) +               &
38                (tv_ma(ks+1)-tv_ma(ks))*TGL_alpha +     &
39                alpha*TGL_tv_ma(ks+1)
40       TGL_pre_ma(ks-1) = 0.
41       tv = tv_ma(ks) * (1.-alpha) + tv_ma(ks+1) * alpha
42       pre_ma(ks-1) = pre
43 
44    else if (pre <= pre_ma(ke)) then
45 
46       ! Above model top
47       layer = ke+1
48       coef1=log(pre/pre_ma(ke))/(pre_ma(ke-1)*           &
49             (log(pre_ma(ke-1)/pre_ma(ke)))**2)
50       coef2=log(pre_ma(ke-1)/pre)/(pre_ma(ke)*           &
51             (log(pre_ma(ke-1)/pre_ma(ke)))**2)
52       TGL_alpha = coef1 * TGL_pre_ma(ke-1) + coef2 * TGL_pre_ma(ke)
53       alpha = log(pre_ma(ke-1)/pre)/log(pre_ma(ke-1)/pre_ma(ke))
54 
55       TGL_tv = (1.-alpha)*TGL_tv_ma(ke-1) +                 &
56                (tv_ma(ke)-tv_ma(ke-1))*TGL_alpha +           &
57                alpha*TGL_tv_ma(ke)
58       TGL_pre_ma(ke+1) = 0.
59       tv = tv_ma(ke-1) * (1.-alpha) + tv_ma(ke) * alpha
60       pre_ma(ke+1) = pre
61 
62    else
63 
64       ! Between model layers
65       do k=ks,ke-1
66          if (pre>=pre_ma(k+1) .and. pre<pre_ma(k)) then
67             layer = k+1
68             coef1=log(pre/pre_ma(k+1))/(pre_ma(k)*   &
69                   (log(pre_ma(k)/pre_ma(k+1)))**2)
70             coef2=log(pre_ma(k)/pre)/(pre_ma(k+1)*   &
71                   (log(pre_ma(k)/pre_ma(k+1)))**2)
72             TGL_alpha = coef1 * TGL_pre_ma(k) + coef2 * TGL_pre_ma(k+1)
73             alpha = log(pre_ma(k)/pre)/log(pre_ma(k)/pre_ma(k+1))
74             TGL_tv = (1.-alpha)*TGL_tv_ma(k) +                 &
75                      (tv_ma(k+1)-tv_ma(k))*TGL_alpha +         &
76                       alpha*TGL_tv_ma(k+1)
77             tv = tv_ma(k) * (1.-alpha) + tv_ma(k+1) * alpha
78             exit
79          end if
80       end do
81 
82    end if
83  
84 end subroutine da_find_layer_tl
85 
86