New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
sbcice_lim.F90 in branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 2370

Last change on this file since 2370 was 2370, checked in by gm, 13 years ago

v3.3beta: ice-ocean stress at kt with VP & EVP (LIM-2 and -3)

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