SICOPOLIS V3.2
 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-2016 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 use sico_vars
43 #if (NETCDF > 1)
44 use netcdf
45 use nc_check
46 #endif
47 
48 implicit none
49 
50 real(dp), intent(in) :: time, dtime, dxi, deta
51 
52 real(dp), intent(out) :: delta_ts, glac_index, dzsl_dtau, z_mar
53 real(dp), intent(inout) :: z_sl
54 
55 ! Further return variables
56 ! (defined as global variables in module sico_variables):
57 !
58 ! accum(j,i), evap(j,i), runoff(j,i), as_perp(j,i), temp_s(j,i)
59 
60 integer(i4b) :: i, j, n
61 integer(i4b) :: i_gr, i_kl
62 integer(i4b) :: nrec_temp_precip
63 integer(i4b), save :: nrec_temp_precip_save = -1
64 real(dp) :: z_sl_old
65 real(dp) :: z_sl_min, t1, t2, t3, t4, t5, t6
66 real(dp) :: time_gr, time_kl
67 real(dp) :: z_sle_present, z_sle_help
68 real(dp), dimension(0:JMAX,0:IMAX,0:12) :: precip
69 real(dp), dimension(0:JMAX,0:IMAX) :: &
70  snowfall, rainfall, melt, melt_star
71 real(dp), dimension(0:JMAX,0:IMAX,12) :: temp_mm
72 real(dp), dimension(0:JMAX,0:IMAX) :: temp_ma, temp_ampl
73 real(dp), dimension(0:JMAX,0:IMAX) :: temp_ma_anom, temp_mj_anom, precip_ma_anom
74 real(dp), dimension(0:IMAX,0:JMAX), save :: temp_ma_anom_tra, temp_mj_anom_tra, &
75  precip_ma_anom_tra
76  ! These arrays are transposed with respect
77  ! to the usual SICOPOLIS convention;
78  ! they have i as first and j as second index
79 real(dp), dimension(12) :: temp_mm_help
80 real(dp) :: temp_jja_help
81 real(dp), dimension(0:JMAX,0:IMAX) :: et
82 real(dp) :: theta_ma, c_ma, gamma_ma, &
83  theta_ma_1, c_ma_1, gamma_ma_1, &
84  theta_ma_2, c_ma_2, gamma_ma_2, &
85  theta_ma_3, c_ma_3, gamma_ma_3, &
86  zs_sep_1, zs_sep_2, &
87  theta_mj, c_mj, gamma_mj
88 real(dp) :: sine_factor
89 real(dp) :: gamma_p, zs_thresh, &
90  alpha_p, beta_p, temp_0, alpha_t, beta_t, &
91  temp_inv, temp_inv_present, &
92  temp_rain, temp_snow, &
93  inv_delta_temp_rain_snow, coeff(0:5), inv_sqrt2_s_stat, &
94  precip_fact, frac_solid
95 real(dp) :: s_stat, beta1, beta2, pmax, mu, lambda_lti, temp_lti
96 logical, dimension(0:JMAX,0:IMAX) :: check_point
97 
98 #if (NETCDF > 1)
99 integer(i4b) :: ncv
100 ! ncv: Variable ID
101 integer(i4b) :: nc3cor(3)
102 ! nc3cor(3): Corner of a 3-d array
103 integer(i4b) :: nc3cnt(3)
104 ! nc3cnt(3): Count of a 3-d array
105 #endif
106 
107 real(dp), parameter :: &
108  inv_twelve = 1.0_dp/12.0_dp, one_third = 1.0_dp/3.0_dp
109 
110 character(len=64), parameter :: thisroutine = 'boundary'
111 
112 !-------- Initialization of intent(out) variables --------
113 
114 z_sl_old = z_sl
115 
116 delta_ts = 0.0_dp
117 glac_index = 0.0_dp
118 z_sl = 0.0_dp
119 dzsl_dtau = 0.0_dp
120 z_mar = 0.0_dp
121 
122 !-------- Surface-temperature deviation from present values --------
123 
124 #if TSURFACE==1
125 delta_ts = delta_ts0
126 ! ! Steady state with prescribed constant
127 ! ! air-temperature deviation
128 #elif TSURFACE==3
129 delta_ts = sine_amplit &
130  *cos(2.0_dp*pi*time/(sine_period*year_sec)) &
131  -sine_amplit
132 ! ! Sinusoidal air-temperature forcing
133 #elif TSURFACE==4
134 
135 ! ------ delta_ts from ice-core record
136 
137 if (time/year_sec < real(grip_time_min,dp)) then
138  delta_ts = griptemp(0)
139 else if (time/year_sec < real(grip_time_max,dp)) then
140 
141  i_kl = floor(((time/year_sec) &
142  -real(grip_time_min,dp))/real(grip_time_stp,dp))
143  i_kl = max(i_kl, 0)
144 
145  i_gr = ceiling(((time/year_sec) &
146  -real(grip_time_min,dp))/real(grip_time_stp,dp))
147  i_gr = min(i_gr, ndata_grip)
148 
149  if (i_kl == i_gr) then
150 
151  delta_ts = griptemp(i_kl)
152 
153  else
154 
155  time_kl = (grip_time_min + i_kl*grip_time_stp) *year_sec
156  time_gr = (grip_time_min + i_gr*grip_time_stp) *year_sec
157 
158  delta_ts = griptemp(i_kl) &
159  +(griptemp(i_gr)-griptemp(i_kl)) &
160  *(time-time_kl)/(time_gr-time_kl)
161  ! linear interpolation of the ice-core data
162 
163  end if
164 
165 else
166  delta_ts = griptemp(ndata_grip)
167 end if
168 
169 delta_ts = delta_ts * grip_temp_fact
170 ! ! modification by constant factor
171 
172 !-------- Glacial index --------
173 
174 #elif TSURFACE==5
175 
176 if (time/year_sec < real(gi_time_min,dp)) then
177  glac_index = glacial_index(0)
178 else if (time/year_sec < real(gi_time_max,dp)) then
179 
180  i_kl = floor(((time/year_sec) &
181  -real(gi_time_min,dp))/real(gi_time_stp,dp))
182  i_kl = max(i_kl, 0)
183 
184  i_gr = ceiling(((time/year_sec) &
185  -real(gi_time_min,dp))/real(gi_time_stp,dp))
186  i_gr = min(i_gr, ndata_gi)
187 
188  if (i_kl == i_gr) then
189 
190  glac_index = glacial_index(i_kl)
191 
192  else
193 
194  time_kl = (gi_time_min + i_kl*gi_time_stp) *year_sec
195  time_gr = (gi_time_min + i_gr*gi_time_stp) *year_sec
196 
197  glac_index = glacial_index(i_kl) &
198  +(glacial_index(i_gr)-glacial_index(i_kl)) &
199  *(time-time_kl)/(time_gr-time_kl)
200  ! linear interpolation of the glacial-index data
201 
202  end if
203 
204 else
205  glac_index = glacial_index(ndata_gi)
206 end if
207 
208 !-------- Reading of surface temperature and precipitation
209 ! directly from data --------
210 
211 #elif ( (TSURFACE==6) && (ACCSURFACE==6) )
212 
213 if (time/year_sec <= temp_precip_time_min) then
214  nrec_temp_precip = 0
215 else if (time/year_sec < temp_precip_time_max) then
216  nrec_temp_precip = nint( ((time/year_sec)-temp_precip_time_min) &
217  / temp_precip_time_stp )
218 else
219  nrec_temp_precip = ndata_temp_precip
220 end if
221 
222 if ( nrec_temp_precip < 0 ) then
223  stop ' boundary: nrec_temp_precip < 0, not allowed!'
224 else if ( nrec_temp_precip > ndata_temp_precip ) then
225  stop ' boundary: nrec_temp_precip > ndata_temp_precip, not allowed!'
226 end if
227 
228 if ( nrec_temp_precip /= nrec_temp_precip_save ) then
229 
230  nc3cor(1) = 1
231  nc3cor(2) = 1
232  nc3cor(3) = nrec_temp_precip + 1
233  nc3cnt(1) = imax + 1
234  nc3cnt(2) = jmax + 1
235  nc3cnt(3) = 1
236 
237  call check( nf90_inq_varid(ncid_temp_precip, 'annualtemp_anom', ncv), &
238  thisroutine )
239  call check( nf90_get_var(ncid_temp_precip, ncv, temp_ma_anom_tra, &
240  start=nc3cor, count=nc3cnt), thisroutine )
241 
242  call check( nf90_inq_varid(ncid_temp_precip, 'januarytemp_anom', ncv), &
243  thisroutine )
244  call check( nf90_get_var(ncid_temp_precip, ncv, temp_mj_anom_tra, &
245  start=nc3cor, count=nc3cnt), thisroutine )
246 
247  call check( nf90_inq_varid(ncid_temp_precip, 'precipitation_anom', ncv), &
248  thisroutine )
249  call check( nf90_get_var(ncid_temp_precip, ncv, precip_ma_anom_tra, &
250  start=nc3cor, count=nc3cnt), thisroutine )
251 
252 end if
253 
254 temp_ma_anom = transpose(temp_ma_anom_tra)
255 temp_mj_anom = transpose(temp_mj_anom_tra)
256 precip_ma_anom = transpose(precip_ma_anom_tra) *(1.0e-03_dp/year_sec)*(rho_w/rho)
257  ! mm/a water equiv. -> m/s ice equiv.
258 
259 nrec_temp_precip_save = nrec_temp_precip
260 
261 #endif
262 
263 !-------- Sea level --------
264 
265 #if SEA_LEVEL==1
266 ! ------ constant sea level
267 z_sl = z_sl0
268 
269 #elif SEA_LEVEL==2
270 ! ------ saw-tooth-shaped palaeoclimatic sea-level forcing
271 
272 z_sl_min = -130.0_dp
273 
274 t1 = -250000.0_dp *year_sec
275 t2 = -140000.0_dp *year_sec
276 t3 = -125000.0_dp *year_sec
277 t4 = -21000.0_dp *year_sec
278 t5 = -8000.0_dp *year_sec
279 t6 = 0.0_dp *year_sec
280 
281 if (time < t1) then
282  z_sl = 0.0_dp
283 else if (time < t2) then
284  z_sl = z_sl_min*(time-t1)/(t2-t1)
285 else if (time < t3) then
286  z_sl = -z_sl_min*(time-t3)/(t3-t2)
287 else if (time < t4) then
288  z_sl = z_sl_min*(time-t3)/(t4-t3)
289 else if (time < t5) then
290  z_sl = -z_sl_min*(time-t5)/(t5-t4)
291 else if (time < t6) then
292  z_sl = 0.0_dp
293 else
294  z_sl = 0.0_dp
295 end if
296 
297 #elif SEA_LEVEL==3
298 
299 ! ------ z_sl from the SPECMAP record
300 
301 if (time/year_sec < real(specmap_time_min,dp)) then
302  z_sl = specmap_zsl(0)
303 else if (time/year_sec < real(specmap_time_max,dp)) then
304 
305  i_kl = floor(((time/year_sec) &
306  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
307  i_kl = max(i_kl, 0)
308 
309  i_gr = ceiling(((time/year_sec) &
310  -real(specmap_time_min,dp))/real(specmap_time_stp,dp))
311  i_gr = min(i_gr, ndata_specmap)
312 
313  if (i_kl == i_gr) then
314 
315  z_sl = specmap_zsl(i_kl)
316 
317  else
318 
319  time_kl = (specmap_time_min + i_kl*specmap_time_stp) *year_sec
320  time_gr = (specmap_time_min + i_gr*specmap_time_stp) *year_sec
321 
322  z_sl = specmap_zsl(i_kl) &
323  +(specmap_zsl(i_gr)-specmap_zsl(i_kl)) &
324  *(time-time_kl)/(time_gr-time_kl)
325  ! linear interpolation of the sea-level data
326 
327  end if
328 
329 else
330  z_sl = specmap_zsl(ndata_specmap)
331 end if
332 
333 #endif
334 
335 ! ------ Time derivative of the sea level
336 
337 if ( z_sl_old > -999999.9_dp ) then
338  dzsl_dtau = (z_sl-z_sl_old)/dtime
339 else ! only dummy value for z_sl_old available
340  dzsl_dtau = 0.0_dp
341 end if
342 
343 ! ------ Minimum bedrock elevation for extent of marine ice
344 
345 #if MARGIN==2
346 
347 #if ( MARINE_ICE_CALVING==2 || MARINE_ICE_CALVING==3 )
348 z_mar = z_mar
349 #elif ( MARINE_ICE_CALVING==4 || MARINE_ICE_CALVING==5 )
350 z_mar = fact_z_mar*z_sl
351 #elif ( MARINE_ICE_CALVING==6 || MARINE_ICE_CALVING==7 )
352 if (z_sl >= -80.0_dp) then
353  z_mar = 2.5_dp*z_sl
354 else
355  z_mar = 10.25_dp*(z_sl+80.0_dp)-200.0_dp
356 end if
357 z_mar = fact_z_mar*z_mar
358 #endif
359 
360 #endif
361 
362 ! ------ Update of the mask according to the sea level
363 
364 ! ---- Check all sea and floating-ice points and their direct
365 ! neighbours
366 
367 do i=0, imax
368 do j=0, jmax
369  check_point(j,i) = .false.
370 end do
371 end do
372 
373 do i=1, imax-1
374 do j=1, jmax-1
375  if (maske(j,i) >= 2) then
376  check_point(j ,i ) = .true.
377  check_point(j ,i+1) = .true.
378  check_point(j ,i-1) = .true.
379  check_point(j+1,i ) = .true.
380  check_point(j-1,i ) = .true.
381  end if
382 end do
383 end do
384 
385 do i=1, imax-1
386 do j=1, jmax-1
387  if (check_point(j,i)) then
388  maske_neu(j,i) = mask_update(z_sl, i, j)
389  end if
390 end do
391 end do
392 
393 ! ---- Assign new values of the mask
394 
395 do i=1, imax-1
396 do j=1, jmax-1
397  if (check_point(j,i)) then
398  maske(j,i) = maske_neu(j,i)
399  end if
400 end do
401 end do
402 
403 !-------- Surface air temperatures --------
404 
405 #if TEMP_PRESENT_PARA == 1 /* Parameterisation by Fortuin and Oerlemans */
406  ! (1990) for the whole ice sheet
407 
408 theta_ma = 34.461_dp
409 gamma_ma = -9.14e-03_dp
410 c_ma = -0.688_dp
411 
412 theta_mj = 16.81_dp
413 gamma_mj = -6.92e-03_dp
414 c_mj = -0.27973_dp
415 
416 #elif TEMP_PRESENT_PARA == 2 /* Parameterisation by Fortuin and Oerlemans */
417  ! (1990), separately for three different
418  ! elevation ranges
419 
420 zs_sep_1 = 200.0_dp
421 zs_sep_2 = 1500.0_dp
422 
423 theta_ma_1 = 49.642_dp
424 gamma_ma_1 = 0.0_dp
425 c_ma_1 = -0.943_dp
426 
427 theta_ma_2 = 36.689_dp
428 gamma_ma_2 = -5.102e-03_dp
429 c_ma_2 = -0.725_dp
430 
431 theta_ma_3 = 7.405_dp
432 gamma_ma_3 = -14.285e-03_dp
433 c_ma_3 = -0.180_dp
434 
435 theta_mj = 16.81_dp
436 gamma_mj = -6.92e-03_dp
437 c_mj = -0.27937_dp
438 
439 #else
440 
441 stop ' boundary: Parameter TEMP_PRESENT_PARA must be either 1 or 2!'
442 
443 #endif
444 
445 do i=0, imax
446 do j=0, jmax
447 
448 ! ------ Present-day mean-annual air temperature
449 
450 #if TEMP_PRESENT_PARA == 1
451  temp_ma_present(j,i) = theta_ma + gamma_ma*zs(j,i) &
452  + c_ma*abs(phi(j,i))*pi_180_inv
453 #elif TEMP_PRESENT_PARA == 2
454  if ( zs(j,i) <= zs_sep_1 ) then
455  temp_ma_present(j,i) = theta_ma_1 + gamma_ma_1*zs(j,i) &
456  + c_ma_1*abs(phi(j,i))*pi_180_inv
457  else if ( zs(j,i) <= zs_sep_2 ) then
458  temp_ma_present(j,i) = theta_ma_2 + gamma_ma_2*zs(j,i) &
459  + c_ma_2*abs(phi(j,i))*pi_180_inv
460  else
461  temp_ma_present(j,i) = theta_ma_3 + gamma_ma_3*zs(j,i) &
462  + c_ma_3*abs(phi(j,i))*pi_180_inv
463  end if
464 #endif
465 
466 ! ------ Present-day mean-January (summer) air temperature
467 
468  temp_mj_present(j,i) = theta_mj + gamma_mj*zs(j,i) &
469  + c_mj*abs(phi(j,i))*pi_180_inv
470 
471 #if (TSURFACE <= 4)
472 
473 ! ------ Correction with deviation delta_ts
474 
475  temp_ma(j,i) = temp_ma_present(j,i) + delta_ts
476  temp_mm(j,i,7) = temp_mj_present(j,i) + delta_ts
477 
478 #elif (TSURFACE == 5)
479 
480 ! ------ Correction with LGM anomaly and glacial index
481 
482  temp_ma(j,i) = temp_ma_present(j,i) + glac_index*temp_ma_lgm_anom(j,i)
483  temp_mm(j,i,7) = temp_mj_present(j,i) + glac_index*temp_mj_lgm_anom(j,i)
484 
485 #elif (TSURFACE == 6)
486 
487 ! ------ Mean-annual and mean-January (summer) surface temperatures
488 ! from data read above --------
489 
490  temp_ma(j,i) = temp_ma_present(j,i) + temp_ma_anom_fact*temp_ma_anom(j,i)
491  temp_mm(j,i,7) = temp_mj_present(j,i) + temp_mj_anom_fact*temp_mj_anom(j,i)
492 
493 #endif
494 
495 ! ------ Amplitude of the annual cycle
496 
497  temp_ampl(j,i) = temp_mm(j,i,7) - temp_ma(j,i)
498 
499  if (temp_ampl(j,i) < eps) then
500  temp_ampl(j,i) = eps ! Correction of amplitude, if required
501  end if
502 
503 end do
504 end do
505 
506 ! ------ Monthly temperatures
507 
508 do n=1, 12 ! month counter
509 
510  sine_factor = sin((real(n,dp)-4.0_dp)*pi/6.0_dp)
511 
512  do i=0, imax
513  do j=0, jmax
514  temp_mm(j,i,n) = temp_ma(j,i) + sine_factor*temp_ampl(j,i)
515  end do
516  end do
517 
518 end do
519 
520 !-------- Accumulation-ablation function as_perp --------
521 
522 #if (ACCSURFACE <= 3)
523 
524 #if (ELEV_DESERT == 1)
525 
526 gamma_p = gamma_p*1.0e-03_dp ! Precipitation lapse rate
527  ! for elevation desertification, in m^(-1)
528 zs_thresh = zs_thresh ! Elevation threshold, in m
529 
530 #endif
531 
532 #elif (ACCSURFACE == 4)
533 
534 alpha_p = 22.47_dp
535 beta_p = 0.046_dp
536 temp_0 = 273.15_dp
537 alpha_t = 0.67_dp
538 beta_t = 88.9_dp
539 
540 #endif
541 
542 #if (SOLID_PRECIP == 1) /* Marsiat (1994) */
543 
544 temp_rain = 7.0_dp ! Threshold monthly mean temperature for
545  ! precipitation = 100% rain, in deg C
546 temp_snow = -10.0_dp ! Threshold monthly mean temperature for &
547  ! precipitation = 100% snow, in deg C
548 
549 inv_delta_temp_rain_snow = 1.0_dp/(temp_rain-temp_snow)
550 
551 #elif (SOLID_PRECIP == 2) /* Bales et al. (2009) */
552 
553 temp_rain = 7.2_dp ! Threshold monthly mean temperature for &
554  ! precipitation = 100% rain, in deg C
555 temp_snow = -11.6_dp ! Threshold monthly mean temperature for &
556  ! precipitation = 100% snow, in deg C
557 
558 coeff(0) = 5.4714e-01_dp ! Coefficients
559 coeff(1) = -9.1603e-02_dp ! of
560 coeff(2) = -3.314e-03_dp ! the
561 coeff(3) = 4.66e-04_dp ! fifth-order
562 coeff(4) = 3.8e-05_dp ! polynomial
563 coeff(5) = 6.0e-07_dp ! fit
564 
565 #elif (SOLID_PRECIP == 3) /* Huybrechts and de Wolde (1999) */
566 
567 temp_rain = 2.0_dp ! Threshold instantaneous temperature for &
568  ! precipitation = 100% rain, in deg C
569 temp_snow = temp_rain ! Threshold instantaneous temperature for &
570  ! precipitation = 100% snow, in deg C
571 
572 s_stat = s_stat_0 ! Standard deviation of the air termperature
573  ! (same parameter as in the PDD model)
574 
575 inv_sqrt2_s_stat = 1.0_dp/(sqrt(2.0_dp)*s_stat)
576 
577 #endif
578 
579 #if (ABLSURFACE==1 || ABLSURFACE==2)
580 
581 s_stat = s_stat_0
582 beta1 = beta1_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
583  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
584 beta2 = beta2_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
585  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
586 pmax = pmax_0
587 mu = mu_0 *(1000.0_dp*86400.0_dp)*(rho/rho_w)
588  ! (d*deg C)/(mm WE) --> (s*deg C)/(m IE)
589 
590 #elif (ABLSURFACE==3)
591 
592 lambda_lti = lambda_lti *(0.001_dp/year_sec)*(rho_w/rho)
593  ! (mm WE)/(a*deg C) --> (m IE)/(s*deg C)
594 temp_lti = temp_lti
595 mu = 0.0_dp ! no superimposed ice considered
596 
597 #endif
598 
599 do i=0, imax
600 do j=0, jmax
601 
602 ! ------ Accumulation
603 
604 #if (ACCSURFACE <= 3)
605 
606 ! ---- Elevation desertification of precipitation
607 
608 #if (ELEV_DESERT == 0)
609 
610  precip_fact = 1.0_dp ! no elevation desertification
611 
612 #elif (ELEV_DESERT == 1)
613 
614  if (zs_ref(j,i) < zs_thresh) then
615  precip_fact &
616  = exp(gamma_p*(max(zs(j,i),zs_thresh)-zs_thresh))
617  else
618  precip_fact &
619  = exp(gamma_p*(max(zs(j,i),zs_thresh)-zs_ref(j,i)))
620  end if
621 
622 #else
623  stop ' boundary: Parameter ELEV_DESERT must be either 0 or 1!'
624 #endif
625 
626  do n=1, 12 ! month counter
627  precip(j,i,n) = precip_present(j,i,n)*precip_fact
628  end do
629 
630 #endif
631 
632 ! ---- Precipitation change related to changing climate
633 
634 #if ACCSURFACE==1
635  precip_fact = accfact
636 #elif ACCSURFACE==2
637  precip_fact = 1.0_dp + gamma_s*delta_ts
638 #elif ACCSURFACE==3
639  precip_fact = exp(gamma_s*delta_ts)
640 #endif
641 
642 #if (ACCSURFACE <= 3)
643 
644  precip(j,i,0) = 0.0_dp ! initialisation value for mean annual precip
645 
646  do n=1, 12 ! month counter
647  precip(j,i,n) = precip(j,i,n)*precip_fact ! monthly precip
648  precip(j,i,0) = precip(j,i,0) + precip(j,i,n)*inv_twelve
649  ! mean annual precip
650  end do
651 
652 #elif (ACCSURFACE == 4)
653 
654  precip(j,i,0) = 0.0_dp ! initialisation value for mean annual precip
655 
656  temp_inv = alpha_t * (temp_ma(j,i)+temp_0) + beta_t ! in K
657  temp_inv_present = alpha_t * (temp_ma_present(j,i)+temp_0) + beta_t ! in K
658 
659  precip_fact = exp(alpha_p*(temp_0/temp_inv_present-temp_0/temp_inv)) &
660  *(temp_inv_present/temp_inv)**2 &
661  *(1.0_dp+beta_p*(temp_inv-temp_inv_present))
662 
663  do n=1, 12 ! month counter
664  precip(j,i,n) = precip_present(j,i,n)*precip_fact ! monthly precip
665  precip(j,i,0) = precip(j,i,0) + precip(j,i,n)*inv_twelve
666  ! mean annual precip
667  end do
668 
669 #elif (ACCSURFACE == 5)
670 
671  precip(j,i,0) = 0.0_dp ! initialisation value for mean annual precip
672 
673  do n=1, 12 ! month counter
674 
675 #if (PRECIP_ANOM_INTERPOL==1)
676  precip_fact = 1.0_dp-glac_index+glac_index*precip_lgm_anom(j,i,n)
677  ! interpolation with a linear function
678 #elif (PRECIP_ANOM_INTERPOL==2)
679  precip_fact = exp(-glac_index*gamma_precip_lgm_anom(j,i,n))
680  ! interpolation with an exponential function
681 #endif
682 
683  precip(j,i,n) = precip_present(j,i,n)*precip_fact ! monthly precip
684  precip(j,i,0) = precip(j,i,0) + precip(j,i,n)*inv_twelve
685  ! mean annual precip
686  end do
687 
688 #elif (ACCSURFACE == 6)
689 
690 ! -- Mean-annual precipitation from data already read above
691 
692  precip(j,i,0) = precip_ma_present(j,i) + precip_anom_fact*precip_ma_anom(j,i)
693  ! mean annual precip
694 
695  do n=1, 12 ! month counter
696  precip(j,i,n) = precip(j,i,0) ! monthly precip, assumed to be equal
697  ! to the mean annual precip
698  end do
699 
700 #endif
701 
702 ! ---- Annual accumulation, snowfall and rainfall rates
703 
704  accum(j,i) = precip(j,i,0)
705 
706  snowfall(j,i) = 0.0_dp ! initialisation value
707 
708  do n=1, 12 ! month counter
709 
710 #if (SOLID_PRECIP == 1) /* Marsiat (1994) */
711 
712  if (temp_mm(j,i,n) >= temp_rain) then
713  frac_solid = 0.0_dp
714  else if (temp_mm(j,i,n) <= temp_snow) then
715  frac_solid = 1.0_dp
716  else
717  frac_solid = (temp_rain-temp_mm(j,i,n))*inv_delta_temp_rain_snow
718  end if
719 
720 #elif (SOLID_PRECIP == 2) /* Bales et al. (2009) */
721 
722  if (temp_mm(j,i,n) >= temp_rain) then
723  frac_solid = 0.0_dp
724  else if (temp_mm(j,i,n) <= temp_snow) then
725  frac_solid = 1.0_dp
726  else
727  frac_solid = coeff(0) + temp_mm(j,i,n) * ( coeff(1) &
728  + temp_mm(j,i,n) * ( coeff(2) &
729  + temp_mm(j,i,n) * ( coeff(3) &
730  + temp_mm(j,i,n) * ( coeff(4) &
731  + temp_mm(j,i,n) * coeff(5) ) ) ) )
732  ! evaluation of 5th-order polynomial by Horner scheme
733  end if
734 
735 #elif (SOLID_PRECIP == 3) /* Huybrechts and de Wolde (1999) */
736 
737  frac_solid = 1.0_dp &
738  - 0.5_dp*erfcc((temp_rain-temp_mm(j,i,n))*inv_sqrt2_s_stat)
739 
740 #endif
741 
742  snowfall(j,i) = snowfall(j,i) + precip(j,i,n)*frac_solid*inv_twelve
743 
744  end do
745 
746  rainfall(j,i) = precip(j,i,0) - snowfall(j,i)
747 
748  if (snowfall(j,i) < 0.0_dp) snowfall(j,i) = 0.0_dp ! correction of
749  if (rainfall(j,i) < 0.0_dp) rainfall(j,i) = 0.0_dp ! negative values
750 
751 ! ------ Ablation
752 
753 ! ---- Runoff
754 
755 #if (ABLSURFACE==1 || ABLSURFACE==2)
756 
757 ! -- Temperature excess ET
758 
759  do n=1, 12 ! month counter
760  temp_mm_help(n) = temp_mm(j,i,n)
761  end do
762 
763  call pdd(temp_mm_help, s_stat, et(j,i))
764 
765 ! -- Formation rate of superimposed ice (melt_star), melt rate (melt)
766 ! and runoff rate (runoff)
767 
768 #if (ABLSURFACE==1)
769 
770  if ((beta1*et(j,i)) <= (pmax*snowfall(j,i))) then
771  melt_star(j,i) = beta1*et(j,i)
772  melt(j,i) = 0.0_dp
773  runoff(j,i) = melt(j,i)+rainfall(j,i)
774  else
775  melt_star(j,i) = pmax*snowfall(j,i)
776  melt(j,i) = beta2*(et(j,i)-melt_star(j,i)/beta1)
777  runoff(j,i) = melt(j,i)+rainfall(j,i)
778  end if
779 
780 #elif (ABLSURFACE==2)
781 
782  if ( rainfall(j,i) <= (pmax*snowfall(j,i)) ) then
783 
784  if ( (rainfall(j,i)+beta1*et(j,i)) <= (pmax*snowfall(j,i)) ) then
785  melt_star(j,i) = rainfall(j,i)+beta1*et(j,i)
786  melt(j,i) = 0.0_dp
787  runoff(j,i) = melt(j,i)
788  else
789  melt_star(j,i) = pmax*snowfall(j,i)
790  melt(j,i) = beta2 &
791  *(et(j,i)-(melt_star(j,i)-rainfall(j,i))/beta1)
792  runoff(j,i) = melt(j,i)
793  end if
794 
795  else
796 
797  melt_star(j,i) = pmax*snowfall(j,i)
798  melt(j,i) = beta2*et(j,i)
799  runoff(j,i) = melt(j,i) + rainfall(j,i)-pmax*snowfall(j,i)
800 
801  end if
802 
803 #endif
804 
805 #elif (ABLSURFACE==3)
806 
807  temp_jja_help = one_third*(temp_mm(j,i,6)+temp_mm(j,i,7)+temp_mm(j,i,8))
808 
809  melt_star(j,i) = 0.0_dp ! no superimposed ice considered
810  melt(j,i) = lambda_lti*max((temp_jja_help-temp_lti), 0.0_dp)
811  runoff(j,i) = melt(j,i) + rainfall(j,i)
812 
813 #endif
814 
815 ! ---- Evaporation
816 
817  evap(j,i) = 0.0_dp
818 
819 ! ------ Accumulation minus ablation
820 
821  as_perp(j,i) = accum(j,i) - evap(j,i) - runoff(j,i)
822 
823 ! ------ Ice-surface temperature (10-m firn temperature) temp_s,
824 ! including empirical firn-warming correction due to
825 ! refreezing meltwater when superimposed ice is formed
826 
827  if (melt_star(j,i) >= melt(j,i)) then
828  temp_s(j,i) = temp_ma(j,i) &
829  +mu*(melt_star(j,i)-melt(j,i))
830  else
831  temp_s(j,i) = temp_ma(j,i)
832  end if
833 
834  if (temp_s(j,i) > -0.001_dp) temp_s(j,i) = -0.001_dp
835  ! Cut-off of positive air temperatures
836 
837 end do
838 end do
839 
840 end subroutine boundary
841 !
NetCDF error capturing.
Definition: nc_check.F90:35
integer(i2b) function mask_update(z_sl, i, j)
Update of the topography mask due to changes of the sea level.
Definition: mask_update.F90:35
Declarations of kind types for SICOPOLIS.
Definition: sico_types.F90:35
subroutine pdd(temp_mm, s_stat, ET)
Computation of the positive degree days (PDD) with statistical temperature fluctuations; based on sem...
Definition: pdd.F90:41
subroutine check(status, ch_calling_routine)
NetCDF error capturing.
Definition: nc_check.F90:45
real(dp) function erfcc(x)
Computation of the complementary error function erfc(x) = 1-erf(x) with a fractional error everywhere...
Definition: erfcc.F90:39
Declarations of global variables for SICOPOLIS (for the ANT domain).
Definition: sico_vars.F90:35
subroutine boundary(time, dtime, dxi, deta, delta_ts, glac_index, z_sl, dzsl_dtau, z_mar)
Computation of the surface temperature (must be less than 0 deg C!) and of the accumulation-ablation ...
Definition: boundary.F90:37
Declarations of global variables for SICOPOLIS.