ROMS
Loading...
Searching...
No Matches
mod_ice Module Reference

Data Types

type  t_ice
 
type  t_ice_avg
 
type  t_ice_lobc
 

Functions/Subroutines

subroutine, public allocate_ice (ng, lbi, ubi, lbj, ubj, ice_kernel)
 
subroutine, public deallocate_ice (ng)
 
subroutine, public initialize_ice (ng, tile, model)
 

Variables

integer idaice
 
integer idaicl
 
integer idhage
 
integer idhice
 
integer idhicl
 
integer idhmel
 
integer idhsno
 
integer idiage
 
integer idiofv
 
integer idiomf
 
integer idiomt
 
integer idisst
 
integer idisxx
 
integer idisxy
 
integer idisyy
 
integer ids0mk
 
integer idt0mk
 
integer idtice
 
integer iduice
 
integer iduicl
 
integer iduier
 
integer idvice
 
integer idvicl
 
integer idvinr
 
integer idwdiv
 
integer idw_ai
 
integer idw_ao
 
integer idw_fr
 
integer idw_io
 
integer idw_ro
 
integer, parameter nices = 19
 
integer, dimension(nicesisice
 
integer, parameter isaice = 1
 
integer, parameter ishice = 2
 
integer, parameter ishmel = 3
 
integer, parameter ishsno = 4
 
integer, parameter isiage = 5
 
integer, parameter isisxx = 6
 
integer, parameter isisxy = 7
 
integer, parameter isisyy = 8
 
integer, parameter istice = 9
 
integer, parameter isuice = 10
 
integer, parameter isvice = 11
 
integer, parameter isenth = 12
 
integer, parameter ishage = 13
 
integer, parameter isuevp = 14
 
integer, parameter isvevp = 15
 
integer, parameter isiphy = 16
 
integer, parameter isino3 = 17
 
integer, parameter isinh4 = 18
 
integer, parameter isilog = 19
 
integer, dimension(nicesibice
 
integer, dimension(4, nicesiceobc
 
integer, parameter nicef = 24
 
integer, dimension(nicefifice
 
integer, parameter icaius = 1
 
integer, parameter icaivs = 2
 
integer, parameter icbvis = 3
 
integer, parameter ichsse = 4
 
integer, parameter iciofv = 5
 
integer, parameter iciomf = 6
 
integer, parameter iciomt = 7
 
integer, parameter iciovs = 8
 
integer, parameter icisst = 9
 
integer, parameter icpgrd = 10
 
integer, parameter icpice = 11
 
integer, parameter icqcon = 12
 
integer, parameter icqrhs = 13
 
integer, parameter icsvis = 14
 
integer, parameter ics0mk = 15
 
integer, parameter ict0mk = 16
 
integer, parameter icuavg = 17
 
integer, parameter icvavg = 18
 
integer, parameter icwdiv = 19
 
integer, parameter icw_ai = 20
 
integer, parameter icw_ao = 21
 
integer, parameter icw_fr = 22
 
integer, parameter icw_io = 23
 
integer, parameter icw_ro = 24
 
logical, dimension(:,:), allocatable licefavg
 
logical, dimension(:,:), allocatable licesavg
 
integer, dimension(:), allocatable ievp
 
integer, dimension(:), allocatable nevp
 
integer, dimension(:), allocatable dtice
 
integer, dimension(:), allocatable dtevp
 
real(r8), dimension(:), allocatable airrho
 
real(r8), dimension(:), allocatable icerho
 
real(r8), dimension(:), allocatable snowdryrho
 
real(r8), dimension(:), allocatable snowwetrho
 
real(r8), dimension(:), allocatable cd_ai
 
real(r8), dimension(:), allocatable cd_io
 
real(r8), dimension(:), allocatable astrength
 
real(r8), dimension(:), allocatable pstar
 
real(r8), dimension(:), allocatable min_ai
 
real(r8), dimension(:), allocatable max_ai
 
real(r8), dimension(:), allocatable min_hi
 
real(r8), dimension(:), allocatable max_hmelt
 
real(r8), dimension(:), allocatable zetamin
 
real(r8), dimension(:), allocatable zetamax
 
real(r8), dimension(:), allocatable stressang
 
real(r8), dimension(:), allocatable ellip_sq
 
real(r8ice_emiss
 
real(r8spec_heat_air
 
real(r8trans_coeff
 
real(r8sublimation
 
type(t_ice), dimension(:), allocatable ice
 
type(t_ice_lobc), dimension(:,:), allocatable ice_lobc
 
type(t_ice_avg), dimension(:,:), allocatable ice_favg
 
type(t_ice_avg), dimension(:,:), allocatable ice_savg
 

Function/Subroutine Documentation

◆ allocate_ice()

subroutine, public mod_ice::allocate_ice ( integer, intent(in) ng,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
logical, intent(in) ice_kernel )

Definition at line 325 of file ice_mod.h.

326!
327!=======================================================================
328! !
329! This routine allocates either the derived-type ice model variables !
330! (ice_kernel=TRUE) or module parameters (ice_kernel=FALSE). !
331! !
332!=======================================================================
333!
334 USE mod_param, ONLY : dmem, lbc, ngrids
335 USE mod_scalars, ONLY : iwest, ieast, isouth, inorth
336!
337! Imported variable declarations.
338!
339 logical, intent(in) :: ice_kernel
340!
341 integer, intent(in) :: ng, LBi, UBi, LBj, UBj
342!
343! Local variable declarations.
344!
345 integer :: i
346!
347 real(r8) :: size2d, Xsize, Ysize
348!
349 size2d=real((ubi-lbi+1)*(ubj-lbj+1),r8)
350 xsize =real(ubi-lbi,r8)
351 ysize =real(ubj-lbj,r8)
352
353!
354!-----------------------------------------------------------------------
355! Allocate ice model parameters.
356!-----------------------------------------------------------------------
357!
358 IF (.not.ice_kernel) THEN
359
360#ifdef AVERAGES
361 IF (.not.allocated(licefavg)) &
362 allocate ( licefavg(nicef,ngrids) )
363
364 IF (.not.allocated(licesavg)) &
365 allocate ( licesavg(nices,ngrids) )
366#endif
367 IF (.not.allocated(ievp)) &
368 & allocate ( ievp(ngrids) )
369
370 IF (.not.allocated(nevp)) &
371 & allocate ( nevp(ngrids) )
372
373 IF (.not.allocated(dtice)) &
374 & allocate ( dtice(ngrids) )
375
376 IF (.not.allocated(dtevp)) &
377 & allocate ( dtevp(ngrids) )
378
379 IF (.not.allocated(airrho)) &
380 & allocate ( airrho(ngrids) )
381
382 IF (.not.allocated(icerho)) &
383 & allocate ( icerho(ngrids) )
384
385 IF (.not.allocated(snowdryrho)) &
386 & allocate ( snowdryrho(ngrids) )
387
388 IF (.not.allocated(snowwetrho)) &
389 & allocate ( snowwetrho(ngrids) )
390
391 IF (.not.allocated(cd_ai)) &
392 & allocate ( cd_ai(ngrids) )
393
394 IF (.not.allocated(cd_io)) &
395 & allocate ( cd_io(ngrids) )
396
397 IF (.not.allocated(astrength)) &
398 & allocate ( astrength(ngrids) )
399
400 IF (.not.allocated(pstar)) &
401 & allocate ( pstar(ngrids) )
402
403 IF (.not.allocated(min_ai)) &
404 & allocate ( min_ai(ngrids) )
405
406 IF (.not.allocated(max_ai)) &
407 & allocate ( max_ai(ngrids) )
408
409 IF (.not.allocated(min_hi)) &
410 & allocate ( min_hi(ngrids) )
411
412 IF (.not.allocated(max_hmelt)) &
413 & allocate ( max_hmelt(ngrids) )
414
415 IF (.not.allocated(zetamin)) &
416 & allocate ( zetamin(ngrids) )
417
418 IF (.not.allocated(zetamax)) &
419 & allocate ( zetamax(ngrids) )
420
421 IF (.not.allocated(stressang)) &
422 & allocate ( stressang(ngrids) )
423
424 IF (.not.allocated(ellip_sq)) &
425 & allocate ( ellip_sq(ngrids) )
426
427 END IF
428!
429!-----------------------------------------------------------------------
430! Allocate derived-type structure ice model kernel variables.
431!-----------------------------------------------------------------------
432!
433 IF (ice_kernel) THEN
434 IF (ng.eq.1) allocate ( ice(ngrids) )
435!
436! Nonlinear ice model state variables (Si) and internal kernel field
437! (Fi) arrays.
438!
439 allocate ( ice(ng) % Si(lbi:ubi,lbj:ubj,2,nices) )
440 dmem(ng)=dmem(ng)+2.0_r8*real(nices,r8)*size2d
441!
442 allocate ( ice(ng) % Fi(lbi:ubi,lbj:ubj,nicef) )
443 dmem(ng)=dmem(ng)+real(nicef,r8)*size2d
444 END IF
445!
446!-----------------------------------------------------------------------
447! Allocate derived-type structure ice model lateral boundary variables.
448!-----------------------------------------------------------------------
449!
450 IF (ice_kernel) THEN
451 IF (ng.eq.1) allocate ( ice_lobc(nices,ngrids) )
452!
453! Nonlinear ice model 'state' lateral boundary arrays.
454!
455 DO i=1,nices
456 IF (lbc(iwest,ibice(i),ng)%acquire) THEN
457 allocate ( ice_lobc(i,ng) % ice_west(lbj:ubj) )
458 dmem(ng)=dmem(ng)+ysize
459
460 allocate ( ice_lobc(i,ng) % iceG_west(lbj:ubj,2) )
461 dmem(ng)=dmem(ng)+2.0_r8*ysize
462 END IF
463!
464 IF (lbc(ieast,ibice(i),ng)%acquire) THEN
465 allocate ( ice_lobc(i,ng) % ice_east(lbj:ubj) )
466 dmem(ng)=dmem(ng)+ysize
467
468 allocate ( ice_lobc(i,ng) % iceG_east(lbj:ubj,2) )
469 dmem(ng)=dmem(ng)+2.0_r8*ysize
470 END IF
471!
472 IF (lbc(isouth,ibice(i),ng)%acquire) THEN
473 allocate ( ice_lobc(i,ng) % ice_south(lbi:ubi) )
474 dmem(ng)=dmem(ng)+xsize
475
476 allocate ( ice_lobc(i,ng) % iceG_south(lbi:ubi,2) )
477 dmem(ng)=dmem(ng)+2.0_r8*xsize
478 END IF
479!
480 IF (lbc(inorth,ibice(i),ng)%acquire) THEN
481 allocate ( ice_lobc(i,ng) % ice_north(lbi:ubi) )
482 dmem(ng)=dmem(ng)+xsize
483
484 allocate ( ice_lobc(i,ng) % iceG_north(lbi:ubi,2) )
485 dmem(ng)=dmem(ng)+2.0_r8*xsize
486 END IF
487 END DO
488 END IF
489
490#ifdef AVERAGES
491!
492!-----------------------------------------------------------------------
493! Allocate derived-type structure ice model state and internal arrays
494! time-averaged variables.
495!-----------------------------------------------------------------------
496!
497 IF (ice_kernel) THEN
498 IF (ng.eq.1) THEN
499 allocate ( ice_favg(nicef,ngrids) )
500 allocate ( ice_savg(nices,ngrids) )
501 END IF
502!
503! Time-averaged internal fields. Only fields with metadata (iFice > 0)
504! are allocated and processed.
505!
506 DO i=1,nicef
507 IF (ifice(i).gt.0) THEN
508 IF (licefavg(i,ng)) THEN
509 allocate ( ice_favg(i,ng) % var(lbi:ubi,lbj:ubj) )
510 dmem(ng)=dmem(ng)+size2d
511 END IF
512 END IF
513 END DO
514!
515! Time-averaged state fields. Only fields with metadata (iSice > 0)
516! are allocated and processed.
517!
518 DO i=1,nices
519 IF (isice(i).gt.0) THEN
520 IF (licesavg(i,ng)) THEN
521 allocate ( ice_savg(i,ng) % var(lbi:ubi,lbj:ubj) )
522 dmem(ng)=dmem(ng)+size2d
523 END IF
524 END IF
525 END DO
526 END IF
527#endif
528!
529 RETURN
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
integer ngrids
Definition mod_param.F:113
integer, parameter iwest
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth

References airrho, astrength, cd_ai, cd_io, mod_param::dmem, dtevp, dtice, ellip_sq, ibice, ice, ice_favg, ice_lobc, ice_savg, icerho, mod_scalars::ieast, ievp, ifice, mod_scalars::inorth, isice, mod_scalars::isouth, mod_scalars::iwest, mod_param::lbc, licefavg, licesavg, max_ai, max_hmelt, min_ai, min_hi, nevp, mod_param::ngrids, nicef, nices, pstar, mod_kinds::r8, snowdryrho, snowwetrho, stressang, zetamax, and zetamin.

Referenced by read_icepar(), and mod_arrays::roms_allocate_arrays().

Here is the caller graph for this function:

◆ deallocate_ice()

subroutine, public mod_ice::deallocate_ice ( integer, intent(in) ng)

Definition at line 532 of file ice_mod.h.

533!
534!=======================================================================
535! !
536! This routine deallocates all variables in the module for all nested !
537! grids. !
538! !
539!=======================================================================
540!
541 USE mod_param, ONLY : ngrids
542#ifdef SUBOBJECT_DEALLOCATION
543 USE mod_param, ONLY : lbc
544 USE mod_scalars, ONLY : iwest, ieast, isouth, inorth
545!
546 USE destroy_mod, ONLY : destroy
547#endif
548!
549! Imported variable declarations.
550!
551 integer, intent(in) :: ng
552!
553! Local variable declarations.
554!
555#ifdef SUBOBJECT_DEALLOCATION
556 integer :: i
557!
558#endif
559 character (len=*), parameter :: MyFile = &
560 & __FILE__//", deallocate_ice"
561
562#ifdef SUBOBJECT_DEALLOCATION
563!
564!-----------------------------------------------------------------------
565! Deallocate each variable in the derived-type T_FORCES structure
566! separately.
567!-----------------------------------------------------------------------
568!
569! Nonlinear ice model state variables 'Si' and internal kernel arrays
570! 'Fi'.
571!
572 IF (.not.destroy(ng, ice(ng)%Si, myfile, &
573 & __line__, 'ICE(ng)%Si')) RETURN
574!
575 IF (.not.destroy(ng, ice(ng)%Fi, myfile, &
576 & __line__, 'ICE(ng)%Fi')) RETURN
577!
578! Nonlinear ice model state lateral boundary arrays.
579!
580 DO i=1,nices
581 IF (lbc(iwest,ibice(i),ng)%acquire) THEN
582 IF (.not.destroy(ng, ice_lobc(i,ng)%ice_west, myfile, &
583 & __line__, 'ICE_LOBC(i,ng)%ice_west')) RETURN
584 IF (.not.destroy(ng, ice_lobc(i,ng)%iceG_west, myfile, &
585 & __line__, 'ICE_LOBC(i,ng)%iceG_west')) RETURN
586 END IF
587!
588 IF (lbc(ieast,ibice(i),ng)%acquire) THEN
589 IF (.not.destroy(ng, ice_lobc(i,ng)%ice_east, myfile, &
590 & __line__, 'ICE_LOBC(i,ng)%ice_west')) RETURN
591 IF (.not.destroy(ng, ice_lobc(i,ng)%iceG_east, myfile, &
592 & __line__, 'ICE_LOBC(i,ng)%iceG_west')) RETURN
593 END IF
594!
595 IF (lbc(isouth,ibice(i),ng)%acquire) THEN
596 IF (.not.destroy(ng, ice_lobc(i,ng)%ice_south, myfile, &
597 & __line__, 'ICE_LOBC(i,ng)%ice_south')) RETURN
598 IF (.not.destroy(ng, ice_lobc(i,ng)%iceG_south, myfile, &
599 & __line__,'ICE_LOBC(i,ng)%iceG_south')) RETURN
600 END IF
601!
602 IF (lbc(inorth,ibice(i),ng)%acquire) THEN
603 IF (.not.destroy(ng, ice_lobc(i,ng)%ice_north, myfile, &
604 & __line__, 'ICE_LOBC(i,ng)%ice_north')) RETURN
605 IF (.not.destroy(ng, ice_lobc(i,ng)%iceG_south, myfile, &
606 & __line__,'ICE_LOBC(i,ng)%iceG_north')) RETURN
607 END IF
608 END DO
609
610# ifdef AVERAGES
611!
612!-----------------------------------------------------------------------
613! Time-averaged ice model state and internal variables.
614!-----------------------------------------------------------------------
615!
616 DO i=1,nicef
617 IF (licefavg(i,ng)) THEN
618 IF (.not.destroy(ng, ice_favg(i,ng)%var, myfile, &
619 & __line__, 'ICE_FAVG(i,ng)%var')) RETURN
620 END IF
621 END DO
622!
623 DO i=1,nices
624 IF (licesavg(i,ng)) THEN
625 IF (.not.destroy(ng, ice_savg(i,ng)%var, myfile, &
626 & __line__, 'ICE_SAVG(i,ng)%var')) RETURN
627 END IF
628 END DO
629# endif
630#endif
631!
632!-----------------------------------------------------------------------
633! Deallocate derived-type ICE and ICE_LOBC structures.
634!-----------------------------------------------------------------------
635!
636 IF (allocated(ice)) deallocate ( ice )
637 IF (allocated(ice_lobc)) deallocate ( ice_lobc )
638#ifdef AVERAGES
639 IF (allocated(ice_favg)) deallocate ( ice_favg )
640 IF (allocated(ice_savg)) deallocate ( ice_savg )
641#endif
642!
643 RETURN

References ibice, ice, ice_favg, ice_lobc, ice_savg, mod_scalars::ieast, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_param::lbc, licefavg, licesavg, mod_param::ngrids, nicef, and nices.

Referenced by mod_arrays::roms_deallocate_arrays().

Here is the caller graph for this function:

◆ initialize_ice()

subroutine, public mod_ice::initialize_ice ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model )

Definition at line 646 of file ice_mod.h.

647!
648!=======================================================================
649! !
650! This routine initialize all variables in the module using first !
651! touch distribution policy. In shared-memory configuration, this !
652! operation actually performs propagation of the "shared arrays" !
653! across the cluster, unless another policy is specified to !
654! override the default. !
655! !
656!=======================================================================
657!
658 USE mod_param, ONLY : bounds, domain, lbc, inlm
659 USE mod_scalars, ONLY : iwest, ieast, isouth, inorth
660!
661! Imported variable declarations.
662!
663 integer, intent(in) :: ng, tile, model
664!
665! Local variable declarations.
666!
667 integer :: i, j, nf, ns
668 integer :: Imin, Imax, Jmin, Jmax
669
670 real(r8), parameter :: IniVal = 0.0_r8
671!
672#include "tile.h"
673!
674! Set array initialization range.
675!
676#ifdef _OPENMP
677 IF (domain(ng)%Western_Edge(tile)) THEN
678 imin=bounds(ng)%LBi(tile)
679 ELSE
680 imin=istr
681 END IF
682 IF (domain(ng)%Eastern_Edge(tile)) THEN
683 imax=bounds(ng)%UBi(tile)
684 ELSE
685 imax=iend
686 END IF
687 IF (domain(ng)%Southern_Edge(tile)) THEN
688 jmin=bounds(ng)%LBj(tile)
689 ELSE
690 jmin=jstr
691 END IF
692 IF (domain(ng)%Northern_Edge(tile)) THEN
693 jmax=bounds(ng)%UBj(tile)
694 ELSE
695 jmax=jend
696 END IF
697#else
698 imin=lbi
699 imax=ubi
700 jmin=lbj
701 jmax=ubj
702#endif
703!
704!-----------------------------------------------------------------------
705! Initialize ice model state (Si) and internal field (Fi) arrays.
706!-----------------------------------------------------------------------
707!
708 DO j=jmin,jmax
709 DO i=imin,imax
710 DO ns=1,nices
711 ice(ng) % Si(i,j,1,ns) = inival
712 ice(ng) % Si(i,j,2,ns) = inival
713 END DO
714!
715 DO nf=1,nicef
716 ice(ng) % Fi(i,j,nf) = inival
717 END DO
718 END DO
719 END DO
720!
721!-----------------------------------------------------------------------
722! Initialize ice model lateral boundary arrays.
723!-----------------------------------------------------------------------
724!
725 IF ((model.eq.0).or.(model.eq.inlm)) THEN
726
727 DO i=1,nices
728 IF (domain(ng)%NorthWest_Test(tile)) THEN
729 IF (lbc(iwest,ibice(i),ng)%acquire) THEN
730 ice_lobc(i,ng) % ice_west = inival
731 ice_lobc(i,ng) % iceG_west = inival
732 END IF
733 END IF
734!
735 IF (domain(ng)%SouthEast_Test(tile)) THEN
736 IF (lbc(ieast,ibice(i),ng)%acquire) THEN
737 ice_lobc(i,ng) % ice_east = inival
738 ice_lobc(i,ng) % iceG_east = inival
739 END IF
740 END IF
741!
742 IF (domain(ng)%SouthWest_Test(tile)) THEN
743 IF (lbc(isouth,ibice(i),ng)%acquire) THEN
744 ice_lobc(i,ng) % ice_south = inival
745 ice_lobc(i,ng) % iceG_south = inival
746 END IF
747 END IF
748!
749 IF (domain(ng)%NorthEast_Test(tile)) THEN
750 IF (lbc(inorth,ibice(i),ng)%acquire) THEN
751 ice_lobc(i,ng) % ice_north = inival
752 ice_lobc(i,ng) % iceG_north = inival
753 END IF
754 END IF
755 END DO
756
757 END IF
758
759#ifdef AVERAGES
760!
761!-----------------------------------------------------------------------
762! Initialize time-averaged variables.
763!-----------------------------------------------------------------------
764!
765! Ice model internal fields.
766!
767 DO nf=1,nicef
768 IF (ifice(nf).gt.0) THEN
769 IF (licefavg(nf,ng)) THEN
770 DO j=jmin,jmax
771 DO i=imin,imax
772 ice_favg(nf,ng) % var(i,j) = inival
773 END DO
774 END DO
775 END IF
776 END IF
777 END DO
778!
779! Time-averaged state fields.
780!
781 DO ns=1,nices
782 IF (isice(ns).gt.0) THEN
783 IF (licesavg(ns,ng)) THEN
784 DO j=jmin,jmax
785 DO i=imin,imax
786 ice_savg(ns,ng) % var(i,j) = inival
787 END DO
788 END DO
789 END IF
790 END IF
791 END DO
792#endif
793!
794 RETURN
integer, parameter inlm
Definition mod_param.F:662
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329

References mod_param::bounds, mod_param::domain, ibice, ice, ice_favg, ice_lobc, ice_savg, mod_scalars::ieast, ifice, mod_param::inlm, mod_scalars::inorth, isice, mod_scalars::isouth, mod_scalars::iwest, mod_param::lbc, licefavg, licesavg, nicef, and nices.

Referenced by mod_arrays::roms_initialize_arrays().

Here is the caller graph for this function:

Variable Documentation

◆ airrho

real(r8), dimension(:), allocatable mod_ice::airrho

Definition at line 222 of file ice_mod.h.

222 real(r8), allocatable :: AirRho(:) ! air density

Referenced by allocate_ice(), and ice_thermo_mod::ice_thermo_tile().

◆ astrength

real(r8), dimension(:), allocatable mod_ice::astrength

Definition at line 235 of file ice_mod.h.

235 real(r8), allocatable :: Astrength(:)

Referenced by allocate_ice().

◆ cd_ai

real(r8), dimension(:), allocatable mod_ice::cd_ai

Definition at line 229 of file ice_mod.h.

229 real(r8), allocatable :: Cd_ai(:) ! air-ice boundary

Referenced by allocate_ice(), and bulk_flux_mod::bulk_flux_tile().

◆ cd_io

real(r8), dimension(:), allocatable mod_ice::cd_io

Definition at line 230 of file ice_mod.h.

230 real(r8), allocatable :: Cd_io(:) ! ice-ocean boundary

Referenced by allocate_ice().

◆ dtevp

integer, dimension(:), allocatable mod_ice::dtevp

Definition at line 218 of file ice_mod.h.

218 integer, allocatable :: dtevp(:) ! elastic

Referenced by allocate_ice().

◆ dtice

integer, dimension(:), allocatable mod_ice::dtice

Definition at line 217 of file ice_mod.h.

217 integer, allocatable :: dtice(:) ! viscous

Referenced by allocate_ice(), ice_advect_mod::ice_mpdata_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ ellip_sq

real(r8), dimension(:), allocatable mod_ice::ellip_sq

Definition at line 263 of file ice_mod.h.

263 real(r8), allocatable :: ellip_sq(:)

Referenced by allocate_ice().

◆ ibice

integer, dimension(nices) mod_ice::ibice

Definition at line 162 of file ice_mod.h.

162 integer :: ibICE(nIceS) ! indices to LBC switch

Referenced by allocate_ice(), deallocate_ice(), ice_advect_mod::ice_advect_tile(), ice_thermo_mod::ice_thermo_tile(), initialize_ice(), and mod_ncparam::initialize_ncparam().

◆ icaius

integer, parameter mod_ice::icaius = 1

Definition at line 171 of file ice_mod.h.

171 integer, parameter :: icAIus = 1 ! surface Air-Ice U-stress

◆ icaivs

integer, parameter mod_ice::icaivs = 2

Definition at line 172 of file ice_mod.h.

172 integer, parameter :: icAIvs = 2 ! surface Air-Ice V-stress

◆ icbvis

integer, parameter mod_ice::icbvis = 3

Definition at line 173 of file ice_mod.h.

173 integer, parameter :: icBvis = 3 ! ice bulk viscosity

◆ ice

type (t_ice), dimension(:), allocatable mod_ice::ice

◆ ice_emiss

real(r8) mod_ice::ice_emiss

Definition at line 267 of file ice_mod.h.

267 real(r8) :: ice_emiss ! ice emissivity (unitless)

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ ice_favg

type (t_ice_avg), dimension(:,:), allocatable mod_ice::ice_favg

Definition at line 319 of file ice_mod.h.

319 TYPE (T_ICE_AVG), allocatable :: ICE_FAVG(:,:) ! [nIceF,Ngrids]

Referenced by allocate_ice(), deallocate_ice(), and initialize_ice().

◆ ice_lobc

type (t_ice_lobc), dimension(:,:), allocatable mod_ice::ice_lobc

Definition at line 303 of file ice_mod.h.

303 TYPE (T_ICE_LOBC), allocatable :: ICE_LOBC(:,:) ! [nIceS,Ngrids]

Referenced by allocate_ice(), deallocate_ice(), and initialize_ice().

◆ ice_savg

type (t_ice_avg), dimension(:,:), allocatable mod_ice::ice_savg

Definition at line 320 of file ice_mod.h.

320 TYPE (T_ICE_AVG), allocatable :: ICE_SAVG(:,:) ! [nIceS,Ngrids]

Referenced by allocate_ice(), deallocate_ice(), and initialize_ice().

◆ iceobc

integer, dimension(4,nices) mod_ice::iceobc

Definition at line 163 of file ice_mod.h.

163 integer :: iceOBC(4,nIceS) ! I/O metadata indices

◆ icerho

real(r8), dimension(:), allocatable mod_ice::icerho

Definition at line 223 of file ice_mod.h.

223 real(r8), allocatable :: IceRho(:) ! ice density

Referenced by allocate_ice(), and ice_thermo_mod::ice_thermo_tile().

◆ ichsse

integer, parameter mod_ice::ichsse = 4

Definition at line 174 of file ice_mod.h.

174 integer, parameter :: icHsse = 4 ! sea surface elevation

◆ iciofv

integer, parameter mod_ice::iciofv = 5

Definition at line 175 of file ice_mod.h.

175 integer, parameter :: icIOfv = 5 ! Ice-Ocean friction velocity

◆ iciomf

integer, parameter mod_ice::iciomf = 6

Definition at line 176 of file ice_mod.h.

176 integer, parameter :: icIOmf = 6 ! Ice-Ocean mass flux

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ iciomt

integer, parameter mod_ice::iciomt = 7

Definition at line 177 of file ice_mod.h.

177 integer, parameter :: icIOmt = 7 ! Ice-Ocean momentum transfer

◆ iciovs

integer, parameter mod_ice::iciovs = 8

Definition at line 178 of file ice_mod.h.

178 integer, parameter :: icIOvs = 8 ! Ice-Ocean velocity shear

◆ icisst

integer, parameter mod_ice::icisst = 9

Definition at line 179 of file ice_mod.h.

179 integer, parameter :: icIsst = 9 ! ice/snow surface temperature

Referenced by bulk_flux_mod::bulk_flux_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ icpgrd

integer, parameter mod_ice::icpgrd = 10

Definition at line 180 of file ice_mod.h.

180 integer, parameter :: icPgrd = 10 ! gridded ice strength

◆ icpice

integer, parameter mod_ice::icpice = 11

Definition at line 181 of file ice_mod.h.

181 integer, parameter :: icPice = 11 ! ice pressure or strength

◆ icqcon

integer, parameter mod_ice::icqcon = 12

Definition at line 182 of file ice_mod.h.

182 integer, parameter :: icQcon = 12 ! ice/snow heat conductivity

Referenced by bulk_flux_mod::bulk_flux_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ icqrhs

integer, parameter mod_ice::icqrhs = 13

Definition at line 183 of file ice_mod.h.

183 integer, parameter :: icQrhs = 13 ! RHS heat flux over ice/snow

Referenced by bulk_flux_mod::bulk_flux_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ ics0mk

integer, parameter mod_ice::ics0mk = 15

Definition at line 185 of file ice_mod.h.

185 integer, parameter :: icS0mk = 15 ! molecular sublayer salinity

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ icsvis

integer, parameter mod_ice::icsvis = 14

Definition at line 184 of file ice_mod.h.

184 integer, parameter :: icSvis = 14 ! ice shear viscosity

◆ ict0mk

integer, parameter mod_ice::ict0mk = 16

Definition at line 186 of file ice_mod.h.

186 integer, parameter :: icT0mk = 16 ! molecular sublayer temperature

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ icuavg

integer, parameter mod_ice::icuavg = 17

Definition at line 187 of file ice_mod.h.

187 integer, parameter :: icUavg = 17 ! average mixed-layer U-velocity

◆ icvavg

integer, parameter mod_ice::icvavg = 18

Definition at line 188 of file ice_mod.h.

188 integer, parameter :: icVavg = 18 ! average mixed-layer V-velocity

◆ icw_ai

integer, parameter mod_ice::icw_ai = 20

Definition at line 190 of file ice_mod.h.

190 integer, parameter :: icW_ai = 20 ! melt/freeze rate at Air/Ice

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ icw_ao

integer, parameter mod_ice::icw_ao = 21

Definition at line 191 of file ice_mod.h.

191 integer, parameter :: icW_ao = 21 ! melt/freeze rate at Air/Ocean

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ icw_fr

integer, parameter mod_ice::icw_fr = 22

Definition at line 192 of file ice_mod.h.

192 integer, parameter :: icW_fr = 22 ! ice accretion rate by frazil

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ icw_io

integer, parameter mod_ice::icw_io = 23

Definition at line 193 of file ice_mod.h.

193 integer, parameter :: icW_io = 23 ! melt/freeze rate at Ice/Ocean

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ icw_ro

integer, parameter mod_ice::icw_ro = 24

Definition at line 194 of file ice_mod.h.

194 integer, parameter :: icW_ro = 24 ! melt/freeze rate runoff

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ icwdiv

integer, parameter mod_ice::icwdiv = 19

Definition at line 189 of file ice_mod.h.

189 integer, parameter :: icWdiv = 19 ! ice divergence rate

Referenced by ice_advect_mod::ice_advect_tile().

◆ idaice

integer mod_ice::idaice

Definition at line 96 of file ice_mod.h.

96 integer :: idAice ! ice concentration

◆ idaicl

integer mod_ice::idaicl

Definition at line 97 of file ice_mod.h.

97 integer :: idAiCL ! ice concentration climatology

◆ idhage

integer mod_ice::idhage

Definition at line 98 of file ice_mod.h.

98 integer :: idHage ! thickness associated with age of ice

◆ idhice

integer mod_ice::idhice

Definition at line 99 of file ice_mod.h.

99 integer :: idHice ! ice thickness

◆ idhicl

integer mod_ice::idhicl

Definition at line 100 of file ice_mod.h.

100 integer :: idHiCL ! ice thickness climatology

◆ idhmel

integer mod_ice::idhmel

Definition at line 101 of file ice_mod.h.

101 integer :: idHmel ! surface meltwater thickness

◆ idhsno

integer mod_ice::idhsno

Definition at line 102 of file ice_mod.h.

102 integer :: idHsno ! snow cover thickness

◆ idiage

integer mod_ice::idiage

Definition at line 103 of file ice_mod.h.

103 integer :: idIage ! ice age

◆ idiofv

integer mod_ice::idiofv

Definition at line 104 of file ice_mod.h.

104 integer :: idIOfv ! ice-ocean friction velocity

◆ idiomf

integer mod_ice::idiomf

Definition at line 105 of file ice_mod.h.

105 integer :: idIOmf ! ice-ocean mass flux

◆ idiomt

integer mod_ice::idiomt

Definition at line 106 of file ice_mod.h.

106 integer :: idIOmt ! ice-ocean momentum transfer coefficient

◆ idisst

integer mod_ice::idisst

Definition at line 107 of file ice_mod.h.

107 integer :: idIsst ! ice/snow surface temperature

◆ idisxx

integer mod_ice::idisxx

Definition at line 108 of file ice_mod.h.

108 integer :: idISxx ! internal ice stress xx-component

◆ idisxy

integer mod_ice::idisxy

Definition at line 109 of file ice_mod.h.

109 integer :: idISxy ! internal ice stress xy-component

◆ idisyy

integer mod_ice::idisyy

Definition at line 110 of file ice_mod.h.

110 integer :: idISyy ! internal ice stress yy-component

◆ ids0mk

integer mod_ice::ids0mk

Definition at line 111 of file ice_mod.h.

111 integer :: idS0mk ! under ice molecular sublayer salinity

◆ idt0mk

integer mod_ice::idt0mk

Definition at line 112 of file ice_mod.h.

112 integer :: idT0mk ! under ice molecular sublayer temperature

◆ idtice

integer mod_ice::idtice

Definition at line 113 of file ice_mod.h.

113 integer :: idTice ! ice interior temperature

◆ iduice

integer mod_ice::iduice

Definition at line 114 of file ice_mod.h.

114 integer :: idUice ! ice U-velocity

◆ iduicl

integer mod_ice::iduicl

Definition at line 115 of file ice_mod.h.

115 integer :: idUiCL ! ice U-velocity climatology

◆ iduier

integer mod_ice::iduier

Definition at line 116 of file ice_mod.h.

116 integer :: idUiER ! ice U-eastward at RHO-points

◆ idvice

integer mod_ice::idvice

Definition at line 117 of file ice_mod.h.

117 integer :: idVice ! ice V-velocity

◆ idvicl

integer mod_ice::idvicl

Definition at line 118 of file ice_mod.h.

118 integer :: idViCL ! ice V-velocity climatology

◆ idvinr

integer mod_ice::idvinr

Definition at line 119 of file ice_mod.h.

119 integer :: idViNR ! ice V-northward at RHO-points

◆ idw_ai

integer mod_ice::idw_ai

Definition at line 121 of file ice_mod.h.

121 integer :: idW_ai ! rate of melt/freeze at air/ice edge

◆ idw_ao

integer mod_ice::idw_ao

Definition at line 122 of file ice_mod.h.

122 integer :: idW_ao ! rate of melt/freeze at air/ocean edge

◆ idw_fr

integer mod_ice::idw_fr

Definition at line 123 of file ice_mod.h.

123 integer :: idW_fr ! rate of ice accretion by frazil growth

◆ idw_io

integer mod_ice::idw_io

Definition at line 124 of file ice_mod.h.

124 integer :: idW_io ! rate of melt/freeze at ice/ocean edge

◆ idw_ro

integer mod_ice::idw_ro

Definition at line 125 of file ice_mod.h.

125 integer :: idW_ro ! rate of melt/freeze runoff into ocean

◆ idwdiv

integer mod_ice::idwdiv

Definition at line 120 of file ice_mod.h.

120 integer :: idWdiv ! rate of ice divergence

◆ ievp

integer, dimension(:), allocatable mod_ice::ievp

Definition at line 212 of file ice_mod.h.

212 integer, allocatable :: iEVP(:) ! timestep counter

Referenced by allocate_ice().

◆ ifice

integer, dimension(nicef) mod_ice::ifice

Definition at line 169 of file ice_mod.h.

169 integer :: iFice(nIceF) ! internal fields I/O indices

Referenced by allocate_ice(), and initialize_ice().

◆ isaice

integer, parameter mod_ice::isaice = 1

Definition at line 137 of file ice_mod.h.

137 integer, parameter :: isAice = 1 ! ice concentration

Referenced by bulk_flux_mod::bulk_flux_tile(), ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ isenth

integer, parameter mod_ice::isenth = 12

Definition at line 148 of file ice_mod.h.

148 integer, parameter :: isEnth = 12 ! ice/brine enthalpy

Referenced by ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ ishage

integer, parameter mod_ice::ishage = 13

Definition at line 149 of file ice_mod.h.

149 integer, parameter :: isHage = 13 ! thickness linked with ice age

Referenced by ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ ishice

integer, parameter mod_ice::ishice = 2

Definition at line 138 of file ice_mod.h.

138 integer, parameter :: isHice = 2 ! ice thickness

Referenced by bulk_flux_mod::bulk_flux_tile(), ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ ishmel

integer, parameter mod_ice::ishmel = 3

Definition at line 139 of file ice_mod.h.

139 integer, parameter :: isHmel = 3 ! meltwater thickness on ice

Referenced by bulk_flux_mod::bulk_flux_tile(), ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ ishsno

integer, parameter mod_ice::ishsno = 4

Definition at line 140 of file ice_mod.h.

140 integer, parameter :: isHsno = 4 ! snow thickness

Referenced by bulk_flux_mod::bulk_flux_tile(), ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ isiage

integer, parameter mod_ice::isiage = 5

Definition at line 141 of file ice_mod.h.

141 integer, parameter :: isIage = 5 ! ice age

Referenced by ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ isice

integer, dimension(nices) mod_ice::isice

◆ isilog

integer, parameter mod_ice::isilog = 19

Definition at line 156 of file ice_mod.h.

156 integer, parameter :: isIlog = 19 ! ice biology log

◆ isinh4

integer, parameter mod_ice::isinh4 = 18

Definition at line 155 of file ice_mod.h.

155 integer, parameter :: isINH4 = 18 ! ice biology, ammonia

Referenced by ice_advect_mod::ice_advect_tile().

◆ isino3

integer, parameter mod_ice::isino3 = 17

Definition at line 154 of file ice_mod.h.

154 integer, parameter :: isINO3 = 17 ! ice biology, nitrate

Referenced by ice_advect_mod::ice_advect_tile().

◆ isiphy

integer, parameter mod_ice::isiphy = 16

Definition at line 153 of file ice_mod.h.

153 integer, parameter :: isIphy = 16 ! ice biology, algae

Referenced by ice_advect_mod::ice_advect_tile().

◆ isisxx

integer, parameter mod_ice::isisxx = 6

Definition at line 142 of file ice_mod.h.

142 integer, parameter :: isISxx = 6 ! internal ice xx-stress

◆ isisxy

integer, parameter mod_ice::isisxy = 7

Definition at line 143 of file ice_mod.h.

143 integer, parameter :: isISxy = 7 ! internal ice xy-stress

◆ isisyy

integer, parameter mod_ice::isisyy = 8

Definition at line 144 of file ice_mod.h.

144 integer, parameter :: isISyy = 8 ! internal ice yy-stress

◆ istice

integer, parameter mod_ice::istice = 9

Definition at line 145 of file ice_mod.h.

145 integer, parameter :: isTice = 9 ! ice interior temperature

Referenced by ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ isuevp

integer, parameter mod_ice::isuevp = 14

Definition at line 150 of file ice_mod.h.

150 integer, parameter :: isUevp = 14 ! EVP ice U-velocity

◆ isuice

integer, parameter mod_ice::isuice = 10

Definition at line 146 of file ice_mod.h.

146 integer, parameter :: isUice = 10 ! ice U-velocity

Referenced by get_state_mod::get_state_nf90(), get_state_mod::get_state_pio(), ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ isvevp

integer, parameter mod_ice::isvevp = 15

Definition at line 151 of file ice_mod.h.

151 integer, parameter :: isVevp = 15 ! EVP ice V-velocity

◆ isvice

integer, parameter mod_ice::isvice = 11

◆ licefavg

logical, dimension(:,:), allocatable mod_ice::licefavg

Definition at line 205 of file ice_mod.h.

205 logical, allocatable :: LiceFavg(:,:) ! internal variables

Referenced by allocate_ice(), deallocate_ice(), and initialize_ice().

◆ licesavg

logical, dimension(:,:), allocatable mod_ice::licesavg

Definition at line 206 of file ice_mod.h.

206 logical, allocatable :: LiceSavg(:,:) ! state variables

Referenced by allocate_ice(), deallocate_ice(), and initialize_ice().

◆ max_ai

real(r8), dimension(:), allocatable mod_ice::max_ai

Definition at line 242 of file ice_mod.h.

242 real(r8), allocatable :: max_ai(:) ! maximun limiter

Referenced by allocate_ice(), and ice_thermo_mod::ice_thermo_tile().

◆ max_hmelt

real(r8), dimension(:), allocatable mod_ice::max_hmelt

Definition at line 250 of file ice_mod.h.

250 real(r8), allocatable :: max_hmelt(:)

Referenced by allocate_ice(), and ice_thermo_mod::ice_thermo_tile().

◆ min_ai

real(r8), dimension(:), allocatable mod_ice::min_ai

Definition at line 241 of file ice_mod.h.

241 real(r8), allocatable :: min_ai(:) ! minimum limiter

Referenced by allocate_ice(), bulk_flux_mod::bulk_flux_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ min_hi

real(r8), dimension(:), allocatable mod_ice::min_hi

Definition at line 246 of file ice_mod.h.

246 real(r8), allocatable :: min_hi(:) ! minimum limiter

Referenced by allocate_ice(), ice_advect_mod::ice_advect_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ nevp

integer, dimension(:), allocatable mod_ice::nevp

Definition at line 213 of file ice_mod.h.

213 integer, allocatable :: nEVP(:) ! number of timesteps

Referenced by allocate_ice().

◆ nicef

integer, parameter mod_ice::nicef = 24

Definition at line 167 of file ice_mod.h.

167 integer, parameter :: nIceF = 24 ! number of ice field variables

Referenced by allocate_ice(), deallocate_ice(), and initialize_ice().

◆ nices

integer parameter mod_ice::nices = 19

Definition at line 130 of file ice_mod.h.

130 integer, parameter :: nIceS = 19 ! number of ice state variables

Referenced by allocate_ice(), checkvars_mod::checkvars::checkvars_nf90(), deallocate_ice(), get_state_mod::get_state_nf90(), get_state_mod::get_state_pio(), initialize_ice(), and mod_ncparam::initialize_ncparam().

◆ pstar

real(r8), dimension(:), allocatable mod_ice::pstar

Definition at line 236 of file ice_mod.h.

236 real(r8), allocatable :: Pstar(:)

Referenced by allocate_ice().

◆ snowdryrho

real(r8), dimension(:), allocatable mod_ice::snowdryrho

Definition at line 224 of file ice_mod.h.

224 real(r8), allocatable :: SnowDryRho(:) ! dry snow density

Referenced by allocate_ice(), bulk_flux_mod::bulk_flux_tile(), and ice_thermo_mod::ice_thermo_tile().

◆ snowwetrho

real(r8), dimension(:), allocatable mod_ice::snowwetrho

Definition at line 225 of file ice_mod.h.

225 real(r8), allocatable :: SnowWetRho(:) ! wet snow density

Referenced by allocate_ice(), and ice_thermo_mod::ice_thermo_tile().

◆ spec_heat_air

real(r8) mod_ice::spec_heat_air

Definition at line 268 of file ice_mod.h.

268 real(r8) :: spec_heat_air ! specific heat of air (J/kg/K)

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ stressang

real(r8), dimension(:), allocatable mod_ice::stressang

Definition at line 259 of file ice_mod.h.

259 real(r8), allocatable :: stressAng(:)

Referenced by allocate_ice().

◆ sublimation

real(r8) mod_ice::sublimation

Definition at line 270 of file ice_mod.h.

270 real(r8) :: sublimation ! latent heat of sublimation (J/kg)

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ trans_coeff

real(r8) mod_ice::trans_coeff

Definition at line 269 of file ice_mod.h.

269 real(r8) :: trans_coeff ! heat transfer coefficient (unitless)

Referenced by ice_thermo_mod::ice_thermo_tile().

◆ zetamax

real(r8), dimension(:), allocatable mod_ice::zetamax

Definition at line 255 of file ice_mod.h.

255 real(r8), allocatable :: zetaMax(:) ! maximum limiter

Referenced by allocate_ice().

◆ zetamin

real(r8), dimension(:), allocatable mod_ice::zetamin

Definition at line 254 of file ice_mod.h.

254 real(r8), allocatable :: zetaMin(:) ! minimum limiter

Referenced by allocate_ice().