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

source: branches/2016/dev_r6711_SIMPLIF_6_aerobulk/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 6723

Last change on this file since 6723 was 6723, checked in by gm, 8 years ago

#1751 - branch SIMPLIF_6_aerobulk: add aerobulk package including NCAR, COARE and ECMWF bulk

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