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

Last change on this file since 1308 was 1270, checked in by rblod, 15 years ago

Light bug in blk_ice, see ticket #289

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