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 trunk/NEMOGCM/NEMO/LIM_SRC_2 – NEMO

source: trunk/NEMOGCM/NEMO/LIM_SRC_2/limdyn_2.F90 @ 3432

Last change on this file since 3432 was 3294, checked in by rblod, 12 years ago

Merge of 3.4beta into the trunk

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