SICOPOLIS V3.0
 All Classes Files Functions Variables Macros
shift_cts_upward.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : s h i f t _ c t s _ u p w a r d
4 !
5 !> @file
6 !!
7 !! Upward shifting of the CTS.
8 !!
9 !! @section Copyright
10 !!
11 !! Copyright 2009-2013 Ralf Greve
12 !!
13 !! @section License
14 !!
15 !! This file is part of SICOPOLIS.
16 !!
17 !! SICOPOLIS is free software: you can redistribute it and/or modify
18 !! it under the terms of the GNU General Public License as published by
19 !! the Free Software Foundation, either version 3 of the License, or
20 !! (at your option) any later version.
21 !!
22 !! SICOPOLIS is distributed in the hope that it will be useful,
23 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
24 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25 !! GNU General Public License for more details.
26 !!
27 !! You should have received a copy of the GNU General Public License
28 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
29 !<
30 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31 
32 !-------------------------------------------------------------------------------
33 !> Upward shifting of the CTS.
34 !<------------------------------------------------------------------------------
35 subroutine shift_cts_upward(at1, at2_1, at2_2, at3_1, at3_2, &
36  at4_1, at4_2, at5, at6, at7, atr1, am1, am2, alb1, &
37  aw1, aw2, aw3, aw4, aw5, aw7, aw8, aw9, aqtld, &
38  ai1, ai2, ai3, mean_accum_inv, dzeta_t, &
39  dtime_temp, dtt_2dxi, dtt_2deta, dtime_temp_inv, &
40  i, j)
41 
42 use sico_types
44 
45 implicit none
46 integer(i4b) :: i, j
47 real(dp) :: zm_shift
48 real(dp) :: at1(0:kcmax), at2_1(0:kcmax), at2_2(0:kcmax), &
49  at3_1(0:kcmax), at3_2(0:kcmax), at4_1(0:kcmax), &
50  at4_2(0:kcmax), at5(0:kcmax), at6(0:kcmax), at7, &
51  ai1(0:kcmax), ai2(0:kcmax), ai3, &
52  atr1, am1, am2, alb1
53 real(dp) :: aw1, aw2, aw3, aw4, aw5, aw7, aw8, aw9, aqtld
54 real(dp) :: mean_accum_inv
55 real(dp) :: dtt_2dxi, dtt_2deta, dtime_temp, dtime_temp_inv, dzeta_t
56 real(dp) :: difftemp_a, difftemp_b, interpol
57 
58 zm_shift = 1.0_dp ! CTS shift in intervals of 1 m
59 
60 !-------- Temperature discrepancy from the computation of the main
61 ! program --------
62 
63 difftemp_a = temp_c_neu(0,j,i)-(-beta*h_c_neu(j,i))
64 if (difftemp_a.le.0.0_dp) return
65 
66 !-------- Shift CTS upward until it is too high --------
67 
68  10 continue
69 
70  zm_neu(j,i) = zm_neu(j,i) + zm_shift
71  if (zm_neu(j,i).ge.zs(j,i)) then
72  zm_neu(j,i) = zm_neu(j,i) - zm_shift
73  return
74  end if
75  h_c_neu(j,i) = h_c_neu(j,i) - zm_shift
76  h_t_neu(j,i) = h_t_neu(j,i) + zm_shift
77 
78  dh_t_dtau(j,i) = (zm_neu(j,i)-zm(j,i))*dtime_temp_inv
79  dzm_dtau(j,i) = dzb_dtau(j,i)+dh_t_dtau(j,i)
80  dh_c_dtau(j,i) = dzs_dtau(j,i)-dzm_dtau(j,i)
81 
82  am_perp(j,i) = am_perp_st(j,i) + dzm_dtau(j,i)
83 
84  call calc_temp3(at1, at2_1, at2_2, at3_1, at3_2, &
85  at4_1, at4_2, at5, at6, at7, atr1, &
86  am1, am2, alb1, &
87  aw1, aw2, aw3, aw4, aw5, aw7, aw8, aw9, aqtld, &
88  ai1, ai2, ai3, mean_accum_inv, dzeta_t, &
89  dtime_temp, dtt_2dxi, dtt_2deta, i, j)
90 
91  difftemp_b = difftemp_a
92  difftemp_a = temp_c_neu(0,j,i)-(-beta*h_c_neu(j,i))
93 
94 if (difftemp_a.gt.0.0_dp) go to 10
95 
96 !-------- Interpolate the CTS position from the last (_a) and the
97 ! last but one (_b) value, weighed with the temperature
98 ! discrepancies at the CTS --------
99 
100 interpol = difftemp_a/(difftemp_b-difftemp_a)*zm_shift
101 
102 zm_neu(j,i) = zm_neu(j,i) + interpol
103 h_c_neu(j,i) = h_c_neu(j,i) - interpol
104 h_t_neu(j,i) = h_t_neu(j,i) + interpol
105 
106 dh_t_dtau(j,i) = (zm_neu(j,i)-zm(j,i))*dtime_temp_inv
107 dzm_dtau(j,i) = dzb_dtau(j,i)+dh_t_dtau(j,i)
108 dh_c_dtau(j,i) = dzs_dtau(j,i)-dzm_dtau(j,i)
109 
110 am_perp(j,i) = am_perp_st(j,i) + dzm_dtau(j,i)
111 
112 call calc_temp3(at1, at2_1, at2_2, at3_1, at3_2, &
113  at4_1, at4_2, at5, at6, at7, atr1, &
114  am1, am2, alb1, &
115  aw1, aw2, aw3, aw4, aw5, aw7, aw8, aw9, aqtld, &
116  ai1, ai2, ai3, mean_accum_inv, dzeta_t, &
117  dtime_temp, dtt_2dxi, dtt_2deta, i, j)
118 
119 end subroutine shift_cts_upward
120 !