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

source: branches/2016/dev_CNRS_2016/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 7366

Last change on this file since 7366 was 7366, checked in by clem, 7 years ago

make the option Cd-Lupkes restartable

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