SICOPOLIS V3.0
 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(i2b) :: mask_update
56 integer(i4b) :: i, j
57 integer(i4b) :: i_gr, i_kl
58 real(dp) :: z_sl_old
59 real(dp) :: z_sl_min, t1, t2, t3, t4, t5, t6
60 real(dp) :: time_gr, time_kl
61 real(dp), dimension(0:JMAX,0:IMAX) :: dist
62 logical, dimension(0:JMAX,0:IMAX) :: check_point
63 
64 !-------- Initialization of intent(out) variables --------
65 
66 z_sl_old = z_sl
67 
68 delta_ts = 0.0_dp
69 glac_index = 0.0_dp
70 z_sl = 0.0_dp
71 dzsl_dtau = 0.0_dp
72 z_mar = 0.0_dp
73 
74 !-------- Surface-temperature deviation from present values --------
75 
76 #if TSURFACE==1
77 delta_ts = delta_ts0
78 ! ! Steady state with prescribed constant
79 ! ! air-temperature deviation
80 #elif TSURFACE==3
81 delta_ts = sine_amplit &
82  *cos(2.0_dp*pi*time/(sine_period*year_sec)) &
83  -sine_amplit
84 ! ! Sinusoidal air-temperature forcing
85 #elif TSURFACE==4
86 
87 ! ------ delta_ts from the GRIP record
88 
89 if (time/year_sec.lt.real(grip_time_min,dp)) then
90  delta_ts = griptemp(0)
91 else if (time/year_sec.lt.real(grip_time_max,dp)) then
92 
93  i_kl = floor(((time/year_sec) &
94  -real(grip_time_min,dp))/real(grip_time_stp,dp))
95  i_kl = max(i_kl, 0)
96 
97  i_gr = ceiling(((time/year_sec) &
98  -real(grip_time_min,dp))/real(grip_time_stp,dp))
99  i_gr = min(i_gr, ndata_grip)
100 
101  if (i_kl.eq.i_gr) then
102 
103  delta_ts = griptemp(i_kl)
104 
105  else
106 
107  time_kl = (grip_time_min + i_kl*grip_time_stp) *year_sec
108  time_gr = (grip_time_min + i_gr*grip_time_stp) *year_sec
109 
110  delta_ts = griptemp(i_kl) &
111  +(griptemp(i_gr)-griptemp(i_kl)) &
112  *(time-time_kl)/(time_gr-time_kl)
113  ! linear interpolation of the ice-core data
114 
115  end if
116 
117 else
118  delta_ts = griptemp(ndata_grip)
119 end if
120 
121 delta_ts = delta_ts * grip_temp_fact
122 ! ! modification by constant factor
123 
124 #elif TSURFACE==5
125 ! Enter here delta_ts scenario:
126 
127 #endif
128 
129 !-------- Sea level --------
130 
131 #if SEA_LEVEL==1
132 ! ------ constant sea level
133 z_sl = z_sl0
134 
135 #elif SEA_LEVEL==2
136 ! ------ saw-tooth-shaped palaeoclimatic sea-level forcing
137 
138 z_sl_min = -130.0_dp
139 
140 t1 = -250000.0_dp *year_sec
141 t2 = -140000.0_dp *year_sec
142 t3 = -125000.0_dp *year_sec
143 t4 = -21000.0_dp *year_sec
144 t5 = -8000.0_dp *year_sec
145 t6 = 0.0_dp *year_sec
146 
147 if (time.lt.t1) then
148  z_sl = 0.0_dp
149 else if (time.lt.t2) then
150  z_sl = z_sl_min*(time-t1)/(t2-t1)
151 else if (time.lt.t3) then
152  z_sl = -z_sl_min*(time-t3)/(t3-t2)
153 else if (time.lt.t4) then
154  z_sl = z_sl_min*(time-t3)/(t4-t3)
155 else if (time.lt.t5) then
156  z_sl = -z_sl_min*(time-t5)/(t5-t4)
157 else if (time.lt.t6) then
158  z_sl = 0.0_dp
159 else
160  z_sl = 0.0_dp
161 end if
162 
163 #elif SEA_LEVEL==3
164 
165 ! ------ z_sl from the SPECMAP record
166 
167 if (time/year_sec.lt.real(specmap_time_min,dp)) then
168  z_sl = specmap_zsl(0)
169 else if (time/year_sec.lt.real(specmap_time_max,dp)) then
170 
171  i_kl = floor(((time/year_sec) &
172  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
173  i_kl = max(i_kl, 0)
174 
175  i_gr = ceiling(((time/year_sec) &
176  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
177  i_gr = min(i_gr, ndata_specmap)
178 
179  if (i_kl.eq.i_gr) then
180 
181  z_sl = specmap_zsl(i_kl)
182 
183  else
184 
185  time_kl = (specmap_time_min + i_kl*specmap_time_stp) *year_sec
186  time_gr = (specmap_time_min + i_gr*specmap_time_stp) *year_sec
187 
188  z_sl = specmap_zsl(i_kl) &
189  +(specmap_zsl(i_gr)-specmap_zsl(i_kl)) &
190  *(time-time_kl)/(time_gr-time_kl)
191  ! linear interpolation of the sea-level data
192 
193  end if
194 
195 else
196  z_sl = specmap_zsl(ndata_specmap)
197 end if
198 
199 #endif
200 
201 ! ------ Time derivative of the sea level
202 
203 if ( z_sl_old > -999999.9_dp ) then
204  dzsl_dtau = (z_sl-z_sl_old)/dtime
205 else ! only dummy value for z_sl_old available
206  dzsl_dtau = 0.0_dp
207 end if
208 
209 ! ------ Minimum bedrock elevation for extent of marine ice
210 
211 #if MARGIN==2
212 
213 #if ( MARINE_ICE_CALVING==2 || MARINE_ICE_CALVING==3 )
214 z_mar = z_mar
215 #elif ( MARINE_ICE_CALVING==4 || MARINE_ICE_CALVING==5 )
216 z_mar = fact_z_mar*z_sl
217 #elif ( MARINE_ICE_CALVING==6 || MARINE_ICE_CALVING==7 )
218 if (z_sl >= -80.0_dp) then
219  z_mar = 2.5_dp*z_sl
220 else
221  z_mar = 10.25_dp*(z_sl+80.0_dp)-200.0_dp
222 end if
223 z_mar = fact_z_mar*z_mar
224 #endif
225 
226 #endif
227 
228 ! ------ Update of the mask according to the sea level
229 
230 ! ---- Check all sea and floating-ice points and their direct
231 ! neighbours
232 
233 do i=0, imax
234 do j=0, jmax
235  check_point(j,i) = .false.
236 end do
237 end do
238 
239 do i=1, imax-1
240 do j=1, jmax-1
241  if (maske(j,i).ge.2) then
242  check_point(j ,i ) = .true.
243  check_point(j ,i+1) = .true.
244  check_point(j ,i-1) = .true.
245  check_point(j+1,i ) = .true.
246  check_point(j-1,i ) = .true.
247  end if
248 end do
249 end do
250 
251 do i=1, imax-1
252 do j=1, jmax-1
253  if (check_point(j,i)) then
254  maske_neu(j,i) = mask_update(z_sl, i, j)
255  end if
256 end do
257 end do
258 
259 ! ---- Assign new values of the mask
260 
261 do i=1, imax-1
262 do j=1, jmax-1
263  if (check_point(j,i)) then
264  maske(j,i) = maske_neu(j,i)
265  end if
266 end do
267 end do
268 
269 !-------- Surface air temperature and accumulation-ablation function --------
270 
271 do i=0, imax
272 do j=0, jmax
273  dist(j,i) = sqrt( (xi(i)-x_hat)**2 + (eta(j)-y_hat)**2 )
274 end do
275 end do
276 
277 do i=0, imax
278 do j=0, jmax
279 
280  temp_s(j,i) = temp_min + s_t*dist(j,i)
281  temp_s(j,i) = temp_s(j,i) + delta_ts
282  ! Correction with temperature deviation delta_ts
283  if (temp_s(j,i).gt. -0.001_dp) temp_s(j,i) = -0.001_dp
284  ! Cut-off of positive air temperatures
285 
286  accum_present(j,i) = b_max
287  accum(j,i) = accum_present(j,i)
288 
289  as_perp(j,i) = s_b*(eld-dist(j,i))
290  if (as_perp(j,i).gt.accum(j,i)) as_perp(j,i) = accum(j,i)
291 
292 end do
293 end do
294 
295 end subroutine boundary
296 !