da_thickness_adj.inc
References to this file elsewhere.
1 subroutine da_thickness_adj(pre_ma,tv_ma,ks,ke,tv1,tv2,layer1,layer2,pre1,pre2, &
2 ADJ_pre_ma,ADJ_tv_ma,ADJ_tv1,ADJ_tv2,ADJ_thk)
3
4 !-----------------------------------------------------------------------
5 ! Purpose: adjoint routine for thickness
6 !-----------------------------------------------------------------------
7
8 implicit none
9
10 integer, intent(in) :: layer1,layer2
11 integer, intent(in) :: ks,ke
12 real, intent(in) :: pre_ma(ks-1:ke+1)
13 real, intent(in) :: tv_ma(ks-1:ke+1)
14 real, intent(inout) :: ADJ_pre_ma(ks-1:ke+1)
15 real, intent(inout) :: ADJ_tv_ma(ks-1:ke+1)
16 real, intent(in) :: tv1,tv2
17 real, intent(inout) :: ADJ_tv1,ADJ_tv2
18 real, intent(in) :: pre1,pre2
19 real, intent(inout) :: ADJ_thk
20
21 integer :: k
22 real :: p_tmp(ks-1:ke+1)
23 real :: ADJ_p_tmp(ks-1:ke+1)
24
25 if (trace_use_dull) call da_trace_entry("da_thickness_adj")
26
27 ! p_tmp and ADJ_p_tmp are temporary (local) variables
28
29 ADJ_p_tmp(:)=0.0
30
31 p_tmp(layer1) = pre1
32 p_tmp(layer2-1) = pre2
33 do k=layer2,layer1-1
34 p_tmp(k) = pre_ma(k)
35 end do
36
37 do k=layer2,layer1-1
38 ADJ_p_tmp(k+1) = ADJ_p_tmp(k+1) - 0.5*gas_constant/gravity * &
39 ADJ_thk*tv_ma(k)/p_tmp(k+1)
40 ADJ_p_tmp(k-1) = ADJ_p_tmp(k-1) + 0.5*gas_constant/gravity * &
41 ADJ_thk*tv_ma(k)/p_tmp(k-1)
42 ADJ_tv_ma(k) = ADJ_tv_ma(k) + 0.5*gas_constant/gravity * &
43 ADJ_thk*log(p_tmp(k-1)/p_tmp(k+1))
44 end do
45
46 do k=layer2,layer1-1
47 ADJ_pre_ma(k) = ADJ_pre_ma(k) + ADJ_p_tmp(k)
48 end do
49
50 ADJ_thk = 0.5 * gas_constant/gravity * ADJ_thk
51 ADJ_pre_ma(layer2) = ADJ_pre_ma(layer2) - ADJ_thk*tv2/pre_ma(layer2)
52 ADJ_tv2 = ADJ_tv2 + ADJ_thk*log(pre2/pre_ma(layer2))
53 ADJ_pre_ma(layer1-1) = ADJ_pre_ma(layer1-1) + &
54 ADJ_thk*tv1/pre_ma(layer1-1)
55 ADJ_tv1 = ADJ_tv1 + ADJ_thk*log(pre_ma(layer1-1)/pre1)
56
57 if (trace_use_dull) call da_trace_exit("da_thickness_adj")
58
59 end subroutine da_thickness_adj
60
61