ROMS
Loading...
Searching...
No Matches
ice_mod.h
Go to the documentation of this file.
1 MODULE mod_ice
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group Paul Budgell !
6! Licensed under a MIT/X style license Katherine Hedstrom !
7! See License_ROMS.md Scott M. Durski !
8!=======================================================================
9! !
10! The sea-ice model was originally written and tested by Paul Budgell !
11! (2005). It has been updated and improved by Kate Hedstrom and Scott !
12! Durski. !
13! !
14! It uses the elastic-viscous-plastic (EVP) rheology of Hunke and !
15! Dukowicz (1997) and Hunke (2001). It has a simple one-layer ice and !
16! snow thermodynamics with a molecular sublayer under the ice (Mellor !
17! and Kantha, 1989). !
18! !
19! Variables: !
20! !
21! The ice model prognostic 'state' variable are declares in compact !
22! form in the derived-type structure as: !
23! !
24! ICE(ng) % state(:,:,TimeLevel,StateIndex) !
25! !
26! where the 'StateIndex' are as follows: !
27! !
28! isAice ice concentration (grid cell fraction) !
29! isEnth enthalpy of the ice/brine system, ice heat content !
30! isHage thickness associated with age of ice (m) !
31! isHice average ice thickness (m), ice mass divided by area !
32! isHmel surface meltwater thickness on ice (m) !
33! isHsno average thickness of snow coverage (m), mass snow !
34! isIage age of ice(s) !
35! isISxx internal ice stress, xx-component (N/m) !
36! isISxy internal ice stress, xy-component (N/m) !
37! isISyy internal ice stress, yy-component (N/m) !
38! isTice ice interior temperature (Celsius) !
39! isUice ice U-velocity (m/s) !
40! isUevp elastic-viscous-plastic ice U-velocity (m/s) !
41! isVice ice V-velocity (m/s) !
42! isVevp elastic-viscous-plastic ice V-velocity (m/s) !
43#if defined ICE_BIO
44! isIphy ice biology, algae !
45! isINO3 ice biology, nitrate !
46! isINH4 ice biology, ammonia !
47! isIlog ice biology log !
48#endif
49! !
50! The ice model internal 'field' arrays are declares in compact form !
51! in the derived-type structure as: !
52! !
53! ICE(ng) % field(:,:,FieldIndex) !
54! !
55! where the 'FieldIndex' are as follows: !
56! !
57! icAIus surface Air-Ice U-stress (N/m2) !
58! icAIvs surface Air-Ice V-stress (N/m2) !
59! icBvis ice bulk viscosity !
60! icHsse sea surface elevation (m) !
61! icPgrd gridded ice strength parameter (unitless) !
62! icPice ice pressure or strength (N/m2) !
63! icQcon gradient heat conductivity over ice/snow (W/m2/K) !
64! icQrhs RHS surface net heat flux over ice/snow (W/m2) !
65! icSvis ice shear viscosity !
66! icS0mk salinity of molecular sublayer under ice (unitless) !
67! icT0mk temperature of molecular sublayer under ice (Celsius) !
68! icIsst temperature at the snow/atmosphere interface (Celsius) !
69! icUavg vertically averaged mixed-layer U-velocity (m/s) !
70! icVavg vertically averaged mixed-layer V-velocity (m/s) !
71! icWdiv rate of ice divergence (m3/s) !
72! icW_ai rate of melt/freeze at Air/Ice edge (m3/s) !
73! icW_ao rate of melt/freeze at Air/Ocean edge (m3/s) !
74! icW_fr rate of ice accretion due to frazil growth (m3/s) !
75! icW_io rate of melt/freeze at Ice/Ocean edge (m3/s) !
76! icW_ro rate of melt/freeze runoff into ocean (m3/s) !
77! icIOmf Ice-Ocean mass flux (m/s) !
78! icIOvs Ice-Ocean velocity shear magnitude (m/s) !
79! icIOfv Ice-Ocean friction velocity (m/s) !
80! icIOmt Ice-Ocean momentum transfer coefficient (m/s) !
81! !
82!=======================================================================
83!
84 USE mod_kinds
85!
86 implicit none
87!
88 PUBLIC :: allocate_ice
89 PUBLIC :: deallocate_ice
90 PUBLIC :: initialize_ice
91!
92!-----------------------------------------------------------------------
93! Ice model identification indices.
94!-----------------------------------------------------------------------
95!
96 integer :: idaice ! ice concentration
97 integer :: idaicl ! ice concentration climatology
98 integer :: idhage ! thickness associated with age of ice
99 integer :: idhice ! ice thickness
100 integer :: idhicl ! ice thickness climatology
101 integer :: idhmel ! surface meltwater thickness
102 integer :: idhsno ! snow cover thickness
103 integer :: idiage ! ice age
104 integer :: idiofv ! ice-ocean friction velocity
105 integer :: idiomf ! ice-ocean mass flux
106 integer :: idiomt ! ice-ocean momentum transfer coefficient
107 integer :: idisst ! ice/snow surface temperature
108 integer :: idisxx ! internal ice stress xx-component
109 integer :: idisxy ! internal ice stress xy-component
110 integer :: idisyy ! internal ice stress yy-component
111 integer :: ids0mk ! under ice molecular sublayer salinity
112 integer :: idt0mk ! under ice molecular sublayer temperature
113 integer :: idtice ! ice interior temperature
114 integer :: iduice ! ice U-velocity
115 integer :: iduicl ! ice U-velocity climatology
116 integer :: iduier ! ice U-eastward at RHO-points
117 integer :: idvice ! ice V-velocity
118 integer :: idvicl ! ice V-velocity climatology
119 integer :: idvinr ! ice V-northward at RHO-points
120 integer :: idwdiv ! rate of ice divergence
121 integer :: idw_ai ! rate of melt/freeze at air/ice edge
122 integer :: idw_ao ! rate of melt/freeze at air/ocean edge
123 integer :: idw_fr ! rate of ice accretion by frazil growth
124 integer :: idw_io ! rate of melt/freeze at ice/ocean edge
125 integer :: idw_ro ! rate of melt/freeze runoff into ocean
126!
127! Ice model state prognostic variables indices.
128!
129#if defined ICE_BIO
130 integer, parameter :: nices = 19 ! number of ice state variables
131#else
132 integer, parameter :: nices = 15 ! number of ice state variables
133#endif
134!
135 integer :: isice(nices) ! state I/O indices
136!
137 integer, parameter :: isaice = 1 ! ice concentration
138 integer, parameter :: ishice = 2 ! ice thickness
139 integer, parameter :: ishmel = 3 ! meltwater thickness on ice
140 integer, parameter :: ishsno = 4 ! snow thickness
141 integer, parameter :: isiage = 5 ! ice age
142 integer, parameter :: isisxx = 6 ! internal ice xx-stress
143 integer, parameter :: isisxy = 7 ! internal ice xy-stress
144 integer, parameter :: isisyy = 8 ! internal ice yy-stress
145 integer, parameter :: istice = 9 ! ice interior temperature
146 integer, parameter :: isuice = 10 ! ice U-velocity
147 integer, parameter :: isvice = 11 ! ice V-velocity
148 integer, parameter :: isenth = 12 ! ice/brine enthalpy
149 integer, parameter :: ishage = 13 ! thickness linked with ice age
150 integer, parameter :: isuevp = 14 ! EVP ice U-velocity
151 integer, parameter :: isvevp = 15 ! EVP ice V-velocity
152#if defined ICE_BIO
153 integer, parameter :: isiphy = 16 ! ice biology, algae
154 integer, parameter :: isino3 = 17 ! ice biology, nitrate
155 integer, parameter :: isinh4 = 18 ! ice biology, ammonia
156 integer, parameter :: isilog = 19 ! ice biology log
157#endif
158
159!
160! Ice model state lateral boundary conditions indices.
161!
162 integer :: ibice(nices) ! indices to LBC switch
163 integer :: iceobc(4,nices) ! I/O metadata indices
164!
165! Ice model internal variables indices.
166!
167 integer, parameter :: nicef = 24 ! number of ice field variables
168!
169 integer :: ifice(nicef) ! internal fields I/O indices
170!
171 integer, parameter :: icaius = 1 ! surface Air-Ice U-stress
172 integer, parameter :: icaivs = 2 ! surface Air-Ice V-stress
173 integer, parameter :: icbvis = 3 ! ice bulk viscosity
174 integer, parameter :: ichsse = 4 ! sea surface elevation
175 integer, parameter :: iciofv = 5 ! Ice-Ocean friction velocity
176 integer, parameter :: iciomf = 6 ! Ice-Ocean mass flux
177 integer, parameter :: iciomt = 7 ! Ice-Ocean momentum transfer
178 integer, parameter :: iciovs = 8 ! Ice-Ocean velocity shear
179 integer, parameter :: icisst = 9 ! ice/snow surface temperature
180 integer, parameter :: icpgrd = 10 ! gridded ice strength
181 integer, parameter :: icpice = 11 ! ice pressure or strength
182 integer, parameter :: icqcon = 12 ! ice/snow heat conductivity
183 integer, parameter :: icqrhs = 13 ! RHS heat flux over ice/snow
184 integer, parameter :: icsvis = 14 ! ice shear viscosity
185 integer, parameter :: ics0mk = 15 ! molecular sublayer salinity
186 integer, parameter :: ict0mk = 16 ! molecular sublayer temperature
187 integer, parameter :: icuavg = 17 ! average mixed-layer U-velocity
188 integer, parameter :: icvavg = 18 ! average mixed-layer V-velocity
189 integer, parameter :: icwdiv = 19 ! ice divergence rate
190 integer, parameter :: icw_ai = 20 ! melt/freeze rate at Air/Ice
191 integer, parameter :: icw_ao = 21 ! melt/freeze rate at Air/Ocean
192 integer, parameter :: icw_fr = 22 ! ice accretion rate by frazil
193 integer, parameter :: icw_io = 23 ! melt/freeze rate at Ice/Ocean
194 integer, parameter :: icw_ro = 24 ! melt/freeze rate runoff
195!
196!-----------------------------------------------------------------------
197! Ice model parameters.
198!-----------------------------------------------------------------------
199
200#ifdef AVERAGES
201!
202! Switches to process time-averaged ice model state and internal
203! variables.
204!
205 logical, allocatable :: licefavg(:,:) ! internal variables
206 logical, allocatable :: licesavg(:,:) ! state variables
207#endif
208!
209! Counter and number of Elastic-Viscous-Plastic (EVP) rheology
210! equations timesteps to resolve elastic dynamics.
211!
212 integer, allocatable :: ievp(:) ! timestep counter
213 integer, allocatable :: nevp(:) ! number of timesteps
214!
215! Ice equations time step (s).
216!
217 integer, allocatable :: dtice(:) ! viscous
218 integer, allocatable :: dtevp(:) ! elastic
219!
220! Density parameters (kg/m3).
221!
222 real(r8), allocatable :: airrho(:) ! air density
223 real(r8), allocatable :: icerho(:) ! ice density
224 real(r8), allocatable :: snowdryrho(:) ! dry snow density
225 real(r8), allocatable :: snowwetrho(:) ! wet snow density
226!
227! Boundary layer bulk drag coefficients (nondimensional).
228!
229 real(r8), allocatable :: cd_ai(:) ! air-ice boundary
230 real(r8), allocatable :: cd_io(:) ! ice-ocean boundary
231!
232! Ice strength exponential weighting coefficient on concentration
233! (nondimensional).
234!
235 real(r8), allocatable :: astrength(:)
236 real(r8), allocatable :: pstar(:)
237!
238! Minimum and maximum average fractional ice concetration coverage
239! (nomdimensional) limiters.
240!
241 real(r8), allocatable :: min_ai(:) ! minimum limiter
242 real(r8), allocatable :: max_ai(:) ! maximun limiter
243!
244! Minimum average ice thickness (m) limiter.
245!
246 real(r8), allocatable :: min_hi(:) ! minimum limiter
247!
248! Maximum surface melt water thickness (m) limiter.
249!
250 real(r8), allocatable :: max_hmelt(:)
251!
252! Minimum and maximum ice shear strength (N/m2) limiters.
253!
254 real(r8), allocatable :: zetamin(:) ! minimum limiter
255 real(r8), allocatable :: zetamax(:) ! maximum limiter
256!
257! Turning angle for ice-water drag (radians).
258!
259 real(r8), allocatable :: stressang(:)
260!
261! Squared ellipticity/eccentricity of the yield curve (nondimensional).
262!
263 real(r8), allocatable :: ellip_sq(:)
264!
265! Model formulation constants.
266!
267 real(r8) :: ice_emiss ! ice emissivity (unitless)
268 real(r8) :: spec_heat_air ! specific heat of air (J/kg/K)
269 real(r8) :: trans_coeff ! heat transfer coefficient (unitless)
270 real(r8) :: sublimation ! latent heat of sublimation (J/kg)
271!
272!-----------------------------------------------------------------------
273! Define derived-type structure ice model state and internal arrays.
274!-----------------------------------------------------------------------
275!
276 TYPE t_ice
277
278 real(r8), pointer :: fi(:,:,:) ! [i,j,1:nIceF]
279 real(r8), pointer :: si(:,:,:,:) ! [i,j,1:2,1:nIceS]
280
281 END TYPE t_ice
282!
283 TYPE (t_ice), allocatable :: ice(:) ! [Ngrids]
284!
285!-----------------------------------------------------------------------
286! Define derived-type structure ice model lateral boundary variables.
287!-----------------------------------------------------------------------
288!
290
291 real(r8), pointer :: ice_west (:)
292 real(r8), pointer :: ice_east (:)
293 real(r8), pointer :: ice_south(:)
294 real(r8), pointer :: ice_north(:)
295!
296 real(r8), pointer :: iceg_west (:,:)
297 real(r8), pointer :: iceg_east (:,:)
298 real(r8), pointer :: iceg_south(:,:)
299 real(r8), pointer :: iceg_north(:,:)
300
301 END TYPE t_ice_lobc
302!
303 TYPE (t_ice_lobc), allocatable :: ice_lobc(:,:) ! [nIceS,Ngrids]
304
305#ifdef AVERAGES
306!
307!-----------------------------------------------------------------------
308! Define derived-type structure ice model state and internal arrays
309! time-averaged variables. Notice that only the requested arrays are
310! allocated and processed.
311!-----------------------------------------------------------------------
312!
314
315 real(r8), pointer :: var(:,:) ! [i,j]
316
317 END TYPE t_ice_avg
318!
319 TYPE (t_ice_avg), allocatable :: ice_favg(:,:) ! [nIceF,Ngrids]
320 TYPE (t_ice_avg), allocatable :: ice_savg(:,:) ! [nIceS,Ngrids]
321#endif
322!
323 CONTAINS
324!
325 SUBROUTINE allocate_ice (ng, LBi, UBi, LBj, UBj, ice_kernel)
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
530 END SUBROUTINE allocate_ice
531!
532 SUBROUTINE deallocate_ice (ng)
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
644 END SUBROUTINE deallocate_ice
645!
646 SUBROUTINE initialize_ice (ng, tile, model)
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
795 END SUBROUTINE initialize_ice
796!
797 END MODULE mod_ice
integer idhicl
Definition ice_mod.h:100
integer, parameter isvice
Definition ice_mod.h:147
integer, dimension(:), allocatable dtice
Definition ice_mod.h:217
real(r8), dimension(:), allocatable pstar
Definition ice_mod.h:236
real(r8), dimension(:), allocatable min_hi
Definition ice_mod.h:246
integer idisxy
Definition ice_mod.h:109
integer, dimension(:), allocatable ievp
Definition ice_mod.h:212
integer, parameter ishage
Definition ice_mod.h:149
integer idvinr
Definition ice_mod.h:119
integer, parameter isino3
Definition ice_mod.h:154
integer idhmel
Definition ice_mod.h:101
integer, parameter iciomf
Definition ice_mod.h:176
integer, parameter icpice
Definition ice_mod.h:181
integer, parameter icw_ro
Definition ice_mod.h:194
integer, parameter icwdiv
Definition ice_mod.h:189
integer, parameter icqcon
Definition ice_mod.h:182
integer idw_ai
Definition ice_mod.h:121
integer, dimension(:), allocatable dtevp
Definition ice_mod.h:218
real(r8), dimension(:), allocatable ellip_sq
Definition ice_mod.h:263
integer idvicl
Definition ice_mod.h:118
integer, dimension(nicef) ifice
Definition ice_mod.h:169
integer, parameter ics0mk
Definition ice_mod.h:185
integer idw_ro
Definition ice_mod.h:125
integer, dimension(4, nices) iceobc
Definition ice_mod.h:163
integer idiomf
Definition ice_mod.h:105
integer idvice
Definition ice_mod.h:117
integer, parameter isenth
Definition ice_mod.h:148
integer, parameter ichsse
Definition ice_mod.h:174
integer, parameter icuavg
Definition ice_mod.h:187
real(r8), dimension(:), allocatable airrho
Definition ice_mod.h:222
integer, parameter icsvis
Definition ice_mod.h:184
integer, parameter icbvis
Definition ice_mod.h:173
integer, parameter icw_ao
Definition ice_mod.h:191
integer, parameter isisxx
Definition ice_mod.h:142
integer ids0mk
Definition ice_mod.h:111
integer idiofv
Definition ice_mod.h:104
real(r8), dimension(:), allocatable icerho
Definition ice_mod.h:223
integer idwdiv
Definition ice_mod.h:120
real(r8) ice_emiss
Definition ice_mod.h:267
real(r8), dimension(:), allocatable zetamax
Definition ice_mod.h:255
integer idt0mk
Definition ice_mod.h:112
type(t_ice_lobc), dimension(:,:), allocatable ice_lobc
Definition ice_mod.h:303
real(r8), dimension(:), allocatable zetamin
Definition ice_mod.h:254
integer, parameter ishsno
Definition ice_mod.h:140
real(r8), dimension(:), allocatable snowwetrho
Definition ice_mod.h:225
type(t_ice), dimension(:), allocatable ice
Definition ice_mod.h:283
integer, dimension(:), allocatable nevp
Definition ice_mod.h:213
real(r8), dimension(:), allocatable snowdryrho
Definition ice_mod.h:224
integer, parameter icvavg
Definition ice_mod.h:188
integer idw_ao
Definition ice_mod.h:122
integer, parameter istice
Definition ice_mod.h:145
integer, parameter iciomt
Definition ice_mod.h:177
subroutine, public deallocate_ice(ng)
Definition ice_mod.h:533
integer, parameter isiphy
Definition ice_mod.h:153
integer, parameter nicef
Definition ice_mod.h:167
real(r8), dimension(:), allocatable max_hmelt
Definition ice_mod.h:250
integer idw_fr
Definition ice_mod.h:123
integer idisst
Definition ice_mod.h:107
integer, dimension(nices) ibice
Definition ice_mod.h:162
integer idiage
Definition ice_mod.h:103
integer iduicl
Definition ice_mod.h:115
integer idiomt
Definition ice_mod.h:106
integer, parameter isisxy
Definition ice_mod.h:143
real(r8), dimension(:), allocatable min_ai
Definition ice_mod.h:241
integer, parameter iciovs
Definition ice_mod.h:178
integer idhice
Definition ice_mod.h:99
integer idisxx
Definition ice_mod.h:108
real(r8) spec_heat_air
Definition ice_mod.h:268
integer, parameter isuevp
Definition ice_mod.h:150
integer, parameter isilog
Definition ice_mod.h:156
integer, parameter ishmel
Definition ice_mod.h:139
real(r8) trans_coeff
Definition ice_mod.h:269
integer, parameter icpgrd
Definition ice_mod.h:180
integer, parameter iciofv
Definition ice_mod.h:175
integer, dimension(nices) isice
Definition ice_mod.h:135
integer idhage
Definition ice_mod.h:98
integer, parameter icaius
Definition ice_mod.h:171
integer, parameter icw_fr
Definition ice_mod.h:192
integer idisyy
Definition ice_mod.h:110
integer, parameter icqrhs
Definition ice_mod.h:183
integer idhsno
Definition ice_mod.h:102
integer, parameter isiage
Definition ice_mod.h:141
integer, parameter isaice
Definition ice_mod.h:137
integer, parameter isvevp
Definition ice_mod.h:151
integer, parameter ict0mk
Definition ice_mod.h:186
real(r8), dimension(:), allocatable max_ai
Definition ice_mod.h:242
integer idaice
Definition ice_mod.h:96
subroutine, public initialize_ice(ng, tile, model)
Definition ice_mod.h:647
real(r8) sublimation
Definition ice_mod.h:270
integer, parameter isinh4
Definition ice_mod.h:155
integer idw_io
Definition ice_mod.h:124
integer, parameter isuice
Definition ice_mod.h:146
type(t_ice_avg), dimension(:,:), allocatable ice_savg
Definition ice_mod.h:320
real(r8), dimension(:), allocatable cd_io
Definition ice_mod.h:230
real(r8), dimension(:), allocatable astrength
Definition ice_mod.h:235
integer idaicl
Definition ice_mod.h:97
subroutine, public allocate_ice(ng, lbi, ubi, lbj, ubj, ice_kernel)
Definition ice_mod.h:326
integer iduice
Definition ice_mod.h:114
type(t_ice_avg), dimension(:,:), allocatable ice_favg
Definition ice_mod.h:319
integer, parameter ishice
Definition ice_mod.h:138
integer, parameter nices
Definition ice_mod.h:130
integer, parameter icw_io
Definition ice_mod.h:193
logical, dimension(:,:), allocatable licesavg
Definition ice_mod.h:206
integer, parameter icisst
Definition ice_mod.h:179
integer, parameter icaivs
Definition ice_mod.h:172
integer, parameter isisyy
Definition ice_mod.h:144
logical, dimension(:,:), allocatable licefavg
Definition ice_mod.h:205
real(r8), dimension(:), allocatable stressang
Definition ice_mod.h:259
integer, parameter icw_ai
Definition ice_mod.h:190
integer idtice
Definition ice_mod.h:113
integer iduier
Definition ice_mod.h:116
real(r8), dimension(:), allocatable cd_ai
Definition ice_mod.h:229
integer, parameter r8
Definition mod_kinds.F:28
integer, parameter inlm
Definition mod_param.F:662
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer ngrids
Definition mod_param.F:113
integer, parameter iwest
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth