My Project
Functions/Subroutines | Variables
ice_itd Module Reference

Functions/Subroutines

subroutine init_itd
 
subroutine aggregate
 
subroutine aggregate_area
 
subroutine bound_aggregate
 
subroutine bound_state
 
subroutine check_state
 
subroutine rebin
 
subroutine reduce_area (hice1_old, hice1)
 
subroutine shift_ice (donor, daice, dvice, hicen)
 
subroutine column_sum (nsum, xin, xout)
 
subroutine column_conservation_check (x1, x2, max_err, fieldid)
 
subroutine zap_small_areas
 

Variables

integer(kind=int_kind) kitd
 
integer(kind=int_kind) kcatbound
 
integer(kind=int_kind), dimension(ncat) ilyr1
 
integer(kind=int_kind), dimension(ncat) ilyrn
 
real(kind=dbl_kind), parameter hi_min = p01
 
real(kind=dbl_kind), dimension(0:ncat) hin_max
 

Function/Subroutine Documentation

◆ aggregate()

subroutine ice_itd::aggregate ( )

Definition at line 234 of file ice_itd.f90.

234 !
235 ! !DESCRIPTION:
236 !
237 ! Aggregate ice state variables over thickness categories.
238 !
239 ! !REVISION HISTORY:
240 !
241 ! authors: C. M. Bitz, UW
242 ! W. H. Lipscomb, LANL
243 !
244 ! !USES:
245 !
246  use ice_domain
247  use ice_flux, only : tf
248  use ice_grid
249 !
250 ! !INPUT/OUTPUT PARAMETERS:
251 !
252 !EOP
253 !
254  integer (kind=int_kind) :: i, j, k, ni
255 
256  integer (kind=int_kind), dimension (1:(ihi-ilo+1)*(jhi-jlo+1)) :: &
257  indxi & ! compressed indices in i/j directions
258  , indxj
259 
260  integer (kind=int_kind) :: &
261  icells & ! number of ocean/ice cells
262  , ij ! combined i/j horizontal index
263 
264 !! ggao 60162008
265  real (kind=dbl_kind) ::eps
266 !! change end
267 
268  !-----------------------------------------------------------------
269  ! Initialize
270  !-----------------------------------------------------------------
271  do j=jlo,jhi
272  do i=ilo,ihi
273  aice0(i,j) = c1i
274  aice(i,j) = c0i
275  vice(i,j) = c0i
276  vsno(i,j) = c0i
277  eice(i,j) = c0i
278  esno(i,j) = c0i
279  tsfc(i,j) = c0i
280  enddo
281  enddo
282 
283  !-----------------------------------------------------------------
284  ! Aggregate
285  !-----------------------------------------------------------------
286 
287  icells = 0
288  do j = jlo, jhi
289  do i = ilo, ihi
290  if (tmask(i,j)) then
291  icells = icells + 1
292  indxi(icells) = i
293  indxj(icells) = j
294  endif ! tmask
295  enddo
296  enddo
297 
298  do ni = 1, ncat
299 
300 !DIR$ CONCURRENT !Cray
301 !cdir nodep !NEC
302 !ocl novrec !Fujitsu
303  do ij = 1, icells
304  i = indxi(ij)
305  j = indxj(ij)
306  aice(i,j) = aice(i,j) + aicen(i,j,ni)
307  vice(i,j) = vice(i,j) + vicen(i,j,ni)
308  vsno(i,j) = vsno(i,j) + vsnon(i,j,ni)
309  esno(i,j) = esno(i,j) + esnon(i,j,ni)
310  tsfc(i,j) = tsfc(i,j) + tsfcn(i,j,ni)*aicen(i,j,ni)
311  enddo ! ij
312 
313  do k = 1, nilyr
314 !DIR$ CONCURRENT !Cray
315 !cdir nodep !NEC
316 !ocl novrec !Fujitsu
317  do ij = 1, icells
318  i = indxi(ij)
319  j = indxj(ij)
320  eice(i,j) = eice(i,j) + eicen(i,j,(ni-1)*nilyr+k)
321  enddo
322  enddo ! k
323 
324  enddo ! n
325 
326  !-----------------------------------------------------------------
327  ! Temperature, open water fraction
328  !-----------------------------------------------------------------
329  do j=jlo,jhi
330  do i=ilo,ihi
331  if (aice(i,j) > c0i) then
332  aice0(i,j) = max(c1i - aice(i,j), c0i)
333 ! Tsfc(i,j) = Tsfc(i,j) / aice(i,j)
334 !! ggao change 0616-2008
335  tsfc(i,j) = tsfc(i,j) /(aice(i,j)+epsilon(eps))
336  else
337  tsfc(i,j) = tf(i,j)
338  endif
339  enddo ! i
340  enddo ! j
341 
real(kind=dbl_kind), dimension(:,:), allocatable, save tf
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, target, save eice
Definition: ice_state.f90:82
integer(kind=int_kind), parameter nilyr
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save esnon
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:), allocatable, target, save esno
Definition: ice_state.f90:82
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:), allocatable, target, save vsno
Definition: ice_state.f90:82
integer(kind=int_kind) ilo
Definition: ice_domain.f90:101
integer(kind=int_kind) jhi
Definition: ice_domain.f90:101
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), dimension(:,:), allocatable, target, save aice
Definition: ice_state.f90:82
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save eicen
Definition: ice_state.f90:108
real(kind=dbl_kind), dimension(:,:), allocatable, target, save aice0
Definition: ice_state.f90:105
integer(kind=int_kind), parameter ncat
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:), allocatable, target, save vice
Definition: ice_state.f90:82
real(kind=dbl_kind), parameter c1i
real(kind=dbl_kind), dimension(:,:), allocatable, target, save tsfc
Definition: ice_state.f90:82
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
Here is the caller graph for this function:

◆ aggregate_area()

subroutine ice_itd::aggregate_area ( )

Definition at line 352 of file ice_itd.f90.

352 !
353 ! !DESCRIPTION:
354 !
355 ! Aggregate ice area (but not other state variables) over thickness categories
356 !
357 ! !REVISION HISTORY:
358 !
359 ! authors: William H. Lipscomb, LANL
360 ! modified Jan 2004 by Clifford Chen, Fujitsu
361 !
362 ! !USES:
363 !
364 ! !INPUT/OUTPUT PARAMETERS:
365 !
366 !EOP
367 !
368 
369  use ice_grid, only : tlat,tlon
370 ! ggao
371  integer (kind=int_kind) :: i, j, ni
372 
373  logical (kind=log_kind) :: &
374  outbound ! = .true. if aggregate ice area out of bounds
375 
376 
377  !-----------------------------------------------------------------
378  ! Initialize
379  !-----------------------------------------------------------------
380  do j=jlo,jhi
381  do i=ilo,ihi
382  aice(i,j) = c0i
383  enddo
384  enddo
385 
386  !-----------------------------------------------------------------
387  ! Aggregate
388  !-----------------------------------------------------------------
389  do ni = 1, ncat
390  do j=jlo,jhi
391  do i=ilo,ihi
392  aice(i,j) = aice(i,j) + aicen(i,j,ni)
393  enddo ! i
394  enddo ! j
395  enddo ! n
396 
397  outbound = .false.
398 
399  do j = jlo, jhi
400  do i = ilo, ihi
401 
402  !-----------------------------------------------------------------
403  ! Bug check
404  !-----------------------------------------------------------------
405  if (aice(i,j) > c1i+puny .or. &
406  aice(i,j) < -puny) then
407  outbound = .true.
408  endif
409 
410  !-----------------------------------------------------------------
411  ! open water fraction
412  !-----------------------------------------------------------------
413  aice0(i,j) = max(c1i - aice(i,j), c0i)
414 
415  enddo ! i
416  enddo ! j
417 
418  if ( outbound ) then ! area out of bounds
419  do j = jlo, jhi
420  do i = ilo, ihi
421 
422  if (aice(i,j) > c1i+puny .or. &
423  aice(i,j) < -puny) then
424 ! write(nu_diag,*) ' '
425 ! write(nu_diag,*) 'aggregate ice area out of bounds'
426 ! write(nu_diag,*) 'my_task, i, j, aice:', &
427 ! my_task, i, j, aice(i,j)
428 
429 !! write(nu_diag,*) 'aicen =', aicen(i,j,:)
430 ! call abort_ice('ice_itd: aggregate_area')
431  endif
432 
433  enddo ! i
434  enddo ! j
435  endif ! outbound
436 
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind), dimension(:,:), allocatable tlat
Definition: ice_grid.f90:122
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
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), dimension(:,:), allocatable, target, save aice
Definition: ice_state.f90:82
real(kind=dbl_kind), dimension(:,:), allocatable, target, save aice0
Definition: ice_state.f90:105
integer(kind=int_kind), parameter ncat
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
real(kind=dbl_kind), parameter c1i
real(kind=dbl_kind), dimension(:,:), allocatable tlon
Definition: ice_grid.f90:122
Here is the caller graph for this function:

◆ bound_aggregate()

subroutine ice_itd::bound_aggregate ( )

Definition at line 447 of file ice_itd.f90.

447 !
448 ! !DESCRIPTION:
449 !
450 ! Get ghost cell values for aggregate ice state variables
451 ! NOTE: This subroutine is called only at initialization. It could be
452 ! eliminated if the aggregate variables were defined only in
453 ! physical grid cells.
454 !
455 ! !REVISION HISTORY:
456 !
457 ! authors: William H. Lipscomb, LANL
458 !
459 ! !USES:
460 !
461 ! !INPUT/OUTPUT PARAMETERS:
462 !
463 !EOP
464 !
465  use ice_grid
466 
467 ! call bound (aice0)
468 ! call bound (aice)
469 ! call bound (vice)
470 ! call bound (vsno)
471 ! call bound (eice)
472 ! call bound (esno)
473 ! call bound (Tsfc)
474 
475 ! need inter processors exchange
476 ! ggao
477 
478 
Here is the caller graph for this function:

◆ bound_state()

subroutine ice_itd::bound_state ( )

Definition at line 489 of file ice_itd.f90.

489 !
490 ! !DESCRIPTION:
491 !
492 ! Get ghost cell values for ice state variables in each thickness category
493 !
494 ! !REVISION HISTORY:
495 !
496 ! authors: William H. Lipscomb, LANL
497 !
498 ! !USES:
499 !
500 ! !INPUT/OUTPUT PARAMETERS:
501 !
502 !EOP
503 !
504  use ice_grid
505 
506 ! call bound_narr (ncat,aicen)
507 ! call bound_narr (ncat,vicen)
508 ! call bound_narr (ncat,vsnon)
509 ! call bound_narr (ntilay,eicen)
510 ! call bound_narr (ncat,esnon)
511 ! call bound_narr (ncat,Tsfcn)
512 ! need inter processors exchange
513 ! ggao
514 
515 

◆ check_state()

subroutine ice_itd::check_state ( )

Definition at line 526 of file ice_itd.f90.

526 !
527 ! !DESCRIPTION:
528 !
529 ! Insist that certain fields are monotone.
530 ! Should not be necessary if all is well,
531 ! but best to keep going. Model will not conserve
532 ! energy and water if fields are zeroed here.
533 !
534 ! !REVISION HISTORY:
535 !
536 ! author: C. M. Bitz, UW
537 !
538 ! !USES:
539 !
540  use ice_flux
541 !
542 ! !INPUT/OUTPUT PARAMETERS:
543 !
544 !EOP
545 !
546  integer (kind=int_kind) :: &
547  i, j & ! horizontal indices
548  , ni & ! thickness category index
549  , k ! ice layer index
550 
551  integer (kind=int_kind), dimension (ilo:ihi,jlo:jhi) :: &
552  zerout ! 0=false, 1=true
553 
554  do ni = 1, ncat
555 
556  do j = jlo, jhi
557  do i = ilo, ihi
558  if (aicen(i,j,ni) < puny .or. vicen(i,j,ni) < puny) then
559  zerout(i,j) = 1
560  else
561  zerout(i,j) = 0
562  endif
563  enddo
564  enddo
565 
566  do k = 1, nilyr
567  do j = jlo, jhi
568  do i = ilo, ihi
569  if (eicen(i,j,ilyr1(ni)+k-1) > -puny) zerout(i,j) = 1
570  enddo ! i
571  enddo ! j
572  enddo ! k
573 
574  do j = jlo, jhi
575  do i = ilo, ihi
576 
577  if (zerout(i,j)==1) then
578  aice0(i,j) = aice0(i,j) + aicen(i,j,ni)
579  aicen(i,j,ni) = c0i
580  vicen(i,j,ni) = c0i
581  vsnon(i,j,ni) = c0i
582  esnon(i,j,ni) = c0i
583  tsfcn(i,j,ni) = tf(i,j)
584  elseif (vsnon(i,j,ni) <= puny) then
585  vsnon(i,j,ni) = c0i
586  esnon(i,j,ni) = c0i
587  endif
588 
589  if (vsnon(i,j,ni) > puny) then
590  if (-esnon(i,j,ni)/vsnon(i,j,ni)-lfresh*rhos < eps04) &
591  esnon(i,j,ni) = -vsnon(i,j,ni)*(lfresh*rhos + eps04)
592  endif
593 
594  enddo ! i
595  enddo ! j
596 
597  do k = 1, nilyr
598  do j = jlo, jhi
599  do i = ilo, ihi
600  if (zerout(i,j)==1) eicen(i,j,ilyr1(ni)+k-1) = c0i
601  enddo
602  enddo
603  enddo ! k
604 
605  enddo ! n
606 
607  do j=jlo,jhi
608  do i=ilo,ihi
609  if (aice0(i,j) < puny) aice0(i,j) = c0i
610  enddo
611  enddo
612 
real(kind=dbl_kind), dimension(:,:), allocatable, save tf
Definition: ice_flux.f90:91
integer(kind=int_kind), parameter nilyr
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save esnon
Definition: ice_state.f90:97
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind), parameter rhos
real(kind=dbl_kind), parameter lfresh
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
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), 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), dimension(:,:,:), allocatable, target, save eicen
Definition: ice_state.f90:108
real(kind=dbl_kind), dimension(:,:), allocatable, target, save aice0
Definition: ice_state.f90:105
integer(kind=int_kind), parameter ncat
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vsnon
Definition: ice_state.f90:97
real(kind=dbl_kind), parameter eps04

◆ column_conservation_check()

subroutine ice_itd::column_conservation_check ( real (kind=dbl_kind), dimension (imt_local,jmt_local), intent(in)  x1,
real (kind=dbl_kind), dimension (imt_local,jmt_local), intent(in)  x2,
real (kind=dbl_kind), intent(in)  max_err,
character (len=char_len), intent(in)  fieldid 
)

Definition at line 1243 of file ice_itd.f90.

1243 !
1244 ! !DESCRIPTION:
1245 !
1246 ! For each physical grid cell, check that initial and final values
1247 ! of a conserved field are equal to within a small value.
1248 !
1249 ! !REVISION HISTORY:
1250 !
1251 ! author: William H. Lipscomb, LANL
1252 !
1253 ! !USES:
1254 !
1255 ! !INPUT/OUTPUT PARAMETERS:
1256 !
1257  real (kind=dbl_kind), intent(in) :: &
1258  x1(imt_local,jmt_local) ! initial field
1259 
1260  real (kind=dbl_kind), intent(in) :: &
1261  x2(imt_local,jmt_local) ! final field
1262 
1263  real (kind=dbl_kind), intent(in) :: &
1264  max_err ! max allowed error
1265 
1266  character (len=char_len), intent(in) :: &
1267  fieldid ! field identifier
1268 !
1269 !EOP
1270 !
1271  integer (kind=int_kind) :: &
1272  i, j ! horizontal indices
1273 
1274  logical (kind=log_kind) :: &
1275  conserv_err ! = .true. if conservation check failed
1276 
1277  conserv_err = .false.
1278 
1279  do j = jlo, jhi
1280  do i = ilo, ihi
1281  if (abs(x2(i,j) - x1(i,j)) > max_err) then
1282  conserv_err = .true.
1283  endif
1284  enddo
1285  enddo
1286 
1287  if ( conserv_err ) then
1288  do j = jlo, jhi
1289  do i = ilo, ihi
1290  if (abs(x2(i,j) - x1(i,j)) > max_err) then
1291  write (nu_diag,*) ' '
1292  write (nu_diag,*) 'Conservation error: ', fieldid
1293  write (nu_diag,*) my_task, ':', i, j
1294  write (nu_diag,*) 'Initial value =', x1(i,j)
1295  write (nu_diag,*) 'Final value =', x2(i,j)
1296  write (nu_diag,*) 'Difference =', x2(i,j) - x1(i,j)
1297 ! call abort_ice ('ice: Conservation error')
1298  endif
1299  enddo
1300  enddo
1301  endif
1302 
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
integer(kind=int_kind) ilo
Definition: ice_domain.f90:101
integer(kind=int_kind) jhi
Definition: ice_domain.f90:101
integer(kind=int_kind), save my_task
Definition: ice_domain.f90:95
integer(kind=int_kind) imt_local
Definition: ice_domain.f90:101
integer(kind=int_kind) jmt_local
Definition: ice_domain.f90:101
Here is the caller graph for this function:

◆ column_sum()

subroutine ice_itd::column_sum ( integer (kind=int_kind), intent(in)  nsum,
real (kind=dbl_kind), dimension (imt_local,jmt_local,nsum), intent(in)  xin,
real (kind=dbl_kind), dimension (imt_local,jmt_local), intent(out)  xout 
)

Definition at line 1190 of file ice_itd.f90.

1190 !
1191 ! !DESCRIPTION:
1192 !
1193 ! For each grid cell, sum field over all ice categories.
1194 !
1195 ! !REVISION HISTORY:
1196 !
1197 ! author: William H. Lipscomb, LANL
1198 !
1199 ! !USES:
1200 !
1201 ! !INPUT/OUTPUT PARAMETERS:
1202 !
1203  integer (kind=int_kind), intent(in) :: &
1204  nsum ! number of categories/layers
1205 
1206  real (kind=dbl_kind), intent(in) :: &
1207  xin(imt_local,jmt_local,nsum) ! input field
1208 
1209 
1210  real (kind=dbl_kind), intent(out) :: &
1211  xout(imt_local,jmt_local) ! output field
1212 !
1213 !EOP
1214 !
1215  integer (kind=int_kind) :: &
1216  i, j & ! horizontal indices
1217  , ni ! category/layer index
1218 
1219  do j = 1, jmt_local
1220  do i = 1, imt_local
1221  xout(i,j) = c0i
1222  enddo
1223  enddo
1224 
1225  do ni = 1, nsum
1226  do j = 1, jmt_local
1227  do i = 1, imt_local
1228  xout(i,j) = xout(i,j) + xin(i,j,ni)
1229  enddo ! i
1230  enddo ! j
1231  enddo ! n
1232 
real(kind=dbl_kind), parameter c0i
integer(kind=int_kind) imt_local
Definition: ice_domain.f90:101
integer(kind=int_kind) jmt_local
Definition: ice_domain.f90:101
Here is the caller graph for this function:

◆ init_itd()

subroutine ice_itd::init_itd ( )

Definition at line 109 of file ice_itd.f90.

109 !
110 ! !DESCRIPTION:
111 !
112 ! Initialize area fraction and thickness boundaries for the itd model
113 !
114 ! !REVISION HISTORY:
115 !
116 ! authors: William H. Lipscomb, LANL
117 ! Elizabeth C. Hunke LANL
118 ! C. M. Bitz UW
119 !
120 ! !USES:
121 !
122 ! !INPUT/OUTPUT PARAMETERS:
123 !
124 !EOP
125 !
126  integer (kind=int_kind) :: &
127  ni ! thickness category index
128 
129  real (kind=dbl_kind) :: &
130  cc1, cc2, cc3 & ! parameters for kcatbound = 0
131  , x1
132 
133  real (kind=dbl_kind) :: &
134  rn &! real(n)
135  , rncat &! real(ncat)
136  , d1 & ! parameters for kcatbound = 1 (m)
137  , d2
138 
139  if (ncat == 1 .and. kitd == 1) then
140  write (nu_diag,*) 'Remapping the ITD is not allowed for ncat=1'
141  write (nu_diag,*) 'Use the delta function ITD option instead'
142 ! call abort_ice ('(init_itd) Linear remapping not allowed')
143  endif
144 
145  rncat = real(ncat, kind=dbl_kind)
146  d1 = 3.0_dbl_kind / rncat
147  d2 = 0.5_dbl_kind / rncat
148 
149  !-----------------------------------------------------------------
150  ! Choose category boundaries based on one of two formulas.
151  !
152  ! The first formula (kcatbound = 0) was used in Lipscomb (2001)
153  ! and in CICE versions 3.0 and 3.1.
154  !
155  ! The second formula is more user-friendly in the sense that it
156  ! is easy to obtain round numbers for category boundaries:
157  !
158  ! H(n) = n * [d1 + d2*(n-1)]
159  !
160  ! Default values are d1 = 300/ncat, d2 = 50/ncat.
161  ! For ncat = 5, boundaries in cm are 60, 140, 240, 360, which are
162  ! close to the standard values given by the first formula.
163  ! For ncat = 10, boundaries in cm are 30, 70, 120, 180, 250, 330,
164  ! 420, 520, 630.
165  !-----------------------------------------------------------------
166 
167  if (kcatbound == 0) then ! original scheme
168 
169  if (kitd == 1) then
170  ! linear remapping itd category limits
171  cc1 = c3i/rncat
172  cc2 = c15*cc1
173  cc3 = c3i
174 
175  hin_max(0) = c0i ! minimum ice thickness, m
176  else
177  ! delta function itd category limits
178  cc1 = max(1.1_dbl_kind/rncat,c1i*hi_min)
179  cc2 = c25*cc1
180  cc3 = 2.25_dbl_kind
181 
182  ! hin_max(0) should not be zero
183  ! use some caution in making it less than 0.10
184  hin_max(0) = hi_min ! minimum ice thickness, m
185  endif ! kitd
186 
187  do ni = 1, ncat
188  x1 = real(ni-1,kind=dbl_kind) / rncat
189  hin_max(ni) = hin_max(ni-1) &
190  + cc1 + cc2*(c1i + tanh(cc3*(x1-c1i)))
191  enddo
192 
193  elseif (kcatbound == 1) then ! new scheme
194 
195  hin_max(0) = c0i
196  do ni = 1, ncat
197  rn = real(ni, kind=dbl_kind)
198  hin_max(ni) = rn * (d1 + (rn-c1i)*d2)
199  enddo
200 
201  endif
202 
203  if (my_task == master_task) then
204  write (nu_diag,*) ' '
205  write (nu_diag,*) 'The thickness categories are:'
206  write (nu_diag,*) ' '
207  write (nu_diag,*) 'hin_max(n-1) < Cat n < hin_max(n)'
208  do ni = 1, ncat
209  write (nu_diag,*) hin_max(ni-1),' < Cat ',ni, ' < ',hin_max(ni)
210  enddo
211  write (nu_diag,*) ' '
212  endif
213 
214  !-----------------------------------------------------------------
215  ! vectors identifying first and last layer in each category
216  !-----------------------------------------------------------------
217  ilyr1(1) = 1 ! if nilyr = 4
218  ilyrn(1) = nilyr ! ilyr1 = { 1,5,9 }
219  do ni = 2,ncat ! ilyrn = { 4,8,12} etc
220  ilyr1(ni) = ilyrn(ni-1) + 1
221  ilyrn(ni) = ilyrn(ni-1) + nilyr
222  enddo
223 
integer(kind=int_kind), parameter nilyr
real(sp), dimension(:), allocatable, target d1
Definition: mod_main.f90:1116
integer(kind=int_kind) kcatbound
Definition: ice_itd.f90:62
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind), parameter c3i
integer(kind=int_kind) kitd
Definition: ice_itd.f90:62
real(kind=dbl_kind), parameter c25
integer(kind=int_kind), parameter ncat
integer(kind=int_kind), save my_task
Definition: ice_domain.f90:95
real(kind=dbl_kind), parameter c1i
real(kind=dbl_kind), parameter c15
integer(kind=int_kind), save master_task
Definition: ice_domain.f90:95

◆ rebin()

subroutine ice_itd::rebin ( )

Definition at line 623 of file ice_itd.f90.

623 !
624 ! !DESCRIPTION:
625 !
626 ! Rebins thicknesses into defined categories
627 !
628 ! !REVISION HISTORY:
629 !
630 ! authors: William H. Lipscomb, LANL
631 ! Elizabeth C. Hunke LANL
632 !
633 ! !USES:
634 !
635  use ice_grid
636 !
637 ! !INPUT/OUTPUT PARAMETERS:
638 !
639 !EOP
640 !
641  integer (kind=int_kind) :: &
642  i,j & ! horizontal indices
643  , ni ! category index
644 
645  logical (kind=log_kind) :: &
646  shiftflag ! = .true. if ice must be shifted
647 
648  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat) :: &
649  hicen ! ice thickness for each cat (m)
650 
651  integer (kind=int_kind), dimension(ilo:ihi,jlo:jhi,ncat) :: &
652  donor ! donor category index
653 
654  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat) :: &
655  daice & ! ice area transferred
656  , dvice ! ice volume transferred
657 
658  !-----------------------------------------------------------------
659  ! Compute ice thickness.
660  !-----------------------------------------------------------------
661  do ni = 1, ncat
662  do j = jlo, jhi
663  do i = ilo, ihi
664  if (aicen(i,j,ni) > puny) then
665  hicen(i,j,ni) = vicen(i,j,ni) / aicen(i,j,ni)
666  else
667  hicen(i,j,ni) = c0i
668  endif
669  enddo ! i
670  enddo ! j
671  enddo ! n
672 
673  !-----------------------------------------------------------------
674  ! make sure thickness of cat 1 is at least hin_max(0)
675  !-----------------------------------------------------------------
676  do j = jlo, jhi
677  do i = ilo, ihi
678  if (aicen(i,j,1) > puny) then
679  if (hicen(i,j,1) <= hin_max(0) .and. hin_max(0) > c0i ) then
680  aicen(i,j,1) = vicen(i,j,1) / hin_max(0)
681  hicen(i,j,1) = hin_max(0)
682  endif
683  endif
684  enddo ! i
685  enddo ! j
686 
687  !-----------------------------------------------------------------
688  ! If a category thickness is not in bounds, shift the
689  ! entire area, volume, and energy to the neighboring category
690  !-----------------------------------------------------------------
691 
692  !-----------------------------------------------------------------
693  ! Initialize shift arrays
694  !-----------------------------------------------------------------
695 
696  donor(:,:,:) = 0
697  daice(:,:,:) = c0i
698  dvice(:,:,:) = c0i
699 
700  !-----------------------------------------------------------------
701  ! Move thin categories up
702  !-----------------------------------------------------------------
703 
704  do ni = 1, ncat-1 ! loop over category boundaries
705 
706  !-----------------------------------------------------------------
707  ! identify thicknesses that are too big
708  !-----------------------------------------------------------------
709  shiftflag = .false.
710  do j = jlo, jhi
711  do i = ilo, ihi
712  if (aicen(i,j,ni) > puny .and. &
713  hicen(i,j,ni) > hin_max(ni)) then
714  shiftflag = .true.
715  donor(i,j,ni) = ni
716  daice(i,j,ni) = aicen(i,j,ni)
717  dvice(i,j,ni) = vicen(i,j,ni)
718  endif
719  enddo ! i
720  enddo ! j
721 
722  if (shiftflag) then
723 
724  !-----------------------------------------------------------------
725  ! shift ice between categories
726  !-----------------------------------------------------------------
727  call shift_ice (donor, daice, dvice, hicen)
728 
729  !-----------------------------------------------------------------
730  ! reset shift parameters
731  !-----------------------------------------------------------------
732  do j = jlo, jhi
733  do i = ilo, ihi
734  donor(i,j,ni) = 0
735  daice(i,j,ni) = c0i
736  dvice(i,j,ni) = c0i
737  enddo
738  enddo
739 
740  endif ! shiftflag
741 
742  enddo ! n
743 
744  !-----------------------------------------------------------------
745  ! Move thick categories down
746  !-----------------------------------------------------------------
747 
748  do ni = ncat-1, 1, -1 ! loop over category boundaries
749 
750  !-----------------------------------------------------------------
751  ! identify thicknesses that are too small
752  !-----------------------------------------------------------------
753  shiftflag = .false.
754  do j = jlo, jhi
755  do i = ilo, ihi
756  if (aicen(i,j,ni+1) > puny .and. &
757  hicen(i,j,ni+1) <= hin_max(ni)) then
758  shiftflag = .true.
759  donor(i,j,ni) = ni+1
760  daice(i,j,ni) = aicen(i,j,ni+1)
761  dvice(i,j,ni) = vicen(i,j,ni+1)
762  endif
763  enddo ! i
764  enddo ! j
765 
766  if (shiftflag) then
767 
768  !-----------------------------------------------------------------
769  ! shift ice between categories
770  !-----------------------------------------------------------------
771  call shift_ice (donor, daice, dvice, hicen)
772 
773  !-----------------------------------------------------------------
774  ! reset shift parameters
775  !-----------------------------------------------------------------
776  do j = jlo, jhi
777  do i = ilo, ihi
778  donor(i,j,ni) = 0
779  daice(i,j,ni) = c0i
780  dvice(i,j,ni) = c0i
781  enddo
782  enddo
783 
784  endif ! shiftflag
785 
786  enddo ! n
787 
788 
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
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), dimension(:,:,:), allocatable, target, save vicen
Definition: ice_state.f90:97
integer(kind=int_kind), parameter ncat
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
Here is the call graph for this function:
Here is the caller graph for this function:

◆ reduce_area()

subroutine ice_itd::reduce_area ( real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi), intent(in)  hice1_old,
real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi), intent(inout)  hice1 
)

Definition at line 799 of file ice_itd.f90.

799 !
800 ! !DESCRIPTION:
801 !
802 ! Reduce area when ice melts for special case of ncat=1
803 !
804 ! Use CSM 1.0-like method of reducing ice area
805 ! when melting occurs: assume only half the ice volume
806 ! change goes to thickness decrease, the other half
807 ! to reduction in ice fraction
808 !
809 ! !REVISION HISTORY:
810 !
811 ! authors: C. M. Bitz, UW
812 ! modified by: Elizabeth C. Hunke, LANL
813 !
814 ! !USES:
815 !
816  use ice_grid
817 !
818 ! !INPUT/OUTPUT PARAMETERS:
819 !
820  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi), &
821  intent(in) :: &
822  hice1_old ! old ice thickness for category 1 (m)
823 
824  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi), &
825  intent(inout) :: &
826  hice1 ! new ice thickness for category 1 (m)
827 !
828 !EOP
829 !
830  integer (kind=int_kind) :: &
831  i, j ! horizontal indices
832 
833  real (kind=dbl_kind) :: &
834  hi0 & ! current hi for ice fraction adjustment
835  , dai0 & ! change in aice for ice fraction adjustment
836  , dhi ! hice1 - hice1_old
837 
838  do j=jlo,jhi
839  do i=ilo,ihi
840  if (tmask(i,j)) then
841 
842  !-----------------------------------------------------------------
843  ! make sure thickness of cat 1 is at least hin_max(0)
844  !-----------------------------------------------------------------
845 
846  if (hice1(i,j) <= hin_max(0) .and. hin_max(0) > c0i ) then
847  aicen(i,j,1) = vicen(i,j,1) / hin_max(0)
848  hice1(i,j) = hin_max(0)
849  endif
850 
851  if (aicen(i,j,1) > c0i) then
852  dhi = hice1(i,j) - hice1_old(i,j)
853  if (dhi < c0i) then
854  hi0 = vicen(i,j,1) / aicen(i,j,1)
855  dai0 = vicen(i,j,1) / (hi0-p5*dhi) &
856  - aicen(i,j,1)
857  aicen(i,j,1) = aicen(i,j,1) + dai0
858  endif
859  endif
860 
861  endif ! tmask
862  enddo ! i
863  enddo ! j
864 
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
integer(kind=int_kind) ilo
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter p5
integer(kind=int_kind) jhi
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vicen
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
logical(kind=log_kind), dimension(:,:), allocatable tmask
Definition: ice_grid.f90:164
Here is the caller graph for this function:

◆ shift_ice()

subroutine ice_itd::shift_ice ( integer (kind=int_kind), dimension(ilo:ihi,jlo:jhi,ncat), intent(in)  donor,
real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat), intent(inout)  daice,
real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat), intent(inout)  dvice,
real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat), intent(inout)  hicen 
)

Definition at line 875 of file ice_itd.f90.

875 !
876 ! !DESCRIPTION:
877 !
878 ! Shift ice across category boundaries, conserving area, volume, and
879 ! energy.
880 !
881 ! !REVISION HISTORY:
882 !
883 ! authors: William H. Lipscomb, LANL
884 ! Elizabeth C. Hunke, LANL
885 !
886 ! !USES:
887 !
888  use ice_flux
889  use ice_work, only: worka
890 !
891 ! !INPUT/OUTPUT PARAMETERS:
892 !
893  ! NOTE: Third index of donor, daice, dvice should be ncat-1,
894  ! except that compilers would have trouble when ncat = 1
895  integer (kind=int_kind), dimension(ilo:ihi,jlo:jhi,ncat), &
896  intent(in) :: &
897  donor ! donor category index
898 
899  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat), &
900  intent(inout) :: &
901  daice & ! ice area transferred across boundary
902  , dvice ! ice volume transferred across boundary
903 
904  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat), &
905  intent(inout) :: &
906  hicen ! ice thickness for each cat (m)
907 !
908 !EOP
909 !
910  integer (kind=int_kind) :: &
911  i, j & ! horizontal indices
912  , ni & ! thickness category index
913  , n2 & ! receiver category
914  , n1 & ! donor category
915  , k ! ice layer index
916 
917  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi,ncat) :: &
918  atsfn
919 
920  real (kind=dbl_kind) :: &
921  dvsnow & ! snow volume transferred
922  , desnow & ! snow energy transferred
923  , deice & ! ice energy transferred
924  , datsf ! aicen*Tsfcn transferred
925 
926  integer (kind=int_kind), dimension (1:(ihi-ilo+1)*(jhi-jlo+1)) :: &
927  indxi & ! compressed indices for i/j directions
928  , indxj
929 
930  integer (kind=int_kind) :: &
931  icells & ! number of cells with ice to transfer
932  , ij ! combined i/j horizontal index
933 
934  logical (kind=log_kind) :: &
935  daice_negative & ! true if daice < -puny
936  , dvice_negative & ! true if dvice < -puny
937  , daice_greater_aicen & ! true if daice > aicen
938  , dvice_greater_vicen ! true if dvice > vicen
939 
940 !! ggao 60162008
941  real (kind=dbl_kind) ::eps
942 !! change end
943 
944 
945  !-----------------------------------------------------------------
946  ! Define a variable equal to aicen*Tsfcn
947  !-----------------------------------------------------------------
948  do ni = 1, ncat
949  do j = jlo,jhi
950  do i = ilo,ihi
951  atsfn(i,j,ni) = aicen(i,j,ni)*tsfcn(i,j,ni)
952  enddo ! i
953  enddo ! j
954  enddo ! n
955 
956  !-----------------------------------------------------------------
957  ! Check for daice or dvice out of range, allowing for roundoff error
958  !-----------------------------------------------------------------
959 
960  do ni = 1, ncat-1
961 
962  daice_negative = .false.
963  dvice_negative = .false.
964  daice_greater_aicen = .false.
965  dvice_greater_vicen = .false.
966 
967  do j = jlo,jhi
968  do i = ilo,ihi
969 
970  if (donor(i,j,ni) > 0) then
971  n1 = donor(i,j,ni)
972 
973  if (daice(i,j,ni) < c0i) then
974  if (daice(i,j,ni) > -puny*aicen(i,j,n1)) then
975  daice(i,j,ni) = c0i ! shift no ice
976  dvice(i,j,ni) = c0i
977  else
978  daice_negative = .true.
979  endif
980  endif
981 
982  if (dvice(i,j,ni) < c0i) then
983  if (dvice(i,j,ni) > -puny*vicen(i,j,n1)) then
984  daice(i,j,ni) = c0i ! shift no ice
985  dvice(i,j,ni) = c0i
986  else
987  dvice_negative = .true.
988  endif
989  endif
990 
991  if (daice(i,j,ni) > aicen(i,j,n1)*(c1i-puny)) then
992  if (daice(i,j,ni) < aicen(i,j,n1)*(c1i+puny)) then
993  daice(i,j,ni) = aicen(i,j,n1)
994  dvice(i,j,ni) = vicen(i,j,n1)
995  else
996  daice_greater_aicen = .true.
997  endif
998  endif
999 
1000  if (dvice(i,j,ni) > vicen(i,j,n1)*(c1i-puny)) then
1001  if (dvice(i,j,ni) < vicen(i,j,n1)*(c1i+puny)) then
1002  daice(i,j,ni) = aicen(i,j,n1)
1003  dvice(i,j,ni) = vicen(i,j,n1)
1004  else
1005  dvice_greater_vicen = .true.
1006  endif
1007  endif
1008 
1009  endif ! donor > 0
1010  enddo ! i
1011  enddo ! j
1012 
1013  !-----------------------------------------------------------------
1014  ! error messages
1015  !-----------------------------------------------------------------
1016 
1017  if (daice_negative) then
1018  do j = jlo,jhi
1019  do i = ilo,ihi
1020  if (donor(i,j,ni) > 0 .and. &
1021  daice(i,j,ni) <= -puny*aicen(i,j,n1)) then
1022  write(nu_diag,*) my_task,':',i,j, &
1023  'ITD Neg daice =',daice(i,j,ni),' boundary',ni
1024 ! call abort_ice ('ice: ITD Neg daice')
1025  endif
1026  enddo
1027  enddo
1028  endif
1029 
1030  if (dvice_negative) then
1031  do j = jlo,jhi
1032  do i = ilo,ihi
1033  if (donor(i,j,ni) > 0 .and. &
1034  dvice(i,j,ni) <= -puny*vicen(i,j,n1)) then
1035  write(nu_diag,*) my_task,':',i,j, &
1036  'ITD Neg dvice =',dvice(i,j,ni),' boundary',ni
1037 ! call abort_ice ('ice: ITD Neg dvice')
1038  endif
1039  enddo
1040  enddo
1041  endif
1042 
1043  if (daice_greater_aicen) then
1044  do j = jlo,jhi
1045  do i = ilo,ihi
1046  if (donor(i,j,ni) > 0) then
1047  n1 = donor(i,j,ni)
1048  if (daice(i,j,ni) >= aicen(i,j,n1)*(c1i+puny)) then
1049  write(nu_diag,*) my_task,':',i,j, &
1050  'ITD daice > aicen, cat',n1
1051  write(nu_diag,*) my_task,':',i,j, &
1052  'daice =', daice(i,j,ni), &
1053  'aicen =', aicen(i,j,n1)
1054 ! call abort_ice ('ice: ITD daice > aicen')
1055  endif
1056  endif
1057  enddo
1058  enddo
1059  endif
1060 
1061  if (dvice_greater_vicen) then
1062  do j = jlo,jhi
1063  do i = ilo,ihi
1064  if (donor(i,j,ni) > 0) then
1065  n1 = donor(i,j,ni)
1066  if (dvice(i,j,ni) >= vicen(i,j,n1)*(c1i+puny)) then
1067  write(nu_diag,*) my_task,':',i,j, &
1068  'ITD dvice > vicen, cat',n1
1069  write(nu_diag,*) my_task,':',i,j, &
1070  'dvice =', dvice(i,j,ni), &
1071  'vicen =', vicen(i,j,n1)
1072 ! call abort_ice ('ice: ITD dvice > vicen')
1073  endif
1074  endif
1075  enddo
1076  enddo
1077  endif
1078 
1079  !-----------------------------------------------------------------
1080  ! transfer volume and energy between categories
1081  !-----------------------------------------------------------------
1082 
1083  icells = 0
1084  do j = jlo, jhi
1085  do i = ilo, ihi
1086  if (daice(i,j,ni) > c0i) then ! daice(n) can be < puny
1087  icells = icells + 1
1088  indxi(icells) = i
1089  indxj(icells) = j
1090  endif ! tmask
1091  enddo
1092  enddo
1093 
1094 !DIR$ CONCURRENT !Cray
1095 !cdir nodep !NEC
1096 !ocl novrec !Fujitsu
1097  do ij = 1, icells
1098  i = indxi(ij)
1099  j = indxj(ij)
1100 
1101  n1 = donor(i,j,ni)
1102  !worka(i,j) = dvice(i,j,ni) / vicen(i,j,n1)
1103 !! ggao change
1104  worka(i,j) = dvice(i,j,ni) / (vicen(i,j,n1)+epsilon(eps))
1105 !! ggao change end 0616-2008
1106 
1107  if (n1 == ni) then
1108  n2 = n1+1
1109  else ! n1 = n+1
1110  n2 = ni
1111  endif
1112 
1113  aicen(i,j,n1) = aicen(i,j,n1) - daice(i,j,ni)
1114  aicen(i,j,n2) = aicen(i,j,n2) + daice(i,j,ni)
1115  vicen(i,j,n1) = vicen(i,j,n1) - dvice(i,j,ni)
1116  vicen(i,j,n2) = vicen(i,j,n2) + dvice(i,j,ni)
1117 
1118  dvsnow = vsnon(i,j,n1) * worka(i,j)
1119  vsnon(i,j,n1) = vsnon(i,j,n1) - dvsnow
1120  vsnon(i,j,n2) = vsnon(i,j,n2) + dvsnow
1121 
1122  datsf = daice(i,j,ni)*tsfcn(i,j,n1)
1123  atsfn(i,j,n1) = atsfn(i,j,n1) - datsf
1124  atsfn(i,j,n2) = atsfn(i,j,n2) + datsf
1125 
1126  desnow = esnon(i,j,n1) * worka(i,j)
1127  esnon(i,j,n1) = esnon(i,j,n1) - desnow
1128  esnon(i,j,n2) = esnon(i,j,n2) + desnow
1129 
1130  enddo ! ij
1131 
1132  do k = 1, nilyr
1133 !DIR$ CONCURRENT !Cray
1134 !cdir nodep !NEC
1135 !ocl novrec !Fujitsu
1136  do ij = 1, icells
1137  i = indxi(ij)
1138  j = indxj(ij)
1139 
1140  n1 = donor(i,j,ni)
1141  if (n1 == ni) then
1142  n2 = n1+1
1143  else ! n1 = n+1
1144  n2 = ni
1145  endif
1146 
1147  deice = eicen(i,j,ilyr1(n1)+k-1) * worka(i,j)
1148  eicen(i,j,ilyr1(n1)+k-1) = &
1149  eicen(i,j,ilyr1(n1)+k-1) - deice
1150  eicen(i,j,ilyr1(n2)+k-1) = &
1151  eicen(i,j,ilyr1(n2)+k-1) + deice
1152  enddo ! ij
1153  enddo ! k
1154 
1155  enddo ! boundaries, 1 to ncat-1
1156 
1157  !-----------------------------------------------------------------
1158  ! Update ice thickness and temperature
1159  !-----------------------------------------------------------------
1160 
1161  do ni = 1, ncat
1162  do j = jlo,jhi
1163  do i = ilo,ihi
1164  if (aicen(i,j,ni) > puny) then
1165 ! hicen(i,j,ni) = vicen(i,j,ni) / aicen(i,j,ni)
1166 ! Tsfcn(i,j,ni) = aTsfn(i,j,ni) / aicen(i,j,ni)
1167 !! ggao change
1168  hicen(i,j,ni) = vicen(i,j,ni) /( aicen(i,j,ni)+epsilon(eps))
1169  tsfcn(i,j,ni) = atsfn(i,j,ni) /( aicen(i,j,ni)+epsilon(eps))
1170 !! ggao change end 06162008
1171 
1172  else
1173  hicen(i,j,ni) = c0i
1174  tsfcn(i,j,ni) = tf(i,j)
1175  endif
1176  enddo ! i
1177  enddo ! j
1178  enddo ! n
1179 
real(kind=dbl_kind), dimension(:,:), allocatable, save tf
Definition: ice_flux.f90:91
integer(kind=int_kind), parameter nilyr
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save esnon
Definition: ice_state.f90:97
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
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), dimension(:,:,:), allocatable, target, save tsfcn
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:), allocatable worka
Definition: ice_work.f90:61
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vicen
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save eicen
Definition: ice_state.f90:108
integer(kind=int_kind), parameter ncat
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
integer(kind=int_kind), save my_task
Definition: ice_domain.f90:95
real(kind=dbl_kind), parameter c1i
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vsnon
Definition: ice_state.f90:97
Here is the caller graph for this function:

◆ zap_small_areas()

subroutine ice_itd::zap_small_areas ( )

Definition at line 1313 of file ice_itd.f90.

1313 !
1314 ! !DESCRIPTION:
1315 !
1316 ! For each ice category in each grid cell, remove ice if the fractional
1317 ! area is less than puny.
1318 !
1319 ! !REVISION HISTORY:
1320 !
1321 ! author: William H. Lipscomb, LANL
1322 ! Nov 2003: Modified by Julie Schramm to conserve volume and energy
1323 ! Sept 2004: Modified by William Lipscomb; replaced normalize_state with
1324 ! additions to local freshwater, salt, and heat fluxes
1325 !
1326 !
1327 ! !USES:
1328 !
1329  use ice_flux
1330 ! use ice_calendar, only: dt
1331  use ice_calendar, only: dtice
1332 
1333 !
1334 ! !INPUT/OUTPUT PARAMETERS:
1335 !
1336 !EOP
1337 !
1338  integer (kind=int_kind) :: &
1339  i,j & ! horizontal indices
1340  , ni & ! ice category index
1341  , k & ! ice layer index
1342  , icells & ! number of cells with ice to zap
1343  , ij ! combined i/j horizontal index
1344 
1345  integer (kind=int_kind), dimension (1:(ihi-ilo+1)*(jhi-jlo+1)) :: &
1346  indxi & ! compressed indices for i/j directions
1347  , indxj
1348 
1349  real (kind=dbl_kind) :: &
1350  xtmp ! temporary variable
1351 
1352  do ni = 1, ncat
1353 
1354  !-----------------------------------------------------------------
1355  ! Count categories to be zapped.
1356  ! Abort model in case of negative area.
1357  !-----------------------------------------------------------------
1358 
1359  icells = 0
1360  do j = jlo, jhi
1361  do i = ilo, ihi
1362  if (aicen(i,j,ni) < -puny) then
1363  write (nu_diag,*) 'Negative ice area: i, j, n:', i, j, ni
1364  write (nu_diag,*) 'aicen =', aicen(i,j,ni)
1365 ! call abort_ice ('zap: negative ice area')
1366  elseif ((aicen(i,j,ni) >= -puny .and. aicen(i,j,ni) < c0i) &
1367  .or. &
1368  (aicen(i,j,ni) > c0i .and. aicen(i,j,ni) <= puny)) then
1369  icells = icells + 1
1370  indxi(icells) = i
1371  indxj(icells) = j
1372  endif
1373  enddo
1374  enddo
1375 
1376  !-----------------------------------------------------------------
1377  ! Zap ice energy and use ocean heat to melt ice
1378  !-----------------------------------------------------------------
1379 
1380  do k = 1, nilyr
1381 !DIR$ CONCURRENT !Cray
1382 !cdir nodep !NEC
1383 !ocl novrec !Fujitsu
1384  do ij = 1, icells
1385  i = indxi(ij)
1386  j = indxj(ij)
1387 
1388 ! xtmp = eicen(i,j,ilyr1(ni)+k-1) / dt ! < 0
1389  xtmp = eicen(i,j,ilyr1(ni)+k-1) / dtice ! < 0
1390  fhnet(i,j) = fhnet(i,j) + xtmp
1391  fhnet_hist(i,j) = fhnet_hist(i,j) + xtmp
1392  eicen(i,j,ilyr1(ni)+k-1) = c0i
1393 
1394  enddo ! ij
1395  enddo ! k
1396 
1397 !DIR$ CONCURRENT !Cray
1398 !cdir nodep !NEC
1399 !ocl novrec !Fujitsu
1400  do ij = 1, icells
1401  i = indxi(ij)
1402  j = indxj(ij)
1403 
1404  !-----------------------------------------------------------------
1405  ! Zap snow energy and use ocean heat to melt snow
1406  !-----------------------------------------------------------------
1407 
1408 ! xtmp = esnon(i,j,ni) / dt ! < 0
1409  xtmp = esnon(i,j,ni) / dtice ! < 0
1410  fhnet(i,j) = fhnet(i,j) + xtmp
1411  fhnet_hist(i,j) = fhnet_hist(i,j) + xtmp
1412  esnon(i,j,ni) = c0i
1413 
1414  !-----------------------------------------------------------------
1415  ! zap ice and snow volume, add water and salt to ocean
1416  !-----------------------------------------------------------------
1417 
1418 ! xtmp = (rhoi*vicen(i,j,ni) + rhos*vsnon(i,j,ni)) / dt
1419  xtmp = (rhoi*vicen(i,j,ni) + rhos*vsnon(i,j,ni)) / dtice
1420  fresh(i,j) = fresh(i,j) + xtmp
1421  fresh_hist(i,j) = fresh_hist(i,j) + xtmp
1422 
1423 ! xtmp = rhoi*vicen(i,j,n)*ice_ref_salinity*p001/ dt
1424  xtmp = rhoi*vicen(i,j,ni)*ice_ref_salinity*p001/ dtice
1425  fsalt(i,j) = fsalt(i,j) + xtmp
1426  fsalt_hist(i,j) = fsalt_hist(i,j) + xtmp
1427 
1428  aice0(i,j) = aice0(i,j) + aicen(i,j,ni)
1429  aicen(i,j,ni) = c0i
1430  vicen(i,j,ni) = c0i
1431  vsnon(i,j,ni) = c0i
1432  tsfcn(i,j,ni) = tf(i,j)
1433 
1434  enddo ! ij
1435  enddo ! n
1436 
real(kind=dbl_kind), dimension(:,:), allocatable, save tf
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save fsalt_hist
Definition: ice_flux.f90:140
integer(kind=int_kind), parameter nilyr
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save esnon
Definition: ice_state.f90:97
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind), parameter rhos
real(kind=dbl_kind), parameter ice_ref_salinity
real(kind=dbl_kind), dimension(:,:), allocatable, save fhnet_hist
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save fresh
Definition: ice_flux.f90:120
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
integer(kind=int_kind) ilo
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter puny
real(kind=dbl_kind) dtice
integer(kind=int_kind) jhi
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save tsfcn
Definition: ice_state.f90:97
real(kind=dbl_kind), parameter p001
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vicen
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:), allocatable, save fresh_hist
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save fsalt
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save eicen
Definition: ice_state.f90:108
real(kind=dbl_kind), dimension(:,:), allocatable, target, save aice0
Definition: ice_state.f90:105
integer(kind=int_kind), parameter ncat
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
real(kind=dbl_kind), parameter rhoi
real(kind=dbl_kind), dimension(:,:), allocatable, save fhnet
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save vsnon
Definition: ice_state.f90:97
Here is the caller graph for this function:

Variable Documentation

◆ hi_min

real (kind=dbl_kind), parameter ice_itd::hi_min = p01

Definition at line 71 of file ice_itd.f90.

71  real (kind=dbl_kind), parameter :: &
72  hi_min = p01 ! minimum ice thickness allowed (m)
real(kind=dbl_kind), parameter p01

◆ hin_max

real (kind=dbl_kind), dimension(0:ncat) ice_itd::hin_max

Definition at line 74 of file ice_itd.f90.

74  real (kind=dbl_kind) :: &
75  hin_max(0:ncat) ! category limits (m)
integer(kind=int_kind), parameter ncat

◆ ilyr1

integer (kind=int_kind), dimension (ncat) ice_itd::ilyr1

Definition at line 62 of file ice_itd.f90.

◆ ilyrn

integer (kind=int_kind), dimension (ncat) ice_itd::ilyrn

Definition at line 62 of file ice_itd.f90.

◆ kcatbound

integer (kind=int_kind) ice_itd::kcatbound

Definition at line 62 of file ice_itd.f90.

◆ kitd

integer (kind=int_kind) ice_itd::kitd

Definition at line 62 of file ice_itd.f90.

62  integer (kind=int_kind) :: &
63  kitd & ! type of itd conversions
64  ! 0 = delta function
65  ! 1 = linear remap
66  , kcatbound & ! 0 = old category boundary formula
67  ! 1 = new formula giving round numbers
68  , ilyr1(ncat) & ! position of the top layer in each cat
69  , ilyrn(ncat) ! position of the bottom layer in each cat
integer(kind=int_kind) kcatbound
Definition: ice_itd.f90:62
integer(kind=int_kind) kitd
Definition: ice_itd.f90:62
integer(kind=int_kind), parameter ncat