SICOPOLIS V3.1
 All Classes Files Functions Variables Macros
boundary.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : b o u n d a r y
4 !
5 !> @file
6 !!
7 !! Computation of the surface temperature (must be less than 0 deg C!)
8 !! and of the accumulation-ablation function.
9 !!
10 !! @section Copyright
11 !!
12 !! Copyright 2009-2013 Ralf Greve
13 !!
14 !! @section License
15 !!
16 !! This file is part of SICOPOLIS.
17 !!
18 !! SICOPOLIS is free software: you can redistribute it and/or modify
19 !! it under the terms of the GNU General Public License as published by
20 !! the Free Software Foundation, either version 3 of the License, or
21 !! (at your option) any later version.
22 !!
23 !! SICOPOLIS is distributed in the hope that it will be useful,
24 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
25 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
26 !! GNU General Public License for more details.
27 !!
28 !! You should have received a copy of the GNU General Public License
29 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
30 !<
31 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 
33 !-------------------------------------------------------------------------------
34 !> Computation of the surface temperature (must be less than 0 deg C!)
35 !! and of the accumulation-ablation function.
36 !<------------------------------------------------------------------------------
37 subroutine boundary(time, dtime, dxi, deta, &
38  delta_ts, glac_index, z_sl, dzsl_dtau, z_mar)
39 
40 use sico_types
42 
43 implicit none
44 
45 real(dp), intent(in) :: time, dtime, dxi, deta
46 
47 real(dp), intent(out) :: delta_ts, glac_index, dzsl_dtau, z_mar
48 real(dp), intent(inout) :: z_sl
49 
50 ! Further return variables
51 ! (defined as global variables in module sico_variables):
52 !
53 ! accum(j,i), as_perp(j,i), temp_s(j,i)
54 
55 integer(i4b) :: i, j
56 integer(i4b) :: i_gr, i_kl
57 real(dp) :: z_sl_old
58 real(dp) :: z_sl_min, t1, t2, t3, t4, t5, t6
59 real(dp) :: time_gr, time_kl
60 real(dp), dimension(0:JMAX,0:IMAX) :: dist
61 logical, dimension(0:JMAX,0:IMAX) :: check_point
62 
63 !-------- Initialization of intent(out) variables --------
64 
65 z_sl_old = z_sl
66 
67 delta_ts = 0.0_dp
68 glac_index = 0.0_dp
69 z_sl = 0.0_dp
70 dzsl_dtau = 0.0_dp
71 z_mar = 0.0_dp
72 
73 !-------- Surface-temperature deviation from present values --------
74 
75 #if TSURFACE==1
76 delta_ts = delta_ts0
77 ! ! Steady state with prescribed constant
78 ! ! air-temperature deviation
79 #elif TSURFACE==3
80 delta_ts = sine_amplit &
81  *cos(2.0_dp*pi*time/(sine_period*year_sec)) &
82  -sine_amplit
83 ! ! Sinusoidal air-temperature forcing
84 #elif TSURFACE==4
85 
86 ! ------ delta_ts from the GRIP record
87 
88 if (time/year_sec.lt.real(grip_time_min,dp)) then
89  delta_ts = griptemp(0)
90 else if (time/year_sec.lt.real(grip_time_max,dp)) then
91 
92  i_kl = floor(((time/year_sec) &
93  -real(grip_time_min,dp))/real(grip_time_stp,dp))
94  i_kl = max(i_kl, 0)
95 
96  i_gr = ceiling(((time/year_sec) &
97  -real(grip_time_min,dp))/real(grip_time_stp,dp))
98  i_gr = min(i_gr, ndata_grip)
99 
100  if (i_kl.eq.i_gr) then
101 
102  delta_ts = griptemp(i_kl)
103 
104  else
105 
106  time_kl = (grip_time_min + i_kl*grip_time_stp) *year_sec
107  time_gr = (grip_time_min + i_gr*grip_time_stp) *year_sec
108 
109  delta_ts = griptemp(i_kl) &
110  +(griptemp(i_gr)-griptemp(i_kl)) &
111  *(time-time_kl)/(time_gr-time_kl)
112  ! linear interpolation of the ice-core data
113 
114  end if
115 
116 else
117  delta_ts = griptemp(ndata_grip)
118 end if
119 
120 delta_ts = delta_ts * grip_temp_fact
121 ! ! modification by constant factor
122 
123 #elif TSURFACE==5
124 ! Enter here delta_ts scenario:
125 
126 #endif
127 
128 !-------- Sea level --------
129 
130 #if SEA_LEVEL==1
131 ! ------ constant sea level
132 z_sl = z_sl0
133 
134 #elif SEA_LEVEL==2
135 ! ------ saw-tooth-shaped palaeoclimatic sea-level forcing
136 
137 z_sl_min = -130.0_dp
138 
139 t1 = -250000.0_dp *year_sec
140 t2 = -140000.0_dp *year_sec
141 t3 = -125000.0_dp *year_sec
142 t4 = -21000.0_dp *year_sec
143 t5 = -8000.0_dp *year_sec
144 t6 = 0.0_dp *year_sec
145 
146 if (time.lt.t1) then
147  z_sl = 0.0_dp
148 else if (time.lt.t2) then
149  z_sl = z_sl_min*(time-t1)/(t2-t1)
150 else if (time.lt.t3) then
151  z_sl = -z_sl_min*(time-t3)/(t3-t2)
152 else if (time.lt.t4) then
153  z_sl = z_sl_min*(time-t3)/(t4-t3)
154 else if (time.lt.t5) then
155  z_sl = -z_sl_min*(time-t5)/(t5-t4)
156 else if (time.lt.t6) then
157  z_sl = 0.0_dp
158 else
159  z_sl = 0.0_dp
160 end if
161 
162 #elif SEA_LEVEL==3
163 
164 ! ------ z_sl from the SPECMAP record
165 
166 if (time/year_sec.lt.real(specmap_time_min,dp)) then
167  z_sl = specmap_zsl(0)
168 else if (time/year_sec.lt.real(specmap_time_max,dp)) then
169 
170  i_kl = floor(((time/year_sec) &
171  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
172  i_kl = max(i_kl, 0)
173 
174  i_gr = ceiling(((time/year_sec) &
175  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
176  i_gr = min(i_gr, ndata_specmap)
177 
178  if (i_kl.eq.i_gr) then
179 
180  z_sl = specmap_zsl(i_kl)
181 
182  else
183 
184  time_kl = (specmap_time_min + i_kl*specmap_time_stp) *year_sec
185  time_gr = (specmap_time_min + i_gr*specmap_time_stp) *year_sec
186 
187  z_sl = specmap_zsl(i_kl) &
188  +(specmap_zsl(i_gr)-specmap_zsl(i_kl)) &
189  *(time-time_kl)/(time_gr-time_kl)
190  ! linear interpolation of the sea-level data
191 
192  end if
193 
194 else
195  z_sl = specmap_zsl(ndata_specmap)
196 end if
197 
198 #endif
199 
200 ! ------ Time derivative of the sea level
201 
202 if ( z_sl_old > -999999.9_dp ) then
203  dzsl_dtau = (z_sl-z_sl_old)/dtime
204 else ! only dummy value for z_sl_old available
205  dzsl_dtau = 0.0_dp
206 end if
207 
208 ! ------ Minimum bedrock elevation for extent of marine ice
209 
210 #if MARGIN==2
211 
212 #if ( MARINE_ICE_CALVING==2 || MARINE_ICE_CALVING==3 )
213 z_mar = z_mar
214 #elif ( MARINE_ICE_CALVING==4 || MARINE_ICE_CALVING==5 )
215 z_mar = fact_z_mar*z_sl
216 #elif ( MARINE_ICE_CALVING==6 || MARINE_ICE_CALVING==7 )
217 if (z_sl >= -80.0_dp) then
218  z_mar = 2.5_dp*z_sl
219 else
220  z_mar = 10.25_dp*(z_sl+80.0_dp)-200.0_dp
221 end if
222 z_mar = fact_z_mar*z_mar
223 #endif
224 
225 #endif
226 
227 ! ------ Update of the mask according to the sea level
228 
229 ! ---- Check all sea and floating-ice points and their direct
230 ! neighbours
231 
232 do i=0, imax
233 do j=0, jmax
234  check_point(j,i) = .false.
235 end do
236 end do
237 
238 do i=1, imax-1
239 do j=1, jmax-1
240  if (maske(j,i).ge.2) then
241  check_point(j ,i ) = .true.
242  check_point(j ,i+1) = .true.
243  check_point(j ,i-1) = .true.
244  check_point(j+1,i ) = .true.
245  check_point(j-1,i ) = .true.
246  end if
247 end do
248 end do
249 
250 do i=1, imax-1
251 do j=1, jmax-1
252  if (check_point(j,i)) then
253  maske_neu(j,i) = mask_update(z_sl, i, j)
254  end if
255 end do
256 end do
257 
258 ! ---- Assign new values of the mask
259 
260 do i=1, imax-1
261 do j=1, jmax-1
262  if (check_point(j,i)) then
263  maske(j,i) = maske_neu(j,i)
264  end if
265 end do
266 end do
267 
268 !-------- Surface air temperature and accumulation-ablation function --------
269 
270 do i=0, imax
271 do j=0, jmax
272  dist(j,i) = sqrt( (xi(i)-x_hat)**2 + (eta(j)-y_hat)**2 )
273 end do
274 end do
275 
276 do i=0, imax
277 do j=0, jmax
278 
279  temp_s(j,i) = temp_min + s_t*dist(j,i)
280  temp_s(j,i) = temp_s(j,i) + delta_ts
281  ! Correction with temperature deviation delta_ts
282  if (temp_s(j,i).gt. -0.001_dp) temp_s(j,i) = -0.001_dp
283  ! Cut-off of positive air temperatures
284 
285  accum_present(j,i) = b_max
286  accum(j,i) = accum_present(j,i)
287 
288  as_perp(j,i) = s_b*(eld-dist(j,i))
289  if (as_perp(j,i).gt.accum(j,i)) as_perp(j,i) = accum(j,i)
290 
291 end do
292 end do
293 
294 end subroutine boundary
295 !