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), evap(j,i), runoff(j,i), as_perp(j,i), temp_s(j,i)
54 
55 integer(i4b) :: i, j, n
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) :: z_sle_present, z_sle_help
61 real(dp), dimension(0:JMAX,0:IMAX,0:12) :: precip
62 real(dp), dimension(0:JMAX,0:IMAX) :: &
63  snowfall, rainfall, melt, melt_star
64 real(dp), dimension(0:JMAX,0:IMAX,12) :: temp_mm
65 real(dp), dimension(0:JMAX,0:IMAX) :: temp_ma
66 real(dp), dimension(12) :: temp_mm_help
67 real(dp) :: temp_jja_help
68 real(dp), dimension(0:JMAX,0:IMAX) :: et
69 real(dp) :: gamma_t, temp_diff
70 real(dp) :: gamma_p, zs_thresh, &
71  temp_rain, temp_snow, &
72  inv_delta_temp_rain_snow, coeff(0:5), inv_sqrt2_s_stat, &
73  precip_fact, frac_solid
74 real(dp) :: s_stat, beta1, beta2, pmax, mu, lambda_lti, temp_lti
75 logical, dimension(0:JMAX,0:IMAX) :: check_point
76 
77 real(dp), parameter :: &
78  inv_twelve = 1.0_dp/12.0_dp, one_third = 1.0_dp/3.0_dp
79 
80 !-------- Initialization of intent(out) variables --------
81 
82 z_sl_old = z_sl
83 
84 delta_ts = 0.0_dp
85 glac_index = 0.0_dp
86 z_sl = 0.0_dp
87 dzsl_dtau = 0.0_dp
88 z_mar = 0.0_dp
89 
90 !-------- Surface-temperature deviation from present values --------
91 
92 #if TSURFACE==1
93 delta_ts = delta_ts0
94 ! ! Steady state with prescribed constant
95 ! ! air-temperature deviation
96 #elif TSURFACE==3
97 delta_ts = sine_amplit &
98  *cos(2.0_dp*pi*time/(sine_period*year_sec)) &
99  -sine_amplit
100 ! ! Sinusoidal air-temperature forcing
101 #elif TSURFACE==4
102 
103 ! ------ delta_ts from ice-core record
104 
105 if (time/year_sec.lt.real(grip_time_min,dp)) then
106  delta_ts = griptemp(0)
107 else if (time/year_sec.lt.real(grip_time_max,dp)) then
108 
109  i_kl = floor(((time/year_sec) &
110  -real(grip_time_min,dp))/real(grip_time_stp,dp))
111  i_kl = max(i_kl, 0)
112 
113  i_gr = ceiling(((time/year_sec) &
114  -real(grip_time_min,dp))/real(grip_time_stp,dp))
115  i_gr = min(i_gr, ndata_grip)
116 
117  if (i_kl.eq.i_gr) then
118 
119  delta_ts = griptemp(i_kl)
120 
121  else
122 
123  time_kl = (grip_time_min + i_kl*grip_time_stp) *year_sec
124  time_gr = (grip_time_min + i_gr*grip_time_stp) *year_sec
125 
126  delta_ts = griptemp(i_kl) &
127  +(griptemp(i_gr)-griptemp(i_kl)) &
128  *(time-time_kl)/(time_gr-time_kl)
129  ! linear interpolation of the ice-core data
130 
131  end if
132 
133 else
134  delta_ts = griptemp(ndata_grip)
135 end if
136 
137 delta_ts = delta_ts * grip_temp_fact
138 ! ! modification by constant factor
139 
140 !-------- Glacial index --------
141 
142 #elif TSURFACE==5
143 
144 if (time/year_sec < real(gi_time_min,dp)) then
145  glac_index = glacial_index(0)
146 else if (time/year_sec < real(gi_time_max,dp)) then
147 
148  i_kl = floor(((time/year_sec) &
149  -real(gi_time_min,dp))/real(gi_time_stp,dp))
150  i_kl = max(i_kl, 0)
151 
152  i_gr = ceiling(((time/year_sec) &
153  -real(gi_time_min,dp))/real(gi_time_stp,dp))
154  i_gr = min(i_gr, ndata_gi)
155 
156  if (i_kl == i_gr) then
157 
158  glac_index = glacial_index(i_kl)
159 
160  else
161 
162  time_kl = (gi_time_min + i_kl*gi_time_stp) *year_sec
163  time_gr = (gi_time_min + i_gr*gi_time_stp) *year_sec
164 
165  glac_index = glacial_index(i_kl) &
166  +(glacial_index(i_gr)-glacial_index(i_kl)) &
167  *(time-time_kl)/(time_gr-time_kl)
168  ! linear interpolation of the glacial-index data
169 
170  end if
171 
172 else
173  glac_index = glacial_index(ndata_gi)
174 end if
175 
176 #endif
177 
178 !-------- Sea level --------
179 
180 #if SEA_LEVEL==1
181 ! ------ constant sea level
182 z_sl = z_sl0
183 
184 #elif SEA_LEVEL==2
185 ! ------ saw-tooth-shaped palaeoclimatic sea-level forcing
186 
187 z_sl_min = -130.0_dp
188 
189 t1 = -250000.0_dp *year_sec
190 t2 = -140000.0_dp *year_sec
191 t3 = -125000.0_dp *year_sec
192 t4 = -21000.0_dp *year_sec
193 t5 = -8000.0_dp *year_sec
194 t6 = 0.0_dp *year_sec
195 
196 if (time.lt.t1) then
197  z_sl = 0.0_dp
198 else if (time.lt.t2) then
199  z_sl = z_sl_min*(time-t1)/(t2-t1)
200 else if (time.lt.t3) then
201  z_sl = -z_sl_min*(time-t3)/(t3-t2)
202 else if (time.lt.t4) then
203  z_sl = z_sl_min*(time-t3)/(t4-t3)
204 else if (time.lt.t5) then
205  z_sl = -z_sl_min*(time-t5)/(t5-t4)
206 else if (time.lt.t6) then
207  z_sl = 0.0_dp
208 else
209  z_sl = 0.0_dp
210 end if
211 
212 #elif SEA_LEVEL==3
213 
214 ! ------ z_sl from the SPECMAP record
215 
216 if (time/year_sec.lt.real(specmap_time_min,dp)) then
217  z_sl = specmap_zsl(0)
218 else if (time/year_sec.lt.real(specmap_time_max,dp)) then
219 
220  i_kl = floor(((time/year_sec) &
221  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
222  i_kl = max(i_kl, 0)
223 
224  i_gr = ceiling(((time/year_sec) &
225  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
226  i_gr = min(i_gr, ndata_specmap)
227 
228  if (i_kl.eq.i_gr) then
229 
230  z_sl = specmap_zsl(i_kl)
231 
232  else
233 
234  time_kl = (specmap_time_min + i_kl*specmap_time_stp) *year_sec
235  time_gr = (specmap_time_min + i_gr*specmap_time_stp) *year_sec
236 
237  z_sl = specmap_zsl(i_kl) &
238  +(specmap_zsl(i_gr)-specmap_zsl(i_kl)) &
239  *(time-time_kl)/(time_gr-time_kl)
240  ! linear interpolation of the sea-level data
241 
242  end if
243 
244 else
245  z_sl = specmap_zsl(ndata_specmap)
246 end if
247 
248 #endif
249 
250 ! ------ Time derivative of the sea level
251 
252 if ( z_sl_old > -999999.9_dp ) then
253  dzsl_dtau = (z_sl-z_sl_old)/dtime
254 else ! only dummy value for z_sl_old available
255  dzsl_dtau = 0.0_dp
256 end if
257 
258 ! ------ Minimum bedrock elevation for extent of marine ice
259 
260 #if MARGIN==2
261 
262 #if ( MARINE_ICE_CALVING==2 || MARINE_ICE_CALVING==3 )
263 z_mar = z_mar
264 #elif ( MARINE_ICE_CALVING==4 || MARINE_ICE_CALVING==5 )
265 z_mar = fact_z_mar*z_sl
266 #elif ( MARINE_ICE_CALVING==6 || MARINE_ICE_CALVING==7 )
267 if (z_sl >= -80.0_dp) then
268  z_mar = 2.5_dp*z_sl
269 else
270  z_mar = 10.25_dp*(z_sl+80.0_dp)-200.0_dp
271 end if
272 z_mar = fact_z_mar*z_mar
273 #endif
274 
275 #endif
276 
277 ! ------ Update of the mask according to the sea level
278 
279 ! ---- Check all sea and floating-ice points and their direct
280 ! neighbours
281 
282 do i=0, imax
283 do j=0, jmax
284  check_point(j,i) = .false.
285 end do
286 end do
287 
288 do i=1, imax-1
289 do j=1, jmax-1
290  if (maske(j,i).ge.2) then
291  check_point(j ,i ) = .true.
292  check_point(j ,i+1) = .true.
293  check_point(j ,i-1) = .true.
294  check_point(j+1,i ) = .true.
295  check_point(j-1,i ) = .true.
296  end if
297 end do
298 end do
299 
300 do i=1, imax-1
301 do j=1, jmax-1
302  if (check_point(j,i)) then
303  maske_neu(j,i) = mask_update(z_sl, i, j)
304  end if
305 end do
306 end do
307 
308 ! ---- Assign new values of the mask
309 
310 do i=1, imax-1
311 do j=1, jmax-1
312  if (check_point(j,i)) then
313  maske(j,i) = maske_neu(j,i)
314  end if
315 end do
316 end do
317 
318 !-------- Surface air temperatures --------
319 
320 gamma_t = -6.5e-03_dp ! atmospheric lapse rate
321 
322 do i=0, imax
323 do j=0, jmax
324 
325 #if (TSURFACE <= 4)
326 
327 ! ------ Correction of present monthly temperatures with elevation changes
328 ! and temperature deviation delta_ts
329 
330  temp_diff = gamma_t*(zs(j,i)-zs_ref(j,i)) + delta_ts
331 
332  do n=1, 12 ! month counter
333  temp_mm(j,i,n) = temp_mm_present(j,i,n) + temp_diff
334  end do
335 
336 #elif (TSURFACE == 5)
337 
338 ! ------ Correction of present monthly temperatures with LGM anomaly and
339 ! glacial index as well as elevation changes
340 
341  temp_diff = gamma_t*(zs(j,i)-zs_ref(j,i))
342 
343  do n=1, 12 ! month counter
344  temp_mm(j,i,n) = temp_mm_present(j,i,n) &
345  + glac_index*temp_mm_lgm_anom(j,i,n) &
346  + temp_diff
347  end do
348 
349 #endif
350 
351 ! ------ Mean annual air temperature
352 
353  temp_ma(j,i) = 0.0_dp ! initialisation value
354 
355  do n=1, 12 ! month counter
356  temp_ma(j,i) = temp_ma(j,i) + temp_mm(j,i,n)*inv_twelve
357  end do
358 
359 end do
360 end do
361 
362 !-------- Accumulation-ablation function as_perp --------
363 
364 #if (ELEV_DESERT == 1)
365 
366 gamma_p = gamma_p*1.0e-03_dp ! Precipitation lapse rate
367  ! for elevation desertification, in m^(-1)
368 zs_thresh = zs_thresh ! Elevation threshold, in m
369 
370 #endif
371 
372 #if (SOLID_PRECIP == 1) /* Marsiat (1994) */
373 
374 temp_rain = 7.0_dp ! Threshold monthly mean temperature for
375  ! precipitation = 100% rain, in deg C
376 temp_snow = -10.0_dp ! Threshold monthly mean temperature for &
377  ! precipitation = 100% snow, in deg C
378 
379 inv_delta_temp_rain_snow = 1.0_dp/(temp_rain-temp_snow)
380 
381 #elif (SOLID_PRECIP == 2) /* Bales et al. (2009) */
382 
383 temp_rain = 7.2_dp ! Threshold monthly mean temperature for &
384  ! precipitation = 100% rain, in deg C
385 temp_snow = -11.6_dp ! Threshold monthly mean temperature for &
386  ! precipitation = 100% snow, in deg C
387 
388 coeff(0) = 5.4714e-01_dp ! Coefficients
389 coeff(1) = -9.1603e-02_dp ! of
390 coeff(2) = -3.314e-03_dp ! the
391 coeff(3) = 4.66e-04_dp ! fifth-order
392 coeff(4) = 3.8e-05_dp ! polynomial
393 coeff(5) = 6.0e-07_dp ! fit
394 
395 #elif (SOLID_PRECIP == 3) /* Huybrechts and de Wolde (1999) */
396 
397 temp_rain = 2.0_dp ! Threshold instantaneous temperature for &
398  ! precipitation = 100% rain, in deg C
399 temp_snow = temp_rain ! Threshold instantaneous temperature for &
400  ! precipitation = 100% snow, in deg C
401 
402 s_stat = s_stat_0 ! Standard deviation of the air termperature
403  ! (same parameter as in the PDD model)
404 
405 inv_sqrt2_s_stat = 1.0_dp/(sqrt(2.0_dp)*s_stat)
406 
407 #endif
408 
409 #if (ABLSURFACE==1 || ABLSURFACE==2)
410 
411 s_stat = s_stat_0
412 beta1 = beta1_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
413  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
414 beta2 = beta2_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
415  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
416 pmax = pmax_0
417 mu = mu_0 *(1000.0_dp*86400.0_dp)*(rho/rho_w)
418  ! (d*deg C)/(mm WE) --> (s*deg C)/(m IE)
419 
420 #elif (ABLSURFACE==3)
421 
422 lambda_lti = lambda_lti *(0.001_dp/year_sec)*(rho_w/rho)
423  ! (mm WE)/(a*deg C) --> (m IE)/(s*deg C)
424 temp_lti = temp_lti
425 mu = 0.0_dp ! no superimposed ice considered
426 
427 #endif
428 
429 do i=0, imax
430 do j=0, jmax
431 
432 ! ------ Accumulation
433 
434 #if (ACCSURFACE <= 3)
435 
436 ! ---- Elevation desertification of precipitation
437 
438 #if (ELEV_DESERT == 0)
439 
440  precip_fact = 1.0_dp ! no elevation desertification
441 
442 #elif (ELEV_DESERT == 1)
443 
444  if (zs_ref(j,i) < zs_thresh) then
445  precip_fact &
446  = exp(gamma_p*(max(zs(j,i),zs_thresh)-zs_thresh))
447  else
448  precip_fact &
449  = exp(gamma_p*(max(zs(j,i),zs_thresh)-zs_ref(j,i)))
450  end if
451 
452 #else
453  stop ' boundary: Parameter ELEV_DESERT must be either 0 or 1!'
454 #endif
455 
456  do n=1, 12 ! month counter
457  precip(j,i,n) = precip_present(j,i,n)*precip_fact
458  end do
459 
460 #endif
461 
462 ! ---- Precipitation change related to changing climate
463 
464 #if ACCSURFACE==1
465  precip_fact = accfact
466 #elif ACCSURFACE==2
467  precip_fact = 1.0_dp + gamma_s*delta_ts
468 #elif ACCSURFACE==3
469  precip_fact = exp(gamma_s*delta_ts)
470 #endif
471 
472 #if (ACCSURFACE <= 3)
473 
474  precip(j,i,0) = 0.0_dp ! initialisation value for mean annual precip
475 
476  do n=1, 12 ! month counter
477  precip(j,i,n) = precip(j,i,n)*precip_fact ! monthly precip
478  precip(j,i,0) = precip(j,i,0) + precip(j,i,n)*inv_twelve
479  ! mean annual precip
480  end do
481 
482 #elif (ACCSURFACE == 5)
483 
484  precip(j,i,0) = 0.0_dp ! initialisation value for mean annual precip
485 
486  do n=1, 12 ! month counter
487 
488 #if (PRECIP_ANOM_INTERPOL==1)
489  precip_fact = 1.0_dp-glac_index+glac_index*precip_lgm_anom(j,i,n)
490  ! interpolation with a linear function
491 #elif (PRECIP_ANOM_INTERPOL==2)
492  precip_fact = exp(-glac_index*gamma_precip_lgm_anom(j,i,n))
493  ! interpolation with an exponential function
494 #endif
495 
496  precip(j,i,n) = precip_present(j,i,n)*precip_fact ! monthly precip
497  precip(j,i,0) = precip(j,i,0) + precip(j,i,n)*inv_twelve
498  ! mean annual precip
499  end do
500 
501 #endif
502 
503 ! ---- Annual accumulation, snowfall and rainfall rates
504 
505  accum(j,i) = precip(j,i,0)
506 
507  snowfall(j,i) = 0.0_dp ! initialisation value
508 
509  do n=1, 12 ! month counter
510 
511 #if (SOLID_PRECIP == 1) /* Marsiat (1994) */
512 
513  if (temp_mm(j,i,n) >= temp_rain) then
514  frac_solid = 0.0_dp
515  else if (temp_mm(j,i,n) <= temp_snow) then
516  frac_solid = 1.0_dp
517  else
518  frac_solid = (temp_rain-temp_mm(j,i,n))*inv_delta_temp_rain_snow
519  end if
520 
521 #elif (SOLID_PRECIP == 2) /* Bales et al. (2009) */
522 
523  if (temp_mm(j,i,n) >= temp_rain) then
524  frac_solid = 0.0_dp
525  else if (temp_mm(j,i,n) <= temp_snow) then
526  frac_solid = 1.0_dp
527  else
528  frac_solid = coeff(0) + temp_mm(j,i,n) * ( coeff(1) &
529  + temp_mm(j,i,n) * ( coeff(2) &
530  + temp_mm(j,i,n) * ( coeff(3) &
531  + temp_mm(j,i,n) * ( coeff(4) &
532  + temp_mm(j,i,n) * coeff(5) ) ) ) )
533  ! evaluation of 5th-order polynomial by Horner scheme
534  end if
535 
536 #elif (SOLID_PRECIP == 3) /* Huybrechts and de Wolde (1999) */
537 
538  frac_solid = 1.0_dp &
539  - 0.5_dp*erfcc((temp_rain-temp_mm(j,i,n))*inv_sqrt2_s_stat)
540 
541 #endif
542 
543  snowfall(j,i) = snowfall(j,i) + precip(j,i,n)*frac_solid*inv_twelve
544 
545  end do
546 
547  rainfall(j,i) = precip(j,i,0) - snowfall(j,i)
548 
549  if (snowfall(j,i) < 0.0_dp) snowfall(j,i) = 0.0_dp ! correction of
550  if (rainfall(j,i) < 0.0_dp) rainfall(j,i) = 0.0_dp ! negative values
551 
552 ! ------ Ablation
553 
554 ! ---- Runoff
555 
556 #if (ABLSURFACE==1 || ABLSURFACE==2)
557 
558 ! -- Temperature excess ET
559 
560  do n=1, 12 ! month counter
561  temp_mm_help(n) = temp_mm(j,i,n)
562  end do
563 
564  call pdd(temp_mm_help, s_stat, et(j,i))
565 
566 ! -- Formation rate of superimposed ice (melt_star), melt rate (melt)
567 ! and runoff rate (runoff)
568 
569 #if (ABLSURFACE==1)
570 
571  if ((beta1*et(j,i)) <= (pmax*snowfall(j,i))) then
572  melt_star(j,i) = beta1*et(j,i)
573  melt(j,i) = 0.0_dp
574  runoff(j,i) = melt(j,i)+rainfall(j,i)
575  else
576  melt_star(j,i) = pmax*snowfall(j,i)
577  melt(j,i) = beta2*(et(j,i)-melt_star(j,i)/beta1)
578  runoff(j,i) = melt(j,i)+rainfall(j,i)
579  end if
580 
581 #elif (ABLSURFACE==2)
582 
583  if ( rainfall(j,i) <= (pmax*snowfall(j,i)) ) then
584 
585  if ( (rainfall(j,i)+beta1*et(j,i)) <= (pmax*snowfall(j,i)) ) then
586  melt_star(j,i) = rainfall(j,i)+beta1*et(j,i)
587  melt(j,i) = 0.0_dp
588  runoff(j,i) = melt(j,i)
589  else
590  melt_star(j,i) = pmax*snowfall(j,i)
591  melt(j,i) = beta2 &
592  *(et(j,i)-(melt_star(j,i)-rainfall(j,i))/beta1)
593  runoff(j,i) = melt(j,i)
594  end if
595 
596  else
597 
598  melt_star(j,i) = pmax*snowfall(j,i)
599  melt(j,i) = beta2*et(j,i)
600  runoff(j,i) = melt(j,i) + rainfall(j,i)-pmax*snowfall(j,i)
601 
602  end if
603 
604 #endif
605 
606 #elif (ABLSURFACE==3)
607 
608  temp_jja_help = one_third*(temp_mm(j,i,6)+temp_mm(j,i,7)+temp_mm(j,i,8))
609 
610  melt_star(j,i) = 0.0_dp ! no superimposed ice considered
611  melt(j,i) = lambda_lti*max((temp_jja_help-temp_lti), 0.0_dp)
612  runoff(j,i) = melt(j,i) + rainfall(j,i)
613 
614 #endif
615 
616 ! ---- Evaporation
617 
618  evap(j,i) = 0.0_dp
619 
620 ! ------ Accumulation minus ablation
621 
622  as_perp(j,i) = accum(j,i) - evap(j,i) - runoff(j,i)
623 
624 ! ------ Ice-surface temperature (10-m firn temperature) temp_s,
625 ! including empirical firn-warming correction due to
626 ! refreezing meltwater when superimposed ice is formed
627 
628  if (melt_star(j,i).ge.melt(j,i)) then
629  temp_s(j,i) = temp_ma(j,i) &
630  +mu*(melt_star(j,i)-melt(j,i))
631  else
632  temp_s(j,i) = temp_ma(j,i)
633  end if
634 
635  if (temp_s(j,i) > -0.001_dp) temp_s(j,i) = -0.001_dp
636  ! Cut-off of positive air temperatures
637 
638 end do
639 end do
640 
641 end subroutine boundary
642 !