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