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 @ 7325

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

Unfinished work toward incorporation of melt ponds

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