My Project
Functions/Subroutines | Variables
ice_albedo Module Reference

Functions/Subroutines

subroutine albedos
 

Variables

real(kind=dbl_kind), parameter albocn = 0.06_dbl_kind
 
real(kind=dbl_kind), parameter awtvdr = 0.29_dbl_kind
 
real(kind=dbl_kind), parameter awtidr = 0.31_dbl_kind
 
real(kind=dbl_kind), parameter awtvdf = 0.24_dbl_kind
 
real(kind=dbl_kind), parameter awtidf = 0.16_dbl_kind
 
real(kind=dbl_kind), parameter snowpatch = 0.02_dbl_kind
 
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alvdrn
 
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alidrn
 
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alvdfn
 
real(kind=dbl_kind), dimension(:,:,:), allocatable, save alidfn
 
real(kind=dbl_kind), dimension(:,:), allocatable, save alvdr
 
real(kind=dbl_kind), dimension(:,:), allocatable, save alidr
 
real(kind=dbl_kind), dimension(:,:), allocatable, save alvdf
 
real(kind=dbl_kind), dimension(:,:), allocatable, save alidf
 
real(kind=dbl_kind) albicev
 
real(kind=dbl_kind) albicei
 
real(kind=dbl_kind) albsnowv
 
real(kind=dbl_kind) albsnowi
 

Function/Subroutine Documentation

◆ albedos()

subroutine ice_albedo::albedos ( )

Definition at line 109 of file ice_albedo.f90.

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 
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind), parameter c4i
real(kind=dbl_kind), parameter puny
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save tsfcn
Definition: ice_state.f90:97
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), parameter c1i
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

Variable Documentation

◆ albicei

real (kind=dbl_kind) ice_albedo::albicei

Definition at line 91 of file ice_albedo.f90.

◆ albicev

real (kind=dbl_kind) ice_albedo::albicev

Definition at line 91 of file ice_albedo.f90.

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

◆ albocn

real (kind=dbl_kind), parameter ice_albedo::albocn = 0.06_dbl_kind

Definition at line 48 of file ice_albedo.f90.

48  real (kind=dbl_kind), parameter :: &
49  albocn = 0.06_dbl_kind ! ocean albedo

◆ albsnowi

real (kind=dbl_kind) ice_albedo::albsnowi

Definition at line 91 of file ice_albedo.f90.

◆ albsnowv

real (kind=dbl_kind) ice_albedo::albsnowv

Definition at line 91 of file ice_albedo.f90.

◆ alidf

real (kind=dbl_kind), dimension(:,:), allocatable, save ice_albedo::alidf

Definition at line 83 of file ice_albedo.f90.

◆ alidfn

real (kind=dbl_kind), dimension(:,:,:), allocatable, save ice_albedo::alidfn

Definition at line 76 of file ice_albedo.f90.

◆ alidr

real (kind=dbl_kind), dimension(:,:), allocatable, save ice_albedo::alidr

Definition at line 83 of file ice_albedo.f90.

◆ alidrn

real (kind=dbl_kind), dimension(:,:,:), allocatable, save ice_albedo::alidrn

Definition at line 76 of file ice_albedo.f90.

◆ alvdf

real (kind=dbl_kind), dimension(:,:), allocatable, save ice_albedo::alvdf

Definition at line 83 of file ice_albedo.f90.

◆ alvdfn

real (kind=dbl_kind), dimension(:,:,:), allocatable, save ice_albedo::alvdfn

Definition at line 76 of file ice_albedo.f90.

◆ alvdr

real (kind=dbl_kind), dimension(:,:), allocatable, save ice_albedo::alvdr

Definition at line 83 of file ice_albedo.f90.

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)

◆ alvdrn

real (kind=dbl_kind), dimension(:,:,:), allocatable, save ice_albedo::alvdrn

Definition at line 76 of file ice_albedo.f90.

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)

◆ awtidf

real (kind=dbl_kind), parameter ice_albedo::awtidf = 0.16_dbl_kind

Definition at line 52 of file ice_albedo.f90.

◆ awtidr

real (kind=dbl_kind), parameter ice_albedo::awtidr = 0.31_dbl_kind

Definition at line 52 of file ice_albedo.f90.

◆ awtvdf

real (kind=dbl_kind), parameter ice_albedo::awtvdf = 0.24_dbl_kind

Definition at line 52 of file ice_albedo.f90.

◆ awtvdr

real (kind=dbl_kind), parameter ice_albedo::awtvdr = 0.29_dbl_kind

Definition at line 52 of file ice_albedo.f90.

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

◆ snowpatch

real (kind=dbl_kind), parameter ice_albedo::snowpatch = 0.02_dbl_kind

Definition at line 59 of file ice_albedo.f90.

59  real (kind=dbl_kind), parameter :: &
60  snowpatch = 0.02_dbl_kind