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 branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/LIM_SRC_3/limdyn.F90 @ 3750

Last change on this file since 3750 was 3625, checked in by acc, 11 years ago

Branch dev_NOC_2012_r3555. #1006. Step 7. Check in code now merged with dev_r3385_NOCS04_HAMF

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