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

source: trunk/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 888

Last change on this file since 888 was 888, checked in by ctlod, 16 years ago

merge dev_001_SBC branche with the trunk to include the New Surface Module package, see ticket: #113

File size: 37.6 KB
RevLine 
[888]1MODULE sbcice_lim
2   !!======================================================================
3   !!                       ***  MODULE  sbcice_lim  ***
4   !! Surface module :  update surface ocean boundary condition over ice
5   !!                   covered area using LIM sea-ice model
6   !! Sea-Ice model  :  LIM 3.0 Sea ice model time-stepping
7   !!======================================================================
8   !! History :  9.0   !  06-12  (M. Vancoppenolle) Original code
9   !!            9.0   !  06-06  (G. Madec)  Surface module from icestp.F90
10   !!----------------------------------------------------------------------
11#if defined key_lim3
12   !!----------------------------------------------------------------------
13   !!   'key_lim3' :                                  LIM 3.0 sea-ice model
14   !!----------------------------------------------------------------------
15   !!----------------------------------------------------------------------
16   !!   sbc_ice_lim  : sea-ice model time-stepping and
17   !!                  update ocean sbc over ice-covered area
18   !!----------------------------------------------------------------------
19   USE oce             ! ocean dynamics and tracers
20   USE dom_oce         ! ocean space and time domain
21   USE par_ice         ! sea-ice parameters
22   USE ice
23   USE iceini
24   USE ice_oce         ! ice variables
25   USE dom_ice
26   USE cpl_oce
27
28   USE sbc_oce         ! Surface boundary condition: ocean fields
29   USE sbc_ice         ! Surface boundary condition: ice   fields
30   USE sbcblk_core     ! Surface boundary condition: CORE bulk
31   USE sbcblk_clio     ! Surface boundary condition: CLIO bulk
32   USE albedo
33
34   USE daymod          ! day calendar
35   USE phycst          ! Define parameters for the routines
36   USE eosbn2          ! equation of state
37   USE limdyn          ! Ice dynamics
38   USE limtrp          ! Ice transport
39   USE limthd          ! Ice thermodynamics
40   USE limitd_th       ! Thermodynamics on ice thickness distribution
41   USE limitd_me       ! Mechanics on ice thickness distribution
42   USE limsbc          ! sea surface boundary condition
43   USE limdia          ! Ice diagnostics
44   USE limwri          ! Ice outputs
45   USE limrst          ! Ice restarts
46   USE limupdate       ! update of global variables
47   USE limvar          ! Ice variables switch
48
49   USE lbclnk
50   USE iom             ! I/O manager library
51   USE in_out_manager  ! I/O manager
52   USE prtctl          ! Print control
53
54   IMPLICIT NONE
55   PRIVATE
56
57   PUBLIC sbc_ice_lim  ! routine called by sbcmod.F90
58   
59   CHARACTER(len=1) ::   cl_grid = 'C'     ! type of grid used in ice dynamics
60
61   !! * Substitutions
62#  include "domzgr_substitute.h90"
63#  include "vectopt_loop_substitute.h90"
64   !!----------------------------------------------------------------------
65   !!   OPA 9.0 , LOCEAN-IPSL (2006)
66   !! $ Id: $
67   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
68   !!----------------------------------------------------------------------
69
70CONTAINS
71
72   SUBROUTINE sbc_ice_lim( kt, kblk )
73      !!---------------------------------------------------------------------
74      !!                  ***  ROUTINE sbc_ice_lim  ***
75      !!                   
76      !! ** Purpose :   update the ocean surface boundary condition via the
77      !!                Louvain la Neuve Sea Ice Model time stepping
78      !!
79      !! ** Method  :   ice model time stepping
80      !!              - call the ice dynamics routine
81      !!              - call the ice advection/diffusion routine
82      !!              - call the ice thermodynamics routine
83      !!              - call the routine that computes mass and
84      !!                heat fluxes at the ice/ocean interface
85      !!              - save the outputs
86      !!              - save the outputs for restart when necessary
87      !!
88      !! ** Action  : - time evolution of the LIM sea-ice model
89      !!              - update all sbc variables below sea-ice:
90      !!                utau, vtau, qns , qsr, emp , emps
91      !!---------------------------------------------------------------------
92      INTEGER, INTENT(in) ::   kt      ! ocean time step
93      INTEGER, INTENT(in) ::   kblk    ! type of bulk (=3 CLIO, =4 CORE)
94      !!
95      INTEGER  ::   ji, jj, jk, jl     ! dummy loop indices
96      INTEGER  ::   indx               ! indexes for ice points
97      INTEGER  ::   numaltests         ! number of alert tests (max 20)
98      INTEGER  ::   alert_id           ! number of the current alert
99      REAL(wp) ::   ztmelts            ! ice layer melting point
100      INTEGER , DIMENSION(20) ::  numal                     ! number of alerts positive
101      CHARACTER (len=30), DIMENSION(20) ::   alname      ! name of alert
102      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   alb_ice_os   ! albedo of the ice under overcast sky
103      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   alb_ice_cs   ! albedo of ice under clear sky
104      !!----------------------------------------------------------------------
105
106      IF( kt == nit000 ) THEN
107         IF(lwp) WRITE(numout,*)
108         IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition' 
109         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   via Louvain la Neuve Ice Model (LIM) time stepping'
110
111         CALL ice_init
112
113         !+++++
114         indx = 12
115         jiindx = 44
116         jjindx = 140
117         WRITE(numout,*) ' The debugging point is : jiindx : ',jiindx,  &
118                                                  ' jjindx : ',jjindx
119         !+++++
120
121      ENDIF
122
123      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN
124         !
125         ! ... mean surface ocean current at ice dynamics point
126         !     C-grid dynamics :  U- & V-points as the ocean
127         DO jj = 2, jpj
128            DO ji = fs_2, jpi
129               u_oce(ji,jj) = ssu_m(ji,jj) * tmu(ji,jj)
130               v_oce(ji,jj) = ssv_m(ji,jj) * tmv(ji,jj)
131            END DO
132         END DO
133         CALL lbc_lnk( u_oce, 'U', -1. )   ! U-point
134         CALL lbc_lnk( v_oce, 'V', -1. )   ! V-point
135
136         ! ... masked sea surface freezing temperature [Kelvin] (set to rt0 over land)
137!old fashion         zsss_m (:,:) = SQRT( sss_m(:,:) )
138!old fashion         zsss2_m(:,:) = sss_m(:,:) * sss_m(:,:)
139!old fashion         zsss3_m(:,:) = zsss_m(:,:) * zsss_m(:,:) * zsss_m(:,:)
140!old fashion
141!old fashion         DO jj = 1, jpj
142!old fashion            DO ji = 1, jpi
143!old fashion               t_bo(ji,jj)  = ABS ( rt0 - 0.0575      * sss_m(ji,jj)      &
144!old fashion                  &                    + 1.710523e-03 * zsss3_m(ji,jj)    &
145!old fashion                  &                    - 2.154996e-04 * zsss2_m(ji,jj) )  &
146!old fashion                            * tms(ji,jj)
147!old fashion            END DO
148!old fashion         END DO
149         t_bo(:,:) = tfreez( sss_m ) +  rt0 
150
151
152         ! ... ice albedo
153         CALL albedo_ice( t_su, ht_i, ht_s, alb_ice_cs, alb_ice_os )
154
155         ! ... Sea-ice surface boundary conditions output from bulk formulae :
156         !     - utaui_ice  ! surface ice stress i-component (I-point)   [N/m2]
157         !     - vtaui_ice  ! surface ice stress j-component (I-point)   [N/m2]
158         !     - qns_ice    ! non solar heat flux over ice   (T-point)   [W/m2]
159         !     - qsr_ice    !     solar heat flux over ice   (T-point)   [W/m2]
160         !     - qla_ice    ! latent    heat flux over ice   (T-point)   [W/m2]
161         !     - dqns_ice   ! non solar heat sensistivity    (T-point)   [W/m2]
162         !     - dqla_ice   ! latent    heat sensistivity    (T-point)   [W/m2]
163         !     - tprecip    ! total precipitation            (T-point)   [Kg/m2/s]
164         !     - sprecip    ! solid precipitation            (T-point)   [Kg/m2/s]
165         !     - fr1_i0     ! 1sr fraction of qsr penetration in ice     [%]
166         !     - fr2_i0     ! 2nd fraction of qsr penetration in ice     [%]
167         !
168         SELECT CASE( kblk )
169         CASE( 3 )           ! CLIO bulk formulation
170            CALL blk_ice_clio( t_su , u_ice , v_ice   , alb_ice_cs, alb_ice_os,            &
171               &                             utaui_ice, vtaui_ice , qns_ice   , qsr_ice,   &
172               &                             qla_ice  , dqns_ice  , dqla_ice  ,            &
173               &                             tprecip  , sprecip   ,                        &
174               &                             fr1_i0   , fr2_i0    , cl_grid  )
175         CASE( 4 )           ! CORE bulk formulation
176            CALL blk_ice_core( t_su , u_ice , v_ice   , alb_ice_cs,                      &
177               &                             utaui_ice, vtaui_ice , qns_ice , qsr_ice,   &
178               &                             qla_ice  , dqns_ice  , dqla_ice,            &
179               &                             tprecip  , sprecip   ,                      &
180               &                             fr1_i0   , fr2_i0    , cl_grid  )
181         END SELECT
182
183         IF(ln_ctl) THEN         ! print mean trends (used for debugging)
184            CALL prt_ctl_info( 'Ice Forcings ' )
185            CALL prt_ctl( tab2d_1=tprecip  ,clinfo1=' sbc_ice_lim: precip   : ' )
186            CALL prt_ctl( tab2d_1=utaui_ice,clinfo1=' sbc_ice_lim: utaui_ice: ', tab2d_2=vtaui_ice, clinfo2=' vtaui_ice: ' )
187            CALL prt_ctl( tab2d_1=sst_m    ,clinfo1=' sbc_ice_lim: sst      : ', tab2d_2=sss_m    , clinfo2=' sss      : ' )
188            CALL prt_ctl( tab2d_1=u_oce    ,clinfo1=' sbc_ice_lim: u_io     : ', tab2d_2=v_oce    , clinfo2=' v_io     : ' )
189            CALL prt_ctl( tab2d_1=frld     ,clinfo1=' sbc_ice_lim: frld   1 : ' )
190            !
191            DO jl = 1, jpl
192               CALL prt_ctl_info('* - category number ', ivar1=jl)
193               CALL prt_ctl(tab3d_1=t_su    , clinfo1=' sbc_ice_lim: t_su      : ', kdim=jl)
194               CALL prt_ctl(tab3d_1=qsr_ice , clinfo1=' sbc_ice_lim: qsr_ice   : ', kdim=jl)
195               CALL prt_ctl(tab3d_1=qns_ice , clinfo1=' sbc_ice_lim: qns_ice   : ', kdim=jl)
196               CALL prt_ctl(tab3d_1=dqns_ice, clinfo1=' sbc_ice_lim: dqns_ice  : ', kdim=jl)
197               CALL prt_ctl(tab3d_1=qla_ice , clinfo1=' sbc_ice_lim: qla_ice   : ', kdim=jl)
198               CALL prt_ctl(tab3d_1=dqla_ice, clinfo1=' sbc_ice_lim: dqla_ice  : ', kdim=jl)
199            END DO
200            !
201         ENDIF
202
203         !------------------------------------------------
204         ! Store old values of ice model global variables
205         !------------------------------------------------
206
207         old_a_i(:,:,:)   = a_i(:,:,:)     ! ice area
208         old_e_i(:,:,:,:) = e_i(:,:,:,:)   ! ice thermal energy
209         old_v_i(:,:,:)   = v_i(:,:,:)     ! ice volume
210         old_v_s(:,:,:)   = v_s(:,:,:)     ! snow volume
211         old_e_s(:,:,:,:) = e_s(:,:,:,:)   ! snow thermal energy
212         old_smv_i(:,:,:) = smv_i(:,:,:)   ! salt content
213         old_oa_i(:,:,:)  = oa_i(:,:,:)    ! areal age content
214
215         d_a_i_thd(:,:,:)   = 0.e0 ; d_a_i_trp(:,:,:)   = 0.e0
216         d_v_i_thd(:,:,:)   = 0.e0 ; d_v_i_trp(:,:,:)   = 0.e0
217         d_e_i_thd(:,:,:,:) = 0.e0 ; d_e_i_trp(:,:,:,:) = 0.e0
218         d_v_s_thd(:,:,:)   = 0.e0 ; d_v_s_trp(:,:,:)   = 0.e0
219         d_e_s_thd(:,:,:,:) = 0.e0 ; d_e_s_trp(:,:,:,:) = 0.e0
220         d_smv_i_thd(:,:,:) = 0.e0 ; d_smv_i_trp(:,:,:) = 0.e0
221         d_oa_i_thd(:,:,:)  = 0.e0 ; d_oa_i_trp(:,:,:)  = 0.e0
222
223         emps(:,:)      = 0.e0     ; fseqv(:,:)     = 0.e0
224         fsbri(:,:)     = 0.e0     ; fsalt_res(:,:) = 0.e0
225         fsalt_rpo(:,:) = 0.e0
226         fhmec(:,:)     = 0.e0     ; fhbri(:,:)     = 0.e0
227         fmmec(:,:)     = 0.e0     ; fheat_res(:,:) = 0.e0
228         fheat_rpo(:,:) = 0.e0     ; focea2D(:,:)   = 0.e0
229         fsup2D(:,:)    = 0.e0
230
231         diag_sni_gr(:,:) = 0.e0   ; diag_lat_gr(:,:) = 0.e0
232         diag_bot_gr(:,:) = 0.e0   ; diag_dyn_gr(:,:) = 0.e0
233         diag_bot_me(:,:) = 0.e0   ; diag_sur_me(:,:) = 0.e0
234
235         ! dynamical invariants
236         delta_i(:,:) = 0.e0
237         divu_i(:,:)  = 0.e0
238         shear_i(:,:) = 0.e0
239
240         !----------------!
241         ! Ice model step !
242         !----------------!
243                         numit = numit + nn_fsbc
244                         CALL lim_rst_opn( kt )          ! Open Ice restart file
245                         !+++++
246                         WRITE(numout,*) ' - Beginning the time step - '
247                         CALL lim_inst_state(jiindx,jjindx,1)
248                         WRITE(numout,*) ' '
249                         !+++++
250         !---------------------|
251         ! Dynamical processes |       
252         !---------------------|
253                         CALL lim_dyn                    ! Ice dynamics    ( rheology/dynamics )
254                         CALL lim_trp                    ! Ice transport   ( Advection/diffusion )
255                         CALL lim_var_agg(1)             ! aggregate categories, requested
256                         CALL lim_var_glo2eqv            ! equivalent variables, requested for rafting
257                         !+++++
258                         WRITE(numout,*) ' - After ice dynamics and transport '
259                         CALL lim_inst_state(jiindx,jjindx,1)
260                         WRITE(numout,*) ' '
261                         WRITE(numout,*)
262                         WRITE(numout,*) ' Mechanical Check ************** '
263                         WRITE(numout,*) ' Check what means ice divergence '
264                         WRITE(numout,*) ' Total ice concentration ', at_i (jiindx,jjindx)
265                         WRITE(numout,*) ' Total lead fraction     ', ato_i(jiindx,jjindx)
266                         WRITE(numout,*) ' Sum of both             ', ato_i(jiindx,jjindx) + at_i(jiindx,jjindx)
267                         WRITE(numout,*) ' Sum of both minus 1     ', ato_i(jiindx,jjindx) + at_i(jiindx,jjindx) - 1.00
268                         !+++++
269                         CALL lim_itd_me                 ! Mechanical redistribution ! (ridging/rafting)
270         !--------------------|
271         ! Ice thermodynamics |
272         !--------------------|
273                         CALL lim_var_glo2eqv            ! equivalent variables
274                         CALL lim_var_agg(1)             ! aggregate ice categories
275                         CALL lim_var_bv                 ! bulk brine volume (diag)
276                         CALL lim_thd                    ! Ice thermodynamics
277                         oa_i(:,:,:) = oa_i(:,:,:)   &
278                            &         + a_i(:,:,:)   &
279                            &         * rdt_ice      &
280                            &         / 86400.00         ! Ice natural aging
281                         CALL lim_var_glo2eqv            ! except info message that follows,
282                         !                               ! this CALL is maybe not necessary
283                         !+++++
284                         WRITE(numout,*) ' - After ice thermodynamics '
285                         CALL lim_inst_state(jiindx,jjindx,1)
286                         !+++++
287                         CALL lim_itd_th                 !  Remap ice categories, lateral accretion  !
288         !-------------------------|
289         ! Global variables update |
290         !-------------------------|
291                         CALL lim_var_agg(1)             ! requested by limupdate
292                         CALL lim_update                 ! Global variables update
293                         CALL lim_var_glo2eqv            ! equivalent variables (outputs)
294                         CALL lim_var_agg(2)             ! aggregate ice thickness categories
295                         !+++++
296                         IF(ln_nicep) THEN
297                            WRITE(numout,*) ' - Final ice state after lim_update '
298                            CALL lim_inst_state(jiindx,jjindx,2)
299                            WRITE(numout,*) ' '
300                         ENDIF
301                         !+++++
302         !--------------------------------------|
303         ! Fluxes of mass and heat to the ocean |
304         !--------------------------------------|
305                         CALL lim_sbc( kt )              ! Ice/Ocean Mass & Heat fluxes
306                         !+++++
307                         WRITE(numout,*) ' - Final ice state after lim_flx    '
308                         CALL lim_inst_state(jiindx,jjindx,3)
309                         WRITE(numout,*) ' '
310                         !+++++
311         !-------------------------|
312         ! Diagnostics and outputs |
313         !-------------------------|
314         IF( MOD( kt+nn_fsbc-1, ninfo ) == 0 .OR. ntmoy == 1 )   &
315          &              CALL lim_dia                    ! Ice Diagnostics
316                         CALL lim_wri      ( 1  )        ! Ice outputs
317         IF( lrst_ice )  CALL lim_rst_write( kt )        ! Ice restart file
318                         CALL lim_var_glo2eqv
319         !
320         ! Re-initialization of forcings
321         qsr     (:,:)   = 0.e0
322         qns     (:,:)   = 0.e0
323         tprecip (:,:)   = 0.e0 
324         sprecip (:,:)   = 0.e0
325         qsr_ice (:,:,:) = 0.e0
326         qns_ice (:,:,:) = 0.e0 
327         dqns_ice(:,:,:) = 0.e0 
328#if ! defined key_coupled
329         qla_ice (:,:,:) = 0.e0
330         dqla_ice(:,:,:) = 0.e0
331         fr1_i0  (:,:)   = 0.e0
332         fr2_i0  (:,:)   = 0.e0
333#else
334         rrunoff (:,:)   = 0.e0
335         calving (:,:)   = 0.e0
336#endif
337         !
338         !-------------------------------|
339         ! Alerts in case of model crash |
340         !-------------------------------|
341
342         numaltests = 10
343         numal(:) = 0
344
345         ! Alert if incompatible volume and concentration
346         alert_id = 2 ! reference number of this alert
347         alname(alert_id) = ' Incompat vol and con         ' ! name of the alert
348
349         DO jl = 1, jpl
350            DO jj = 1, jpj
351               DO ji = 1, jpi
352                     IF ((v_i(ji,jj,jl).NE.0.0).AND.(a_i(ji,jj,jl).EQ.0.0)) THEN
353                        WRITE(numout,*) ' ALERTE 2 '
354                        WRITE(numout,*) ' Incompatible volume and concentration '
355                        WRITE(numout,*) ' at_i     ', at_i(ji,jj)
356                        WRITE(numout,*) ' Point - category', ji, jj, jl
357                        WRITE(numout,*)
358                        WRITE(numout,*) ' a_i *** a_i_old ', a_i(ji,jj,jl), old_a_i(ji,jj,jl)
359                        WRITE(numout,*) ' v_i *** v_i_old ', v_i(ji,jj,jl), old_v_i(ji,jj,jl)
360                        WRITE(numout,*) ' d_a_i_thd/trp   ', d_a_i_thd(ji,jj,jl), d_a_i_trp(ji,jj,jl)
361                        WRITE(numout,*) ' d_v_i_thd/trp   ', d_v_i_thd(ji,jj,jl), d_v_i_trp(ji,jj,jl)
362                        numal(alert_id) = numal(alert_id) + 1
363                     ENDIF
364               END DO
365            END DO
366         END DO
367
368         ! Alerte if very thick ice
369         alert_id = 3 ! reference number of this alert
370         alname(alert_id) = ' Very thick ice               ' ! name of the alert
371         jl = jpl 
372         DO jj = 1, jpj
373            DO ji = 1, jpi
374               IF (ht_i(ji,jj,jl) .GT. 50.0 ) THEN
375                  WRITE(numout,*) ' ALERTE 3 '
376                  WRITE(numout,*) ' Very thick ice '
377                  CALL lim_inst_state(ji,jj,2)
378                  WRITE(numout,*) ' '
379                  numal(alert_id) = numal(alert_id) + 1
380               ENDIF
381            END DO
382         END DO
383
384         ! Alert if very fast ice
385         alert_id = 4 ! reference number of this alert
386         alname(alert_id) = ' Very fast ice               ' ! name of the alert
387         DO jj = 1, jpj
388            DO ji = 1, jpi
389               IF ( ( ( ABS( u_ice(ji,jj) ) .GT. 0.50) .OR. &
390                      ( ABS( v_ice(ji,jj) ) .GT. 0.50) ) .AND. &
391                    ( at_i(ji,jj) .GT. 0.0 ) ) THEN
392                  WRITE(numout,*) ' ALERTE 4 '
393                  WRITE(numout,*) ' Very fast ice '
394                  CALL lim_inst_state(ji,jj,1)
395                  WRITE(numout,*) ' ice strength             : ', strength(ji,jj)
396                  WRITE(numout,*) ' oceanic stress utau      : ', utau(ji,jj) 
397                  WRITE(numout,*) ' oceanic stress vtau      : ', vtau(ji,jj)
398                  WRITE(numout,*) ' sea-ice stress utaui_ice : ', utaui_ice(ji,jj) 
399                  WRITE(numout,*) ' sea-ice stress vtaui_ice : ', vtaui_ice(ji,jj)
400                  WRITE(numout,*) ' oceanic speed u          : ', u_oce(ji,jj)
401                  WRITE(numout,*) ' oceanic speed v          : ', v_oce(ji,jj)
402                  WRITE(numout,*) ' sst                      : ', sst_m(ji,jj)
403                  WRITE(numout,*) ' sss                      : ', sss_m(ji,jj)
404                  WRITE(numout,*) 
405                  numal(alert_id) = numal(alert_id) + 1
406               ENDIF
407            END DO
408         END DO
409
410         ! Alert if there is ice on continents
411         alert_id = 6 ! reference number of this alert
412         alname(alert_id) = ' Ice on continents           ' ! name of the alert
413         DO jj = 1, jpj
414            DO ji = 1, jpi
415               IF ( ( tms(ji,jj) .LE. 0.0 ) .AND. ( at_i(ji,jj) .GT. 0.0 ) ) THEN
416                  WRITE(numout,*) ' ALERTE 6 '
417                  WRITE(numout,*) ' Ice on continents '
418                  CALL lim_inst_state(ji,jj,1)
419                  WRITE(numout,*) ' masks s, u, v        : ', tms(ji,jj), &
420                                                              tmu(ji,jj), &
421                                                              tmv(ji,jj) 
422                  WRITE(numout,*) ' sst                  : ', sst_m(ji,jj)
423                  WRITE(numout,*) ' sss                  : ', sss_m(ji,jj)
424                  WRITE(numout,*) ' at_i(ji,jj)          : ', at_i(ji,jj)
425                  WRITE(numout,*) ' v_ice(ji,jj)         : ', v_ice(ji,jj)
426                  WRITE(numout,*) ' v_ice(ji,jj-1)       : ', v_ice(ji,jj-1)
427                  WRITE(numout,*) ' u_ice(ji-1,jj)       : ', u_ice(ji-1,jj)
428                  WRITE(numout,*) ' u_ice(ji,jj)         : ', v_ice(ji,jj)
429
430                  numal(alert_id) = numal(alert_id) + 1
431
432               ENDIF
433            END DO
434         END DO
435
436         ! Alert if very fresh ice
437         alert_id = 7 ! reference number of this alert
438         alname(alert_id) = ' Very fresh ice               ' ! name of the alert
439         DO jl = 1, jpl
440         DO jj = 1, jpj
441            DO ji = 1, jpi
442               IF ( ( ( ABS( sm_i(ji,jj,jl) ) .LT. 0.50) .OR. &
443                      ( ABS( sm_i(ji,jj,jl) ) .LT. 0.50) ) .AND. &
444                    ( a_i(ji,jj,jl) .GT. 0.0 ) ) THEN
445!                 WRITE(numout,*) ' ALERTE 7 '
446!                 WRITE(numout,*) ' Very fresh ice '
447!                 CALL lim_inst_state(ji,jj,1)
448!                 WRITE(numout,*) ' sst                  : ', sst_m(ji,jj)
449!                 WRITE(numout,*) ' sss                  : ', sss_m(ji,jj)
450!                 WRITE(numout,*) ' s_i_newice           : ', s_i_newice(ji,jj,1:jpl)
451!                 WRITE(numout,*)
452                  numal(alert_id) = numal(alert_id) + 1
453               ENDIF
454            END DO
455         END DO
456         END DO
457
458         ! Alert if too old ice
459         alert_id = 9 ! reference number of this alert
460         alname(alert_id) = ' Very old   ice               ' ! name of the alert
461         DO jl = 1, jpl
462         DO jj = 1, jpj
463            DO ji = 1, jpi
464               IF ( ( ( ABS( o_i(ji,jj,jl) ) .GT. rdt_ice ) .OR. &
465                      ( ABS( o_i(ji,jj,jl) ) .LT. 0.00) ) .AND. &
466                    ( a_i(ji,jj,jl) .GT. 0.0 ) ) THEN
467                  WRITE(numout,*) ' ALERTE 9 '
468                  WRITE(numout,*) ' Wrong ice age '
469                  CALL lim_inst_state(ji,jj,1)
470                  WRITE(numout,*) 
471                  numal(alert_id) = numal(alert_id) + 1
472               ENDIF
473            END DO
474         END DO
475
476         END DO
477
478         ! Alert on salt flux
479         alert_id = 5 ! reference number of this alert
480         alname(alert_id) = ' High salt flux               ' ! name of the alert
481         DO jj = 1, jpj
482            DO ji = 1, jpi
483               IF (ABS(emps(ji,jj)).gt.1.0e-2) THEN
484               WRITE(numout,*) ' ALERTE 5 '
485               WRITE(numout,*) ' High salt flux '
486               CALL lim_inst_state(ji,jj,3)
487               WRITE(numout,*) ' '
488               DO jl = 1, jpl
489                  WRITE(numout,*) ' Category no: ', jl
490                  WRITE(numout,*) ' a_i        : ', a_i(ji,jj,jl)     , &
491                                  ' old_a_i    : ', old_a_i(ji,jj,jl)   
492                  WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ji,jj,jl) , &
493                                  ' d_a_i_thd  : ', d_a_i_thd(ji,jj,jl) 
494                  WRITE(numout,*) ' v_i        : ', v_i(ji,jj,jl)    , &
495                                  ' old_v_i    : ', old_v_i(ji,jj,jl)   
496                  WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ji,jj,jl) , &
497                                  ' d_v_i_thd  : ', d_v_i_thd(ji,jj,jl) 
498                  WRITE(numout,*) ' '
499                  END DO
500                  numal(alert_id) = numal(alert_id) + 1
501               ENDIF
502            END DO
503         END DO
504
505         ! Alert if qns very big
506         alert_id = 8 ! reference number of this alert
507         alname(alert_id) = ' qns very big             ' ! name of the alert
508         DO jj = 1, jpj
509            DO ji = 1, jpi
510               IF ( ( ABS( qns(ji,jj) ) .GT. 1500.0 ) .AND. &
511                    ( at_i(ji,jj) .GT. 0.0 ) )  THEN
512
513                  WRITE(numout,*) ' ALERTE 8 '
514                  WRITE(numout,*) ' ji, jj    : ', ji, jj
515                  WRITE(numout,*) ' qns       : ', qns(ji,jj)
516                  WRITE(numout,*) ' Very high non-solar heat flux '
517                  WRITE(numout,*) ' sst       : ', sst_m(ji,jj)
518                  WRITE(numout,*) ' sss       : ', sss_m(ji,jj)
519                  WRITE(numout,*) ' qcmif     : ', qcmif(jiindx,jjindx)
520                  WRITE(numout,*) ' qldif     : ', qldif(jiindx,jjindx)
521                  WRITE(numout,*) ' qcmif     : ', qcmif(jiindx,jjindx) / rdt_ice
522                  WRITE(numout,*) ' qldif     : ', qldif(jiindx,jjindx) / rdt_ice
523                  WRITE(numout,*) ' qfvbq     : ', qfvbq(jiindx,jjindx)
524                  WRITE(numout,*) ' qdtcn     : ', qdtcn(jiindx,jjindx)
525                  WRITE(numout,*) ' qfvbq / dt: ', qfvbq(jiindx,jjindx) / rdt_ice
526                  WRITE(numout,*) ' qdtcn / dt: ', qdtcn(jiindx,jjindx) / rdt_ice
527                  WRITE(numout,*) ' fdtcn     : ', fdtcn(jiindx,jjindx) 
528                  WRITE(numout,*) ' fhmec     : ', fhmec(jiindx,jjindx) 
529                  WRITE(numout,*) ' fheat_rpo : ', fheat_rpo(jiindx,jjindx) 
530                  WRITE(numout,*) ' fheat_res : ', fheat_res(jiindx,jjindx) 
531                  WRITE(numout,*) ' fhbri     : ', fhbri(jiindx,jjindx) 
532
533                  CALL lim_inst_state(ji,jj,2)
534                  numal(alert_id) = numal(alert_id) + 1
535
536               ENDIF
537            END DO
538         END DO
539         !+++++
540 
541         ! Alert if very warm ice
542         alert_id = 10 ! reference number of this alert
543         alname(alert_id) = ' Very warm ice                ' ! name of the alert
544         numal(alert_id) = 0
545         DO jl = 1, jpl
546         DO jk = 1, nlay_i
547         DO jj = 1, jpj
548            DO ji = 1, jpi
549               ztmelts    =  -tmut*s_i(ji,jj,jk,jl) + rtt
550               IF ( ( t_i(ji,jj,jk,jl) .GE. ztmelts) .AND. &
551                    ( v_i(ji,jj,jl)   .GT. 1.0e-6)   .AND. &
552                    ( a_i(ji,jj,jl) .GT. 0.0     )  ) THEN
553                  WRITE(numout,*) ' ALERTE 10 '
554                  WRITE(numout,*) ' ji, jj, jk, jl : ', ji, jj, jk, jl
555                  WRITE(numout,*) ' Very warm ice '
556                  WRITE(numout,*) ' t_i : ', t_i(ji,jj,jk,jl)
557                  WRITE(numout,*) ' e_i : ', e_i(ji,jj,jk,jl)
558                  WRITE(numout,*) ' s_i : ', s_i(ji,jj,jk,jl)
559                  WRITE(numout,*) ' ztmelts : ', ztmelts
560                  numal(alert_id) = numal(alert_id) + 1
561               ENDIF
562            END DO
563         END DO
564         END DO
565         END DO
566
567         alert_id = 1 ! reference number of this alert
568         alname(alert_id) = ' Il n''y a pas d''alerte 1       ' ! name of the alert
569         WRITE(numout,*)
570         WRITE(numout,*) ' All alerts at the end of ice model '
571         DO alert_id = 1, numaltests
572            WRITE(numout,*) alert_id, alname(alert_id)//' : ', numal(alert_id), ' times ! '
573         END DO
574         !
575      ENDIF  ! End sea-ice coupling
576      !
577   END SUBROUTINE sbc_ice_lim
578
579
580   SUBROUTINE lim_inst_state( ki, kj, kn )
581      !!-----------------------------------------------------------------------
582      !!                   ***  ROUTINE lim_inst_state ***
583      !!                 
584      !! ** Purpose :   Writes global ice state on the (i,j) point
585      !!                in ocean.ouput
586      !!                3 possibilities exist
587      !!                n = 1 -> simple ice state
588      !!                n = 2 -> exhaustive state
589      !!                n = 3 -> ice/ocean salt fluxes
590      !!
591      !! ** input   :   point coordinates (i,j)
592      !!                n : number of the option
593      !!-------------------------------------------------------------------
594      INTEGER, INTENT( in ) ::   ki, kj, kn     ! ocean time-step index
595      INTEGER :: jl
596      !!-------------------------------------------------------------------
597
598      !----------------
599      !  Simple state
600      !----------------
601
602      IF ( kn .EQ. 1 ) THEN
603         WRITE(numout,*) ' lim_inst_state - Point : ',ki,kj
604         WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
605         WRITE(numout,*) ' Simple state '
606         WRITE(numout,*) ' masks s,u,v   : ', tms(ki,kj), tmu(ki,kj), tmv(ki,kj)
607         WRITE(numout,*) ' lat - long    : ', gphit(ki,kj), glamt(ki,kj)
608         WRITE(numout,*) ' Time step     : ', numit
609         WRITE(numout,*) ' - Ice drift   '
610         WRITE(numout,*) '   ~~~~~~~~~~~ '
611         WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ki-1,kj)
612         WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ki,kj)
613         WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ki,kj-1)
614         WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ki,kj)
615         WRITE(numout,*) ' strength      : ', strength(ki,kj)
616         WRITE(numout,*)
617         WRITE(numout,*) ' - Cell values '
618         WRITE(numout,*) '   ~~~~~~~~~~~ '
619         WRITE(numout,*) ' cell area     : ', area(ki,kj)
620         WRITE(numout,*) ' at_i          : ', at_i(ki,kj)       
621         WRITE(numout,*) ' vt_i          : ', vt_i(ki,kj)       
622         WRITE(numout,*) ' vt_s          : ', vt_s(ki,kj)       
623         DO jl = 1, jpl
624            WRITE(numout,*) ' - Category (', jl,')'
625            WRITE(numout,*) ' a_i           : ', a_i(ki,kj,jl)
626            WRITE(numout,*) ' ht_i          : ', ht_i(ki,kj,jl)
627            WRITE(numout,*) ' ht_s          : ', ht_s(ki,kj,jl)
628            WRITE(numout,*) ' v_i           : ', v_i(ki,kj,jl)
629            WRITE(numout,*) ' v_s           : ', v_s(ki,kj,jl)
630            WRITE(numout,*) ' e_s           : ', e_s(ki,kj,1,jl)/1.0e9
631            WRITE(numout,*) ' e_i           : ', e_i(ki,kj,1:nlay_i,jl)/1.0e9
632            WRITE(numout,*) ' t_su          : ', t_su(ki,kj,jl)
633            WRITE(numout,*) ' t_snow        : ', t_s(ki,kj,1,jl)
634            WRITE(numout,*) ' t_i           : ', t_i(ki,kj,1:nlay_i,jl)
635            WRITE(numout,*) ' sm_i          : ', sm_i(ki,kj,jl)
636            WRITE(numout,*) ' smv_i         : ', smv_i(ki,kj,jl)
637            WRITE(numout,*)
638            WRITE(numout,*) ' Pathological case : ', patho_case(ki,kj,jl)
639         END DO
640      ENDIF
641
642      !--------------------
643      !  Exhaustive state
644      !--------------------
645
646      IF ( kn .EQ. 2 ) THEN
647         WRITE(numout,*) ' lim_inst_state - Point : ',ki,kj
648         WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
649         WRITE(numout,*) ' Exhaustive state '
650         WRITE(numout,*) ' lat - long ', gphit(ki,kj), glamt(ki,kj)
651         WRITE(numout,*) ' Time step ', numit
652         WRITE(numout,*) 
653         WRITE(numout,*) ' - Cell values '
654         WRITE(numout,*) '   ~~~~~~~~~~~ '
655         WRITE(numout,*) ' cell area     : ', area(ki,kj)
656         WRITE(numout,*) ' at_i          : ', at_i(ki,kj)       
657         WRITE(numout,*) ' vt_i          : ', vt_i(ki,kj)       
658         WRITE(numout,*) ' vt_s          : ', vt_s(ki,kj)       
659         WRITE(numout,*) ' sdvt          : ', sdvt(ki,kj)       
660         WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ki-1,kj)
661         WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ki,kj)
662         WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ki,kj-1)
663         WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ki,kj)
664         WRITE(numout,*) ' strength      : ', strength(ki,kj)
665         WRITE(numout,*) ' d_u_ice_dyn   : ', d_u_ice_dyn(ki,kj), ' d_v_ice_dyn   : ', d_v_ice_dyn(ki,kj)
666         WRITE(numout,*) ' old_u_ice     : ', old_u_ice(ki,kj)  , ' old_v_ice     : ', old_v_ice(ki,kj) 
667         WRITE(numout,*)
668
669         DO jl = 1, jpl
670              WRITE(numout,*) ' - Category (',jl,')'
671              WRITE(numout,*) '   ~~~~~~~~         ' 
672              WRITE(numout,*) ' ht_i       : ', ht_i(ki,kj,jl)             , ' ht_s       : ', ht_s(ki,kj,jl)
673              WRITE(numout,*) ' t_i        : ', t_i(ki,kj,1:nlay_i,jl)
674              WRITE(numout,*) ' t_su       : ', t_su(ki,kj,jl)             , ' t_s        : ', t_s(ki,kj,1,jl)
675              WRITE(numout,*) ' sm_i       : ', sm_i(ki,kj,jl)             , ' o_i        : ', o_i(ki,kj,jl)
676              WRITE(numout,*) ' a_i        : ', a_i(ki,kj,jl)              , ' old_a_i    : ', old_a_i(ki,kj,jl)   
677              WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ki,kj,jl)        , ' d_a_i_thd  : ', d_a_i_thd(ki,kj,jl) 
678              WRITE(numout,*) ' v_i        : ', v_i(ki,kj,jl)              , ' old_v_i    : ', old_v_i(ki,kj,jl)   
679              WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ki,kj,jl)        , ' d_v_i_thd  : ', d_v_i_thd(ki,kj,jl) 
680              WRITE(numout,*) ' v_s        : ', v_s(ki,kj,jl)              , ' old_v_s    : ', old_v_s(ki,kj,jl) 
681              WRITE(numout,*) ' d_v_s_trp  : ', d_v_s_trp(ki,kj,jl)        , ' d_v_s_thd  : ', d_v_s_thd(ki,kj,jl)
682              WRITE(numout,*) ' e_i1       : ', e_i(ki,kj,1,jl)/1.0e9      , ' old_ei1    : ', old_e_i(ki,kj,1,jl)/1.0e9 
683              WRITE(numout,*) ' de_i1_trp  : ', d_e_i_trp(ki,kj,1,jl)/1.0e9, ' de_i1_thd  : ', d_e_i_thd(ki,kj,1,jl)/1.0e9
684              WRITE(numout,*) ' e_i2       : ', e_i(ki,kj,2,jl)/1.0e9      , ' old_ei2    : ', old_e_i(ki,kj,2,jl)/1.0e9 
685              WRITE(numout,*) ' de_i2_trp  : ', d_e_i_trp(ki,kj,2,jl)/1.0e9, ' de_i2_thd  : ', d_e_i_thd(ki,kj,2,jl)/1.0e9
686              WRITE(numout,*) ' e_snow     : ', e_s(ki,kj,1,jl)            , ' old_e_snow : ', old_e_s(ki,kj,1,jl) 
687              WRITE(numout,*) ' d_e_s_trp  : ', d_e_s_trp(ki,kj,1,jl)      , ' d_e_s_thd  : ', d_e_s_thd(ki,kj,1,jl)
688              WRITE(numout,*) ' smv_i      : ', smv_i(ki,kj,jl)            , ' old_smv_i  : ', old_smv_i(ki,kj,jl)   
689              WRITE(numout,*) ' d_smv_i_trp: ', d_smv_i_trp(ki,kj,jl)      , ' d_smv_i_thd: ', d_smv_i_thd(ki,kj,jl) 
690              WRITE(numout,*) ' oa_i       : ', oa_i(ki,kj,jl)             , ' old_oa_i   : ', old_oa_i(ki,kj,jl)
691              WRITE(numout,*) ' d_oa_i_trp : ', d_oa_i_trp(ki,kj,jl)       , ' d_oa_i_thd : ', d_oa_i_thd(ki,kj,jl)
692              WRITE(numout,*) ' Path. case : ', patho_case(ki,kj,jl)
693        END DO !jl
694
695        WRITE(numout,*)
696        WRITE(numout,*) ' - Heat / FW fluxes '
697        WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ '
698!       WRITE(numout,*) ' fsbri      : ', fsbri(ki,kj)
699!       WRITE(numout,*) ' fseqv      : ', fseqv(ki,kj)
700!       WRITE(numout,*) ' fsalt_res  : ', fsalt_res(ki,kj)
701        WRITE(numout,*) ' fmmec      : ', fmmec(ki,kj)
702        WRITE(numout,*) ' fhmec      : ', fhmec(ki,kj)
703        WRITE(numout,*) ' fhbri      : ', fhbri(ki,kj)
704        WRITE(numout,*) ' fheat_rpo  : ', fheat_rpo(ki,kj)
705        WRITE(numout,*) 
706        WRITE(numout,*) ' sst        : ', sst_m(ki,kj) 
707        WRITE(numout,*) ' sss        : ', sss_m(ki,kj) 
708        WRITE(numout,*) 
709        WRITE(numout,*) ' - Stresses '
710        WRITE(numout,*) '   ~~~~~~~~ '
711        WRITE(numout,*) ' utaui_ice  : ', utaui_ice(ki,kj) 
712        WRITE(numout,*) ' vtaui_ice  : ', vtaui_ice(ki,kj)
713        WRITE(numout,*) ' utau       : ', utau(ki,kj) 
714        WRITE(numout,*) ' vtau       : ', vtau(ki,kj)
715        WRITE(numout,*) ' oc. vel. u : ', u_oce(ki,kj)
716        WRITE(numout,*) ' oc. vel. v : ', v_oce(ki,kj)
717     ENDIF
718
719     !---------------------
720     ! Salt / heat fluxes
721     !---------------------
722
723     IF ( kn .EQ. 3 ) THEN
724        WRITE(numout,*) ' lim_inst_state - Point : ',ki,kj
725        WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
726        WRITE(numout,*) ' - Salt / Heat Fluxes '
727        WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ '
728        WRITE(numout,*) ' lat - long ', gphit(ki,kj), glamt(ki,kj)
729        WRITE(numout,*) ' Time step ', numit
730        WRITE(numout,*)
731        WRITE(numout,*) ' - Heat fluxes at bottom interface ***'
732        WRITE(numout,*) ' qsr        : ', qsr(ki,kj)
733        WRITE(numout,*) ' qns        : ', qns(ki,kj)
734        WRITE(numout,*)
735        WRITE(numout,*) ' - Salt fluxes at bottom interface ***'
736        WRITE(numout,*) ' emps       : ', emps(ki,kj)
737        WRITE(numout,*) ' emp        : ', emp(ki,kj)
738        WRITE(numout,*) ' fsbri      : ', fsbri(ki,kj)
739        WRITE(numout,*) ' fseqv      : ', fseqv(ki,kj)
740        WRITE(numout,*) ' fsalt_res  : ', fsalt_res(ki,kj)
741        WRITE(numout,*) ' fsalt_rpo  : ', fsalt_rpo(ki,kj)
742        WRITE(numout,*) ' - Heat fluxes at bottom interface ***'
743        WRITE(numout,*) ' fheat_res  : ', fheat_res(ki,kj)
744        WRITE(numout,*)
745        WRITE(numout,*) ' - Momentum fluxes '
746        WRITE(numout,*) ' utau      : ', utau(ki,kj) 
747        WRITE(numout,*) ' vtau      : ', vtau(ki,kj)
748      ENDIF
749
750   END SUBROUTINE lim_inst_state
751
752#else
753   !!----------------------------------------------------------------------
754   !!   Default option           Dummy module      NO LIM 3.0 sea-ice model
755   !!----------------------------------------------------------------------
756CONTAINS
757   SUBROUTINE sbc_ice_lim ( kt, kblk )     ! Dummy routine
758      WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk
759   END SUBROUTINE sbc_ice_lim
760#endif
761
762   !!======================================================================
763END MODULE sbcice_lim
Note: See TracBrowser for help on using the repository browser.