SICOPOLIS V3.0
 All Classes Files Functions Variables Macros
calc_enhance.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : c a l c _ e n h a n c e
4 !
5 !> @file
6 !!
7 !! Computation of the flow enhancement factor.
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 the flow enhancement factor.
34 !<------------------------------------------------------------------------------
35 subroutine calc_enhance(time)
36 
37 use sico_types
39 
40 implicit none
41 
42 real(dp), intent(in) :: time
43 
44 integer(i4b) :: i, j, kc, kt
45 real(dp) :: age_trans, date_trans1, date_trans2, date_trans3
46 real(dp) :: enh_mult
47 
48 #if ENHMOD==1
49 
50 enh_t = enh_fact
51 enh_c = enh_fact
52 
53 #elif ENHMOD==2
54 
55 age_trans = age_trans_0*year_sec
56 
57 do i=0, imax
58 do j=0, jmax
59 
60  do kt=0, ktmax
61  if (age_t(kt,j,i) < age_trans) then
62  enh_t(kt,j,i) = enh_intg ! Holocene ice
63  else
64  enh_t(kt,j,i) = enh_fact ! Pleistocene ice
65  end if
66  end do
67 
68  do kc=0, kcmax
69  if (age_c(kc,j,i) < age_trans) then
70  enh_c(kc,j,i) = enh_intg ! Holocene ice
71  else
72  enh_c(kc,j,i) = enh_fact ! Pleistocene ice
73  end if
74  end do
75 
76 end do
77 end do
78 
79 #elif ENHMOD==3
80 
81 date_trans1 = date_trans1_0*year_sec
82 date_trans2 = date_trans2_0*year_sec
83 date_trans3 = date_trans3_0*year_sec
84 
85 do i=0, imax
86 do j=0, jmax
87 
88  do kt=0, ktmax
89  if ( (time-age_t(kt,j,i)) < date_trans1 ) then
90  enh_t(kt,j,i) = enh_fact ! pre-Eemian ice
91  else if ( ((time-age_t(kt,j,i)) >= date_trans1).and. &
92  ((time-age_t(kt,j,i)) < date_trans2) ) then
93  enh_t(kt,j,i) = enh_intg ! Eemian ice
94  else if ( ((time-age_t(kt,j,i)) >= date_trans2).and. &
95  ((time-age_t(kt,j,i)) < date_trans3) ) then
96  enh_t(kt,j,i) = enh_fact ! Weichselian ice
97  else
98  enh_t(kt,j,i) = enh_intg ! Holocene ice
99  end if
100  end do
101 
102  do kc=0, kcmax
103  if ( (time-age_c(kc,j,i)) < date_trans1 ) then
104  enh_c(kc,j,i) = enh_fact ! pre-Eemian ice
105  else if ( ((time-age_c(kc,j,i)) >= date_trans1).and. &
106  ((time-age_c(kc,j,i)) < date_trans2) ) then
107  enh_c(kc,j,i) = enh_intg ! Eemian ice
108  else if ( ((time-age_c(kc,j,i)) >= date_trans2).and. &
109  ((time-age_c(kc,j,i)) < date_trans3) ) then
110  enh_c(kc,j,i) = enh_fact ! Weichselian ice
111  else
112  enh_c(kc,j,i) = enh_intg ! Holocene ice
113  end if
114  end do
115 
116 end do
117 end do
118 
119 #endif
120 
121 !-------- Modification for floating ice --------
122 
123 #if MARGIN==3
124 
125 do i=0, imax
126 do j=0, jmax
127 
128  if ( maske(j,i)==3_i2b ) then ! floating ice
129 
130  do kt=0, ktmax
131  enh_t(kt,j,i) = enh_shelf
132  end do
133 
134  do kc=0, kcmax
135  enh_c(kc,j,i) = enh_shelf
136  end do
137 
138  end if
139 
140 end do
141 end do
142 
143 #endif
144 
145 !-------- Modification due to dust content (Martian polar caps only) --------
146 
147 #if ( defined(NMARS) \
148  || defined(smars) ) /* martian ice sheet */
149 
150 #if (FLOW_LAW==1)
151  enh_mult = exp(-3.0_dp*2.0_dp*frac_dust)
152 #elif (FLOW_LAW==2)
153  enh_mult = exp(-1.8_dp*2.0_dp*frac_dust)
154 #elif (FLOW_LAW==3)
155  enh_mult = exp(-4.0_dp*2.0_dp*frac_dust)
156 #elif (FLOW_LAW==4)
157  stop ' calc_enhance: FLOW_LAW==4 has not been implemented yet!'
158 #endif
159 
160 enh_t = enh_t * enh_mult
161 enh_c = enh_c * enh_mult
162 
163 #endif
164 
165 end subroutine calc_enhance
166 !