SICOPOLIS V3.1
 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 use sico_sle_solvers, only : tri_sle
40 
41 implicit none
42 
43 integer(i4b), intent(in) :: i, j
44 real(dp), intent(in) :: atr1, alb1
45 
46 integer(i4b) :: kc, kt, kr
47 real(dp) :: ctr1, clb1
48 real(dp) :: lgs_a0(0:kcmax+ktmax+krmax+imax+jmax), &
49  lgs_a1(0:kcmax+ktmax+krmax+imax+jmax), &
50  lgs_a2(0:kcmax+ktmax+krmax+imax+jmax), &
51  lgs_x(0:kcmax+ktmax+krmax+imax+jmax), &
52  lgs_b(0:kcmax+ktmax+krmax+imax+jmax)
53 
54 !-------- Abbreviations --------
55 
56 ctr1 = atr1
57 clb1 = alb1*q_geo(j,i)
58 
59 !-------- Set up the equations for the bedrock temperature --------
60 
61 kr=0
62 lgs_a1(kr) = 1.0_dp
63 lgs_a2(kr) = -1.0_dp
64 lgs_b(kr) = clb1
65 
66 #if Q_LITHO==1
67 ! (coupled heat-conducting bedrock)
68 
69 do kr=1, krmax-1
70  lgs_a0(kr) = - ctr1
71  lgs_a1(kr) = 1.0_dp + 2.0_dp*ctr1
72  lgs_a2(kr) = - ctr1
73  lgs_b(kr) = temp_r(kr,j,i)
74 end do
75 
76 #elif Q_LITHO==0
77 ! (no coupled heat-conducting bedrock)
78 
79 do kr=1, krmax-1
80  lgs_a0(kr) = 1.0_dp
81  lgs_a1(kr) = 0.0_dp
82  lgs_a2(kr) = -1.0_dp
83  lgs_b(kr) = 2.0_dp*clb1
84 end do
85 
86 #endif
87 
88 kr=krmax
89 lgs_a0(kr) = 0.0_dp
90 lgs_a1(kr) = 1.0_dp
91 lgs_b(kr) = temp_s(j,i)
92 
93 !-------- Solve system of linear equations --------
94 
95 call tri_sle(lgs_a0, lgs_a1, lgs_a2, lgs_x, lgs_b, krmax)
96 
97 !-------- Assign the result --------
98 
99 do kr=0, krmax
100  temp_r_neu(kr,j,i) = lgs_x(kr)
101 end do
102 
103 !-------- Water content and age in the non-existing temperate layer --------
104 
105 do kt=0, ktmax
106  omega_t_neu(kt,j,i) = 0.0_dp
107  age_t_neu(kt,j,i) = 0.0_dp
108 end do
109 
110 !-------- Temperature and age in the non-existing cold layer --------
111 
112 do kc=0, kcmax
113  temp_c_neu(kc,j,i) = temp_s(j,i)
114  age_c_neu(kc,j,i) = 0.0_dp
115 end do
116 
117 end subroutine calc_temp_r
118 !