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.
limthd_2.F90 in branches/2015/dev_r5204_CNRS_PISCES_dcy/NEMOGCM/NEMO/LIM_SRC_2 – NEMO

source: branches/2015/dev_r5204_CNRS_PISCES_dcy/NEMOGCM/NEMO/LIM_SRC_2/limthd_2.F90 @ 5207

Last change on this file since 5207 was 5207, checked in by cetlod, 9 years ago

dev_r5204_CNRS_PISCES_dcy:minor bugs corrections

  • Property svn:keywords set to Id
File size: 35.8 KB
Line 
1MODULE limthd_2
2   !!======================================================================
3   !!                  ***  MODULE limthd_2   ***
4   !!              LIM thermo ice model : ice thermodynamic
5   !!======================================================================
6   !! History :  1.0  ! 2000-01 (LIM)
7   !!            2.0  ! 2002-07 (C. Ethe, G. Madec) F90
8   !!            2.0  ! 2003-08 (C. Ethe)  add lim_thd_init
9   !!             -   ! 2008-2008  (A. Caubel, G. Madec, E. Maisonnave, S. Masson ) generic coupled interface
10   !!---------------------------------------------------------------------
11#if defined key_lim2
12   !!----------------------------------------------------------------------
13   !!   'key_lim2' :                                  LIM 2.0 sea-ice model
14   !!----------------------------------------------------------------------
15   !!   lim_thd_2       : thermodynamic of sea ice
16   !!   lim_thd_init_2  : initialisation of sea-ice thermodynamic
17   !!----------------------------------------------------------------------
18   USE phycst           ! physical constants
19   USE dom_oce          ! ocean space and time domain variables
20   USE domvvl
21   USE lbclnk
22   USE in_out_manager   ! I/O manager
23   USE lib_mpp
24   USE wrk_nemo         ! work arrays
25   USE iom              ! IOM library
26   USE ice_2            ! LIM sea-ice variables
27   USE sbc_oce          !
28   USE sbc_ice          !
29   USE thd_ice_2        ! LIM thermodynamic sea-ice variables
30   USE dom_ice_2        ! LIM sea-ice domain
31   USE limthd_zdf_2
32   USE limthd_lac_2
33   USE limtab_2
34   USE prtctl           ! Print control
35   USE lib_fortran      ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
36   
37   IMPLICIT NONE
38   PRIVATE
39
40   PUBLIC   lim_thd_2  ! called by lim_step
41
42   REAL(wp) ::   epsi20 = 1.e-20   ! constant values
43   REAL(wp) ::   epsi16 = 1.e-16   !
44   REAL(wp) ::   epsi04 = 1.e-04   !
45   REAL(wp) ::   rzero  = 0.e0     !
46   REAL(wp) ::   rone   = 1.e0     !
47
48   !! * Substitutions
49#  include "domzgr_substitute.h90"
50#  include "vectopt_loop_substitute.h90"
51   !!-------- -------------------------------------------------------------
52   !! NEMO/LIM2 3.3 , UCL - NEMO Consortium (2010)
53   !! $Id$
54   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
55   !!----------------------------------------------------------------------
56CONTAINS
57
58   SUBROUTINE lim_thd_2( kt )
59      !!-------------------------------------------------------------------
60      !!                ***  ROUTINE lim_thd_2  ***       
61      !! 
62      !! ** Purpose : This routine manages the ice thermodynamic.
63      !!         
64      !! ** Action : - Initialisation of some variables
65      !!             - Some preliminary computation (oceanic heat flux
66      !!               at the ice base, snow acc.,heat budget of the leads)
67      !!             - selection of the icy points and put them in an array
68      !!             - call lim_vert_ther for vert ice thermodynamic
69      !!             - back to the geographic grid
70      !!             - selection of points for lateral accretion
71      !!             - call lim_lat_acc  for the ice accretion
72      !!             - back to the geographic grid
73      !!
74      !! References :   Goosse et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90
75      !!---------------------------------------------------------------------
76      INTEGER, INTENT(in) ::   kt     ! number of iteration
77      !!
78      INTEGER  ::   ji, jj               ! dummy loop indices
79      INTEGER  ::   nbpb                 ! nb of icy pts for thermo. cal.
80      INTEGER  ::   nbpac                ! nb of pts for lateral accretion
81      CHARACTER (len=22) :: charout
82      REAL(wp) ::   zfric_umin = 5e-03   ! lower bound for the friction velocity
83      REAL(wp) ::   zfric_umax = 2e-02   ! upper bound for the friction velocity
84      REAL(wp) ::   zinda                ! switch for test. the val. of concen.
85      REAL(wp) ::   zindb, zindg         ! switches for test. the val of arg
86      REAL(wp) ::   zfricp               ! temporary scalar
87      REAL(wp) ::   za , zh, zthsnice    !
88      REAL(wp) ::   zfric_u              ! friction velocity
89      REAL(wp) ::   zfntlat, zpareff     ! test. the val. of lead heat budget
90
91      REAL(wp) ::   zuice_m, zvice_m     ! Sea-ice velocities at U & V-points
92      REAL(wp) ::   zhice_u, zhice_v     ! Sea-ice volume at U & V-points
93      REAL(wp) ::   ztr_fram             ! Sea-ice transport through Fram strait
94      REAL(wp) ::   zrhoij, zrhoijm1     ! temporary scalars
95      REAL(wp) ::   zztmp                ! temporary scalars within a loop
96      REAL(wp), POINTER, DIMENSION(:,:)     ::   ztmp      ! 2D workspace
97      REAL(wp), POINTER, DIMENSION(:,:)     ::   zqlbsbq   ! link with lead energy budget qldif
98      REAL(wp), POINTER, DIMENSION(:,:)     ::   zlicegr   ! link with lateral ice growth
99!!$      REAL(wp), DIMENSION(:,:) ::   firic         ! IR flux over the ice            (outputs only)
100!!$      REAL(wp), DIMENSION(:,:) ::   fcsic         ! Sensible heat flux over the ice (outputs only)
101!!$      REAL(wp), DIMENSION(:,:) ::   fleic         ! Latent heat flux over the ice   (outputs only)
102!!$      REAL(wp), DIMENSION(:,:) ::   qlatic        ! latent flux                     (outputs only)
103      REAL(wp), POINTER, DIMENSION(:,:) ::   zdvosif       ! Variation of volume at surface                (outputs only)
104      REAL(wp), POINTER, DIMENSION(:,:) ::   zdvobif       ! Variation of ice volume at the bottom ice     (outputs only)
105      REAL(wp), POINTER, DIMENSION(:,:) ::   zdvolif       ! Total variation of ice volume                 (outputs only)
106      REAL(wp), POINTER, DIMENSION(:,:) ::   zdvonif       ! Surface accretion Snow to Ice transformation  (outputs only)
107      REAL(wp), POINTER, DIMENSION(:,:) ::   zdvomif       ! Bottom variation of ice volume due to melting (outputs only)
108      REAL(wp), POINTER, DIMENSION(:,:) ::   zu_imasstr    ! Sea-ice transport along i-axis at U-point     (outputs only)
109      REAL(wp), POINTER, DIMENSION(:,:) ::   zv_imasstr    ! Sea-ice transport along j-axis at V-point     (outputs only)
110      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zmsk        ! 3D workspace
111      !!-------------------------------------------------------------------
112
113      CALL wrk_alloc( jpi, jpj, ztmp, zqlbsbq, zlicegr, zdvosif, zdvobif, zdvolif, zdvonif, zdvomif, zu_imasstr, zv_imasstr )
114      CALL wrk_alloc( jpi, jpj, jpk, zmsk )
115
116      IF( kt == nit000 ) THEN
117         CALL lim_thd_init_2  ! Initialization (first time-step only)
118         IF( l_trcdm2dc ) ALLOCATE( fstric_mean(jpi,jpj), fstbif_mean_1d(jpij), qsr_ice_mean_1d(jpij) )
119      ENDIF
120   
121      !-------------------------------------------!
122      !   Initilization of diagnostic variables   !
123      !-------------------------------------------!
124     
125!!gm needed?  yes at least for some of these arrays
126      zdvosif(:,:) = 0.e0   ! variation of ice volume at surface
127      zdvobif(:,:) = 0.e0   ! variation of ice volume at bottom
128      zdvolif(:,:) = 0.e0   ! total variation of ice volume
129      zdvonif(:,:) = 0.e0   ! transformation of snow to sea-ice volume
130      zlicegr(:,:) = 0.e0   ! lateral variation of ice volume
131      zdvomif(:,:) = 0.e0   ! variation of ice volume at bottom due to melting only
132      ztr_fram     = 0.e0   ! sea-ice transport through Fram strait
133      fstric (:,:) = 0.e0   ! part of solar radiation absorbing inside the ice
134      fscmbq (:,:) = 0.e0   ! linked with fstric
135      ffltbif(:,:) = 0.e0   ! linked with fstric
136      qfvbq  (:,:) = 0.e0   ! linked with fstric
137      rdm_snw(:,:) = 0.e0   ! variation of snow mass over 1 time step
138      rdq_snw(:,:) = 0.e0   ! heat content associated with rdm_snw
139      rdm_ice(:,:) = 0.e0   ! variation of ice mass over 1 time step
140      rdq_ice(:,:) = 0.e0   ! heat content associated with rdm_ice
141      zmsk (:,:,:) = 0.e0
142      !
143      IF( l_trcdm2dc  ) fstric_mean(:,:) = 0.e0   ! part of solar radiation absorbing inside the ice
144
145      ! set to zero snow thickness smaller than epsi04
146      DO jj = 1, jpj
147         DO ji = 1, jpi
148            hsnif(ji,jj)  = hsnif(ji,jj) *  MAX( rzero, SIGN( rone , hsnif(ji,jj) - epsi04 ) )
149         END DO
150      END DO
151!!gm better coded (do not use SIGN...)
152!     WHERE( hsnif(:,:) < epsi04 )   hsnif(:,:) = 0.e0
153!!gm
154
155      IF(ln_ctl)   CALL prt_ctl( tab2d_1=hsnif, clinfo1=' lim_thd: hsnif   : ' )
156     
157      !-----------------------------------!
158      !   Treatment of particular cases   !
159      !-----------------------------------!
160     
161      DO jj = 1, jpj
162         DO ji = 1, jpi
163            !  snow is transformed into ice if the original ice cover disappears.
164            zindg         = tms(ji,jj) *  MAX( rzero , SIGN( rone , -hicif(ji,jj) ) )
165            hicif(ji,jj)  = hicif(ji,jj) + zindg * rhosn * hsnif(ji,jj) / rau0
166            hsnif(ji,jj)  = ( rone - zindg ) * hsnif(ji,jj) + zindg * hicif(ji,jj) * ( rau0 - rhoic ) / rhosn
167            dmgwi(ji,jj)  = zindg * (1.0 - frld(ji,jj)) * rhoic * hicif(ji,jj)   ! snow/ice mass
168           
169            !  the lead fraction, frld, must be little than or equal to amax (ice ridging).
170            zthsnice      = hsnif(ji,jj) + hicif(ji,jj)
171            zindb         = tms(ji,jj) * ( 1.0 - MAX( rzero , SIGN( rone , - zthsnice ) ) ) 
172            za            = zindb * MIN( rone, ( 1.0 - frld(ji,jj) ) * uscomi )
173            hsnif (ji,jj) = hsnif(ji,jj)  * za
174            hicif (ji,jj) = hicif(ji,jj)  * za
175            qstoif(ji,jj) = qstoif(ji,jj) * za
176            frld  (ji,jj) = 1.0 - zindb * ( 1.0 - frld(ji,jj) ) / MAX( za, epsi20 )
177           
178            !  the in situ ice thickness, hicif, must be equal to or greater than hiclim.
179            zh            = MAX( rone , zindb * hiclim  / MAX( hicif(ji,jj), epsi20 ) )
180            hsnif (ji,jj) = hsnif(ji,jj)  * zh
181            hicif (ji,jj) = hicif(ji,jj)  * zh
182            qstoif(ji,jj) = qstoif(ji,jj) * zh
183            frld  (ji,jj) = ( frld(ji,jj) + ( zh - 1.0 ) ) / zh
184         END DO
185      END DO
186
187      IF(ln_ctl) THEN
188         CALL prt_ctl( tab2d_1=hicif , clinfo1=' lim_thd: hicif   : ' )
189         CALL prt_ctl( tab2d_1=hsnif , clinfo1=' lim_thd: hsnif   : ' )
190         CALL prt_ctl( tab2d_1=dmgwi , clinfo1=' lim_thd: dmgwi   : ' )
191         CALL prt_ctl( tab2d_1=qstoif, clinfo1=' lim_thd: qstoif  : ' )
192         CALL prt_ctl( tab2d_1=frld  , clinfo1=' lim_thd: frld    : ' )
193      ENDIF
194
195     
196      !-------------------------------!
197      !   Thermodynamics of sea ice   !
198      !-------------------------------!
199     
200      !      Partial computation of forcing for the thermodynamic sea ice model.
201      !--------------------------------------------------------------------------
202
203      !CDIR NOVERRCHK
204      DO jj = 1, jpj
205         !CDIR NOVERRCHK
206         DO ji = 1, jpi
207            zthsnice       = hsnif(ji,jj) + hicif(ji,jj)
208            zindb          = tms(ji,jj) * ( 1.0 - MAX( rzero , SIGN( rone , - zthsnice ) ) ) 
209            pfrld(ji,jj)   = frld(ji,jj)
210            zfricp         = 1.0 - frld(ji,jj)
211            zinda          = 1.0 - MAX( rzero , SIGN( rone , - zfricp ) )
212           
213            !  solar irradiance transmission at the mixed layer bottom and used in the lead heat budget
214            thcm(ji,jj)    = 0.e0 
215           
216            !  net downward heat flux from the ice to the ocean, expressed as a function of ocean
217            !  temperature and turbulent mixing (McPhee, 1992)
218            zfric_u        = MAX ( MIN( SQRT( ust2s(ji,jj) ) , zfric_umax ) , zfric_umin )  ! friction velocity
219            fdtcn(ji,jj)  = zindb * rau0 * rcp * 0.006  * zfric_u * ( sst_m(ji,jj) + rt0 - tfu(ji,jj) ) 
220            qdtcn(ji,jj)  = zindb * fdtcn(ji,jj) * frld(ji,jj) * rdt_ice
221                       
222            !  partial computation of the lead energy budget (qldif)
223            IF( lk_cpl ) THEN
224               qldif(ji,jj)   = tms(ji,jj) * rdt_ice                                                  &
225                  &    * (   ( qsr_tot(ji,jj) - qsr_ice(ji,jj,1) * zfricp ) * ( 1.0 - thcm(ji,jj) )   &
226                  &        + ( qns_tot(ji,jj) - qns_ice(ji,jj,1) * zfricp )                           &
227                  &        + frld(ji,jj) * ( fdtcn(ji,jj) + ( 1.0 - zindb ) * fsbbq(ji,jj) )   )
228            ELSE
229               qldif(ji,jj)   = tms(ji,jj) * rdt_ice * frld(ji,jj)                    &
230                  &                        * (  qsr(ji,jj) * ( 1.0 - thcm(ji,jj) )    &
231                  &                           + qns(ji,jj)  +  fdtcn(ji,jj)           &
232                  &                           + ( 1.0 - zindb ) * fsbbq(ji,jj)      )
233            ENDIF
234            !  parlat : percentage of energy used for lateral ablation (0.0)
235            zfntlat        = 1.0 - MAX( rzero , SIGN( rone ,  - qldif(ji,jj) ) )
236            zpareff        = 1.0 + ( parlat - 1.0 ) * zinda * zfntlat
237            zqlbsbq(ji,jj) = qldif(ji,jj) * ( 1.0 - zpareff ) / MAX( (1.0 - frld(ji,jj)) * rdt_ice , epsi16 )
238            qldif  (ji,jj) = zpareff *  qldif(ji,jj)
239            qdtcn  (ji,jj) = zpareff * qdtcn(ji,jj)
240           
241            !  energy needed to bring ocean surface layer until its freezing
242            qcmif  (ji,jj) =  rau0 * rcp * fse3t_m(ji,jj) * ( tfu(ji,jj) - sst_m(ji,jj) - rt0 ) * ( 1 - zinda )
243           
244            !  calculate oceanic heat flux.
245            fbif   (ji,jj) = zindb * (  fsbbq(ji,jj) / MAX( (1.0 - frld(ji,jj)) , epsi20 ) + fdtcn(ji,jj) )
246           
247            ! computation of the thermodynamic ice production (only needed for output)
248            hicifp(ji,jj) = hicif(ji,jj) * ( 1.0 - frld(ji,jj) )
249         END DO
250      END DO
251     
252      !         Select icy points and fulfill arrays for the vectorial grid.
253      !----------------------------------------------------------------------
254      nbpb = 0
255      DO jj = 1, jpj
256         DO ji = 1, jpi
257            IF ( frld(ji,jj) < 1.0 ) THEN     
258               nbpb      = nbpb + 1
259               npb(nbpb) = (jj - 1) * jpi + ji
260            ENDIF
261         END DO
262      END DO
263
264      IF(ln_ctl) THEN
265         CALL prt_ctl(tab2d_1=pfrld, clinfo1=' lim_thd: pfrld   : ', tab2d_2=thcm   , clinfo2='  thcm    : ')
266         CALL prt_ctl(tab2d_1=fdtcn, clinfo1=' lim_thd: fdtcn   : ', tab2d_2=qdtcn  , clinfo2='  qdtcn   : ')
267         CALL prt_ctl(tab2d_1=qldif, clinfo1=' lim_thd: qldif   : ', tab2d_2=zqlbsbq, clinfo2='  zqlbsbq : ')
268         CALL prt_ctl(tab2d_1=qcmif, clinfo1=' lim_thd: qcmif   : ', tab2d_2=fbif   , clinfo2='  fbif    : ')
269         zmsk(:,:,1) = tms(:,:)
270         CALL prt_ctl(tab2d_1=qcmif , clinfo1=' lim_thd: qcmif  : ', mask1=zmsk)
271         CALL prt_ctl(tab2d_1=hicifp, clinfo1=' lim_thd: hicifp : ')
272         WRITE(charout, FMT="('lim_thd: nbpb = ',I4)") nbpb
273         CALL prt_ctl_info(charout)
274      ENDIF
275     
276     
277      ! If there is no ice, do nothing. Otherwise, compute Top and Bottom accretion/ablation
278      !------------------------------------------------------------------------------------
279
280      IF( nbpb > 0 ) THEN
281         !   
282         !  put the variable in a 1-D array for thermodynamics process
283         CALL tab_2d_1d_2( nbpb, frld_1d    (1:nbpb)     , frld           , jpi, jpj, npb(1:nbpb) )
284         CALL tab_2d_1d_2( nbpb, h_ice_1d   (1:nbpb)     , hicif          , jpi, jpj, npb(1:nbpb) )
285         CALL tab_2d_1d_2( nbpb, h_snow_1d  (1:nbpb)     , hsnif          , jpi, jpj, npb(1:nbpb) )
286         CALL tab_2d_1d_2( nbpb, sist_1d    (1:nbpb)     , sist           , jpi, jpj, npb(1:nbpb) )
287         CALL tab_2d_1d_2( nbpb, tbif_1d    (1:nbpb , 1 ), tbif(:,:,1)    , jpi, jpj, npb(1:nbpb) )
288         CALL tab_2d_1d_2( nbpb, tbif_1d    (1:nbpb , 2 ), tbif(:,:,2)    , jpi, jpj, npb(1:nbpb) )
289         CALL tab_2d_1d_2( nbpb, tbif_1d    (1:nbpb , 3 ), tbif(:,:,3)    , jpi, jpj, npb(1:nbpb) )
290         CALL tab_2d_1d_2( nbpb, qsr_ice_1d (1:nbpb)     , qsr_ice(:,:,1) , jpi, jpj, npb(1:nbpb) )
291         IF( l_trcdm2dc ) &
292         & CALL tab_2d_1d_2( nbpb, qsr_ice_mean_1d (1:nbpb), qsr_ice_mean(:,:,1), jpi, jpj, npb(1:nbpb) )
293         CALL tab_2d_1d_2( nbpb, fr1_i0_1d  (1:nbpb)     , fr1_i0         , jpi, jpj, npb(1:nbpb) )
294         CALL tab_2d_1d_2( nbpb, fr2_i0_1d  (1:nbpb)     , fr2_i0         , jpi, jpj, npb(1:nbpb) )
295         CALL tab_2d_1d_2( nbpb,  qns_ice_1d(1:nbpb)     ,  qns_ice(:,:,1), jpi, jpj, npb(1:nbpb) )
296         CALL tab_2d_1d_2( nbpb, dqns_ice_1d(1:nbpb)     , dqns_ice(:,:,1), jpi, jpj, npb(1:nbpb) )
297         IF( .NOT. lk_cpl ) THEN
298            CALL tab_2d_1d_2( nbpb, qla_ice_1d (1:nbpb)     ,  qla_ice(:,:,1), jpi, jpj, npb(1:nbpb) )
299            CALL tab_2d_1d_2( nbpb, dqla_ice_1d(1:nbpb)     , dqla_ice(:,:,1), jpi, jpj, npb(1:nbpb) )
300         ENDIF
301         CALL tab_2d_1d_2( nbpb, tfu_1d     (1:nbpb)     , tfu        , jpi, jpj, npb(1:nbpb) )
302         CALL tab_2d_1d_2( nbpb, sprecip_1d (1:nbpb)     , sprecip    , jpi, jpj, npb(1:nbpb) ) 
303         CALL tab_2d_1d_2( nbpb, fbif_1d    (1:nbpb)     , fbif       , jpi, jpj, npb(1:nbpb) )
304         CALL tab_2d_1d_2( nbpb, thcm_1d    (1:nbpb)     , thcm       , jpi, jpj, npb(1:nbpb) )
305         CALL tab_2d_1d_2( nbpb, qldif_1d   (1:nbpb)     , qldif      , jpi, jpj, npb(1:nbpb) )
306         CALL tab_2d_1d_2( nbpb, qstbif_1d  (1:nbpb)     , qstoif     , jpi, jpj, npb(1:nbpb) )
307         CALL tab_2d_1d_2( nbpb, rdm_ice_1d (1:nbpb)     , rdm_ice    , jpi, jpj, npb(1:nbpb) )
308         CALL tab_2d_1d_2( nbpb, rdq_ice_1d (1:nbpb)     , rdq_ice    , jpi, jpj, npb(1:nbpb) )
309         CALL tab_2d_1d_2( nbpb, dmgwi_1d   (1:nbpb)     , dmgwi      , jpi, jpj, npb(1:nbpb) )
310         CALL tab_2d_1d_2( nbpb, rdm_snw_1d (1:nbpb)     , rdm_snw    , jpi, jpj, npb(1:nbpb) )
311         CALL tab_2d_1d_2( nbpb, rdq_snw_1d (1:nbpb)     , rdq_snw    , jpi, jpj, npb(1:nbpb) )
312         CALL tab_2d_1d_2( nbpb, qlbbq_1d   (1:nbpb)     , zqlbsbq    , jpi, jpj, npb(1:nbpb) )
313         !
314         CALL lim_thd_zdf_2( 1, nbpb )       !  compute ice growth
315         !
316         !  back to the geographic grid.
317         CALL tab_1d_2d_2( nbpb, frld       , npb, frld_1d   (1:nbpb)     , jpi, jpj )
318         CALL tab_1d_2d_2( nbpb, hicif      , npb, h_ice_1d  (1:nbpb)     , jpi, jpj )
319         CALL tab_1d_2d_2( nbpb, hsnif      , npb, h_snow_1d (1:nbpb)     , jpi, jpj )
320         CALL tab_1d_2d_2( nbpb, sist       , npb, sist_1d   (1:nbpb)     , jpi, jpj )
321         CALL tab_1d_2d_2( nbpb, tbif(:,:,1), npb, tbif_1d   (1:nbpb , 1 ), jpi, jpj )   
322         CALL tab_1d_2d_2( nbpb, tbif(:,:,2), npb, tbif_1d   (1:nbpb , 2 ), jpi, jpj )   
323         CALL tab_1d_2d_2( nbpb, tbif(:,:,3), npb, tbif_1d   (1:nbpb , 3 ), jpi, jpj )   
324         CALL tab_1d_2d_2( nbpb, fscmbq     , npb, fscbq_1d  (1:nbpb)     , jpi, jpj )
325         CALL tab_1d_2d_2( nbpb, ffltbif    , npb, fltbif_1d (1:nbpb)     , jpi, jpj )
326         CALL tab_1d_2d_2( nbpb, fstric     , npb, fstbif_1d (1:nbpb)     , jpi, jpj )
327         CALL tab_1d_2d_2( nbpb, qldif      , npb, qldif_1d  (1:nbpb)     , jpi, jpj )
328         CALL tab_1d_2d_2( nbpb, qfvbq      , npb, qfvbq_1d  (1:nbpb)     , jpi, jpj )
329         CALL tab_1d_2d_2( nbpb, qstoif     , npb, qstbif_1d (1:nbpb)     , jpi, jpj )
330         CALL tab_1d_2d_2( nbpb, rdm_ice    , npb, rdm_ice_1d(1:nbpb)     , jpi, jpj )
331         CALL tab_1d_2d_2( nbpb, rdq_ice    , npb, rdq_ice_1d(1:nbpb)     , jpi, jpj )
332         CALL tab_1d_2d_2( nbpb, dmgwi      , npb, dmgwi_1d  (1:nbpb)     , jpi, jpj )
333         CALL tab_1d_2d_2( nbpb, rdm_snw    , npb, rdm_snw_1d(1:nbpb)     , jpi, jpj )
334         CALL tab_1d_2d_2( nbpb, rdq_snw    , npb, rdq_snw_1d(1:nbpb)     , jpi, jpj )
335         CALL tab_1d_2d_2( nbpb, zdvosif    , npb, dvsbq_1d  (1:nbpb)     , jpi, jpj )
336         CALL tab_1d_2d_2( nbpb, zdvobif    , npb, dvbbq_1d  (1:nbpb)     , jpi, jpj )
337         CALL tab_1d_2d_2( nbpb, zdvomif    , npb, rdvomif_1d(1:nbpb)     , jpi, jpj )
338         CALL tab_1d_2d_2( nbpb, zdvolif    , npb, dvlbq_1d  (1:nbpb)     , jpi, jpj )
339         CALL tab_1d_2d_2( nbpb, zdvonif    , npb, dvnbq_1d  (1:nbpb)     , jpi, jpj ) 
340         CALL tab_1d_2d_2( nbpb, qsr_ice(:,:,1), npb, qsr_ice_1d(1:nbpb)  , jpi, jpj )
341         CALL tab_1d_2d_2( nbpb, qns_ice(:,:,1), npb, qns_ice_1d(1:nbpb)  , jpi, jpj )
342         IF( l_trcdm2dc ) THEN
343            CALL tab_1d_2d_2( nbpb, fstric_mean        , npb, fstbif_mean_1d (1:nbpb), jpi, jpj )
344            CALL tab_1d_2d_2( nbpb, qsr_ice_mean(:,:,1), npb, qsr_ice_mean_1d(1:nbpb), jpi, jpj )
345         ENDIF
346         IF( .NOT. lk_cpl )   CALL tab_1d_2d_2( nbpb, qla_ice(:,:,1), npb, qla_ice_1d(1:nbpb), jpi, jpj )
347         !
348      ENDIF
349
350      ! Up-date sea ice thickness
351      !--------------------------
352      DO jj = 1, jpj
353         DO ji = 1, jpi
354            phicif(ji,jj) = hicif(ji,jj) 
355            hicif(ji,jj)  = hicif(ji,jj) *  ( rone -  MAX( rzero, SIGN( rone, - ( 1.0 - frld(ji,jj) ) ) ) )
356         END DO
357      END DO
358
359     
360      ! Tricky trick : add 2 to frld in the Southern Hemisphere
361      !--------------------------------------------------------
362      IF( fcor(1,1) < 0.e0 ) THEN
363         DO jj = 1, njeqm1
364            DO ji = 1, jpi
365               frld(ji,jj) = frld(ji,jj) + 2.0
366            END DO
367         END DO
368      ENDIF
369
370      CALL lbc_lnk( frld , 'T', 1. )     
371     
372      ! Select points for lateral accretion (this occurs when heat exchange
373      ! between ice and ocean is negative; ocean losing heat)
374      !-----------------------------------------------------------------
375      nbpac = 0
376      DO jj = 1, jpj
377         DO ji = 1, jpi
378!i yes!     IF ( ( qcmif(ji,jj) - qldif(ji,jj) ) > 0.e0 ) THEN
379            IF ( tms(ji,jj) * ( qcmif(ji,jj) - qldif(ji,jj) ) > 0.e0 ) THEN
380               nbpac = nbpac + 1
381               npac( nbpac ) = (jj - 1) * jpi + ji
382            ENDIF
383         END DO
384      END DO
385     
386      IF(ln_ctl) THEN
387         CALL prt_ctl(tab2d_1=phicif, clinfo1=' lim_thd: phicif  : ', tab2d_2=hicif, clinfo2=' hicif : ')
388         WRITE(charout, FMT="('lim_thd: nbpac = ',I4)") nbpac
389         CALL prt_ctl_info(charout)
390      ENDIF
391
392
393      ! If ocean gains heat do nothing ; otherwise, one performs lateral accretion
394      !--------------------------------------------------------------------------------
395      IF( nbpac > 0 ) THEN
396         !
397         zlicegr(:,:) = rdm_ice(:,:)      ! to output the lateral sea-ice growth
398         !...Put the variable in a 1-D array for lateral accretion
399         CALL tab_2d_1d_2( nbpac, frld_1d   (1:nbpac)     , frld       , jpi, jpj, npac(1:nbpac) )
400         CALL tab_2d_1d_2( nbpac, h_snow_1d (1:nbpac)     , hsnif      , jpi, jpj, npac(1:nbpac) )
401         CALL tab_2d_1d_2( nbpac, h_ice_1d  (1:nbpac)     , hicif      , jpi, jpj, npac(1:nbpac) )
402         CALL tab_2d_1d_2( nbpac, tbif_1d   (1:nbpac , 1 ), tbif(:,:,1), jpi, jpj, npac(1:nbpac) )   
403         CALL tab_2d_1d_2( nbpac, tbif_1d   (1:nbpac , 2 ), tbif(:,:,2), jpi, jpj, npac(1:nbpac) )   
404         CALL tab_2d_1d_2( nbpac, tbif_1d   (1:nbpac , 3 ), tbif(:,:,3), jpi, jpj, npac(1:nbpac) )   
405         CALL tab_2d_1d_2( nbpac, qldif_1d  (1:nbpac)     , qldif      , jpi, jpj, npac(1:nbpac) )
406         CALL tab_2d_1d_2( nbpac, qcmif_1d  (1:nbpac)     , qcmif      , jpi, jpj, npac(1:nbpac) )
407         CALL tab_2d_1d_2( nbpac, qstbif_1d (1:nbpac)     , qstoif     , jpi, jpj, npac(1:nbpac) )
408         CALL tab_2d_1d_2( nbpac, rdm_ice_1d(1:nbpac)     , rdm_ice    , jpi, jpj, npac(1:nbpac) )
409         CALL tab_2d_1d_2( nbpac, rdq_ice_1d(1:nbpac)     , rdq_ice    , jpi, jpj, npac(1:nbpac) )
410         CALL tab_2d_1d_2( nbpac, dvlbq_1d  (1:nbpac)     , zdvolif    , jpi, jpj, npac(1:nbpac) )
411         CALL tab_2d_1d_2( nbpac, tfu_1d    (1:nbpac)     , tfu        , jpi, jpj, npac(1:nbpac) )
412         !
413         CALL lim_thd_lac_2( 1 , nbpac )         ! lateral accretion routine.
414         !
415         !   back to the geographic grid
416         CALL tab_1d_2d_2( nbpac, frld       , npac(1:nbpac), frld_1d   (1:nbpac)     , jpi, jpj )
417         CALL tab_1d_2d_2( nbpac, hsnif      , npac(1:nbpac), h_snow_1d (1:nbpac)     , jpi, jpj )
418         CALL tab_1d_2d_2( nbpac, hicif      , npac(1:nbpac), h_ice_1d  (1:nbpac)     , jpi, jpj )
419         CALL tab_1d_2d_2( nbpac, tbif(:,:,1), npac(1:nbpac), tbif_1d   (1:nbpac , 1 ), jpi, jpj )
420         CALL tab_1d_2d_2( nbpac, tbif(:,:,2), npac(1:nbpac), tbif_1d   (1:nbpac , 2 ), jpi, jpj )
421         CALL tab_1d_2d_2( nbpac, tbif(:,:,3), npac(1:nbpac), tbif_1d   (1:nbpac , 3 ), jpi, jpj )
422         CALL tab_1d_2d_2( nbpac, qstoif     , npac(1:nbpac), qstbif_1d (1:nbpac)     , jpi, jpj )
423         CALL tab_1d_2d_2( nbpac, rdm_ice    , npac(1:nbpac), rdm_ice_1d(1:nbpac)     , jpi, jpj )
424         CALL tab_1d_2d_2( nbpac, rdq_ice    , npac(1:nbpac), rdq_ice_1d(1:nbpac)     , jpi, jpj )
425         CALL tab_1d_2d_2( nbpac, zdvolif    , npac(1:nbpac), dvlbq_1d  (1:nbpac)     , jpi, jpj )
426         !
427      ENDIF
428       
429       
430      ! Recover frld values between 0 and 1 in the Southern Hemisphere (tricky trick)
431      ! Update daily thermodynamic ice production.   
432      !------------------------------------------------------------------------------
433      DO jj = 1, jpj
434         DO ji = 1, jpi
435            frld  (ji,jj) = MIN( frld(ji,jj), ABS( frld(ji,jj) - 2.0 ) )
436            fr_i  (ji,jj) = 1.0 - frld(ji,jj) 
437            hicifp(ji,jj) = hicif(ji,jj) * fr_i(ji,jj) - hicifp(ji,jj)
438         END DO
439      END DO
440
441      ! Outputs
442      !--------------------------------------------------------------------------------
443      ztmp(:,:) = 1. - pfrld(:,:)                                ! fraction of ice after the dynamic, before the thermodynamic
444      IF( iom_use('ist_cea'    ) )   CALL iom_put( 'ist_cea', (sist(:,:) - rt0) * ztmp(:,:) )   ! Ice surface temperature   [Celius]
445      IF( iom_use('qsr_ai_cea' ) )   CALL iom_put( 'qsr_ai_cea', qsr_ice(:,:,1) * ztmp(:,:) )   ! Solar flux over the ice     [W/m2]
446      IF( iom_use('qns_ai_cea' ) )   CALL iom_put( 'qns_ai_cea', qns_ice(:,:,1) * ztmp(:,:) )   ! Non-solar flux over the ice [W/m2]
447      IF( iom_use('qla_ai_cea' ) .AND. .NOT. lk_cpl ) &
448         &                           CALL iom_put( 'qla_ai_cea', qla_ice(:,:,1) * ztmp(:,:) )   ! Latent flux over the ice [W/m2]
449      !
450      IF( iom_use('snowthic_cea'))   CALL iom_put( 'snowthic_cea', hsnif  (:,:) * fr_i(:,:) )   ! Snow thickness           [m]
451      IF( iom_use('icethic_cea' ))   CALL iom_put( 'icethic_cea' , hicif  (:,:) * fr_i(:,:) )   ! Ice thickness            [m]
452      zztmp = 1.0 / rdt_ice
453      IF( iom_use('iceprod_cea') )   CALL iom_put( 'iceprod_cea' , hicifp (:,:) * zztmp     )   ! Ice produced             [m/s]
454      IF( iom_use('iiceconc'   ) )   CALL iom_put( 'iiceconc'    , fr_i(:,:)                )   ! Ice concentration        [-]
455      IF( iom_use('snowmel_cea') )   CALL iom_put( 'snowmel_cea' , rdm_snw(:,:) * zztmp     )   ! Snow melt                [kg/m2/s]
456      zztmp = rhoic / rdt_ice
457      IF( iom_use('sntoice_cea') )   CALL iom_put( 'sntoice_cea' , zdvonif(:,:) * zztmp     ) ! Snow to Ice transformation [kg/m2/s]
458      IF( iom_use('ticemel_cea') )   CALL iom_put( 'ticemel_cea' , zdvosif(:,:) * zztmp     )   ! Melt at Sea Ice top      [kg/m2/s]
459      IF( iom_use('bicemel_cea') )   CALL iom_put( 'bicemel_cea' , zdvomif(:,:) * zztmp     )   ! Melt at Sea Ice bottom   [kg/m2/s]
460      IF( iom_use('licepro_cea') ) THEN
461         zlicegr(:,:) = MAX( 0.e0, rdm_ice(:,:)-zlicegr(:,:) )
462                                     CALL iom_put( 'licepro_cea' , zlicegr(:,:) * zztmp     )   ! Lateral sea ice growth   [kg/m2/s]
463      ENDIF
464      !
465      ! Compute the Eastward & Northward sea-ice transport
466      IF( iom_use('u_imasstr') ) THEN
467         zztmp = 0.25 * rhoic
468         DO jj = 1, jpjm1 
469            DO ji = 1, jpim1   ! NO vector opt.
470               ! Ice velocities, volume & transport at U-points
471               zuice_m = u_ice(ji+1,jj+1) + u_ice(ji+1,jj )
472               zhice_u = hicif(ji,jj)*e2t(ji,jj)*fr_i(ji,jj) + hicif(ji+1,jj  )*e2t(ji+1,jj  )*fr_i(ji+1,jj  )
473               zu_imasstr(ji,jj) = zztmp * zhice_u * zuice_m 
474            END DO
475         END DO
476         CALL lbc_lnk( zu_imasstr, 'U', -1. )
477         CALL iom_put( 'u_imasstr',  zu_imasstr(:,:) )   ! Ice transport along i-axis at U-point [kg/s]
478      ENDIF
479      IF( iom_use('v_imasstr') ) THEN
480         zztmp = 0.25 * rhoic
481         DO jj = 1, jpjm1 
482            DO ji = 1, jpim1   ! NO vector opt.
483               ! Ice velocities, volume & transport at V-points
484               zvice_m = v_ice(ji+1,jj+1) + v_ice(ji ,jj+1)
485               zhice_v = hicif(ji,jj)*e1t(ji,jj)*fr_i(ji,jj) + hicif(ji  ,jj+1)*e1t(ji  ,jj+1)*fr_i(ji  ,jj+1)
486               zv_imasstr(ji,jj) = zztmp * zhice_v * zvice_m 
487            END DO
488         END DO
489         CALL lbc_lnk( zv_imasstr, 'V', -1. )
490         CALL iom_put( 'v_imasstr',  zv_imasstr(:,:) )   ! Ice transport along j-axis at V-point [kg/s]
491      ENDIF
492
493      !! Fram Strait sea-ice transport (sea-ice + snow)  (in ORCA2 = 5 points)
494      IF( iom_use('fram_trans') .and. cp_cfg == "orca" .AND. jp_cfg == 2 ) THEN    ! ORCA R2 configuration
495         DO jj = mj0(137), mj1(137) ! B grid
496            IF( mj0(jj-1) >= nldj ) THEN
497               DO ji = MAX(mi0(134),nldi), MIN(mi1(138),nlei)
498                  zrhoij    = e1t(ji,jj  ) * fr_i(ji,jj  ) * ( rhoic*hicif(ji,jj  ) + rhosn*hsnif(ji,jj  ) ) 
499                  zrhoijm1  = e1t(ji,jj-1) * fr_i(ji,jj-1) * ( rhoic*hicif(ji,jj-1) + rhosn*hsnif(ji,jj-1) ) 
500                  ztr_fram  = ztr_fram - 0.25 * ( v_ice(ji,jj)+ v_ice(ji+1,jj) ) * ( zrhoij + zrhoijm1 )
501               END DO
502            ENDIF
503         END DO
504         IF( lk_mpp )   CALL mpp_sum( ztr_fram )
505         CALL iom_put( 'fram_trans', ztr_fram )   ! Ice transport through Fram strait     [kg/s]
506      ENDIF
507
508      IF( iom_use('ice_pres') .OR. iom_use('ist_ipa') .OR. iom_use('uice_ipa') .OR. iom_use('vice_ipa') ) THEN
509!! ce     ztmp(:,:) = 1. - AINT( frld(:,:), wp )                        ! return 1 as soon as there is ice
510!! ce     A big warning because the model crashes on IDRIS/IBM SP6 with xlf 13.1.0.3, see ticket #761
511!! ce     We Unroll the loop and everything works fine     
512         DO jj = 1, jpj
513            DO ji = 1, jpi
514               ztmp(ji,jj) = 1. - AINT( frld(ji,jj), wp )                ! return 1 as soon as there is ice
515            END DO
516         END DO
517         !
518         IF( iom_use('ice_pres') ) CALL iom_put( 'ice_pres', ztmp                            )   ! Ice presence                 [-]
519         IF( iom_use('ist_ipa' ) ) CALL iom_put( 'ist_ipa' , ( sist(:,:) - rt0 ) * ztmp(:,:) )   ! Ice surface temperature [Celius]
520         IF( iom_use('uice_ipa') ) CALL iom_put( 'uice_ipa', u_ice(:,:) * ztmp(:,:) ) ! Ice velocity along i-axis at I-point  [m/s]
521         IF( iom_use('vice_ipa') ) CALL iom_put( 'vice_ipa', v_ice(:,:) * ztmp(:,:) ) ! Ice velocity along j-axis at I-point  [m/s]
522      ENDIF
523
524      IF(ln_ctl) THEN
525         CALL prt_ctl_info(' lim_thd  end  ')
526         CALL prt_ctl( tab2d_1=hicif      , clinfo1=' lim_thd: hicif   : ', tab2d_2=hsnif , clinfo2=' hsnif  : ' )
527         CALL prt_ctl( tab2d_1=frld       , clinfo1=' lim_thd: frld    : ', tab2d_2=hicifp, clinfo2=' hicifp : ' )
528         CALL prt_ctl( tab2d_1=phicif     , clinfo1=' lim_thd: phicif  : ', tab2d_2=pfrld , clinfo2=' pfrld  : ' )
529         CALL prt_ctl( tab2d_1=sist       , clinfo1=' lim_thd: sist    : ' )
530         CALL prt_ctl( tab2d_1=tbif(:,:,1), clinfo1=' lim_thd: tbif 1  : ' )
531         CALL prt_ctl( tab2d_1=tbif(:,:,2), clinfo1=' lim_thd: tbif 2  : ' )
532         CALL prt_ctl( tab2d_1=tbif(:,:,3), clinfo1=' lim_thd: tbif 3  : ' )
533         CALL prt_ctl( tab2d_1=fdtcn      , clinfo1=' lim_thd: fdtcn   : ', tab2d_2=qdtcn , clinfo2=' qdtcn  : ' )
534         CALL prt_ctl( tab2d_1=qstoif     , clinfo1=' lim_thd: qstoif  : ', tab2d_2=fsbbq , clinfo2=' fsbbq  : ' )
535      ENDIF
536       !
537      CALL wrk_dealloc( jpi, jpj, ztmp, zqlbsbq, zlicegr, zdvosif, zdvobif, zdvolif, zdvonif, zdvomif, zu_imasstr, zv_imasstr )
538      CALL wrk_dealloc( jpi, jpj, jpk, zmsk )
539      !
540    END SUBROUTINE lim_thd_2
541
542
543    SUBROUTINE lim_thd_init_2
544      !!-------------------------------------------------------------------
545      !!                   ***  ROUTINE lim_thd_init_2 ***
546      !!                 
547      !! ** Purpose :   Physical constants and parameters linked to the ice
548      !!      thermodynamics
549      !!
550      !! ** Method  :   Read the namicethd namelist and check the ice-thermo
551      !!       parameter values called at the first timestep (nit000)
552      !!
553      !! ** input   :   Namelist namicether
554      !!-------------------------------------------------------------------
555      INTEGER  ::   ios                 ! Local integer output status for namelist read
556      NAMELIST/namicethd/ hmelt , hiccrit, hicmin, hiclim, amax  ,        &
557         &                swiqst, sbeta  , parlat, hakspl, hibspl, exld,  &
558         &                hakdif, hnzst  , thth  , parsub, alphs
559      !!-------------------------------------------------------------------
560                   
561      REWIND( numnam_ice_ref )              ! Namelist namicethd in reference namelist : Ice thermodynamics
562      READ  ( numnam_ice_ref, namicethd, IOSTAT = ios, ERR = 901)
563901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicethd in reference namelist', lwp )
564
565      REWIND( numnam_ice_cfg )              ! Namelist namicethd in configuration namelist : Ice thermodynamics
566      READ  ( numnam_ice_cfg, namicethd, IOSTAT = ios, ERR = 902 )
567902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicethd in configuration namelist', lwp )
568      IF(lwm) WRITE ( numoni, namicethd )
569
570      IF( lk_cpl .AND. parsub /= 0.0 )   CALL ctl_stop( 'In coupled mode, use parsub = 0. or send dqla' )
571      !
572      IF(lwp) THEN                          ! control print
573         WRITE(numout,*)
574         WRITE(numout,*)'lim_thd_init_2: ice parameters for ice thermodynamic computation '
575         WRITE(numout,*)'~~~~~~~~~~~~~~'
576         WRITE(numout,*)'       maximum melting at the bottom                           hmelt        = ', hmelt
577         WRITE(numout,*)'       ice thick. for lateral accretion in NH (SH)             hiccrit(1/2) = ', hiccrit
578         WRITE(numout,*)'       ice thick. corr. to max. energy stored in brine pocket  hicmin       = ', hicmin 
579         WRITE(numout,*)'       minimum ice thickness                                   hiclim       = ', hiclim 
580         WRITE(numout,*)'       maximum lead fraction                                   amax         = ', amax
581         WRITE(numout,*)'       energy stored in brine pocket (=1) or not (=0)          swiqst       = ', swiqst 
582         WRITE(numout,*)'       numerical carac. of the scheme for diffusion in ice '
583         WRITE(numout,*)'       Cranck-Nicholson (=0.5), implicit (=1), explicit (=0)   sbeta        = ', sbeta
584         WRITE(numout,*)'       percentage of energy used for lateral ablation          parlat       = ', parlat
585         WRITE(numout,*)'       slope of distr. for Hakkinen-Mellor lateral melting     hakspl       = ', hakspl 
586         WRITE(numout,*)'       slope of distribution for Hibler lateral melting        hibspl       = ', hibspl
587         WRITE(numout,*)'       exponent for leads-closure rate                         exld         = ', exld
588         WRITE(numout,*)'       coefficient for diffusions of ice and snow              hakdif       = ', hakdif
589         WRITE(numout,*)'       threshold thick. for comp. of eq. thermal conductivity  zhth         = ', thth 
590         WRITE(numout,*)'       thickness of the surf. layer in temp. computation       hnzst        = ', hnzst
591         WRITE(numout,*)'       switch for snow sublimation  (=1) or not (=0)           parsub       = ', parsub 
592         WRITE(numout,*)'       coefficient for snow density when snow ice formation    alphs        = ', alphs
593      ENDIF
594      !         
595      uscomi = 1.0 / ( 1.0 - amax )   ! inverse of minimum lead fraction
596      rcdsn = hakdif * rcdsn 
597      rcdic = hakdif * rcdic
598      !
599      IF( hsndif > 100.e0 .OR. hicdif > 100.e0 ) THEN
600         cnscg = 0.e0
601      ELSE
602         cnscg = rcpsn / rcpic   ! ratio  rcpsn/rcpic
603      ENDIF
604      !
605   END SUBROUTINE lim_thd_init_2
606
607#else
608   !!----------------------------------------------------------------------
609   !!   Default option          Dummy module       NO LIM 2.0 sea-ice model
610   !!----------------------------------------------------------------------
611CONTAINS
612   SUBROUTINE lim_thd_2         ! Dummy routine
613   END SUBROUTINE lim_thd_2
614#endif
615
616   !!======================================================================
617END MODULE limthd_2
Note: See TracBrowser for help on using the repository browser.