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