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