Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_physics / da_find_layer.inc
bloba22d366d1c36aa45bf17d2d41779512477acfe1d
1 subroutine da_find_layer(layer,tv,pre,pre_ma,tv_ma,ks,ke)
3    !-----------------------------------------------------------------------
4    ! Purpose: find the vertical location in the Tv profile given
5    ! a specific pressure and vertically interpolate Tv to that height.
6    ! pre_ma,tv_ma give vertical profile of virtual temperature
7    ! pre is a given pressure, alpha is the percentage of pre in the layer.
8    ! layer,tv are calculated vertical layer and interpolated virtual temp.
9    !-----------------------------------------------------------------------
11    implicit none
13    integer, intent(in)    :: ks, ke
14    integer, intent(inout) :: layer
15    real,    intent(inout) :: pre_ma(ks-1:ke+1)
16    real,    intent(in)    :: tv_ma(ks-1:ke+1)
17    real,    intent(in)    :: pre
18    real,    intent(inout) :: tv
20    integer :: k
21    real    :: alpha
23    if (trace_use_frequent) call da_trace_entry("da_find_layer")
25    if (pre >= pre_ma(ks)) then
26       ! Below model bottom
27       layer = ks
28       alpha = log(pre_ma(ks)/pre)/log(pre_ma(ks)/pre_ma(ks+1))
29       tv = tv_ma(ks) * (1.0-alpha) + tv_ma(ks+1) * alpha
30       pre_ma(ks-1)=pre
31    else if (pre <= pre_ma(ke)) then
32       ! Above model top
33       layer = ke+1
34       alpha = log(pre_ma(ke-1)/pre)/log(pre_ma(ke-1)/pre_ma(ke))
35       tv = tv_ma(ke-1) * (1.0-alpha) + tv_ma(ke) * alpha
36       pre_ma(ke+1) = pre
37    else
38       ! Between model layers 
39       do k=ks,ke-1
40          if (pre>=pre_ma(k+1) .and. pre<pre_ma(k)) then
41             layer = k+1
42             alpha = log(pre_ma(k)/pre)/log(pre_ma(k)/pre_ma(k+1))
43             tv = tv_ma(k) * (1.0-alpha) + tv_ma(k+1) * alpha
44             exit
45          end if
46       end do
47    end if
49    if (trace_use_frequent) call da_trace_exit("da_find_layer")
51 end subroutine da_find_layer