New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
limistate.F90 in branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/LIM_SRC_3/limistate.F90 @ 13354

Last change on this file since 13354 was 7993, checked in by frrh, 7 years ago

Merge in missing revisions 6428:2477 inclusive and 6482 from nemo_v3_6_STABLE
branch. In ptic, this includes the fix for restartability of runoff fields in coupled
models. Evolution of coupled models will therefor be affected.

These changes donot affect evolution of the current stand-alone NEMO-CICE GO6
standard configuration.

Work and testing documented in Met Office GMED ticket 320.

File size: 23.3 KB
RevLine 
[825]1MODULE limistate
2   !!======================================================================
3   !!                     ***  MODULE  limistate  ***
4   !!              Initialisation of diagnostics ice variables
5   !!======================================================================
[2528]6   !! History :  2.0  ! 2004-01 (C. Ethe, G. Madec)  Original code
[4161]7   !!            4.0  ! 2011-02 (G. Madec) dynamical allocation
[4990]8   !!             -   ! 2014    (C. Rousset) add N/S initializations
[2528]9   !!----------------------------------------------------------------------
[825]10#if defined key_lim3
11   !!----------------------------------------------------------------------
[834]12   !!   'key_lim3' :                                    LIM3 sea-ice model
[825]13   !!----------------------------------------------------------------------
14   !!   lim_istate      :  Initialisation of diagnostics ice variables
15   !!   lim_istate_init :  initialization of ice state and namelist read
16   !!----------------------------------------------------------------------
[2528]17   USE phycst           ! physical constant
18   USE oce              ! dynamics and tracers variables
19   USE dom_oce          ! ocean domain
20   USE sbc_oce          ! Surface boundary condition: ocean fields
[4161]21   USE sbc_ice          ! Surface boundary condition: ice fields
[2528]22   USE eosbn2           ! equation of state
23   USE ice              ! sea-ice variables
[4161]24   USE par_oce          ! ocean parameters
[2528]25   USE dom_ice          ! sea-ice domain
[7993]26   USE limvar           ! lim_var_salprof
[2528]27   USE in_out_manager   ! I/O manager
[2715]28   USE lib_mpp          ! MPP library
[4161]29   USE lib_fortran      ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
[3294]30   USE wrk_nemo         ! work arrays
[825]31
32   IMPLICIT NONE
33   PRIVATE
34
[2528]35   PUBLIC   lim_istate      ! routine called by lim_init.F90
[825]36
[4166]37   !                          !!** init namelist (namiceini) **
[5123]38   REAL(wp) ::   rn_thres_sst   ! threshold water temperature for initial sea ice
39   REAL(wp) ::   rn_hts_ini_n   ! initial snow thickness in the north
40   REAL(wp) ::   rn_hts_ini_s   ! initial snow thickness in the south
41   REAL(wp) ::   rn_hti_ini_n   ! initial ice thickness in the north
42   REAL(wp) ::   rn_hti_ini_s   ! initial ice thickness in the south
43   REAL(wp) ::   rn_ati_ini_n   ! initial leads area in the north
44   REAL(wp) ::   rn_ati_ini_s   ! initial leads area in the south
45   REAL(wp) ::   rn_smi_ini_n   ! initial salinity
46   REAL(wp) ::   rn_smi_ini_s   ! initial salinity
47   REAL(wp) ::   rn_tmi_ini_n   ! initial temperature
48   REAL(wp) ::   rn_tmi_ini_s   ! initial temperature
[4166]49
[5123]50   LOGICAL  ::  ln_iceini    ! initialization or not
[825]51   !!----------------------------------------------------------------------
[4161]52   !!   LIM 3.0,  UCL-LOCEAN-IPSL (2008)
[1156]53   !! $Id$
[4161]54   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
[825]55   !!----------------------------------------------------------------------
56CONTAINS
57
58   SUBROUTINE lim_istate
59      !!-------------------------------------------------------------------
60      !!                    ***  ROUTINE lim_istate  ***
61      !!
62      !! ** Purpose :   defined the sea-ice initial state
63      !!
[4161]64      !! ** Method  :   
65      !!                This routine will put some ice where ocean
66      !!                is at the freezing point, then fill in ice
67      !!                state variables using prescribed initial
68      !!                values in the namelist           
69      !!
70      !! ** Steps   :   
71      !!                1) Read namelist
72      !!                2) Basal temperature; ice and hemisphere masks
73      !!                3) Fill in the ice thickness distribution using gaussian
74      !!                4) Fill in space-dependent arrays for state variables
75      !!                5) Diagnostic arrays
76      !!                6) Lateral boundary conditions
77      !!
[4333]78      !! ** Notes   : o_i, t_su, t_s, t_i, s_i must be filled everywhere, even
[4990]79      !!              where there is no ice (clem: I do not know why, is it mandatory?)
[4333]80      !!
[4161]81      !! History :
82      !!   2.0  !  01-04  (C. Ethe, G. Madec)  Original code
83      !!   3.0  !  2007   (M. Vancoppenolle)   Rewrite for ice cats
84      !!   4.0  !  09-11  (M. Vancoppenolle)   Enhanced version for ice cats
85      !!--------------------------------------------------------------------
86
87      !! * Local variables
88      INTEGER    :: ji, jj, jk, jl             ! dummy loop indices
[5123]89      REAL(wp)   :: ztmelts, zdh
[4161]90      INTEGER    :: i_hemis, i_fill, jl0 
91      REAL(wp)   :: ztest_1, ztest_2, ztest_3, ztest_4, ztests, zsigma, zarg, zA, zV, zA_cons, zV_cons, zconv
[4688]92      REAL(wp), POINTER, DIMENSION(:)     :: zht_i_ini, zat_i_ini, zvt_i_ini, zht_s_ini, zsm_i_ini, ztm_i_ini
93      REAL(wp), POINTER, DIMENSION(:,:)   :: zh_i_ini, za_i_ini, zv_i_ini
94      REAL(wp), POINTER, DIMENSION(:,:)   :: zswitch    ! ice indicator
[4161]95      INTEGER,  POINTER, DIMENSION(:,:)   :: zhemis   ! hemispheric index
[825]96      !--------------------------------------------------------------------
[921]97
[4688]98      CALL wrk_alloc( jpi, jpj, zswitch )
[4161]99      CALL wrk_alloc( jpi, jpj, zhemis )
[4688]100      CALL wrk_alloc( jpl,   2, zh_i_ini,  za_i_ini,  zv_i_ini )
101      CALL wrk_alloc(   2,      zht_i_ini, zat_i_ini, zvt_i_ini, zht_s_ini, zsm_i_ini, ztm_i_ini )
[2715]102
[4161]103      IF(lwp) WRITE(numout,*)
104      IF(lwp) WRITE(numout,*) 'lim_istate : Ice initialization '
105      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ '
106
[825]107      !--------------------------------------------------------------------
[4161]108      ! 1) Read namelist
[825]109      !--------------------------------------------------------------------
110
111      CALL lim_istate_init     !  reading the initials parameters of the ice
112
[4688]113      ! surface temperature
114      DO jl = 1, jpl ! loop over categories
[5123]115         t_su  (:,:,jl) = rt0 * tmask(:,:,1)
116         tn_ice(:,:,jl) = rt0 * tmask(:,:,1)
[4688]117      END DO
118
[4990]119      ! basal temperature (considered at freezing point)
[6498]120      CALL eos_fzp( sss_m(:,:), t_bo(:,:) )
121      t_bo(:,:) = ( t_bo(:,:) + rt0 ) * tmask(:,:,1) 
[4990]122
[5123]123      IF( ln_iceini ) THEN
[4688]124
[825]125      !--------------------------------------------------------------------
[4161]126      ! 2) Basal temperature, ice mask and hemispheric index
[825]127      !--------------------------------------------------------------------
[4990]128
129      DO jj = 1, jpj                                       ! ice if sst <= t-freez + ttest
[1037]130         DO ji = 1, jpi
[5407]131            IF( ( sst_m(ji,jj)  - ( t_bo(ji,jj) - rt0 ) ) * tmask(ji,jj,1) >= rn_thres_sst ) THEN
[5123]132               zswitch(ji,jj) = 0._wp * tmask(ji,jj,1)    ! no ice
[4765]133            ELSE                                                                                   
[5123]134               zswitch(ji,jj) = 1._wp * tmask(ji,jj,1)    !    ice
[1037]135            ENDIF
136         END DO
137      END DO
138
139
[4161]140      ! Hemispheric index
[825]141      DO jj = 1, jpj
142         DO ji = 1, jpi
[2528]143            IF( fcor(ji,jj) >= 0._wp ) THEN   
[4161]144               zhemis(ji,jj) = 1 ! Northern hemisphere
145            ELSE
146               zhemis(ji,jj) = 2 ! Southern hemisphere
147            ENDIF
148         END DO
149      END DO
[825]150
[4161]151      !--------------------------------------------------------------------
152      ! 3) Initialization of sea ice state variables
153      !--------------------------------------------------------------------
[825]154
[4161]155      !-----------------------------
156      ! 3.1) Hemisphere-dependent arrays
157      !-----------------------------
[4688]158      ! assign initial thickness, concentration, snow depth and salinity to an hemisphere-dependent array
[5123]159      zht_i_ini(1) = rn_hti_ini_n ; zht_i_ini(2) = rn_hti_ini_s  ! ice thickness
160      zht_s_ini(1) = rn_hts_ini_n ; zht_s_ini(2) = rn_hts_ini_s  ! snow depth
161      zat_i_ini(1) = rn_ati_ini_n ; zat_i_ini(2) = rn_ati_ini_s  ! ice concentration
162      zsm_i_ini(1) = rn_smi_ini_n ; zsm_i_ini(2) = rn_smi_ini_s  ! bulk ice salinity
163      ztm_i_ini(1) = rn_tmi_ini_n ; ztm_i_ini(2) = rn_tmi_ini_s  ! temperature (ice and snow)
[825]164
[4688]165      zvt_i_ini(:) = zht_i_ini(:) * zat_i_ini(:)   ! ice volume
166
[4161]167      !---------------------------------------------------------------------
168      ! 3.2) Distribute ice concentration and thickness into the categories
169      !---------------------------------------------------------------------
170      ! a gaussian distribution for ice concentration is used
171      ! then we check whether the distribution fullfills
172      ! volume and area conservation, positivity and ice categories bounds
173      DO i_hemis = 1, 2
[825]174
[4161]175      ztest_1 = 0 ; ztest_2 = 0 ; ztest_3 = 0 ; ztest_4 = 0
[825]176
[4161]177      ! note for the great nemo engineers:
178      ! only very few of the WRITE statements are necessary for the reference version
179      ! they were one day useful, but now i personally doubt of their
180      ! potential for bringing anything useful
[825]181
[4161]182      DO i_fill = jpl, 1, -1
183         IF ( ( ztest_1 + ztest_2 + ztest_3 + ztest_4 ) .NE. 4 ) THEN
184            !----------------------------
185            ! fill the i_fill categories
186            !----------------------------
187            ! *** 1 category to fill
188            IF ( i_fill .EQ. 1 ) THEN
[4688]189               zh_i_ini(1,i_hemis)       = zht_i_ini(i_hemis)
190               za_i_ini(1,i_hemis)       = zat_i_ini(i_hemis)
191               zh_i_ini(2:jpl,i_hemis)   = 0._wp
192               za_i_ini(2:jpl,i_hemis)   = 0._wp
[4161]193            ELSE
[825]194
[4688]195               ! *** >1 categores to fill
196               !--- Ice thicknesses in the i_fill - 1 first categories
[4161]197               DO jl = 1, i_fill - 1
[5123]198                  zh_i_ini(jl,i_hemis) = hi_mean(jl)
[4161]199               END DO
[4688]200               
201               !--- jl0: most likely index where cc will be maximum
[825]202               DO jl = 1, jpl
[5123]203                  IF ( ( zht_i_ini(i_hemis) >  hi_max(jl-1) ) .AND. &
204                     & ( zht_i_ini(i_hemis) <= hi_max(jl)   ) ) THEN
[4161]205                     jl0 = jl
206                  ENDIF
207               END DO
208               jl0 = MIN(jl0, i_fill)
[4688]209               
210               !--- Concentrations
[4161]211               za_i_ini(jl0,i_hemis)      = zat_i_ini(i_hemis) / SQRT(REAL(jpl))
212               DO jl = 1, i_fill - 1
213                  IF ( jl .NE. jl0 ) THEN
[4688]214                     zsigma               = 0.5 * zht_i_ini(i_hemis)
215                     zarg                 = ( zh_i_ini(jl,i_hemis) - zht_i_ini(i_hemis) ) / zsigma
[4161]216                     za_i_ini(jl,i_hemis) = za_i_ini(jl0,i_hemis) * EXP(-zarg**2)
217                  ENDIF
[4688]218               END DO
219               
[4161]220               zA = 0. ! sum of the areas in the jpl categories
221               DO jl = 1, i_fill - 1
222                 zA = zA + za_i_ini(jl,i_hemis)
223               END DO
224               za_i_ini(i_fill,i_hemis)   = zat_i_ini(i_hemis) - zA ! ice conc in the last category
225               IF ( i_fill .LT. jpl ) za_i_ini(i_fill+1:jpl, i_hemis) = 0._wp
226         
[4688]227               !--- Ice thickness in the last category
[4161]228               zV = 0. ! sum of the volumes of the N-1 categories
229               DO jl = 1, i_fill - 1
[4688]230                  zV = zV + za_i_ini(jl,i_hemis)*zh_i_ini(jl,i_hemis)
[4161]231               END DO
[4688]232               zh_i_ini(i_fill,i_hemis) = ( zvt_i_ini(i_hemis) - zV ) / za_i_ini(i_fill,i_hemis) 
233               IF ( i_fill .LT. jpl ) zh_i_ini(i_fill+1:jpl, i_hemis) = 0._wp
[825]234
[4688]235               !--- volumes
236               zv_i_ini(:,i_hemis) = za_i_ini(:,i_hemis) * zh_i_ini(:,i_hemis)
[4161]237               IF ( i_fill .LT. jpl ) zv_i_ini(i_fill+1:jpl, i_hemis) = 0._wp
[921]238
[4161]239            ENDIF ! i_fill
[825]240
[4161]241            !---------------------
242            ! Compatibility tests
243            !---------------------
244            ! Test 1: area conservation
245            zA_cons = SUM(za_i_ini(:,i_hemis)) ; zconv = ABS(zat_i_ini(i_hemis) - zA_cons )
246            IF ( zconv .LT. 1.0e-6 ) THEN
247               ztest_1 = 1
248            ELSE
249               ztest_1 = 0
250            ENDIF
[825]251
[4161]252            ! Test 2: volume conservation
253            zV_cons = SUM(zv_i_ini(:,i_hemis))
254            zconv = ABS(zvt_i_ini(i_hemis) - zV_cons)
[825]255
[4161]256            IF ( zconv .LT. 1.0e-6 ) THEN
257               ztest_2 = 1
258            ELSE
259               ztest_2 = 0
260            ENDIF
[825]261
[4161]262            ! Test 3: thickness of the last category is in-bounds ?
[5123]263            IF ( zh_i_ini(i_fill, i_hemis) > hi_max(i_fill-1) ) THEN
[4161]264               ztest_3 = 1
265            ELSE
266               ztest_3 = 0
267            ENDIF
[825]268
[4161]269            ! Test 4: positivity of ice concentrations
270            ztest_4 = 1
271            DO jl = 1, jpl
272               IF ( za_i_ini(jl,i_hemis) .LT. 0._wp ) THEN
273                  ztest_4 = 0
[825]274               ENDIF
[4161]275            END DO
[825]276
[4161]277         ENDIF ! ztest_1 + ztest_2 + ztest_3 + ztest_4
278 
279         ztests = ztest_1 + ztest_2 + ztest_3 + ztest_4
[825]280
[4161]281      END DO ! i_fill
[825]282
[4161]283      IF(lwp) THEN
[4908]284         WRITE(numout,*) ' ztests : ', ztests
[4161]285         IF ( ztests .NE. 4 ) THEN
286            WRITE(numout,*)
[4908]287            WRITE(numout,*) ' !!!! ALERT                  !!! '
288            WRITE(numout,*) ' !!!! Something is wrong in the LIM3 initialization procedure '
[4161]289            WRITE(numout,*)
[4908]290            WRITE(numout,*) ' *** ztests is not equal to 4 '
291            WRITE(numout,*) ' *** ztest_i (i=1,4) = ', ztest_1, ztest_2, ztest_3, ztest_4
292            WRITE(numout,*) ' zat_i_ini : ', zat_i_ini(i_hemis)
293            WRITE(numout,*) ' zht_i_ini : ', zht_i_ini(i_hemis)
[4161]294         ENDIF ! ztests .NE. 4
295      ENDIF
296     
297      END DO ! i_hemis
[825]298
[4161]299      !---------------------------------------------------------------------
300      ! 3.3) Space-dependent arrays for ice state variables
301      !---------------------------------------------------------------------
[825]302
[4333]303      ! Ice concentration, thickness and volume, ice salinity, ice age, surface temperature
[4161]304      DO jl = 1, jpl ! loop over categories
305         DO jj = 1, jpj
306            DO ji = 1, jpi
[4688]307               a_i(ji,jj,jl)   = zswitch(ji,jj) * za_i_ini (jl,zhemis(ji,jj))  ! concentration
[5202]308               ht_i(ji,jj,jl)  = zswitch(ji,jj) * zh_i_ini(jl,zhemis(ji,jj))   ! ice thickness
[4688]309               ht_s(ji,jj,jl)  = ht_i(ji,jj,jl) * ( zht_s_ini( zhemis(ji,jj) ) / zht_i_ini( zhemis(ji,jj) ) )  ! snow depth
[5202]310               sm_i(ji,jj,jl)  = zswitch(ji,jj) * zsm_i_ini(zhemis(ji,jj))     ! salinity
311               o_i(ji,jj,jl)   = zswitch(ji,jj) * 1._wp                        ! age (1 day)
[5123]312               t_su(ji,jj,jl)  = zswitch(ji,jj) * ztm_i_ini(zhemis(ji,jj)) + ( 1._wp - zswitch(ji,jj) ) * rt0 ! surf temp
[4333]313
314               ! This case below should not be used if (ht_s/ht_i) is ok in namelist
315               ! In case snow load is in excess that would lead to transformation from snow to ice
316               ! Then, transfer the snow excess into the ice (different from limthd_dh)
317               zdh = MAX( 0._wp, ( rhosn * ht_s(ji,jj,jl) + ( rhoic - rau0 ) * ht_i(ji,jj,jl) ) * r1_rau0 ) 
318               ! recompute ht_i, ht_s avoiding out of bounds values
319               ht_i(ji,jj,jl) = MIN( hi_max(jl), ht_i(ji,jj,jl) + zdh )
[5123]320               ht_s(ji,jj,jl) = MAX( 0._wp, ht_s(ji,jj,jl) - zdh * rhoic * r1_rhosn )
[4333]321
322               ! ice volume, salt content, age content
[4161]323               v_i(ji,jj,jl)   = ht_i(ji,jj,jl) * a_i(ji,jj,jl)              ! ice volume
324               v_s(ji,jj,jl)   = ht_s(ji,jj,jl) * a_i(ji,jj,jl)              ! snow volume
325               smv_i(ji,jj,jl) = MIN( sm_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) ! salt content
326               oa_i(ji,jj,jl)  = o_i(ji,jj,jl) * a_i(ji,jj,jl)               ! age content
[5202]327            END DO
328         END DO
329      END DO
[825]330
[7993]331      ! for constant salinity in time
332      IF( nn_icesal == 1 .OR. nn_icesal == 3 )  THEN
333         CALL lim_var_salprof
334         smv_i = sm_i * v_i
335      ENDIF
336
[4161]337      ! Snow temperature and heat content
338      DO jk = 1, nlay_s
339         DO jl = 1, jpl ! loop over categories
340            DO jj = 1, jpj
341               DO ji = 1, jpi
[5123]342                   t_s(ji,jj,jk,jl) = zswitch(ji,jj) * ztm_i_ini(zhemis(ji,jj)) + ( 1._wp - zswitch(ji,jj) ) * rt0
[4161]343                   ! Snow energy of melting
[5123]344                   e_s(ji,jj,jk,jl) = zswitch(ji,jj) * rhosn * ( cpic * ( rt0 - t_s(ji,jj,jk,jl) ) + lfus )
345
346                   ! Mutliply by volume, and divide by number of layers to get heat content in J/m2
347                   e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s
[5202]348               END DO
349            END DO
350         END DO
351      END DO
[825]352
[4161]353      ! Ice salinity, temperature and heat content
354      DO jk = 1, nlay_i
355         DO jl = 1, jpl ! loop over categories
356            DO jj = 1, jpj
357               DO ji = 1, jpi
[5123]358                   t_i(ji,jj,jk,jl) = zswitch(ji,jj) * ztm_i_ini(zhemis(ji,jj)) + ( 1._wp - zswitch(ji,jj) ) * rt0 
359                   s_i(ji,jj,jk,jl) = zswitch(ji,jj) * zsm_i_ini(zhemis(ji,jj)) !+ ( 1._wp - zswitch(ji,jj) ) * rn_simin
360                   ztmelts          = - tmut * s_i(ji,jj,jk,jl) + rt0 !Melting temperature in K
[825]361
[4161]362                   ! heat content per unit volume
[4688]363                   e_i(ji,jj,jk,jl) = zswitch(ji,jj) * rhoic * (   cpic    * ( ztmelts - t_i(ji,jj,jk,jl) ) &
[5123]364                      +   lfus    * ( 1._wp - (ztmelts-rt0) / MIN((t_i(ji,jj,jk,jl)-rt0),-epsi20) ) &
365                      -   rcp     * ( ztmelts - rt0 ) )
[825]366
[5123]367                   ! Mutliply by ice volume, and divide by number of layers to get heat content in J/m2
368                   e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * v_i(ji,jj,jl) * r1_nlay_i
[5202]369               END DO
370            END DO
371         END DO
372      END DO
[825]373
[4688]374      tn_ice (:,:,:) = t_su (:,:,:)
375
376      ELSE 
[5123]377         ! if ln_iceini=false
[4688]378         a_i  (:,:,:) = 0._wp
379         v_i  (:,:,:) = 0._wp
380         v_s  (:,:,:) = 0._wp
381         smv_i(:,:,:) = 0._wp
382         oa_i (:,:,:) = 0._wp
383         ht_i (:,:,:) = 0._wp
384         ht_s (:,:,:) = 0._wp
385         sm_i (:,:,:) = 0._wp
386         o_i  (:,:,:) = 0._wp
387
388         e_i(:,:,:,:) = 0._wp
389         e_s(:,:,:,:) = 0._wp
390
391         DO jl = 1, jpl
392            DO jk = 1, nlay_i
[5123]393               t_i(:,:,jk,jl) = rt0 * tmask(:,:,1)
[4688]394            END DO
395            DO jk = 1, nlay_s
[5123]396               t_s(:,:,jk,jl) = rt0 * tmask(:,:,1)
[4688]397            END DO
398         END DO
399     
[5123]400      ENDIF ! ln_iceini
[4688]401     
402      at_i (:,:) = 0.0_wp
403      DO jl = 1, jpl
404         at_i (:,:) = at_i (:,:) + a_i (:,:,jl)
405      END DO
406      !
[825]407      !--------------------------------------------------------------------
[4161]408      ! 4) Global ice variables for output diagnostics                    |
[825]409      !--------------------------------------------------------------------
[4161]410      u_ice (:,:)     = 0._wp
411      v_ice (:,:)     = 0._wp
412      stress1_i(:,:)  = 0._wp
413      stress2_i(:,:)  = 0._wp
414      stress12_i(:,:) = 0._wp
[825]415
416      !--------------------------------------------------------------------
[4161]417      ! 5) Moments for advection
[825]418      !--------------------------------------------------------------------
419
[4161]420      sxopw (:,:) = 0._wp 
421      syopw (:,:) = 0._wp 
422      sxxopw(:,:) = 0._wp 
423      syyopw(:,:) = 0._wp 
424      sxyopw(:,:) = 0._wp
[3610]425
[4161]426      sxice (:,:,:)  = 0._wp   ;   sxsn (:,:,:)  = 0._wp   ;   sxa  (:,:,:)  = 0._wp
427      syice (:,:,:)  = 0._wp   ;   sysn (:,:,:)  = 0._wp   ;   sya  (:,:,:)  = 0._wp
428      sxxice(:,:,:)  = 0._wp   ;   sxxsn(:,:,:)  = 0._wp   ;   sxxa (:,:,:)  = 0._wp
429      syyice(:,:,:)  = 0._wp   ;   syysn(:,:,:)  = 0._wp   ;   syya (:,:,:)  = 0._wp
430      sxyice(:,:,:)  = 0._wp   ;   sxysn(:,:,:)  = 0._wp   ;   sxya (:,:,:)  = 0._wp
[825]431
[4161]432      sxc0  (:,:,:)  = 0._wp   ;   sxe  (:,:,:,:)= 0._wp   
433      syc0  (:,:,:)  = 0._wp   ;   sye  (:,:,:,:)= 0._wp   
434      sxxc0 (:,:,:)  = 0._wp   ;   sxxe (:,:,:,:)= 0._wp   
435      syyc0 (:,:,:)  = 0._wp   ;   syye (:,:,:,:)= 0._wp   
436      sxyc0 (:,:,:)  = 0._wp   ;   sxye (:,:,:,:)= 0._wp   
[825]437
[4161]438      sxsal  (:,:,:)  = 0._wp
439      sysal  (:,:,:)  = 0._wp
440      sxxsal (:,:,:)  = 0._wp
441      syysal (:,:,:)  = 0._wp
442      sxysal (:,:,:)  = 0._wp
[825]443
[4335]444      sxage  (:,:,:)  = 0._wp
445      syage  (:,:,:)  = 0._wp
446      sxxage (:,:,:)  = 0._wp
447      syyage (:,:,:)  = 0._wp
448      sxyage (:,:,:)  = 0._wp
449
[825]450
[4688]451      CALL wrk_dealloc( jpi, jpj, zswitch )
[4161]452      CALL wrk_dealloc( jpi, jpj, zhemis )
[4688]453      CALL wrk_dealloc( jpl,   2, zh_i_ini,  za_i_ini,  zv_i_ini )
454      CALL wrk_dealloc(   2,      zht_i_ini, zat_i_ini, zvt_i_ini, zht_s_ini, zsm_i_ini, ztm_i_ini )
[4161]455
[825]456   END SUBROUTINE lim_istate
457
458   SUBROUTINE lim_istate_init
459      !!-------------------------------------------------------------------
460      !!                   ***  ROUTINE lim_istate_init  ***
461      !!       
462      !! ** Purpose : Definition of initial state of the ice
463      !!
[4161]464      !! ** Method : Read the namiceini namelist and check the parameter
465      !!       values called at the first timestep (nit000)
[825]466      !!
[4161]467      !! ** input :
468      !!        Namelist namiceini
469      !!
470      !! history :
471      !!  8.5  ! 03-08 (C. Ethe) original code
472      !!  8.5  ! 07-11 (M. Vancoppenolle) rewritten initialization
[825]473      !!-----------------------------------------------------------------------------
[5123]474      NAMELIST/namiceini/ ln_iceini, rn_thres_sst, rn_hts_ini_n, rn_hts_ini_s, rn_hti_ini_n, rn_hti_ini_s,  &
475         &                                      rn_ati_ini_n, rn_ati_ini_s, rn_smi_ini_n, rn_smi_ini_s, rn_tmi_ini_n, rn_tmi_ini_s
[4147]476      INTEGER :: ios                 ! Local integer output status for namelist read
[825]477      !!-----------------------------------------------------------------------------
[2528]478      !
[4147]479      REWIND( numnam_ice_ref )              ! Namelist namiceini in reference namelist : Ice initial state
480      READ  ( numnam_ice_ref, namiceini, IOSTAT = ios, ERR = 901)
481901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceini in reference namelist', lwp )
482
483      REWIND( numnam_ice_cfg )              ! Namelist namiceini in configuration namelist : Ice initial state
484      READ  ( numnam_ice_cfg, namiceini, IOSTAT = ios, ERR = 902 )
485902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceini in configuration namelist', lwp )
[4624]486      IF(lwm) WRITE ( numoni, namiceini )
[4161]487
488      ! Define the initial parameters
489      ! -------------------------
490
491      IF(lwp) THEN
[825]492         WRITE(numout,*)
493         WRITE(numout,*) 'lim_istate_init : ice parameters inititialisation '
494         WRITE(numout,*) '~~~~~~~~~~~~~~~'
[5123]495         WRITE(numout,*) '   initialization with ice (T) or not (F)       ln_iceini     = ', ln_iceini
496         WRITE(numout,*) '   threshold water temp. for initial sea-ice    rn_thres_sst  = ', rn_thres_sst
497         WRITE(numout,*) '   initial snow thickness in the north          rn_hts_ini_n  = ', rn_hts_ini_n
498         WRITE(numout,*) '   initial snow thickness in the south          rn_hts_ini_s  = ', rn_hts_ini_s 
499         WRITE(numout,*) '   initial ice thickness  in the north          rn_hti_ini_n  = ', rn_hti_ini_n
500         WRITE(numout,*) '   initial ice thickness  in the south          rn_hti_ini_s  = ', rn_hti_ini_s
501         WRITE(numout,*) '   initial ice concentr.  in the north          rn_ati_ini_n  = ', rn_ati_ini_n
502         WRITE(numout,*) '   initial ice concentr.  in the north          rn_ati_ini_s  = ', rn_ati_ini_s
503         WRITE(numout,*) '   initial  ice salinity  in the north          rn_smi_ini_n  = ', rn_smi_ini_n
504         WRITE(numout,*) '   initial  ice salinity  in the south          rn_smi_ini_s  = ', rn_smi_ini_s
505         WRITE(numout,*) '   initial  ice/snw temp  in the north          rn_tmi_ini_n  = ', rn_tmi_ini_n
506         WRITE(numout,*) '   initial  ice/snw temp  in the south          rn_tmi_ini_s  = ', rn_tmi_ini_s
[825]507      ENDIF
[4161]508
[825]509   END SUBROUTINE lim_istate_init
510
511#else
512   !!----------------------------------------------------------------------
513   !!   Default option :         Empty module          NO LIM sea-ice model
514   !!----------------------------------------------------------------------
515CONTAINS
516   SUBROUTINE lim_istate          ! Empty routine
517   END SUBROUTINE lim_istate
518#endif
519
520   !!======================================================================
521END MODULE limistate
Note: See TracBrowser for help on using the repository browser.