My Project
ice_albedo.f90
Go to the documentation of this file.
1 
2 
3 
4 
5 
6 
7 
8 
9 
10 
11 
12 !/===========================================================================/
13 ! CVS VERSION INFORMATION
14 ! $Id$
15 ! $Name$
16 ! $Revision$
17 !/===========================================================================/
18 
19 !=======================================================================
20 !BOP
21 !
22 ! !MODULE: ice_albedo - snow and ice albedo parameterization
23 !
24 ! !DESCRIPTION:
25 !
26 ! The albedo parameterization
27 !
28 ! !REVISION HISTORY:
29 !
30 ! authors: Bruce P. Briegleb, NCAR
31 ! Elizabeth C. Hunke, LANL
32 !
33 ! Vectorized by Clifford Chen (Fujitsu) and William H. Lipscomb (LANL)
34 !
35 ! !INTERFACE:
36 !
37  module ice_albedo
38 !
39 ! !USES:
40 !
41  use ice_kinds_mod
42  use ice_domain
43 !
44 !EOP
45 !
46  implicit none
47 
48  real (kind=dbl_kind), parameter :: &
49  albocn = 0.06_dbl_kind ! ocean albedo
50 
51  ! weights for albedos match those for isccp shortwave forcing
52  real (kind=dbl_kind), parameter :: & ! currently used only
53  awtvdr = 0.29_dbl_kind &! visible, direct ! for history and
54  , awtidr = 0.31_dbl_kind &! near IR, direct ! diagnostics
55  , awtvdf = 0.24_dbl_kind &! visible, diffuse
56  , awtidf = 0.16_dbl_kind ! near IR, diffuse
57 
58  ! parameter for fractional snow area
59  real (kind=dbl_kind), parameter :: &
60  snowpatch = 0.02_dbl_kind
61 
62 ! ! albedos for ice in each category
63 ! real (kind=dbl_kind) :: &
64 ! alvdrn (ilo:ihi,jlo:jhi,ncat) &! visible, direct (fraction)
65 ! , alidrn (ilo:ihi,jlo:jhi,ncat) &! near-ir, direct (fraction)
66 ! , alvdfn (ilo:ihi,jlo:jhi,ncat) &! visible, diffuse (fraction)
67 ! , alidfn (ilo:ihi,jlo:jhi,ncat) ! near-ir, diffuse (fraction)
68 
69  ! albedos aggregated over categories
70 ! real (kind=dbl_kind) :: &
71 ! alvdr (ilo:ihi,jlo:jhi) &! visible, direct (fraction)
72 ! , alidr (ilo:ihi,jlo:jhi) &! near-ir, direct (fraction)
73 ! , alvdf (ilo:ihi,jlo:jhi) &! visible, diffuse (fraction)
74 ! , alidf (ilo:ihi,jlo:jhi) ! near-ir, diffuse (fraction)
75 
76  real (kind=dbl_kind),dimension(:,:,:),allocatable,save :: &
77  alvdrn &! visible, direct (fraction)
78  , alidrn &! near-ir, direct (fraction)
79  , alvdfn &! visible, diffuse (fraction)
80  , alidfn ! near-ir, diffuse (fraction)
81 
82  ! albedos aggregated over categories
83  real (kind=dbl_kind),dimension(:,:),allocatable,save :: &
84  alvdr &! visible, direct (fraction)
85  , alidr &! near-ir, direct (fraction)
86  , alvdf &! visible, diffuse (fraction)
87  , alidf ! near-ir, diffuse (fraction)
88 
89 
90  ! baseline albedos for thick cases
91  real (kind=dbl_kind) :: &
92  albicev &! visible ice albedo for h > ahmax
93  , albicei &! near-ir ice albedo for h > ahmax
94  , albsnowv &! cold snow albedo, visible
95  , albsnowi ! cold snow albedo, near IR
96 
97 !=======================================================================
98 
99  contains
100 
101 !=======================================================================
102 !BOP
103 !
104 ! !IROUTINE: albedos - compute snow/ice albedos and aggregate
105 !
106 ! !INTERFACE:
107 !
108  subroutine albedos
109 !
110 ! !DESCRIPTION:
111 !
112 ! Compute albedos and aggregate them \! note: ice albedo is zero if no ice present
113 !
114 ! !REVISION HISTORY:
115 !
116 ! authors: Bruce P. Briegleb, NCAR
117 ! Elizabeth C. Hunke, LANL
118 !
119 ! !USES:
120 !
121  use ice_constants
122  use ice_grid
123  use ice_state
124 !
125 ! !INPUT/OUTPUT PARAMETERS:
126 !
127 !EOP
128 !
129  real (kind=dbl_kind), parameter :: &
130  ahmax = 0.5_dbl_kind &! thickness above which ice albedo constant (m)
131  , dt_mlt = 1._dbl_kind &! change in temp to give dalb_mlt albedo change
132  , dalb_mlt = -0.075_dbl_kind &! albedo change per dT_mlt change
133  ! in temp for ice
134  , dalb_mltv = -0.100_dbl_kind &! albedo vis change per dT_mlt change
135  ! in temp for snow
136  , dalb_mlti = -0.150_dbl_kind ! albedo nir change per dT_mlt change
137  ! in temp for snow
138 
139  integer (kind=int_kind) :: i, j, ni
140 
141  real (kind=dbl_kind) :: &
142  hi &! ice thickness (m)
143  , hs &! snow thickness (m)
144  , albo &! effective ocean albedo, function of ice thickness
145  , asnow &! snow-covered area fraction
146  , asnwv &! snow albedo, visible
147  , asnwi &! snow albedo, near IR
148  , fh &! piecewise linear function of thickness
149  , ft &! piecewise linear function of surface temperature
150  , dts &! difference of Tsfc and Timelt
151  , fhtan ! factor used in albedo dependence on ice thickness
152 
153  integer (kind=int_kind) :: &
154  icells &! number of ice/ocean grid cells
155  , ij ! horizontal index, combines i and j loops
156 
157  integer (kind=int_kind), dimension (1:(ihi-ilo+1)*(jhi-jlo+1)) :: &
158  indxi &! compressed indices for ice/ocean cells
159  , indxj
160 
161  fhtan = atan(ahmax*c4i)
162 
163  icells = 0
164  do j = jlo, jhi
165  do i = ilo, ihi
166  if (tmask(i,j)) then
167  icells = icells + 1
168  indxi(icells) = i
169  indxj(icells) = j
170  endif ! tmask
171  enddo
172  enddo
173 
174  !-----------------------------------------------------------------
175  ! albedo for each thickness category
176  !-----------------------------------------------------------------
177 
178  do ni = 1, ncat
179 !cdir$ ivdep !Cray
180 !cdir nodep !NEC
181 !ocl novrec !Fujitsu
182  do ij = 1, icells
183  i = indxi(ij)
184  j = indxj(ij)
185 
186  if (aicen(i,j,ni) > puny) then
187  hi = vicen(i,j,ni) / aicen(i,j,ni)
188  hs = vsnon(i,j,ni) / aicen(i,j,ni)
189 
190  ! bare ice, thickness dependence
191  fh = min(atan(hi*c4i)/fhtan,c1i)
192  albo = albocn*(c1i-fh)
193  alvdfn(i,j,ni) = albicev*fh + albo
194  alidfn(i,j,ni) = albicei*fh + albo
195 
196  ! bare ice, temperature dependence
197  dts = timelt - tsfcn(i,j,ni)
198  ft = min(dts/dt_mlt-c1i,c0i)
199  alvdfn(i,j,ni) = alvdfn(i,j,ni) - dalb_mlt*ft
200  alidfn(i,j,ni) = alidfn(i,j,ni) - dalb_mlt*ft
201 
202  ! avoid negative albedos for thin, bare, melting ice
203  alvdfn(i,j,ni) = max(alvdfn(i,j,ni), albocn)
204  alidfn(i,j,ni) = max(alidfn(i,j,ni), albocn)
205 
206  if( hs > puny ) then
207 
208  ! fractional area of snow on ice (thickness dependent)
209  asnow = hs / ( hs + snowpatch )
210  asnwv = albsnowv
211  asnwi = albsnowi
212 
213  ! snow on ice, temperature dependence
214  asnwv = asnwv - dalb_mltv*ft
215  asnwi = asnwi - dalb_mlti*ft
216 
217  ! combine ice and snow albedos
218  alvdfn(i,j,ni) = alvdfn(i,j,ni)*(c1i-asnow) + &
219  asnwv*asnow
220  alidfn(i,j,ni) = alidfn(i,j,ni)*(c1i-asnow) + &
221  asnwi*asnow
222  endif ! hs > puny
223 
224  alvdrn(i,j,ni) = alvdfn(i,j,ni)
225  alidrn(i,j,ni) = alidfn(i,j,ni)
226 
227  else ! no ice
228  alvdfn(i,j,ni) = albocn
229  alidfn(i,j,ni) = albocn
230  alvdrn(i,j,ni) = albocn
231  alidrn(i,j,ni) = albocn
232  endif ! aicen > puny
233  enddo ! ij
234  enddo ! ncat
235 
236  !-----------------------------------------------------------------
237  ! aggregate
238  !-----------------------------------------------------------------
239 
240  do j = jlo, jhi
241  do i = ilo, ihi
242  alvdf(i,j) = c0i
243  alidf(i,j) = c0i
244  alvdr(i,j) = c0i
245  alidr(i,j) = c0i
246  enddo
247  enddo
248 
249  do ni = 1, ncat
250 !cdir$ ivdep !Cray
251 !cdir nodep !NEC
252 !ocl novrec !Fujitsu
253  do ij = 1, icells
254  i = indxi(ij)
255  j = indxj(ij)
256  if (aice(i,j) > puny) then
257  alvdf(i,j) = alvdf(i,j) + alvdfn(i,j,ni)*aicen(i,j,ni)
258  alidf(i,j) = alidf(i,j) + alidfn(i,j,ni)*aicen(i,j,ni)
259  alvdr(i,j) = alvdr(i,j) + alvdrn(i,j,ni)*aicen(i,j,ni)
260  alidr(i,j) = alidr(i,j) + alidrn(i,j,ni)*aicen(i,j,ni)
261  endif
262  enddo ! ij
263  enddo ! ncat
264 
265 !cdir$ ivdep !Cray
266 !cdir nodep !NEC
267 !ocl novrec !Fujitsu
268  do ij = 1, icells
269  i = indxi(ij)
270  j = indxj(ij)
271  if (aice(i,j) > puny) then
272  alvdf(i,j) = alvdf(i,j) / aice(i,j)
273  alidf(i,j) = alidf(i,j) / aice(i,j)
274  alvdr(i,j) = alvdr(i,j) / aice(i,j)
275  alidr(i,j) = alidr(i,j) / aice(i,j)
276  endif ! aicen > puny
277  enddo ! ij
278 
279  end subroutine albedos
280 
281 !=======================================================================
282 
283  end module ice_albedo
284 
285 !=======================================================================
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alidrn
Definition: ice_albedo.f90:76
integer, parameter dbl_kind
real(kind=dbl_kind), dimension(:,:), allocatable, save alidf
Definition: ice_albedo.f90:83
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind) albicei
Definition: ice_albedo.f90:91
real(kind=dbl_kind), parameter c4i
real(kind=dbl_kind), dimension(:,:), allocatable, save alvdr
Definition: ice_albedo.f90:83
real(kind=dbl_kind), parameter snowpatch
Definition: ice_albedo.f90:59
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:), allocatable, save alvdf
Definition: ice_albedo.f90:83
real(kind=dbl_kind), parameter albocn
Definition: ice_albedo.f90:48
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alidfn
Definition: ice_albedo.f90:76
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alvdrn
Definition: ice_albedo.f90:76
integer(kind=int_kind) ilo
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter puny
integer(kind=int_kind) jhi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter awtvdf
Definition: ice_albedo.f90:52
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save tsfcn
Definition: ice_state.f90:97
subroutine albedos
Definition: ice_albedo.f90:109
real(kind=dbl_kind) albsnowv
Definition: ice_albedo.f90:91
real(kind=dbl_kind), parameter awtidr
Definition: ice_albedo.f90:52
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vicen
Definition: ice_state.f90:97
real(kind=dbl_kind), parameter timelt
real(kind=dbl_kind), dimension(:,:), allocatable, target, save aice
Definition: ice_state.f90:82
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:), allocatable, save alidr
Definition: ice_albedo.f90:83
real(kind=dbl_kind), parameter c1i
real(kind=dbl_kind), parameter awtvdr
Definition: ice_albedo.f90:52
real(kind=dbl_kind) albsnowi
Definition: ice_albedo.f90:91
real(kind=dbl_kind) albicev
Definition: ice_albedo.f90:91
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alvdfn
Definition: ice_albedo.f90:76
logical(kind=log_kind), dimension(:,:), allocatable tmask
Definition: ice_grid.f90:164
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vsnon
Definition: ice_state.f90:97
real(kind=dbl_kind), parameter awtidf
Definition: ice_albedo.f90:52