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

source: trunk/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 7698

Last change on this file since 7698 was 7698, checked in by mocavero, 7 years ago

update trunk with OpenMP parallelization

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