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