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/2016/dev_r7012_ROBUST5_CNRS/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2016/dev_r7012_ROBUST5_CNRS/NEMOGCM/NEMO/LIM_SRC_3/limistate.F90 @ 7177

Last change on this file since 7177 was 7050, checked in by cetlod, 8 years ago

new top interface : update namelists for PISCES

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