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, 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 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, 'julytemp_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 /* Parameterization by Ritz et al. (1997) */
405 
406 theta_ma = 49.13_dp
407 gamma_ma = -7.992e-03_dp
408 c_ma = -0.7576_dp
409 kappa_ma = 0.0_dp
410 
411 theta_mj = 30.38_dp
412 gamma_mj = -6.277e-03_dp
413 c_mj = -0.3262_dp
414 kappa_mj = 0.0_dp
415 
416 #elif TEMP_PRESENT_PARA == 2 /* Parameterization by Fausto et al. (2009) */
417 
418 theta_ma = 41.83_dp
419 gamma_ma = -6.309e-03_dp
420 c_ma = -0.7189_dp
421 kappa_ma = -0.0672_dp
422 
423 theta_mj = 14.70_dp
424 gamma_mj = -5.426e-03_dp
425 c_mj = -0.1585_dp
426 kappa_mj = -0.0518_dp
427 
428 #else
429 
430 stop ' boundary: Parameter TEMP_PRESENT_PARA must be either 1 or 2!'
431 
432 #endif
433 
434 do i=0, imax
435 do j=0, jmax
436 
437 ! ------ Present-day mean-annual air temperature
438 
439  temp_ma_present(j,i) = theta_ma &
440  + gamma_ma*zs(j,i) &
441  + c_ma*phi(j,i)*pi_180_inv &
442  + kappa_ma*(modulo(lambda(j,i)+pi,2.0_dp*pi)-pi)*pi_180_inv
443  ! west longitudes counted negatively
444 
445 ! ------ Present-day mean-July (summer) air temperature
446 
447  temp_mj_present(j,i) = theta_mj &
448  + gamma_mj*zs(j,i) &
449  + c_mj*phi(j,i)*pi_180_inv &
450  + kappa_mj*(modulo(lambda(j,i)+pi,2.0_dp*pi)-pi)*pi_180_inv
451  ! west longitudes counted negatively
452 
453 #if (TSURFACE <= 4)
454 
455 ! ------ Correction with deviation delta_ts
456 
457  temp_ma(j,i) = temp_ma_present(j,i) + delta_ts
458  temp_mm(j,i,7) = temp_mj_present(j,i) + delta_ts
459 
460 #elif (TSURFACE == 5)
461 
462 ! ------ Correction with LGM anomaly and glacial index
463 
464  temp_ma(j,i) = temp_ma_present(j,i) + glac_index*temp_ma_lgm_anom(j,i)
465  temp_mm(j,i,7) = temp_mj_present(j,i) + glac_index*temp_mj_lgm_anom(j,i)
466 
467 #elif (TSURFACE == 6)
468 
469 ! ------ Mean-annual and mean-July (summer) surface temperatures
470 ! from data read above --------
471 
472  temp_ma(j,i) = temp_ma_present(j,i) + temp_ma_anom_fact*temp_ma_anom(j,i)
473  temp_mm(j,i,7) = temp_mj_present(j,i) + temp_mj_anom_fact*temp_mj_anom(j,i)
474 
475 #endif
476 
477 ! ------ Amplitude of the annual cycle
478 
479  temp_ampl(j,i) = temp_mm(j,i,7) - temp_ma(j,i)
480 
481  if (temp_ampl(j,i) < eps) then
482  temp_ampl(j,i) = eps ! Correction of amplitude, if required
483  end if
484 
485 end do
486 end do
487 
488 ! ------ Monthly temperatures
489 
490 do n=1, 12 ! month counter
491 
492  sine_factor = sin((real(n,dp)-4.0_dp)*pi/6.0_dp)
493 
494  do i=0, imax
495  do j=0, jmax
496  temp_mm(j,i,n) = temp_ma(j,i) + sine_factor*temp_ampl(j,i)
497  end do
498  end do
499 
500 end do
501 
502 !-------- Accumulation-ablation function as_perp --------
503 
504 #if (ELEV_DESERT == 1)
505 
506 gamma_p = gamma_p*1.0e-03_dp ! Precipitation lapse rate
507  ! for elevation desertification, in m^(-1)
508 zs_thresh = zs_thresh ! Elevation threshold, in m
509 
510 #endif
511 
512 #if (SOLID_PRECIP == 1) /* Marsiat (1994) */
513 
514 temp_rain = 7.0_dp ! Threshold monthly mean temperature for
515  ! precipitation = 100% rain, in deg C
516 temp_snow = -10.0_dp ! Threshold monthly mean temperature for &
517  ! precipitation = 100% snow, in deg C
518 
519 inv_delta_temp_rain_snow = 1.0_dp/(temp_rain-temp_snow)
520 
521 #elif (SOLID_PRECIP == 2) /* Bales et al. (2009) */
522 
523 temp_rain = 7.2_dp ! Threshold monthly mean temperature for &
524  ! precipitation = 100% rain, in deg C
525 temp_snow = -11.6_dp ! Threshold monthly mean temperature for &
526  ! precipitation = 100% snow, in deg C
527 
528 coeff(0) = 5.4714e-01_dp ! Coefficients
529 coeff(1) = -9.1603e-02_dp ! of
530 coeff(2) = -3.314e-03_dp ! the
531 coeff(3) = 4.66e-04_dp ! fifth-order
532 coeff(4) = 3.8e-05_dp ! polynomial
533 coeff(5) = 6.0e-07_dp ! fit
534 
535 #elif (SOLID_PRECIP == 3) /* Huybrechts and de Wolde (1999) */
536 
537 temp_rain = 2.0_dp ! Threshold instantaneous temperature for &
538  ! precipitation = 100% rain, in deg C
539 temp_snow = temp_rain ! Threshold instantaneous temperature for &
540  ! precipitation = 100% snow, in deg C
541 
542 s_stat = s_stat_0 ! Standard deviation of the air termperature
543  ! (same parameter as in the PDD model)
544 
545 inv_sqrt2_s_stat = 1.0_dp/(sqrt(2.0_dp)*s_stat)
546 
547 #endif
548 
549 #if (ABLSURFACE==1 || ABLSURFACE==2)
550 
551 s_stat = s_stat_0
552 
553 phi_sep = phi_sep_0*pi_180 ! separates different domains for computation of
554  ! degree-day factors beta1 and beta2
555  ! deg N --> rad
556 temp_lt = -1.0_dp
557 temp_ht = 10.0_dp
558 inv_delta_temp_ht_lt = 1.0_dp/(temp_ht-temp_lt)
559 
560 beta1_lt = beta1_lt_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
561  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
562 beta1_ht = beta1_ht_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
563  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
564 beta2_lt = beta2_lt_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
565  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
566 beta2_ht = beta2_ht_0 *(0.001_dp/86400.0_dp)*(rho_w/rho)
567  ! (mm WE)/(d*deg C) --> (m IE)/(s*deg C)
568 pmax = pmax_0
569 mu = mu_0 *(1000.0_dp*86400.0_dp)*(rho/rho_w)
570  ! (d*deg C)/(mm WE) --> (s*deg C)/(m IE)
571 #if (PDD_MODIFIER==2)
572 
573 lon_w_e_sep = lon_w_e_sep *pi_180 ! deg N -> rad
574 
575 if (lon_w_e_sep < 0.0_dp) then
576  lon_w_e_sep = lon_w_e_sep + 2.0_dp*pi
577 else if (lon_w_e_sep >= (2.0_dp*pi)) then
578  lon_w_e_sep = lon_w_e_sep - 2.0_dp*pi
579 end if
580 
581 pdd_mod_lat = 0.0_dp ! initialisation
582 pdd_mod_fac_w = 0.0_dp ! initialisation
583 pdd_mod_fac_e = 0.0_dp ! initialisation
584 
585 pdd_mod_lat( 1) = pdd_mod_lat_01 *pi_180 ! deg N -> rad
586 pdd_mod_fac_w( 1) = pdd_mod_fac_w_01
587 pdd_mod_fac_e( 1) = pdd_mod_fac_e_01
588 pdd_mod_lat( 2) = pdd_mod_lat_02 *pi_180 ! deg N -> rad
589 pdd_mod_fac_w( 2) = pdd_mod_fac_w_02
590 pdd_mod_fac_e( 2) = pdd_mod_fac_e_02
591 pdd_mod_lat( 3) = pdd_mod_lat_03 *pi_180 ! deg N -> rad
592 pdd_mod_fac_w( 3) = pdd_mod_fac_w_03
593 pdd_mod_fac_e( 3) = pdd_mod_fac_e_03
594 pdd_mod_lat( 4) = pdd_mod_lat_04 *pi_180 ! deg N -> rad
595 pdd_mod_fac_w( 4) = pdd_mod_fac_w_04
596 pdd_mod_fac_e( 4) = pdd_mod_fac_e_04
597 pdd_mod_lat( 5) = pdd_mod_lat_05 *pi_180 ! deg N -> rad
598 pdd_mod_fac_w( 5) = pdd_mod_fac_w_05
599 pdd_mod_fac_e( 5) = pdd_mod_fac_e_05
600 pdd_mod_lat( 6) = pdd_mod_lat_06 *pi_180 ! deg N -> rad
601 pdd_mod_fac_w( 6) = pdd_mod_fac_w_06
602 pdd_mod_fac_e( 6) = pdd_mod_fac_e_06
603 pdd_mod_lat( 7) = pdd_mod_lat_07 *pi_180 ! deg N -> rad
604 pdd_mod_fac_w( 7) = pdd_mod_fac_w_07
605 pdd_mod_fac_e( 7) = pdd_mod_fac_e_07
606 pdd_mod_lat( 8) = pdd_mod_lat_08 *pi_180 ! deg N -> rad
607 pdd_mod_fac_w( 8) = pdd_mod_fac_w_08
608 pdd_mod_fac_e( 8) = pdd_mod_fac_e_08
609 pdd_mod_lat( 9) = pdd_mod_lat_09 *pi_180 ! deg N -> rad
610 pdd_mod_fac_w( 9) = pdd_mod_fac_w_09
611 pdd_mod_fac_e( 9) = pdd_mod_fac_e_09
612 pdd_mod_lat( 10) = pdd_mod_lat_10 *pi_180 ! deg N -> rad
613 pdd_mod_fac_w(10) = pdd_mod_fac_w_10
614 pdd_mod_fac_e(10) = pdd_mod_fac_e_10
615 
616 delta_pdd_mod_lat_inv = 0.0_dp ! initialisation
617 
618 do n=1, n_pdd_mod-1
619  delta_pdd_mod_lat_inv(n) = 1.0_dp/(pdd_mod_lat(n+1)-pdd_mod_lat(n))
620 end do
621 
622 #endif
623 
624 #elif (ABLSURFACE==3)
625 
626 lambda_lti = lambda_lti *(0.001_dp/year_sec)*(rho_w/rho)
627  ! (mm WE)/(a*deg C) --> (m IE)/(s*deg C)
628 temp_lti = temp_lti
629 mu = 0.0_dp ! no superimposed ice considered
630 
631 #endif
632 
633 do i=0, imax
634 do j=0, jmax
635 
636 #if (ABLSURFACE==1 || ABLSURFACE==2)
637 
638  if (phi(j,i) <= phi_sep) then
639 
640  beta1 = beta1_ht
641  beta2 = beta2_ht
642 
643  else
644 
645  if (temp_mm(j,i,7) >= temp_ht) then
646  beta1 = beta1_ht
647  beta2 = beta2_ht
648  else if (temp_mm(j,i,7) <= temp_lt) then
649  beta1 = beta1_lt
650  beta2 = beta2_lt
651  else
652  beta1 = beta1_lt &
653  + (beta1_ht-beta1_lt) &
654  *inv_delta_temp_ht_lt*(temp_mm(j,i,7)-temp_lt)
655  beta2 = beta2_ht &
656  + (beta2_lt-beta2_ht) &
657  *(inv_delta_temp_ht_lt*(temp_ht-temp_mm(j,i,7)))**3
658  end if
659 
660  end if
661 
662 #if (PDD_MODIFIER==2)
663 
664  if ( lambda(j,i) <= lon_w_e_sep ) then ! West Greenland
665  pdd_mod_fac = pdd_mod_fac_w
666  else ! lambda(j,i) > lon_w_e_sep, East Greenland
667  pdd_mod_fac = pdd_mod_fac_e
668  end if
669 
670  if (phi(j,i) <= pdd_mod_lat(1)) then
671 
672  beta1 = beta1 * pdd_mod_fac(1)
673  beta2 = beta2 * pdd_mod_fac(1)
674 
675  else if (phi(j,i) >= pdd_mod_lat(n_pdd_mod)) then
676 
677  beta1 = beta1 * pdd_mod_fac(n_pdd_mod)
678  beta2 = beta2 * pdd_mod_fac(n_pdd_mod)
679 
680  else
681 
682  do n=1, n_pdd_mod-1
683 
684  if ( (phi(j,i) >= pdd_mod_lat(n)) &
685  .and. &
686  (phi(j,i) <= pdd_mod_lat(n+1)) ) then
687 
688  pdd_mod_fac_interpol = pdd_mod_fac(n) &
689  + (pdd_mod_fac(n+1)-pdd_mod_fac(n)) &
690  *(phi(j,i)-pdd_mod_lat(n)) &
691  *delta_pdd_mod_lat_inv(n)
692 
693  beta1 = beta1 * pdd_mod_fac_interpol
694  beta2 = beta2 * pdd_mod_fac_interpol
695 
696  exit
697 
698  end if
699 
700  end do
701 
702  end if
703 
704 #endif
705 
706 #endif
707 
708 ! ------ Accumulation
709 
710 #if (ACCSURFACE <= 3)
711 
712 ! ---- Elevation desertification of precipitation
713 
714 #if (ELEV_DESERT == 0)
715 
716  precip_fact = 1.0_dp ! no elevation desertification
717 
718 #elif (ELEV_DESERT == 1)
719 
720  if (zs_ref(j,i) < zs_thresh) then
721  precip_fact &
722  = exp(gamma_p*(max(zs(j,i),zs_thresh)-zs_thresh))
723  else
724  precip_fact &
725  = exp(gamma_p*(max(zs(j,i),zs_thresh)-zs_ref(j,i)))
726  end if
727 
728 #else
729  stop ' boundary: Parameter ELEV_DESERT must be either 0 or 1!'
730 #endif
731 
732  do n=1, 12 ! month counter
733  precip(j,i,n) = precip_present(j,i,n)*precip_fact
734  end do
735 
736 #endif
737 
738 ! ---- Precipitation change related to changing climate
739 
740 #if ACCSURFACE==1
741  precip_fact = accfact
742 #elif ACCSURFACE==2
743  precip_fact = 1.0_dp + gamma_s*delta_ts
744 #elif ACCSURFACE==3
745  precip_fact = exp(gamma_s*delta_ts)
746 #endif
747 
748 #if (ACCSURFACE <= 3)
749 
750  precip(j,i,0) = 0.0_dp ! initialisation value for mean annual precip
751 
752  do n=1, 12 ! month counter
753  precip(j,i,n) = precip(j,i,n)*precip_fact ! monthly precip
754  precip(j,i,0) = precip(j,i,0) + precip(j,i,n)*inv_twelve
755  ! mean annual precip
756  end do
757 
758 #elif (ACCSURFACE == 5)
759 
760  precip(j,i,0) = 0.0_dp ! initialisation value for mean annual precip
761 
762  do n=1, 12 ! month counter
763 
764 #if (PRECIP_ANOM_INTERPOL==1)
765  precip_fact = 1.0_dp-glac_index+glac_index*precip_lgm_anom(j,i,n)
766  ! interpolation with a linear function
767 #elif (PRECIP_ANOM_INTERPOL==2)
768  precip_fact = exp(-glac_index*gamma_precip_lgm_anom(j,i,n))
769  ! interpolation with an exponential function
770 #endif
771 
772  precip(j,i,n) = precip_present(j,i,n)*precip_fact ! monthly precip
773  precip(j,i,0) = precip(j,i,0) + precip(j,i,n)*inv_twelve
774  ! mean annual precip
775  end do
776 
777 #elif (ACCSURFACE == 6)
778 
779 ! -- Mean-annual precipitation from data already read above
780 
781  precip(j,i,0) = precip_ma_present(j,i) + precip_anom_fact*precip_ma_anom(j,i)
782  ! mean annual precip
783 
784  do n=1, 12 ! month counter
785  precip(j,i,n) = precip(j,i,0) ! monthly precip, assumed to be equal
786  ! to the mean annual precip
787  end do
788 
789 #endif
790 
791 ! ---- Annual accumulation, snowfall and rainfall rates
792 
793  accum(j,i) = precip(j,i,0)
794 
795  snowfall(j,i) = 0.0_dp ! initialisation value
796 
797  do n=1, 12 ! month counter
798 
799 #if (SOLID_PRECIP == 1) /* Marsiat (1994) */
800 
801  if (temp_mm(j,i,n) >= temp_rain) then
802  frac_solid = 0.0_dp
803  else if (temp_mm(j,i,n) <= temp_snow) then
804  frac_solid = 1.0_dp
805  else
806  frac_solid = (temp_rain-temp_mm(j,i,n))*inv_delta_temp_rain_snow
807  end if
808 
809 #elif (SOLID_PRECIP == 2) /* Bales et al. (2009) */
810 
811  if (temp_mm(j,i,n) >= temp_rain) then
812  frac_solid = 0.0_dp
813  else if (temp_mm(j,i,n) <= temp_snow) then
814  frac_solid = 1.0_dp
815  else
816  frac_solid = coeff(0) + temp_mm(j,i,n) * ( coeff(1) &
817  + temp_mm(j,i,n) * ( coeff(2) &
818  + temp_mm(j,i,n) * ( coeff(3) &
819  + temp_mm(j,i,n) * ( coeff(4) &
820  + temp_mm(j,i,n) * coeff(5) ) ) ) )
821  ! evaluation of 5th-order polynomial by Horner scheme
822  end if
823 
824 #elif (SOLID_PRECIP == 3) /* Huybrechts and de Wolde (1999) */
825 
826  frac_solid = 1.0_dp &
827  - 0.5_dp*erfcc((temp_rain-temp_mm(j,i,n))*inv_sqrt2_s_stat)
828 
829 #endif
830 
831  snowfall(j,i) = snowfall(j,i) + precip(j,i,n)*frac_solid*inv_twelve
832 
833  end do
834 
835  rainfall(j,i) = precip(j,i,0) - snowfall(j,i)
836 
837  if (snowfall(j,i) < 0.0_dp) snowfall(j,i) = 0.0_dp ! correction of
838  if (rainfall(j,i) < 0.0_dp) rainfall(j,i) = 0.0_dp ! negative values
839 
840 ! ------ Ablation
841 
842 ! ---- Runoff
843 
844 #if (ABLSURFACE==1 || ABLSURFACE==2)
845 
846 ! -- Temperature excess ET
847 
848  do n=1, 12 ! month counter
849  temp_mm_help(n) = temp_mm(j,i,n)
850  end do
851 
852  call pdd(temp_mm_help, s_stat, et(j,i))
853 
854 ! -- Formation rate of superimposed ice (melt_star), melt rate (melt)
855 ! and runoff rate (runoff)
856 
857 #if (ABLSURFACE==1)
858 
859  if ((beta1*et(j,i)) <= (pmax*snowfall(j,i))) then
860  melt_star(j,i) = beta1*et(j,i)
861  melt(j,i) = 0.0_dp
862  runoff(j,i) = melt(j,i)+rainfall(j,i)
863  else
864  melt_star(j,i) = pmax*snowfall(j,i)
865  melt(j,i) = beta2*(et(j,i)-melt_star(j,i)/beta1)
866  runoff(j,i) = melt(j,i)+rainfall(j,i)
867  end if
868 
869 #elif (ABLSURFACE==2)
870 
871  if ( rainfall(j,i) <= (pmax*snowfall(j,i)) ) then
872 
873  if ( (rainfall(j,i)+beta1*et(j,i)) <= (pmax*snowfall(j,i)) ) then
874  melt_star(j,i) = rainfall(j,i)+beta1*et(j,i)
875  melt(j,i) = 0.0_dp
876  runoff(j,i) = melt(j,i)
877  else
878  melt_star(j,i) = pmax*snowfall(j,i)
879  melt(j,i) = beta2 &
880  *(et(j,i)-(melt_star(j,i)-rainfall(j,i))/beta1)
881  runoff(j,i) = melt(j,i)
882  end if
883 
884  else
885 
886  melt_star(j,i) = pmax*snowfall(j,i)
887  melt(j,i) = beta2*et(j,i)
888  runoff(j,i) = melt(j,i) + rainfall(j,i)-pmax*snowfall(j,i)
889 
890  end if
891 
892 #endif
893 
894 #elif (ABLSURFACE==3)
895 
896  temp_jja_help = one_third*(temp_mm(j,i,6)+temp_mm(j,i,7)+temp_mm(j,i,8))
897 
898  melt_star(j,i) = 0.0_dp ! no superimposed ice considered
899  melt(j,i) = lambda_lti*max((temp_jja_help-temp_lti), 0.0_dp)
900  runoff(j,i) = melt(j,i) + rainfall(j,i)
901 
902 #endif
903 
904 ! ---- Evaporation
905 
906  evap(j,i) = 0.0_dp
907 
908 ! ------ Accumulation minus ablation
909 
910  as_perp(j,i) = accum(j,i) - evap(j,i) - runoff(j,i)
911 
912 ! ------ Ice-surface temperature (10-m firn temperature) temp_s,
913 ! including empirical firn-warming correction due to
914 ! refreezing meltwater when superimposed ice is formed
915 
916  if (melt_star(j,i) >= melt(j,i)) then
917  temp_s(j,i) = temp_ma(j,i) &
918  +mu*(melt_star(j,i)-melt(j,i))
919  else
920  temp_s(j,i) = temp_ma(j,i)
921  end if
922 
923  if (temp_s(j,i) > -0.001_dp) temp_s(j,i) = -0.001_dp
924  ! Cut-off of positive air temperatures
925 
926 end do
927 end do
928 
929 end subroutine boundary
930 !