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_2.F90 in branches/nemo_v3_3_beta/NEMOGCM/NEMO/LIM_SRC_2 – NEMO

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/LIM_SRC_2/limdyn_2.F90 @ 2332

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

correct small bugs in LIM2, see ticket #746

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