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.
sbcice_lim.F90 in branches/2016/dev_r6859_LIM3_meltponds/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2016/dev_r6859_LIM3_meltponds/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 8099

Last change on this file since 8099 was 8099, checked in by vancop, 7 years ago

impact of melt ponds on albedo, first run

  • Property svn:keywords set to Id
File size: 36.0 KB
Line 
1MODULE sbcice_lim
2   !!======================================================================
3   !!                       ***  MODULE  sbcice_lim  ***
4   !! Surface module :  update the ocean surface boundary condition over ice
5   !!       &           covered area using LIM sea-ice model
6   !! Sea-Ice model  :  LIM-3 Sea ice model time-stepping
7   !!=====================================================================
8   !! History :  2.0  ! 2006-12  (M. Vancoppenolle) Original code
9   !!            3.0  ! 2008-02  (C. Talandier)  Surface module from icestp.F90
10   !!             -   ! 2008-04  (G. Madec)  sltyle and lim_ctl routine
11   !!            3.3  ! 2010-11  (G. Madec) ice-ocean stress always computed at each ocean time-step
12   !!            3.4  ! 2011-01  (A Porter)  dynamical allocation
13   !!             -   ! 2012-10  (C. Rousset)  add lim_diahsb
14   !!            3.6  ! 2014-07  (M. Vancoppenolle, G. Madec, O. Marti) revise coupled interface
15   !!----------------------------------------------------------------------
16#if defined key_lim3
17   !!----------------------------------------------------------------------
18   !!   'key_lim3' :                                  LIM 3.0 sea-ice model
19   !!----------------------------------------------------------------------
20   !!   sbc_ice_lim  : sea-ice model time-stepping and update ocean sbc over ice-covered area
21   !!----------------------------------------------------------------------
22   USE oce             ! ocean dynamics and tracers
23   USE dom_oce         ! ocean space and time domain
24   USE ice             ! LIM-3: ice variables
25   USE thd_ice         ! LIM-3: thermodynamical variables
26
27   USE sbc_oce         ! Surface boundary condition: ocean fields
28   USE sbc_ice         ! Surface boundary condition: ice   fields
29   USE sbcblk_core     ! Surface boundary condition: CORE bulk
30   USE sbcblk_clio     ! Surface boundary condition: CLIO bulk
31   USE sbccpl          ! Surface boundary condition: coupled interface
32   USE sbcana          ! Surface boundary condition: analytic formulation
33   USE albedo          ! ocean & ice albedo
34
35   USE phycst          ! Define parameters for the routines
36   USE eosbn2          ! equation of state
37   USE limdyn          ! Ice dynamics
38   USE limtrp          ! Ice transport
39   USE limhdf          ! Ice horizontal diffusion
40   USE limthd          ! Ice thermodynamics
41   USE limitd_me       ! Mechanics on ice thickness distribution
42   USE limsbc          ! sea surface boundary condition
43   USE limdiahsb       ! Ice budget diagnostics
44   USE limwri          ! Ice outputs
45   USE limrst          ! Ice restarts
46   USE limupdate1      ! update of global variables
47   USE limupdate2      ! update of global variables
48   USE limvar          ! Ice variables switch
49
50   ! MV MP 2016
51   USE limmp
52   ! END MV MP 2016
53
54   USE limistate       ! LIM initial state
55   USE limthd_sal      ! LIM ice thermodynamics: salinity
56
57   USE c1d             ! 1D vertical configuration
58   USE lbclnk          ! lateral boundary condition - MPP link
59   USE lib_mpp         ! MPP library
60   USE wrk_nemo        ! work arrays
61   USE timing          ! Timing
62   USE iom             ! I/O manager library
63   USE in_out_manager  ! I/O manager
64   USE prtctl          ! Print control
65   USE lib_fortran     !
66   USE limctl
67
68#if defined key_bdy 
69   USE bdyice_lim       ! unstructured open boundary data  (bdy_ice_lim routine)
70#endif
71# if defined key_agrif
72   USE agrif_ice
73   USE agrif_lim3_update
74   USE agrif_lim3_interp
75# endif
76
77   IMPLICIT NONE
78   PRIVATE
79
80   PUBLIC sbc_ice_lim  ! routine called by sbcmod.F90
81   PUBLIC sbc_lim_init ! routine called by sbcmod.F90
82   
83   !! * Substitutions
84#  include "domzgr_substitute.h90"
85#  include "vectopt_loop_substitute.h90"
86   !!----------------------------------------------------------------------
87   !! NEMO/OPA 4.0 , UCL NEMO Consortium (2011)
88   !! $Id$
89   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
90   !!----------------------------------------------------------------------
91CONTAINS
92
93   !!======================================================================
94
95   SUBROUTINE sbc_ice_lim( kt, kblk )
96      !!---------------------------------------------------------------------
97      !!                  ***  ROUTINE sbc_ice_lim  ***
98      !!                   
99      !! ** Purpose :   update the ocean surface boundary condition via the
100      !!                Louvain la Neuve Sea Ice Model time stepping
101      !!
102      !! ** Method  :   ice model time stepping
103      !!              - call the ice dynamics routine
104      !!              - call the ice advection/diffusion routine
105      !!              - call the ice thermodynamics routine
106      !!              - call the routine that computes mass and
107      !!                heat fluxes at the ice/ocean interface
108      !!              - save the outputs
109      !!              - save the outputs for restart when necessary
110      !!
111      !! ** Action  : - time evolution of the LIM sea-ice model
112      !!              - update all sbc variables below sea-ice:
113      !!                utau, vtau, taum, wndm, qns , qsr, emp , sfx
114      !!---------------------------------------------------------------------
115      INTEGER, INTENT(in) ::   kt      ! ocean time step
116      INTEGER, INTENT(in) ::   kblk    ! type of bulk (=1 ANALYTIC, =3 CLIO, =4 CORE, =5 COUPLED)
117      !!
118      INTEGER  ::   jl                 ! dummy loop index
119      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zalb_os, zalb_cs  ! ice albedo under overcast/clear sky
120      REAL(wp), POINTER, DIMENSION(:,:  )   ::   zutau_ice, zvtau_ice 
121      !!----------------------------------------------------------------------
122
123      IF( nn_timing == 1 )  CALL timing_start('sbc_ice_lim')
124
125      ! clem: it is important to initialize agrif_lim3 variables here and not in sbc_lim_init
126# if defined key_agrif
127      IF( kt == nit000 ) THEN
128         IF( .NOT. Agrif_Root() )   CALL Agrif_InitValues_cont_lim3
129      ENDIF
130# endif
131
132      !-----------------------!
133      ! --- Ice time step --- !
134      !-----------------------!
135      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN
136
137# if defined key_agrif
138         IF( .NOT. Agrif_Root() )  lim_nbstep = MOD( lim_nbstep, Agrif_irhot() * Agrif_Parent(nn_fsbc) / nn_fsbc ) + 1
139# endif
140
141         ! mean surface ocean current at ice velocity point (C-grid dynamics :  U- & V-points as the ocean)
142         u_oce(:,:) = ssu_m(:,:) * umask(:,:,1)
143         v_oce(:,:) = ssv_m(:,:) * vmask(:,:,1)
144         
145         ! masked sea surface freezing temperature [Kelvin] (set to rt0 over land)
146         CALL eos_fzp( sss_m(:,:) , t_bo(:,:) )
147         t_bo(:,:) = ( t_bo(:,:) + rt0 ) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) )
148         
149         ! Mask sea ice surface temperature (set to rt0 over land)
150         DO jl = 1, jpl
151            t_su(:,:,jl) = t_su(:,:,jl) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) )
152         END DO     
153         !
154         !------------------------------------------------!                                           
155         ! --- Dynamical coupling with the atmosphere --- !                                           
156         !------------------------------------------------!
157         ! It provides the following fields:
158         ! utau_ice, vtau_ice : surface ice stress (U- & V-points)   [N/m2]
159         !-----------------------------------------------------------------
160         SELECT CASE( kblk )
161            CASE( jp_ana     )   ;    CALL ana_ice_tau                              ! analytic formulation           
162            CASE( jp_clio    )   ;    CALL blk_ice_clio_tau                         ! CLIO bulk formulation           
163            CASE( jp_core    )   ;    CALL blk_ice_core_tau                         ! CORE bulk formulation
164            CASE( jp_purecpl )   ;    CALL sbc_cpl_ice_tau( utau_ice , vtau_ice )   ! Coupled   formulation
165         END SELECT
166         
167         IF( ln_mixcpl) THEN                                                       ! Case of a mixed Bulk/Coupled formulation
168            CALL wrk_alloc( jpi,jpj    , zutau_ice, zvtau_ice)
169                                      CALL sbc_cpl_ice_tau( zutau_ice , zvtau_ice )
170            utau_ice(:,:) = utau_ice(:,:) * xcplmask(:,:,0) + zutau_ice(:,:) * ( 1. - xcplmask(:,:,0) )
171            vtau_ice(:,:) = vtau_ice(:,:) * xcplmask(:,:,0) + zvtau_ice(:,:) * ( 1. - xcplmask(:,:,0) )
172            CALL wrk_dealloc( jpi,jpj  , zutau_ice, zvtau_ice)
173         ENDIF
174
175         !-------------------------------------------------------!
176         ! --- ice dynamics and transport (except in 1D case) ---!
177         !-------------------------------------------------------!
178         numit = numit + nn_fsbc                  ! Ice model time step
179         !                                                   
180                                      CALL sbc_lim_bef         ! Store previous ice values
181                                      CALL sbc_lim_diag0       ! set diag of mass, heat and salt fluxes to 0
182                                      CALL lim_rst_opn( kt )   ! Open Ice restart file
183         !
184         ! --- zap this if no ice dynamics --- !
185         IF( .NOT. lk_c1d .AND. ln_limdyn ) THEN
186            !
187            IF( nn_limdyn /= 0 ) THEN                          ! -- Ice dynamics
188                                      CALL lim_dyn( kt )       !     rheology 
189            ELSE
190               u_ice(:,:) = rn_uice * umask(:,:,1)             !     or prescribed velocity
191               v_ice(:,:) = rn_vice * vmask(:,:,1)
192            ENDIF
193                                      CALL lim_trp( kt )       ! -- Ice transport (Advection/diffusion)
194            IF( nn_limdyn == 2 .AND. nn_monocat /= 2 )  &      ! -- Mechanical redistribution (ridging/rafting)
195               &                      CALL lim_itd_me         
196            IF( nn_limdyn == 2 )      CALL lim_update1( kt )   ! -- Corrections
197            !
198         ENDIF
199
200         ! ---
201#if defined key_agrif
202         IF( .NOT. Agrif_Root() )     CALL agrif_interp_lim3('T')
203#endif
204#if defined key_bdy
205         IF( ln_limthd )              CALL bdy_ice_lim( kt )   ! -- bdy ice thermo
206#endif
207
208         ! previous lead fraction and ice volume for flux calculations
209                                      CALL sbc_lim_bef                       
210                                      CALL lim_var_glo2eqv     ! ht_i and ht_s for ice albedo calculation
211                                      CALL lim_var_agg(1)      ! at_i for coupling (via pfrld)
212         !
213         pfrld(:,:)   = 1._wp - at_i(:,:)
214         phicif(:,:)  = vt_i(:,:)
215         
216         !------------------------------------------------------!                                           
217         ! --- Thermodynamical coupling with the atmosphere --- !                                           
218         !------------------------------------------------------!
219         ! It provides the following fields:
220         ! qsr_ice , qns_ice  : solar & non solar heat flux over ice   (T-point)         [W/m2]
221         ! qla_ice            : latent heat flux over ice              (T-point)         [W/m2]
222         ! dqns_ice, dqla_ice : non solar & latent heat sensistivity   (T-point)         [W/m2]
223         ! tprecip , sprecip  : total & solid precipitation            (T-point)         [Kg/m2/s]
224         ! fr1_i0  , fr2_i0   : 1sr & 2nd fraction of qsr penetration in ice             [%]
225         !----------------------------------------------------------------------------------------
226         CALL wrk_alloc( jpi,jpj,jpl, zalb_os, zalb_cs )
227         
228                                      !CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os ) ! cloud-sky and overcast-sky ice albedos
229                                      CALL albedo_ice( t_su, ht_i, a_i, ht_s, a_ip, v_ip, zalb_cs, zalb_os, ln_pnd_rad ) ! cloud-sky and overcast-sky ice albedos
230
231         SELECT CASE( kblk )
232
233            CASE( jp_ana )                                        ! analytic formulation
234                                      CALL ana_ice_flx
235               
236            CASE( jp_clio )                                       ! CLIO bulk formulation
237               ! In CLIO the cloud fraction is read in the climatology and the all-sky albedo
238               ! (alb_ice) is computed within the bulk routine
239                                      CALL blk_ice_clio_flx( t_su, zalb_cs, zalb_os, alb_ice )
240               IF( ln_mixcpl      )   CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su )
241               IF( nn_limflx /= 2 )   CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx )
242               
243            CASE( jp_core )                                       ! CORE bulk formulation
244               ! albedo depends on cloud fraction because of non-linear spectral effects
245               alb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:)
246                                      CALL blk_ice_core_flx( t_su, alb_ice )
247               IF( ln_mixcpl      )   CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su )
248               IF( nn_limflx /= 2 )   CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx )
249               
250            CASE ( jp_purecpl )                                    ! Coupled formulation
251               ! albedo depends on cloud fraction because of non-linear spectral effects
252               alb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:)
253                                      CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su )
254               IF( nn_limflx == 2 )   CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx )
255
256         END SELECT
257
258         CALL wrk_dealloc( jpi,jpj,jpl, zalb_os, zalb_cs )
259
260         !----------------------------!
261         ! --- ice thermodynamics --- !
262         !----------------------------!
263         ! --- zap this if no ice thermo --- !
264         IF( ln_limthd )              CALL lim_thd( kt )        ! -- Ice thermodynamics     
265
266         ! MV MP 2016
267         IF ( ln_pnd )                CALL lim_mp( kt )         ! -- Melt ponds
268         ! END MV MP 2016
269
270         IF( ln_limthd )              CALL lim_update2( kt )    ! -- Corrections
271         ! ---
272# if defined key_agrif
273         IF( .NOT. Agrif_Root() )     CALL agrif_update_lim3( kt )
274# endif
275                                      CALL lim_var_glo2eqv      ! necessary calls (at least for coupling)
276                                      CALL lim_var_agg( 2 )     ! necessary calls (at least for coupling)
277                                      !
278# if defined key_agrif
279!!         IF( .NOT. Agrif_Root() )     CALL Agrif_ChildGrid_To_ParentGrid()  ! clem: should be called at the update frequency only (cf agrif_lim3_update)
280# endif
281                                      CALL lim_sbc_flx( kt )    ! -- Update surface ocean mass, heat and salt fluxes
282# if defined key_agrif
283!!         IF( .NOT. Agrif_Root() )     CALL Agrif_ParentGrid_To_ChildGrid()  ! clem: should be called at the update frequency only (cf agrif_lim3_update)
284# endif
285         IF( ln_limdiahsb )           CALL lim_diahsb( kt )     ! -- Diagnostics and outputs
286         !
287                                      CALL lim_wri( 1 )         ! -- Ice outputs
288         !
289         IF( kt == nit000 .AND. ln_rstart )   &
290            &                         CALL iom_close( numrir )  ! close input ice restart file
291         !
292         IF( lrst_ice )               CALL lim_rst_write( kt )  ! -- Ice restart file
293         !
294         IF( ln_limctl )              CALL lim_ctl( kt )        ! alerts in case of model crash
295         !
296      ENDIF   ! End sea-ice time step only
297
298      !-------------------------!
299      ! --- Ocean time step --- !
300      !-------------------------!
301      ! Update surface ocean stresses (only in ice-dynamic case) otherwise the atm.-ocean stresses are used everywhere
302      !    using before instantaneous surf. currents
303      IF( ln_limdyn )                 CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) )
304!!gm   remark, the ocean-ice stress is not saved in ice diag call above .....  find a solution!!!
305      !
306      IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim')
307      !
308   END SUBROUTINE sbc_ice_lim
309   
310
311   SUBROUTINE sbc_lim_init
312      !!----------------------------------------------------------------------
313      !!                  ***  ROUTINE sbc_lim_init  ***
314      !!
315      !! ** purpose :   Allocate all the dynamic arrays of the LIM-3 modules
316      !!----------------------------------------------------------------------
317      INTEGER :: ierr
318      INTEGER :: ji, jj
319      !!----------------------------------------------------------------------
320      IF(lwp) WRITE(numout,*)
321      IF(lwp) WRITE(numout,*) 'sbc_lim_init : update ocean surface boudary condition' 
322      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   via Louvain la Neuve Ice Model (LIM-3) time stepping'
323      !
324                                       ! Open the reference and configuration namelist files and namelist output file
325      CALL ctl_opn( numnam_ice_ref, 'namelist_ice_ref',    'OLD',     'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) 
326      CALL ctl_opn( numnam_ice_cfg, 'namelist_ice_cfg',    'OLD',     'FORMATTED', 'SEQUENTIAL', -1, numout, lwp )
327      IF(lwm) CALL ctl_opn( numoni, 'output.namelist.ice', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, 1 )
328
329      CALL ice_run                     ! set some ice run parameters
330      !
331      !                                ! Allocate the ice arrays
332      ierr =        ice_alloc        ()      ! ice variables
333      ierr = ierr + sbc_ice_alloc    ()      ! surface forcing
334      ierr = ierr + thd_ice_alloc    ()      ! thermodynamics
335      IF( ln_limdyn )   ierr = ierr + lim_itd_me_alloc ()      ! ice thickness distribution - mechanics
336      !
337      IF( lk_mpp    )   CALL mpp_sum( ierr )
338      IF( ierr /= 0 )   CALL ctl_stop('STOP', 'sbc_lim_init : unable to allocate ice arrays')
339      !
340      CALL lim_dyn_init                ! set ice dynamics parameters
341      !
342      CALL lim_itd_init                ! ice thickness distribution initialization
343      !
344      CALL lim_hdf_init                ! set ice horizontal diffusion computation parameters
345      !
346      CALL lim_thd_init                ! set ice thermodynics parameters
347      !
348      CALL lim_thd_sal_init            ! set ice salinity parameters
349       
350      ! MV MP 2016
351      CALL lim_mp_init                 ! set melt ponds parameters
352      ! END MV MP 2016
353
354      IF( ln_limdyn )   CALL lim_itd_me_init             ! ice thickness distribution initialization for mecanical deformation
355      !                                ! Initial sea-ice state
356      IF( .NOT. ln_rstart ) THEN              ! start from rest: sea-ice deduced from sst
357         numit = 0
358         numit = nit000 - 1
359         CALL lim_istate
360      ELSE                                    ! start from a restart file
361         CALL lim_rst_read
362         numit = nit000 - 1
363      ENDIF
364      CALL lim_var_agg(2)
365      CALL lim_var_glo2eqv
366      !
367      CALL lim_sbc_init                 ! ice surface boundary condition   
368      !
369      IF( ln_limdiahsb) CALL lim_diahsb_init  ! initialization for diags
370      !
371      fr_i(:,:)     = at_i(:,:)         ! initialisation of sea-ice fraction
372      tn_ice(:,:,:) = t_su(:,:,:)       ! initialisation of surface temp for coupled simu
373      !
374      DO jj = 1, jpj
375         DO ji = 1, jpi
376            IF( gphit(ji,jj) > 0._wp ) THEN  ;  rn_amax_2d(ji,jj) = rn_amax_n  ! NH
377            ELSE                             ;  rn_amax_2d(ji,jj) = rn_amax_s  ! SH
378            ENDIF
379        ENDDO
380      ENDDO 
381      !
382      nstart = numit  + nn_fsbc     
383      nitrun = nitend - nit000 + 1 
384      nlast  = numit  + nitrun 
385      !
386      IF( nstock == 0 )   nstock = nlast + 1
387      !
388      !
389   END SUBROUTINE sbc_lim_init
390
391
392   SUBROUTINE ice_run
393      !!-------------------------------------------------------------------
394      !!                  ***  ROUTINE ice_run ***
395      !!                 
396      !! ** Purpose :   Definition some run parameter for ice model
397      !!
398      !! ** Method  :   Read the namicerun namelist and check the parameter
399      !!              values called at the first timestep (nit000)
400      !!
401      !! ** input   :   Namelist namicerun
402      !!-------------------------------------------------------------------
403      INTEGER  ::   ios                 ! Local integer output status for namelist read
404      NAMELIST/namicerun/ jpl, nlay_i, nlay_s, rn_amax_n, rn_amax_s, cn_icerst_in, cn_icerst_indir,   &
405         &                cn_icerst_out, cn_icerst_outdir, ln_limthd, ln_limdyn, nn_limdyn, rn_uice, rn_vice 
406      NAMELIST/namicediag/ ln_limdiachk, ln_limdiahsb, ln_limctl, iiceprt, jiceprt 
407      !!-------------------------------------------------------------------
408      !                   
409      REWIND( numnam_ice_ref )              ! Namelist namicerun in reference namelist : Parameters for ice
410      READ  ( numnam_ice_ref, namicerun, IOSTAT = ios, ERR = 901)
411901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in reference namelist', lwp )
412
413      REWIND( numnam_ice_cfg )              ! Namelist namicerun in configuration namelist : Parameters for ice
414      READ  ( numnam_ice_cfg, namicerun, IOSTAT = ios, ERR = 902 )
415902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in configuration namelist', lwp )
416      IF(lwm) WRITE ( numoni, namicerun )
417      !
418      REWIND( numnam_ice_ref )              ! Namelist namicediag in reference namelist : Parameters for ice
419      READ  ( numnam_ice_ref, namicediag, IOSTAT = ios, ERR = 903)
420903   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicediag in reference namelist', lwp )
421
422      REWIND( numnam_ice_cfg )              ! Namelist namicediag in configuration namelist : Parameters for ice
423      READ  ( numnam_ice_cfg, namicediag, IOSTAT = ios, ERR = 904 )
424904   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicediag in configuration namelist', lwp )
425      IF(lwm) WRITE ( numoni, namicediag )
426      !
427      IF(lwp) THEN                        ! control print
428         WRITE(numout,*)
429         WRITE(numout,*) 'ice_run : ice share~d parameters for dynamics/advection/thermo of sea-ice'
430         WRITE(numout,*) ' ~~~~~~'
431         WRITE(numout,*) '   number of ice  categories                               = ', jpl
432         WRITE(numout,*) '   number of ice  layers                                   = ', nlay_i
433         WRITE(numout,*) '   number of snow layers                                   = ', nlay_s
434         WRITE(numout,*) '   maximum ice concentration for NH                        = ', rn_amax_n 
435         WRITE(numout,*) '   maximum ice concentration for SH                        = ', rn_amax_s
436         WRITE(numout,*) '   Ice thermodynamics (T) or not (F)            ln_limthd  = ', ln_limthd
437         WRITE(numout,*) '   Ice dynamics       (T) or not (F)            ln_limdyn  = ', ln_limdyn
438         WRITE(numout,*) '     (ln_limdyn=T) Ice dynamics switch          nn_limdyn  = ', nn_limdyn
439         WRITE(numout,*) '       2: total'
440         WRITE(numout,*) '       1: advection only (no diffusion, no ridging/rafting)'
441         WRITE(numout,*) '       0: advection only (as 1 + prescribed velocity, bypass rheology)'
442         WRITE(numout,*) '     (ln_limdyn=T) prescribed u-vel (case nn_limdyn=0)     = ', rn_uice
443         WRITE(numout,*) '     (ln_limdyn=T) prescribed v-vel (case nn_limdyn=0)     = ', rn_vice
444         WRITE(numout,*)
445         WRITE(numout,*) '...and ice diagnostics'
446         WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~'
447         WRITE(numout,*) '   Diagnose online heat/mass/salt budget     ln_limdiachk  = ', ln_limdiachk
448         WRITE(numout,*) '   Output          heat/mass/salt budget     ln_limdiahsb  = ', ln_limdiahsb
449         WRITE(numout,*) '   control prints in ocean.out for (i,j)=(iiceprt,jiceprt) = ', ln_limctl
450         WRITE(numout,*) '   i-index for control prints (ln_limctl=true)             = ', iiceprt
451         WRITE(numout,*) '   j-index for control prints (ln_limctl=true)             = ', jiceprt
452      ENDIF
453      !
454      ! sea-ice timestep and inverse
455      rdt_ice   = REAL(nn_fsbc) * rdt 
456      r1_rdtice = 1._wp / rdt_ice 
457
458      ! inverse of nlay_i and nlay_s
459      r1_nlay_i = 1._wp / REAL( nlay_i, wp )
460      r1_nlay_s = 1._wp / REAL( nlay_s, wp )
461      !
462#if defined key_bdy
463      IF( lwp .AND. ln_limdiachk )  CALL ctl_warn('online conservation check activated but it does not work with BDY')
464#endif
465      !
466      IF( lwp ) WRITE(numout,*) '   ice timestep rdt_ice  = ', rdt_ice
467      !
468   END SUBROUTINE ice_run
469
470
471   SUBROUTINE lim_itd_init
472      !!------------------------------------------------------------------
473      !!                ***  ROUTINE lim_itd_init ***
474      !!
475      !! ** Purpose :   Initializes the ice thickness distribution
476      !! ** Method  :   ...
477      !! ** input   :   Namelist namiceitd
478      !!-------------------------------------------------------------------
479      INTEGER  ::   ios                 ! Local integer output status for namelist read
480      NAMELIST/namiceitd/ nn_catbnd, rn_himean
481      !
482      INTEGER  ::   jl                   ! dummy loop index
483      REAL(wp) ::   zc1, zc2, zc3, zx1   ! local scalars
484      REAL(wp) ::   zhmax, znum, zden, zalpha !
485      !!------------------------------------------------------------------
486      !
487      REWIND( numnam_ice_ref )              ! Namelist namiceitd in reference namelist : Parameters for ice
488      READ  ( numnam_ice_ref, namiceitd, IOSTAT = ios, ERR = 905)
489905   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in reference namelist', lwp )
490
491      REWIND( numnam_ice_cfg )              ! Namelist namiceitd in configuration namelist : Parameters for ice
492      READ  ( numnam_ice_cfg, namiceitd, IOSTAT = ios, ERR = 906 )
493906   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in configuration namelist', lwp )
494      IF(lwm) WRITE ( numoni, namiceitd )
495      !
496      IF(lwp) THEN                        ! control print
497         WRITE(numout,*)
498         WRITE(numout,*) 'lim_itd_init : Initialization of ice cat distribution '
499         WRITE(numout,*) '~~~~~~~~~~~~'
500         WRITE(numout,*) '   shape of ice categories distribution                          nn_catbnd = ', nn_catbnd
501         WRITE(numout,*) '   mean ice thickness in the domain (only active if nn_catbnd=2) rn_himean = ', rn_himean
502      ENDIF
503
504      !----------------------------------
505      !- Thickness categories boundaries
506      !----------------------------------
507      hi_max(:) = 0._wp
508
509      SELECT CASE ( nn_catbnd  )       
510                                   !----------------------
511         CASE (1)                  ! tanh function (CICE)
512                                   !----------------------
513         zc1 =  3._wp / REAL( jpl, wp )
514         zc2 = 10._wp * zc1
515         zc3 =  3._wp
516
517         DO jl = 1, jpl
518            zx1 = REAL( jl-1, wp ) / REAL( jpl, wp )
519            hi_max(jl) = hi_max(jl-1) + zc1 + zc2 * (1._wp + TANH( zc3 * (zx1 - 1._wp ) ) )
520         END DO
521
522                                   !----------------------
523         CASE (2)                  ! h^(-alpha) function
524                                   !----------------------
525         zalpha = 0.05             ! exponent of the transform function
526
527         zhmax  = 3.*rn_himean
528
529         DO jl = 1, jpl 
530            znum = jpl * ( zhmax+1 )**zalpha
531            zden = ( jpl - jl ) * ( zhmax+1 )**zalpha + jl
532            hi_max(jl) = ( znum / zden )**(1./zalpha) - 1
533         END DO
534
535      END SELECT
536
537      DO jl = 1, jpl
538         hi_mean(jl) = ( hi_max(jl) + hi_max(jl-1) ) * 0.5_wp
539      END DO
540
541      ! Set hi_max(jpl) to a big value to ensure that all ice is thinner than hi_max(jpl)
542      hi_max(jpl) = 99._wp
543
544      IF(lwp) WRITE(numout,*) ' Thickness category boundaries '
545      IF(lwp) WRITE(numout,*) ' hi_max ', hi_max(0:jpl)
546      !
547   END SUBROUTINE lim_itd_init
548
549   
550   SUBROUTINE ice_lim_flx( ptn_ice, palb_ice, pqns_ice, pqsr_ice, pdqn_ice, pevap_ice, pdevap_ice, k_limflx )
551      !!---------------------------------------------------------------------
552      !!                  ***  ROUTINE ice_lim_flx  ***
553      !!                   
554      !! ** Purpose :   update the ice surface boundary condition by averaging and / or
555      !!                redistributing fluxes on ice categories                   
556      !!
557      !! ** Method  :   average then redistribute
558      !!
559      !! ** Action  :   
560      !!---------------------------------------------------------------------
561      INTEGER                   , INTENT(in   ) ::   k_limflx   ! =-1 do nothing; =0 average ;
562                                                                ! =1 average and redistribute ; =2 redistribute
563      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   ptn_ice    ! ice surface temperature
564      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   palb_ice   ! ice albedo
565      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqns_ice   ! non solar flux
566      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqsr_ice   ! net solar flux
567      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdqn_ice   ! non solar flux sensitivity
568      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pevap_ice  ! sublimation
569      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdevap_ice ! sublimation sensitivity
570      !
571      INTEGER  ::   jl      ! dummy loop index
572      !
573      REAL(wp), POINTER, DIMENSION(:,:) :: zalb_m    ! Mean albedo over all categories
574      REAL(wp), POINTER, DIMENSION(:,:) :: ztem_m    ! Mean temperature over all categories
575      !
576      REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_m   ! Mean solar heat flux over all categories
577      REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_m   ! Mean non solar heat flux over all categories
578      REAL(wp), POINTER, DIMENSION(:,:) :: z_evap_m  ! Mean sublimation over all categories
579      REAL(wp), POINTER, DIMENSION(:,:) :: z_dqn_m   ! Mean d(qns)/dT over all categories
580      REAL(wp), POINTER, DIMENSION(:,:) :: z_devap_m ! Mean d(evap)/dT over all categories
581      !!----------------------------------------------------------------------
582
583      IF( nn_timing == 1 )  CALL timing_start('ice_lim_flx')
584      !
585      !
586      SELECT CASE( k_limflx )                              !==  averaged on all ice categories  ==!
587      CASE( 0 , 1 )
588         CALL wrk_alloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m)
589         !
590         z_qns_m  (:,:) = fice_ice_ave ( pqns_ice (:,:,:) )
591         z_qsr_m  (:,:) = fice_ice_ave ( pqsr_ice (:,:,:) )
592         z_dqn_m  (:,:) = fice_ice_ave ( pdqn_ice (:,:,:) )
593         z_evap_m (:,:) = fice_ice_ave ( pevap_ice (:,:,:) )
594         z_devap_m(:,:) = fice_ice_ave ( pdevap_ice (:,:,:) )
595         DO jl = 1, jpl
596            pdqn_ice  (:,:,jl) = z_dqn_m(:,:)
597            pdevap_ice(:,:,jl) = z_devap_m(:,:)
598         END DO
599         !
600         DO jl = 1, jpl
601            pqns_ice (:,:,jl) = z_qns_m(:,:)
602            pqsr_ice (:,:,jl) = z_qsr_m(:,:)
603            pevap_ice(:,:,jl) = z_evap_m(:,:)
604         END DO
605         !
606         CALL wrk_dealloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m)
607      END SELECT
608
609      SELECT CASE( k_limflx )                              !==  redistribution on all ice categories  ==!
610      CASE( 1 , 2 )
611         CALL wrk_alloc( jpi,jpj, zalb_m, ztem_m )
612         !
613         zalb_m(:,:) = fice_ice_ave ( palb_ice (:,:,:) ) 
614         ztem_m(:,:) = fice_ice_ave ( ptn_ice  (:,:,:) ) 
615         DO jl = 1, jpl
616            pqns_ice (:,:,jl) = pqns_ice (:,:,jl) + pdqn_ice  (:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) )
617            pevap_ice(:,:,jl) = pevap_ice(:,:,jl) + pdevap_ice(:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) )
618            pqsr_ice (:,:,jl) = pqsr_ice (:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) ) 
619         END DO
620         !
621         CALL wrk_dealloc( jpi,jpj, zalb_m, ztem_m )
622      END SELECT
623      !
624      IF( nn_timing == 1 )  CALL timing_stop('ice_lim_flx')
625      !
626   END SUBROUTINE ice_lim_flx
627
628   SUBROUTINE sbc_lim_bef
629      !!----------------------------------------------------------------------
630      !!                  ***  ROUTINE sbc_lim_bef  ***
631      !!
632      !! ** purpose :  store ice variables at "before" time step
633      !!----------------------------------------------------------------------
634      a_i_b  (:,:,:)   = a_i  (:,:,:)     ! ice area
635      e_i_b  (:,:,:,:) = e_i  (:,:,:,:)   ! ice thermal energy
636      v_i_b  (:,:,:)   = v_i  (:,:,:)     ! ice volume
637      v_s_b  (:,:,:)   = v_s  (:,:,:)     ! snow volume
638      e_s_b  (:,:,:,:) = e_s  (:,:,:,:)   ! snow thermal energy
639      smv_i_b(:,:,:)   = smv_i(:,:,:)     ! salt content
640      oa_i_b (:,:,:)   = oa_i (:,:,:)     ! areal age content
641      u_ice_b(:,:)     = u_ice(:,:)
642      v_ice_b(:,:)     = v_ice(:,:)
643      at_i_b (:,:)     = SUM( a_i_b(:,:,:), dim=3 )
644     
645   END SUBROUTINE sbc_lim_bef
646
647   SUBROUTINE sbc_lim_diag0
648      !!----------------------------------------------------------------------
649      !!                  ***  ROUTINE sbc_lim_diag0  ***
650      !!
651      !! ** purpose :  set ice-ocean and ice-atm. fluxes to zeros at the beggining
652      !!               of the time step
653      !!----------------------------------------------------------------------
654      sfx    (:,:) = 0._wp   ;
655      sfx_bri(:,:) = 0._wp   ;   sfx_lam(:,:) = 0._wp
656      sfx_sni(:,:) = 0._wp   ;   sfx_opw(:,:) = 0._wp
657      sfx_bog(:,:) = 0._wp   ;   sfx_dyn(:,:) = 0._wp
658      sfx_bom(:,:) = 0._wp   ;   sfx_sum(:,:) = 0._wp
659      sfx_res(:,:) = 0._wp   ;   sfx_sub(:,:) = 0._wp
660     
661      wfx_snw(:,:) = 0._wp   ;   wfx_ice(:,:) = 0._wp
662      wfx_sni(:,:) = 0._wp   ;   wfx_opw(:,:) = 0._wp
663      wfx_bog(:,:) = 0._wp   ;   wfx_dyn(:,:) = 0._wp
664      wfx_bom(:,:) = 0._wp   ;   wfx_sum(:,:) = 0._wp
665      wfx_res(:,:) = 0._wp   ;   wfx_sub(:,:) = 0._wp
666      wfx_spr(:,:) = 0._wp   ;   wfx_lam(:,:) = 0._wp 
667
668      ! MV MP 2016
669      wfx_pnd(:,:) = 0._wp
670      ! END MV MP 2016
671     
672      hfx_thd(:,:) = 0._wp   ;   
673      hfx_snw(:,:) = 0._wp   ;   hfx_opw(:,:) = 0._wp
674      hfx_bog(:,:) = 0._wp   ;   hfx_dyn(:,:) = 0._wp
675      hfx_bom(:,:) = 0._wp   ;   hfx_sum(:,:) = 0._wp
676      hfx_res(:,:) = 0._wp   ;   hfx_sub(:,:) = 0._wp
677      hfx_spr(:,:) = 0._wp   ;   hfx_dif(:,:) = 0._wp 
678      hfx_err(:,:) = 0._wp   ;   hfx_err_rem(:,:) = 0._wp
679      hfx_err_dif(:,:) = 0._wp
680      wfx_err_sub(:,:) = 0._wp
681     
682      afx_tot(:,:) = 0._wp   ;
683      afx_dyn(:,:) = 0._wp   ;   afx_thd(:,:) = 0._wp
684
685      diag_heat(:,:) = 0._wp ;   diag_smvi(:,:) = 0._wp ;
686      diag_vice(:,:) = 0._wp ;   diag_vsnw(:,:) = 0._wp ;
687
688      tau_icebfr(:,:) = 0._wp; ! landfast ice param only (clem: important to keep the init here)
689     
690   END SUBROUTINE sbc_lim_diag0
691
692     
693   FUNCTION fice_cell_ave ( ptab )
694      !!--------------------------------------------------------------------------
695      !! * Compute average over categories, for grid cell (ice covered and free ocean)
696      !!--------------------------------------------------------------------------
697      REAL (wp), DIMENSION (jpi,jpj) :: fice_cell_ave
698      REAL (wp), DIMENSION (jpi,jpj,jpl), INTENT (in) :: ptab
699      INTEGER :: jl ! Dummy loop index
700     
701      fice_cell_ave (:,:) = 0.0_wp
702      DO jl = 1, jpl
703         fice_cell_ave (:,:) = fice_cell_ave (:,:) + a_i (:,:,jl) * ptab (:,:,jl)
704      END DO
705     
706   END FUNCTION fice_cell_ave
707   
708   
709   FUNCTION fice_ice_ave ( ptab )
710      !!--------------------------------------------------------------------------
711      !! * Compute average over categories, for ice covered part of grid cell
712      !!--------------------------------------------------------------------------
713      REAL (kind=wp), DIMENSION (jpi,jpj) :: fice_ice_ave
714      REAL (kind=wp), DIMENSION (jpi,jpj,jpl), INTENT(in) :: ptab
715
716      fice_ice_ave (:,:) = 0.0_wp
717      WHERE ( at_i (:,:) > 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:)
718
719   END FUNCTION fice_ice_ave
720
721
722#else
723   !!----------------------------------------------------------------------
724   !!   Default option           Dummy module      NO LIM 3.0 sea-ice model
725   !!----------------------------------------------------------------------
726CONTAINS
727   SUBROUTINE sbc_ice_lim ( kt, kblk )     ! Dummy routine
728      WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk
729   END SUBROUTINE sbc_ice_lim
730   SUBROUTINE sbc_lim_init                 ! Dummy routine
731   END SUBROUTINE sbc_lim_init
732#endif
733
734   !!======================================================================
735END MODULE sbcice_lim
Note: See TracBrowser for help on using the repository browser.