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

source: branches/2013/dev_r4028_CNRS_LIM3/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 4042

Last change on this file since 4042 was 4042, checked in by clem, 11 years ago

new LIM3

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