/[MITgcm]/MITgcm/pkg/thsice/thsice_albedo.F
ViewVC logotype

Annotation of /MITgcm/pkg/thsice/thsice_albedo.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (hide annotations) (download)
Wed Dec 31 17:44:32 2003 UTC (20 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint52l_pre, checkpoint52e_pre, hrcube4, checkpoint52j_post, checkpoint52e_post, hrcube_1, checkpoint52l_post, checkpoint52k_post, checkpoint52f_post, hrcube5, checkpoint52i_post, checkpoint52j_pre, checkpoint52i_pre, checkpoint52h_pre, checkpoint52f_pre, hrcube_2, hrcube_3
change how Albedo is computed over sea-ice:
 - separate albedo S/R (thsice_albedo.F) from thsice_therm.F file.
 - use GISS formulation for snow-covered area ; turn on snow-aging ;
 - bare-ice albedo as in LANL code ;
 - addd diagnostic of surface albedo.

1 jmc 1.1 C $Header: $
2     C $Name: $
3    
4     #include "THSICE_OPTIONS.h"
5    
6     SUBROUTINE THSICE_ALBEDO(hi,hs,Tsf,age,albedo)
7     C *==========================================================*
8     C | S/R THSICE_ALBEDO
9     C *==========================================================*
10     C | Compute surface albedo
11     C *==========================================================*
12     IMPLICIT NONE
13    
14     C == Global data ==
15     #include "THSICE_PARAMS.h"
16    
17     C == Routine Arguments ==
18     _RL hi ! ice height
19     _RL hs ! snow height
20     _RL Tsf ! surface temperature
21     _RL age ! snow age
22     _RL albedo ! surface albedo
23    
24     #ifdef ALLOW_THSICE
25     C == Local variables ==
26     C frsnow :: fractional snow cover
27     C albsno :: albedo of snow
28     C albice :: albedo of ice
29     C albNewSnow :: albedo of new (fresh) snow
30     _RL frsnow
31     _RL albsno, albice
32     _RL albNewSnow
33    
34     C-- Albedo of Bare Sea-Ice
35     albice = albIceMax + (albIceMin-albIceMax)*exp(-hi/hAlbIce)
36    
37     C-- LANL albedo calculation
38     c frsnow = 0.
39     c if (hs .gt. 0.) frsnow = 1.
40     c if (Tsf .lt. 0.) then
41     c albedo = frsnow*albColdSnow + (1.-frsnow)*albice
42     c else
43     c albedo = frsnow*albWarmSnow + (1.-frsnow)*albice
44     c endif
45     C-end LANL albedo calculation
46    
47     C-- GISS model albedo calculation
48     c albice = 0.7 _d 0
49    
50     C- New snow: (linear) transition between -10.oC and 0.oC
51     C from cold/dry snow albedo to warm/wet snow albedo
52     albNewSnow = albColdSnow
53     & + (albWarmSnow - albColdSnow)
54     & *MAX( 0. _d 0, 1. _d 0 + MIN(Tsf/10. _d 0, 0. _d 0) )
55     C- albedo of snow is function of snow-age (make age units into days):
56     albsno = albOldSnow
57     & +(albNewSnow-albOldSnow)*exp(-0.2 _d 0*age/86400. _d 0)
58     C- layer of snow over the ice:
59     albedo = albsno + (albice-albsno)*exp(-hs/hAlbSnow)
60    
61     if (albedo.gt.1. _d 0 .or. albedo.lt. .2 _d 0) then
62     print*,'QQ - albedo problem', albedo, age, hs, albsno
63     stop
64     endif
65    
66     #endif /* ALLOW_THSICE */
67    
68     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
69    
70     RETURN
71     END

  ViewVC Help
Powered by ViewVC 1.1.22