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/2015/dev_r5218_CNRS17_coupling/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2015/dev_r5218_CNRS17_coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 5220

Last change on this file since 5220 was 5220, checked in by smasson, 9 years ago

dev_r5218_CNRS17_coupling: first update

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