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.
limdyn.F90 in trunk/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: trunk/NEMOGCM/NEMO/LIM_SRC_3/limdyn.F90 @ 2528

Last change on this file since 2528 was 2528, checked in by rblod, 13 years ago

Update NEMOGCM from branch nemo_v3_3_beta

  • Property svn:keywords set to Id
File size: 13.5 KB
Line 
1MODULE limdyn
2   !!======================================================================
3   !!                     ***  MODULE  limdyn  ***
4   !!   Sea-Ice dynamics : 
5   !!======================================================================
6   !! history :  1.0  ! 2002-08 (C. Ethe, G. Madec)  original VP code
7   !!            3.0  ! 2007-03 (MA Morales Maqueda, S. Bouillon, M. Vancoppenolle)  LIM3: EVP-Cgrid
8   !!----------------------------------------------------------------------
9#if defined key_lim3
10   !!----------------------------------------------------------------------
11   !!   'key_lim3' :                                 LIM3 sea-ice model
12   !!----------------------------------------------------------------------
13   !!    lim_dyn      : computes ice velocities
14   !!    lim_dyn_init : initialization and namelist read
15   !!----------------------------------------------------------------------
16   USE phycst
17   USE in_out_manager  ! I/O manager
18   USE dom_ice
19   USE dom_oce         ! ocean space and time domain
20   USE ice
21   USE par_ice
22   USE sbc_oce         ! Surface boundary condition: ocean fields
23   USE sbc_ice         ! Surface boundary condition: ice fields
24   USE limrhg          ! ice rheology
25   USE lbclnk
26   USE lib_mpp
27   USE prtctl          ! Print control
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC   lim_dyn   ! routine called by ice_step
33
34   !! * Substitutions
35#  include "vectopt_loop_substitute.h90"
36   !!----------------------------------------------------------------------
37   !! NEMO/LIM3 3.3 , UCL - NEMO Consortium (2010)
38   !! $Id$
39   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
40   !!----------------------------------------------------------------------
41CONTAINS
42
43   SUBROUTINE lim_dyn( kt )
44      !!-------------------------------------------------------------------
45      !!               ***  ROUTINE lim_dyn  ***
46      !!               
47      !! ** Purpose :   compute ice velocity and ocean-ice stress
48      !!               
49      !! ** Method  :
50      !!
51      !! ** Action  : - Initialisation
52      !!              - Call of the dynamic routine for each hemisphere
53      !!              - computation of the stress at the ocean surface         
54      !!              - treatment of the case if no ice dynamic
55      !!------------------------------------------------------------------------------------
56      INTEGER, INTENT(in) ::   kt     ! number of iteration
57      !!
58      INTEGER  ::   ji, jj, jl, ja    ! dummy loop indices
59      INTEGER  ::   i_j1, i_jpj       ! Starting/ending j-indices for rheology
60      REAL(wp) ::   zcoef             ! local scalar
61      REAL(wp), DIMENSION(jpj)     ::   zind           ! i-averaged indicator of sea-ice
62      REAL(wp), DIMENSION(jpj)     ::   zmsk           ! i-averaged of tmask
63      REAL(wp), DIMENSION(jpi,jpj) ::   zu_io, zv_io   ! ice-ocean velocity
64      !!---------------------------------------------------------------------
65
66      IF( kt == nit000 .AND. lwp ) THEN
67         WRITE(numout,*) ' lim_dyn : Ice dynamics '
68         WRITE(numout,*) ' ~~~~~~~ '
69      ENDIF
70
71      IF( numit == nstart  )   CALL lim_dyn_init   ! Initialization (first time-step only)
72
73      IF ( ln_limdyn ) THEN
74
75         old_u_ice(:,:) = u_ice(:,:) * tmu(:,:)
76         old_v_ice(:,:) = v_ice(:,:) * tmv(:,:)
77
78         ! Rheology (ice dynamics)
79         ! ========
80
81         !  Define the j-limits where ice rheology is computed
82         ! ---------------------------------------------------
83
84         IF( lk_mpp .OR. lk_mpp_rep ) THEN                    ! mpp: compute over the whole domain
85            i_j1 = 1
86            i_jpj = jpj
87            IF(ln_ctl) CALL prt_ctl_info( 'lim_dyn  :    i_j1 = ', ivar1=i_j1, clinfo2=' ij_jpj = ', ivar2=i_jpj )
88            CALL lim_rhg( i_j1, i_jpj )
89         ELSE                                 ! optimization of the computational area
90
91            DO jj = 1, jpj
92               zind(jj) = SUM( 1.0 - at_i (:,jj  ) )   ! = FLOAT(jpj) if ocean everywhere on a j-line
93               zmsk(jj) = SUM( tmask(:,jj,1) )   ! = 0          if land  everywhere on a j-line
94            END DO
95
96            IF( l_jeq ) THEN                     ! local domain include both hemisphere
97               !                                 ! Rheology is computed in each hemisphere
98               !                                 ! only over the ice cover latitude strip
99               ! Northern hemisphere
100               i_j1  = njeq
101               i_jpj = jpj
102               DO WHILE ( i_j1 <= jpj .AND. zind(i_j1) == FLOAT(jpi) .AND. zmsk(i_j1) /=0 )
103                  i_j1 = i_j1 + 1
104               END DO
105               i_j1 = MAX( 1, i_j1-2 )
106               IF(ln_ctl) CALL prt_ctl_info( 'lim_dyn  : NH  i_j1 = ', ivar1=i_j1, clinfo2=' ij_jpj = ', ivar2=i_jpj )
107               CALL lim_rhg( i_j1, i_jpj )
108
109               ! Southern hemisphere
110               i_j1  =  1
111               i_jpj = njeq
112               DO WHILE ( i_jpj >= 1 .AND. zind(i_jpj) == FLOAT(jpi) .AND. zmsk(i_jpj) /=0 )
113                  i_jpj = i_jpj - 1
114               END DO
115               i_jpj = MIN( jpj, i_jpj+1 )
116               IF(ln_ctl) CALL prt_ctl_info( 'lim_dyn  : SH  i_j1 = ', ivar1=i_j1, clinfo2=' ij_jpj = ', ivar2=i_jpj )
117
118       CALL lim_rhg( i_j1, i_jpj )
119
120    ELSE                                 ! local domain extends over one hemisphere only
121       !                                 ! Rheology is computed only over the ice cover
122       !                                 ! latitude strip
123       i_j1  = 1
124               DO WHILE ( i_j1 <= jpj .AND. zind(i_j1) == FLOAT(jpi) .AND. zmsk(i_j1) /=0 )
125                  i_j1 = i_j1 + 1
126               END DO
127               i_j1 = MAX( 1, i_j1-2 )
128
129               i_jpj  = jpj
130               DO WHILE ( i_jpj >= 1  .AND. zind(i_jpj) == FLOAT(jpi) .AND. zmsk(i_jpj) /=0 )
131                  i_jpj = i_jpj - 1
132               END DO
133               i_jpj = MIN( jpj, i_jpj+1)
134
135               IF(ln_ctl) CALL prt_ctl_info( 'lim_dyn  : one hemisphere:  i_j1 = ', ivar1=i_j1, clinfo2=' ij_jpj = ', ivar2=i_jpj )
136
137               CALL lim_rhg( i_j1, i_jpj )
138
139            ENDIF
140
141         ENDIF
142
143         ! computation of friction velocity
144         ! --------------------------------
145         ! ice-ocean velocity at U & V-points (u_ice v_ice at U- & V-points ; ssu_m, ssv_m at U- & V-points)
146         zu_io(:,:) = u_ice(:,:) - ssu_m(:,:)
147         zv_io(:,:) = v_ice(:,:) - ssv_m(:,:)
148         ! frictional velocity at T-point
149         zcoef = 0.5 * cw
150         DO jj = 2, jpjm1 
151            DO ji = fs_2, fs_jpim1   ! vector opt.
152               ust2s(ji,jj) = zcoef * (  zu_io(ji,jj) * zu_io(ji,jj) + zu_io(ji-1,jj) * zu_io(ji-1,jj)   &
153                  &                    + zv_io(ji,jj) * zv_io(ji,jj) + zv_io(ji,jj-1) * zv_io(ji,jj-1)   ) * tms(ji,jj)
154            END DO
155         END DO
156         !
157      ELSE      ! no ice dynamics : transmit directly the atmospheric stress to the ocean
158         !
159         zcoef = SQRT( 0.5 ) / rau0
160         DO jj = 2, jpjm1
161            DO ji = fs_2, fs_jpim1   ! vector opt.
162               ust2s(ji,jj) = zcoef * SQRT(  utau(ji,jj) * utau(ji,jj) + utau(ji-1,jj) * utau(ji-1,jj)   &
163                  &                        + vtau(ji,jj) * vtau(ji,jj) + vtau(ji,jj-1) * vtau(ji,jj-1)   ) * tms(ji,jj)
164            END DO
165         END DO
166         !
167      ENDIF
168      CALL lbc_lnk( ust2s, 'T',  1. )   ! T-point
169
170      IF(ln_ctl) THEN   ! Control print
171         CALL prt_ctl_info(' ')
172         CALL prt_ctl_info(' - Cell values : ')
173         CALL prt_ctl_info('   ~~~~~~~~~~~~~ ')
174         CALL prt_ctl(tab2d_1=ust2s     , clinfo1=' lim_dyn  : ust2s     :')
175         CALL prt_ctl(tab2d_1=divu_i    , clinfo1=' lim_dyn  : divu_i    :')
176         CALL prt_ctl(tab2d_1=delta_i   , clinfo1=' lim_dyn  : delta_i   :')
177         CALL prt_ctl(tab2d_1=strength  , clinfo1=' lim_dyn  : strength  :')
178         CALL prt_ctl(tab2d_1=area      , clinfo1=' lim_dyn  : cell area :')
179         CALL prt_ctl(tab2d_1=at_i      , clinfo1=' lim_dyn  : at_i      :')
180         CALL prt_ctl(tab2d_1=vt_i      , clinfo1=' lim_dyn  : vt_i      :')
181         CALL prt_ctl(tab2d_1=vt_s      , clinfo1=' lim_dyn  : vt_s      :')
182         CALL prt_ctl(tab2d_1=stress1_i , clinfo1=' lim_dyn  : stress1_i :')
183         CALL prt_ctl(tab2d_1=stress2_i , clinfo1=' lim_dyn  : stress2_i :')
184         CALL prt_ctl(tab2d_1=stress12_i, clinfo1=' lim_dyn  : stress12_i:')
185         DO jl = 1, jpl
186            CALL prt_ctl_info(' ')
187            CALL prt_ctl_info(' - Category : ', ivar1=jl)
188            CALL prt_ctl_info('   ~~~~~~~~~~')
189            CALL prt_ctl(tab2d_1=a_i   (:,:,jl)   , clinfo1= ' lim_dyn  : a_i      : ')
190            CALL prt_ctl(tab2d_1=ht_i  (:,:,jl)   , clinfo1= ' lim_dyn  : ht_i     : ')
191            CALL prt_ctl(tab2d_1=ht_s  (:,:,jl)   , clinfo1= ' lim_dyn  : ht_s     : ')
192            CALL prt_ctl(tab2d_1=v_i   (:,:,jl)   , clinfo1= ' lim_dyn  : v_i      : ')
193            CALL prt_ctl(tab2d_1=v_s   (:,:,jl)   , clinfo1= ' lim_dyn  : v_s      : ')
194            CALL prt_ctl(tab2d_1=e_s   (:,:,1,jl) , clinfo1= ' lim_dyn  : e_s      : ')
195            CALL prt_ctl(tab2d_1=t_su  (:,:,jl)   , clinfo1= ' lim_dyn  : t_su     : ')
196            CALL prt_ctl(tab2d_1=t_s   (:,:,1,jl) , clinfo1= ' lim_dyn  : t_snow   : ')
197            CALL prt_ctl(tab2d_1=sm_i  (:,:,jl)   , clinfo1= ' lim_dyn  : sm_i     : ')
198            CALL prt_ctl(tab2d_1=smv_i (:,:,jl)   , clinfo1= ' lim_dyn  : smv_i    : ')
199            DO ja = 1, nlay_i
200               CALL prt_ctl_info(' ')
201               CALL prt_ctl_info(' - Layer : ', ivar1=ja)
202               CALL prt_ctl_info('   ~~~~~~~')
203               CALL prt_ctl(tab2d_1=t_i(:,:,ja,jl) , clinfo1= ' lim_dyn  : t_i      : ')
204               CALL prt_ctl(tab2d_1=e_i(:,:,ja,jl) , clinfo1= ' lim_dyn  : e_i      : ')
205            END DO
206         END DO
207      ENDIF
208      !
209   END SUBROUTINE lim_dyn
210
211
212   SUBROUTINE lim_dyn_init
213      !!-------------------------------------------------------------------
214      !!                  ***  ROUTINE lim_dyn_init  ***
215      !!
216      !! ** Purpose : Physical constants and parameters linked to the ice
217      !!      dynamics
218      !!
219      !! ** Method  :  Read the namicedyn namelist and check the ice-dynamic
220      !!       parameter values called at the first timestep (nit000)
221      !!
222      !! ** input   :   Namelist namicedyn
223      !!-------------------------------------------------------------------
224      NAMELIST/namicedyn/ epsd, alpha,     &
225         &                dm, nbiter, nbitdr, om, resl, cw, angvg, pstar,   &
226         &                c_rhg, etamn, creepl, ecc, ahi0, &
227         &                nevp, telast, alphaevp
228      !!-------------------------------------------------------------------
229
230      REWIND( numnam_ice )                ! Read Namelist namicedyn
231      READ  ( numnam_ice  , namicedyn )
232     
233      IF(lwp) THEN                        ! control print
234         WRITE(numout,*)
235         WRITE(numout,*) 'lim_dyn_init : ice parameters for ice dynamics '
236         WRITE(numout,*) '~~~~~~~~~~~~'
237         WRITE(numout,*) '   tolerance parameter                              epsd   = ', epsd
238         WRITE(numout,*) '   coefficient for semi-implicit coriolis           alpha  = ', alpha
239         WRITE(numout,*) '   diffusion constant for dynamics                  dm     = ', dm
240         WRITE(numout,*) '   number of sub-time steps for relaxation          nbiter = ', nbiter
241         WRITE(numout,*) '   maximum number of iterations for relaxation      nbitdr = ', nbitdr
242         WRITE(numout,*) '   relaxation constant                              om     = ', om
243         WRITE(numout,*) '   maximum value for the residual of relaxation     resl   = ', resl
244         WRITE(numout,*) '   drag coefficient for oceanic stress              cw     = ', cw
245         WRITE(numout,*) '   turning angle for oceanic stress                 angvg  = ', angvg
246         WRITE(numout,*) '   first bulk-rheology parameter                    pstar  = ', pstar
247         WRITE(numout,*) '   second bulk-rhelogy parameter                    c_rhg  = ', c_rhg
248         WRITE(numout,*) '   minimun value for viscosity                      etamn  = ', etamn
249         WRITE(numout,*) '   creep limit                                      creepl = ', creepl
250         WRITE(numout,*) '   eccentricity of the elliptical yield curve       ecc    = ', ecc
251         WRITE(numout,*) '   horizontal diffusivity coeff. for sea-ice        ahi0   = ', ahi0
252         WRITE(numout,*) '   number of iterations for subcycling              nevp   = ', nevp
253         WRITE(numout,*) '   timescale for elastic waves                      telast = ', telast
254         WRITE(numout,*) '   coefficient for the solution of int. stresses  alphaevp = ', alphaevp
255      ENDIF
256      !
257      IF( angvg /= 0._wp ) THEN
258         CALL ctl_warn( 'lim_dyn_init: turning angle for oceanic stress not properly coded for EVP ',   &
259            &           '(see limsbc module). We force  angvg = 0._wp'  )
260         angvg = 0._wp
261      ENDIF
262     
263      usecc2 = 1._wp / ( ecc * ecc )
264      rhoco  = rau0  * cw
265      angvg  = angvg * rad
266      sangvg = SIN( angvg )
267      cangvg = COS( angvg )
268      pstarh = pstar * 0.5_wp
269
270      !  Diffusion coefficients.
271      ahiu(:,:) = ahi0 * umask(:,:,1)
272      ahiv(:,:) = ahi0 * vmask(:,:,1)
273
274   END SUBROUTINE lim_dyn_init
275
276#else
277   !!----------------------------------------------------------------------
278   !!   Default option          Empty module           NO LIM sea-ice model
279   !!----------------------------------------------------------------------
280CONTAINS
281   SUBROUTINE lim_dyn         ! Empty routine
282   END SUBROUTINE lim_dyn
283#endif 
284
285   !!======================================================================
286END MODULE limdyn
Note: See TracBrowser for help on using the repository browser.