SICOPOLIS V3.0
 All Classes Files Functions Variables Macros
calc_temp_r.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : c a l c _ t e m p _ r
4 !
5 !> @file
6 !!
7 !! Computation of temperature for an ice-free column.
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 !> Computation of temperature for an ice-free column.
34 !<------------------------------------------------------------------------------
35 subroutine calc_temp_r(atr1, alb1, i, j)
36 
37 use sico_types
39 
40 implicit none
41 integer(i4b) :: i, j, kc, kt, kr
42 real(dp) :: atr1, alb1, ctr1, clb1
43 real(dp) :: lgs_a0(0:kcmax+ktmax+krmax+imax+jmax), &
44  lgs_a1(0:kcmax+ktmax+krmax+imax+jmax), &
45  lgs_a2(0:kcmax+ktmax+krmax+imax+jmax), &
46  lgs_x(0:kcmax+ktmax+krmax+imax+jmax), &
47  lgs_b(0:kcmax+ktmax+krmax+imax+jmax)
48 
49 !-------- Abbreviations --------
50 
51 ctr1 = atr1
52 clb1 = alb1*q_geo(j,i)
53 
54 !-------- Set up the equations for the bedrock temperature --------
55 
56 kr=0
57 lgs_a1(kr) = 1.0_dp
58 lgs_a2(kr) = -1.0_dp
59 lgs_b(kr) = clb1
60 
61 #if Q_LITHO==1
62 ! (coupled heat-conducting bedrock)
63 
64 do kr=1, krmax-1
65  lgs_a0(kr) = - ctr1
66  lgs_a1(kr) = 1.0_dp + 2.0_dp*ctr1
67  lgs_a2(kr) = - ctr1
68  lgs_b(kr) = temp_r(kr,j,i)
69 end do
70 
71 #elif Q_LITHO==0
72 ! (no coupled heat-conducting bedrock)
73 
74 do kr=1, krmax-1
75  lgs_a0(kr) = 1.0_dp
76  lgs_a1(kr) = 0.0_dp
77  lgs_a2(kr) = -1.0_dp
78  lgs_b(kr) = 2.0_dp*clb1
79 end do
80 
81 #endif
82 
83 kr=krmax
84 lgs_a0(kr) = 0.0_dp
85 lgs_a1(kr) = 1.0_dp
86 lgs_b(kr) = temp_s(j,i)
87 
88 !-------- Solve system of linear equations --------
89 
90 call tri_sle(lgs_a0, lgs_a1, lgs_a2, lgs_x, lgs_b, krmax)
91 
92 !-------- Assign the result --------
93 
94 do kr=0, krmax
95  temp_r_neu(kr,j,i) = lgs_x(kr)
96 end do
97 
98 !-------- Water content and age in the non-existing temperate layer --------
99 
100 do kt=0, ktmax
101  omega_t_neu(kt,j,i) = 0.0_dp
102  age_t_neu(kt,j,i) = 0.0_dp
103 end do
104 
105 !-------- Temperature and age in the non-existing cold layer --------
106 
107 do kc=0, kcmax
108  temp_c_neu(kc,j,i) = temp_s(j,i)
109  age_c_neu(kc,j,i) = 0.0_dp
110 end do
111 
112 end subroutine calc_temp_r
113 !