123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836 |
- MODULE limrhg
- !!======================================================================
- !! *** MODULE limrhg ***
- !! Ice rheology : sea ice rheology
- !!======================================================================
- !! History : - ! 2007-03 (M.A. Morales Maqueda, S. Bouillon) Original code
- !! 3.0 ! 2008-03 (M. Vancoppenolle) LIM3
- !! - ! 2008-11 (M. Vancoppenolle, S. Bouillon, Y. Aksenov) add surface tilt in ice rheolohy
- !! 3.3 ! 2009-05 (G.Garric) addition of the lim2_evp cas
- !! 3.4 ! 2011-01 (A. Porter) dynamical allocation
- !! 3.5 ! 2012-08 (R. Benshila) AGRIF
- !! 3.6 ! 2016-06 (C. Rousset) Rewriting (conserves energy)
- !!----------------------------------------------------------------------
- #if defined key_lim3 || ( defined key_lim2 && ! defined key_lim2_vp )
- !!----------------------------------------------------------------------
- !! 'key_lim3' OR LIM-3 sea-ice model
- !! 'key_lim2' AND NOT 'key_lim2_vp' EVP LIM-2 sea-ice model
- !!----------------------------------------------------------------------
- !! lim_rhg : computes ice velocities
- !!----------------------------------------------------------------------
- USE phycst ! Physical constant
- USE oce , ONLY : snwice_mass, snwice_mass_b
- USE par_oce ! Ocean parameters
- USE dom_oce ! Ocean domain
- USE sbc_oce ! Surface boundary condition: ocean fields
- USE sbc_ice ! Surface boundary condition: ice fields
- #if defined key_lim3
- USE ice ! LIM-3: ice variables
- USE dom_ice ! LIM-3: ice domain
- USE limitd_me ! LIM-3:
- #else
- USE ice_2 ! LIM-2: ice variables
- USE dom_ice_2 ! LIM-2: ice domain
- #endif
- USE lbclnk ! Lateral Boundary Condition / MPP link
- USE lib_mpp ! MPP library
- USE wrk_nemo ! work arrays
- USE in_out_manager ! I/O manager
- USE prtctl ! Print control
- USE iom
- USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
- #if defined key_agrif && defined key_lim2
- USE agrif_lim2_interp
- #endif
- #if defined key_bdy
- USE bdyice_lim
- #endif
- IMPLICIT NONE
- PRIVATE
- PUBLIC lim_rhg ! routine called by lim_dyn (or lim_dyn_2)
- !! * Substitutions
- # include "vectopt_loop_substitute.h90"
- !!----------------------------------------------------------------------
- !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011)
- !! $Id: limrhg.F90 8285 2017-07-06 06:40:51Z vancop $
- !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
- !!----------------------------------------------------------------------
- CONTAINS
- SUBROUTINE lim_rhg( k_j1, k_jpj )
- !!-------------------------------------------------------------------
- !! *** SUBROUTINE lim_rhg ***
- !! EVP-C-grid
- !!
- !! ** purpose : determines sea ice drift from wind stress, ice-ocean
- !! stress and sea-surface slope. Ice-ice interaction is described by
- !! a non-linear elasto-viscous-plastic (EVP) law including shear
- !! strength and a bulk rheology (Hunke and Dukowicz, 2002).
- !!
- !! The points in the C-grid look like this, dear reader
- !!
- !! (ji,jj)
- !! |
- !! |
- !! (ji-1,jj) | (ji,jj)
- !! ---------
- !! | |
- !! | (ji,jj) |------(ji,jj)
- !! | |
- !! ---------
- !! (ji-1,jj-1) (ji,jj-1)
- !!
- !! ** Inputs : - wind forcing (stress), oceanic currents
- !! ice total volume (vt_i) per unit area
- !! snow total volume (vt_s) per unit area
- !!
- !! ** Action : - compute u_ice, v_ice : the components of the
- !! sea-ice velocity vector
- !! - compute delta_i, shear_i, divu_i, which are inputs
- !! of the ice thickness distribution
- !!
- !! ** Steps : 1) Compute ice snow mass, ice strength
- !! 2) Compute wind, oceanic stresses, mass terms and
- !! coriolis terms of the momentum equation
- !! 3) Solve the momentum equation (iterative procedure)
- !! 4) Recompute invariants of the strain rate tensor
- !! which are inputs of the ITD, store stress
- !! for the next time step
- !! 5) Control prints of residual (convergence)
- !! and charge ellipse.
- !! The user should make sure that the parameters
- !! nn_nevp, elastic time scale and rn_creepl maintain stress state
- !! on the charge ellipse for plastic flow
- !! e.g. in the Canadian Archipelago
- !!
- !! ** Notes : Boundary condition for ice is chosen no-slip
- !! but can be adjusted with param rn_shlat
- !!
- !! References : Hunke and Dukowicz, JPO97
- !! Bouillon et al., Ocean Modelling 2009
- !!-------------------------------------------------------------------
- INTEGER, INTENT(in) :: k_j1 ! southern j-index for ice computation
- INTEGER, INTENT(in) :: k_jpj ! northern j-index for ice computation
- !!
- INTEGER :: ji, jj ! dummy loop indices
- INTEGER :: jter ! local integers
- CHARACTER (len=50) :: charout
- REAL(wp) :: zdtevp, z1_dtevp ! time step for subcycling
- REAL(wp) :: ecc2, z1_ecc2 ! square of yield ellipse eccenticity
- REAL(wp) :: zbeta, zalph1, z1_alph1, zalph2, z1_alph2 ! alpha and beta from Bouillon 2009 and 2013
- REAL(wp) :: zm1, zm2, zm3, zmassU, zmassV ! ice/snow mass
- REAL(wp) :: zdelta, zp_delf, zds2, zdt, zdt2, zdiv, zdiv2 ! temporary scalars
- REAL(wp) :: zTauO, zTauE ! temporary scalars
- REAL(wp) :: zsig1, zsig2 ! internal ice stress
- REAL(wp) :: zresm ! Maximal error on ice velocity
- REAL(wp) :: zintb, zintn ! dummy argument
- REAL(wp) :: zfac_x, zfac_y
-
- REAL(wp), POINTER, DIMENSION(:,:) :: zpresh ! temporary array for ice strength
- REAL(wp), POINTER, DIMENSION(:,:) :: z1_e1t0, z1_e2t0 ! scale factors
- REAL(wp), POINTER, DIMENSION(:,:) :: zp_delt ! P/delta at T points
- !
- REAL(wp), POINTER, DIMENSION(:,:) :: zaU , zaV ! ice fraction on U/V points
- REAL(wp), POINTER, DIMENSION(:,:) :: zmU_t, zmV_t ! ice/snow mass/dt on U/V points
- REAL(wp), POINTER, DIMENSION(:,:) :: zmf ! coriolis parameter at T points
- REAL(wp), POINTER, DIMENSION(:,:) :: zTauU_ia , ztauV_ia ! ice-atm. stress at U-V points
- REAL(wp), POINTER, DIMENSION(:,:) :: zspgU , zspgV ! surface pressure gradient at U/V points
- REAL(wp), POINTER, DIMENSION(:,:) :: v_oceU, u_oceV, v_iceU, u_iceV ! ocean/ice u/v component on V/U points
- REAL(wp), POINTER, DIMENSION(:,:) :: zfU , zfV ! internal stresses
-
- REAL(wp), POINTER, DIMENSION(:,:) :: zds ! shear
- REAL(wp), POINTER, DIMENSION(:,:) :: zs1, zs2, zs12 ! stress tensor components
- REAL(wp), POINTER, DIMENSION(:,:) :: zu_ice, zv_ice, zresr ! check convergence
- REAL(wp), POINTER, DIMENSION(:,:) :: zpice ! array used for the calculation of ice surface slope:
- ! ocean surface (ssh_m) if ice is not embedded
- ! ice top surface if ice is embedded
- REAL(wp), POINTER, DIMENSION(:,:) :: zCorx, zCory ! Coriolis stress array
- REAL(wp), POINTER, DIMENSION(:,:) :: ztaux_oi, ztauy_oi ! Ocean-to-ice stress array
- REAL(wp), POINTER, DIMENSION(:,:) :: zswitchU, zswitchV ! dummy arrays
- REAL(wp), POINTER, DIMENSION(:,:) :: zmaskU, zmaskV ! mask for ice presence
- REAL(wp), POINTER, DIMENSION(:,:) :: zfmask, zwf ! mask at F points for the ice
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_xmtrp_ice ! X-component of ice mass transport (kg/s)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_ymtrp_ice ! Y-component of ice mass transport (kg/s)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_xmtrp_snw ! X-component of snow mass transport (kg/s)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_ymtrp_snw ! Y-component of snow mass transport (kg/s)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_xatrp ! X-component of area transport (m2/s)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_yatrp ! Y-component of area transport (m2/s)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_utau_oi ! X-direction ocean-ice stress
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_vtau_oi ! Y-direction ocean-ice stress
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_dssh_dx ! X-direction sea-surface tilt term (N/m2)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_dssh_dy ! X-direction sea-surface tilt term (N/m2)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_corstrx ! X-direction coriolis stress (N/m2)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_corstry ! Y-direction coriolis stress (N/m2)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_intstrx ! X-direction internal stress (N/m2)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_intstry ! Y-direction internal stress (N/m2)
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_sig1 ! Average normal stress in sea ice
- REAL(wp), POINTER, DIMENSION(:,:) :: zdiag_sig2 ! Maximum shear stress in sea ice
- REAL(wp), POINTER, DIMENSION(:,:) :: zswi, zmiss ! Switch & missing value array
- REAL(wp), PARAMETER :: zepsi = 1.0e-20_wp ! tolerance parameter
- REAL(wp), PARAMETER :: zmmin = 1._wp ! ice mass (kg/m2) below which ice velocity equals ocean velocity
- REAL(wp), PARAMETER :: zshlat = 2._wp ! boundary condition for sea-ice velocity (2=no slip ; 0=free slip)
- REAL(wp), PARAMETER :: zmiss_val = 1.0e+20 ! missing value for outputs
- !!-------------------------------------------------------------------
- CALL wrk_alloc( jpi,jpj, zpresh, z1_e1t0, z1_e2t0, zp_delt )
- CALL wrk_alloc( jpi,jpj, zaU, zaV, zmU_t, zmV_t, zmf, zTauU_ia, ztauV_ia )
- CALL wrk_alloc( jpi,jpj, zspgU, zspgV, v_oceU, u_oceV, v_iceU, u_iceV, zfU, zfV )
- CALL wrk_alloc( jpi,jpj, zds, zs1, zs2, zs12, zu_ice, zv_ice, zresr, zpice )
- CALL wrk_alloc( jpi,jpj, zswitchU, zswitchV, zmaskU, zmaskV, zfmask, zwf )
- CALL wrk_alloc( jpi,jpj, zCorx, zCory)
- CALL wrk_alloc( jpi,jpj, ztaux_oi, ztauy_oi)
- CALL wrk_alloc( jpi,jpj, zdiag_xmtrp_ice, zdiag_ymtrp_ice )
- CALL wrk_alloc( jpi,jpj, zdiag_xmtrp_snw, zdiag_ymtrp_snw )
- CALL wrk_alloc( jpi,jpj, zdiag_xatrp , zdiag_yatrp )
- CALL wrk_alloc( jpi,jpj, zdiag_utau_oi , zdiag_vtau_oi )
- CALL wrk_alloc( jpi,jpj, zdiag_dssh_dx , zdiag_dssh_dy )
- CALL wrk_alloc( jpi,jpj, zdiag_corstrx , zdiag_corstry )
- CALL wrk_alloc( jpi,jpj, zdiag_intstrx , zdiag_intstry )
- CALL wrk_alloc( jpi,jpj, zdiag_sig1 , zdiag_sig2 )
- CALL wrk_alloc( jpi,jpj, zswi , zmiss )
- #if defined key_lim2 && ! defined key_lim2_vp
- # if defined key_agrif
- USE ice_2, vt_s => hsnm
- USE ice_2, vt_i => hicm
- # else
- vt_s => hsnm
- vt_i => hicm
- # endif
- at_i(:,:) = 1. - frld(:,:)
- #endif
- #if defined key_agrif && defined key_lim2
- CALL agrif_rhg_lim2_load ! First interpolation of coarse values
- #endif
- !
- !------------------------------------------------------------------------------!
- ! 0) mask at F points for the ice (on the whole domain, not only k_j1,k_jpj)
- !------------------------------------------------------------------------------!
- ! ocean/land mask
- DO jj = 1, jpjm1
- DO ji = 1, jpim1 ! NO vector opt.
- zfmask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1)
- END DO
- END DO
- CALL lbc_lnk( zfmask, 'F', 1._wp )
- ! Lateral boundary conditions on velocity (modify zfmask)
- zwf(:,:) = zfmask(:,:)
- DO jj = 2, jpjm1
- DO ji = fs_2, fs_jpim1 ! vector opt.
- IF( zfmask(ji,jj) == 0._wp ) THEN
- zfmask(ji,jj) = zshlat * MIN( 1._wp , MAX( zwf(ji+1,jj), zwf(ji,jj+1), zwf(ji-1,jj), zwf(ji,jj-1) ) )
- ENDIF
- END DO
- END DO
- DO jj = 2, jpjm1
- IF( zfmask(1,jj) == 0._wp ) THEN
- zfmask(1 ,jj) = zshlat * MIN( 1._wp , MAX( zwf(2,jj), zwf(1,jj+1), zwf(1,jj-1) ) )
- ENDIF
- IF( zfmask(jpi,jj) == 0._wp ) THEN
- zfmask(jpi,jj) = zshlat * MIN( 1._wp , MAX( zwf(jpi,jj+1), zwf(jpim1,jj), zwf(jpi,jj-1) ) )
- ENDIF
- END DO
- DO ji = 2, jpim1
- IF( zfmask(ji,1) == 0._wp ) THEN
- zfmask(ji,1 ) = zshlat * MIN( 1._wp , MAX( zwf(ji+1,1), zwf(ji,2), zwf(ji-1,1) ) )
- ENDIF
- IF( zfmask(ji,jpj) == 0._wp ) THEN
- zfmask(ji,jpj) = zshlat * MIN( 1._wp , MAX( zwf(ji+1,jpj), zwf(ji-1,jpj), zwf(ji,jpjm1) ) )
- ENDIF
- END DO
- CALL lbc_lnk( zfmask, 'F', 1._wp )
- !------------------------------------------------------------------------------!
- ! 1) define some variables and initialize arrays
- !------------------------------------------------------------------------------!
- ! ecc2: square of yield ellipse eccenticrity
- ecc2 = rn_ecc * rn_ecc
- z1_ecc2 = 1._wp / ecc2
- ! Time step for subcycling
- zdtevp = rdt_ice / REAL( nn_nevp )
- z1_dtevp = 1._wp / zdtevp
- ! alpha parameters (Bouillon 2009)
- #if defined key_lim3
- zalph1 = ( 2._wp * rn_relast * rdt_ice ) * z1_dtevp
- #else
- zalph1 = ( 2._wp * telast ) * z1_dtevp
- #endif
- zalph2 = zalph1 * z1_ecc2
- z1_alph1 = 1._wp / ( zalph1 + 1._wp )
- z1_alph2 = 1._wp / ( zalph2 + 1._wp )
- ! Initialise stress tensor
- zs1 (:,:) = stress1_i (:,:)
- zs2 (:,:) = stress2_i (:,:)
- zs12(:,:) = stress12_i(:,:)
- ! Ice strength
- #if defined key_lim3
- CALL lim_itd_me_icestrength( nn_icestr )
- zpresh(:,:) = tmask(:,:,1) * strength(:,:)
- #else
- zpresh(:,:) = tmask(:,:,1) * pstar * vt_i(:,:) * EXP( -c_rhg * (1. - at_i(:,:) ) )
- #endif
- ! scale factors
- DO jj = k_j1+1, k_jpj-1
- DO ji = fs_2, fs_jpim1
- z1_e1t0(ji,jj) = 1._wp / ( e1t(ji+1,jj ) + e1t(ji,jj ) )
- z1_e2t0(ji,jj) = 1._wp / ( e2t(ji ,jj+1) + e2t(ji,jj ) )
- END DO
- END DO
-
- !
- !------------------------------------------------------------------------------!
- ! 2) Wind / ocean stress, mass terms, coriolis terms
- !------------------------------------------------------------------------------!
- IF( nn_ice_embd == 2 ) THEN !== embedded sea ice: compute representative ice top surface ==!
- !
- ! average interpolation coeff as used in dynspg = (1/nn_fsbc) * {SUM[n/nn_fsbc], n=0,nn_fsbc-1}
- ! = (1/nn_fsbc)^2 * {SUM[n], n=0,nn_fsbc-1}
- zintn = REAL( nn_fsbc - 1 ) / REAL( nn_fsbc ) * 0.5_wp
- !
- ! average interpolation coeff as used in dynspg = (1/nn_fsbc) * {SUM[1-n/nn_fsbc], n=0,nn_fsbc-1}
- ! = (1/nn_fsbc)^2 * (nn_fsbc^2 - {SUM[n], n=0,nn_fsbc-1})
- zintb = REAL( nn_fsbc + 1 ) / REAL( nn_fsbc ) * 0.5_wp
- !
- zpice(:,:) = ssh_m(:,:) + ( zintn * snwice_mass(:,:) + zintb * snwice_mass_b(:,:) ) * r1_rau0
- !
- ELSE !== non-embedded sea ice: use ocean surface for slope calculation ==!
- zpice(:,:) = ssh_m(:,:)
- ENDIF
- DO jj = k_j1+1, k_jpj-1
- DO ji = fs_2, fs_jpim1
- ! ice fraction at U-V points
- zaU(ji,jj) = 0.5_wp * ( at_i(ji,jj) * e12t(ji,jj) + at_i(ji+1,jj) * e12t(ji+1,jj) ) * r1_e12u(ji,jj) * umask(ji,jj,1)
- zaV(ji,jj) = 0.5_wp * ( at_i(ji,jj) * e12t(ji,jj) + at_i(ji,jj+1) * e12t(ji,jj+1) ) * r1_e12v(ji,jj) * vmask(ji,jj,1)
- ! Ice/snow mass at U-V points
- zm1 = ( rhosn * vt_s(ji ,jj ) + rhoic * vt_i(ji ,jj ) )
- zm2 = ( rhosn * vt_s(ji+1,jj ) + rhoic * vt_i(ji+1,jj ) )
- zm3 = ( rhosn * vt_s(ji ,jj+1) + rhoic * vt_i(ji ,jj+1) )
- zmassU = 0.5_wp * ( zm1 * e12t(ji,jj) + zm2 * e12t(ji+1,jj) ) * r1_e12u(ji,jj) * umask(ji,jj,1)
- zmassV = 0.5_wp * ( zm1 * e12t(ji,jj) + zm3 * e12t(ji,jj+1) ) * r1_e12v(ji,jj) * vmask(ji,jj,1)
- ! Ocean currents at U-V points
- v_oceU(ji,jj) = 0.5_wp * ( ( v_oce(ji ,jj) + v_oce(ji ,jj-1) ) * e1t(ji+1,jj) &
- & + ( v_oce(ji+1,jj) + v_oce(ji+1,jj-1) ) * e1t(ji ,jj) ) * z1_e1t0(ji,jj) * umask(ji,jj,1)
-
- u_oceV(ji,jj) = 0.5_wp * ( ( u_oce(ji,jj ) + u_oce(ji-1,jj ) ) * e2t(ji,jj+1) &
- & + ( u_oce(ji,jj+1) + u_oce(ji-1,jj+1) ) * e2t(ji,jj ) ) * z1_e2t0(ji,jj) * vmask(ji,jj,1)
- ! Coriolis at T points (m*f)
- zmf(ji,jj) = zm1 * fcor(ji,jj)
- ! m/dt
- zmU_t(ji,jj) = zmassU * z1_dtevp
- zmV_t(ji,jj) = zmassV * z1_dtevp
- ! Drag ice-atm.
- zTauU_ia(ji,jj) = zaU(ji,jj) * utau_ice(ji,jj)
- zTauV_ia(ji,jj) = zaV(ji,jj) * vtau_ice(ji,jj)
- ! Surface pressure gradient (- m*g*GRAD(ssh)) at U-V points
- zspgU(ji,jj) = - zmassU * grav * ( zpice(ji+1,jj) - zpice(ji,jj) ) * r1_e1u(ji,jj)
- zspgV(ji,jj) = - zmassV * grav * ( zpice(ji,jj+1) - zpice(ji,jj) ) * r1_e2v(ji,jj)
- ! masks
- zmaskU(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassU ) ) ! 0 if no ice
- zmaskV(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassV ) ) ! 0 if no ice
- ! switches
- zswitchU(ji,jj) = MAX( 0._wp, SIGN( 1._wp, zmassU - zmmin ) ) ! 0 if ice mass < zmmin
- zswitchV(ji,jj) = MAX( 0._wp, SIGN( 1._wp, zmassV - zmmin ) ) ! 0 if ice mass < zmmin
- END DO
- END DO
- CALL lbc_lnk( zmf, 'T', 1. )
- !
- !------------------------------------------------------------------------------!
- ! 3) Solution of the momentum equation, iterative procedure
- !------------------------------------------------------------------------------!
- !
- ! !----------------------!
- DO jter = 1 , nn_nevp ! loop over jter !
- ! !----------------------!
- IF(ln_ctl) THEN ! Convergence test
- DO jj = k_j1, k_jpj-1
- zu_ice(:,jj) = u_ice(:,jj) ! velocity at previous time step
- zv_ice(:,jj) = v_ice(:,jj)
- END DO
- ENDIF
- ! --- divergence, tension & shear (Appendix B of Hunke & Dukowicz, 2002) --- !
- DO jj = k_j1, k_jpj-1 ! loops start at 1 since there is no boundary condition (lbc_lnk) at i=1 and j=1 for F points
- DO ji = 1, jpim1
- ! shear at F points
- zds(ji,jj) = ( ( u_ice(ji,jj+1) * r1_e1u(ji,jj+1) - u_ice(ji,jj) * r1_e1u(ji,jj) ) * e1f(ji,jj) * e1f(ji,jj) &
- & + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj) &
- & ) * r1_e12f(ji,jj) * zfmask(ji,jj)
- END DO
- END DO
- CALL lbc_lnk( zds, 'F', 1. )
- DO jj = k_j1+1, k_jpj-1
- DO ji = 2, jpim1 ! no vector loop
- ! shear**2 at T points (doc eq. A16)
- zds2 = ( zds(ji,jj ) * zds(ji,jj ) * e12f(ji,jj ) + zds(ji-1,jj ) * zds(ji-1,jj ) * e12f(ji-1,jj ) &
- & + zds(ji,jj-1) * zds(ji,jj-1) * e12f(ji,jj-1) + zds(ji-1,jj-1) * zds(ji-1,jj-1) * e12f(ji-1,jj-1) &
- & ) * 0.25_wp * r1_e12t(ji,jj)
-
- ! divergence at T points
- zdiv = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj) &
- & + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1) &
- & ) * r1_e12t(ji,jj)
- zdiv2 = zdiv * zdiv
-
- ! tension at T points
- zdt = ( ( u_ice(ji,jj) * r1_e2u(ji,jj) - u_ice(ji-1,jj) * r1_e2u(ji-1,jj) ) * e2t(ji,jj) * e2t(ji,jj) &
- & - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj) &
- & ) * r1_e12t(ji,jj)
- zdt2 = zdt * zdt
-
- ! delta at T points
- zdelta = SQRT( zdiv2 + ( zdt2 + zds2 ) * usecc2 )
- ! P/delta at T points
- zp_delt(ji,jj) = zpresh(ji,jj) / ( zdelta + rn_creepl )
-
- ! stress at T points
- zs1(ji,jj) = ( zs1(ji,jj) * zalph1 + zp_delt(ji,jj) * ( zdiv - zdelta ) ) * z1_alph1
- zs2(ji,jj) = ( zs2(ji,jj) * zalph2 + zp_delt(ji,jj) * ( zdt * z1_ecc2 ) ) * z1_alph2
-
- END DO
- END DO
- CALL lbc_lnk( zp_delt, 'T', 1. )
- DO jj = k_j1, k_jpj-1
- DO ji = 1, jpim1
- ! P/delta at F points
- zp_delf = 0.25_wp * ( zp_delt(ji,jj) + zp_delt(ji+1,jj) + zp_delt(ji,jj+1) + zp_delt(ji+1,jj+1) )
-
- ! stress at F points
- zs12(ji,jj)= ( zs12(ji,jj) * zalph2 + zp_delf * ( zds(ji,jj) * z1_ecc2 ) * 0.5_wp ) * z1_alph2
- END DO
- END DO
- CALL lbc_lnk_multi( zs1, 'T', 1., zs2, 'T', 1., zs12, 'F', 1. )
-
- ! --- Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) --- !
- DO jj = k_j1+1, k_jpj-1
- DO ji = fs_2, fs_jpim1
- ! U points
- zfU(ji,jj) = 0.5_wp * ( ( zs1(ji+1,jj) - zs1(ji,jj) ) * e2u(ji,jj) &
- & + ( zs2(ji+1,jj) * e2t(ji+1,jj) * e2t(ji+1,jj) - zs2(ji,jj) * e2t(ji,jj) * e2t(ji,jj) &
- & ) * r1_e2u(ji,jj) &
- & + ( zs12(ji,jj) * e1f(ji,jj) * e1f(ji,jj) - zs12(ji,jj-1) * e1f(ji,jj-1) * e1f(ji,jj-1) &
- & ) * 2._wp * r1_e1u(ji,jj) &
- & ) * r1_e12u(ji,jj)
- ! V points
- zfV(ji,jj) = 0.5_wp * ( ( zs1(ji,jj+1) - zs1(ji,jj) ) * e1v(ji,jj) &
- & - ( zs2(ji,jj+1) * e1t(ji,jj+1) * e1t(ji,jj+1) - zs2(ji,jj) * e1t(ji,jj) * e1t(ji,jj) &
- & ) * r1_e1v(ji,jj) &
- & + ( zs12(ji,jj) * e2f(ji,jj) * e2f(ji,jj) - zs12(ji-1,jj) * e2f(ji-1,jj) * e2f(ji-1,jj) &
- & ) * 2._wp * r1_e2v(ji,jj) &
- & ) * r1_e12v(ji,jj)
- ! u_ice at V point
- u_iceV(ji,jj) = 0.5_wp * ( ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * e2t(ji,jj+1) &
- & + ( u_ice(ji,jj+1) + u_ice(ji-1,jj+1) ) * e2t(ji,jj ) ) * z1_e2t0(ji,jj) * vmask(ji,jj,1)
-
- ! v_ice at U point
- v_iceU(ji,jj) = 0.5_wp * ( ( v_ice(ji ,jj) + v_ice(ji ,jj-1) ) * e1t(ji+1,jj) &
- & + ( v_ice(ji+1,jj) + v_ice(ji+1,jj-1) ) * e1t(ji ,jj) ) * z1_e1t0(ji,jj) * umask(ji,jj,1)
- END DO
- END DO
- !
- ! --- Computation of ice velocity --- !
- ! Bouillon et al. 2013 (eq 47-48) => unstable unless alpha, beta are chosen wisely and large nn_nevp
- ! Bouillon et al. 2009 (eq 34-35) => stable
- IF( MOD(jter,2) .EQ. 0 ) THEN ! even iterations
-
- DO jj = k_j1+1, k_jpj-1
- DO ji = fs_2, fs_jpim1
- ! tau_io/(v_oce - v_ice)
- zTauO = zaV(ji,jj) * rhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) ) &
- & + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
- ! Ocean-to-Ice stress
- ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
- ! Coriolis at V-points (energy conserving formulation)
- zCory(ji,jj) = - 0.25_wp * r1_e2v(ji,jj) * &
- & ( zmf(ji,jj ) * ( e2u(ji,jj ) * u_ice(ji,jj ) + e2u(ji-1,jj ) * u_ice(ji-1,jj ) ) &
- & + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
- ! Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
- zTauE = zfV(ji,jj) + zTauV_ia(ji,jj) + zCory(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
-
- ! ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
- v_ice(ji,jj) = ( ( zmV_t(ji,jj) * v_ice(ji,jj) + zTauE + zTauO * v_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
- & ) / MAX( zepsi, zmV_t(ji,jj) + zTauO ) * zswitchV(ji,jj) & ! m/dt + tau_io(only ice part)
- & + v_oce(ji,jj) * ( 1._wp - zswitchV(ji,jj) ) & ! v_ice = v_oce if mass < zmmin
- & ) * zmaskV(ji,jj)
- END DO
- END DO
- CALL lbc_lnk( v_ice, 'V', -1. )
- #if defined key_agrif && defined key_lim2
- CALL agrif_rhg_lim2( jter, nn_nevp, 'V' )
- #endif
- #if defined key_bdy
- CALL bdy_ice_lim_dyn( 'V' )
- #endif
- DO jj = k_j1+1, k_jpj-1
- DO ji = fs_2, fs_jpim1
-
- ! tau_io/(u_oce - u_ice)
- zTauO = zaU(ji,jj) * rhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) ) &
- & + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
- ! Ocean-to-Ice stress
- ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
- ! Coriolis at U-points (energy conserving formulation)
- zCorx(ji,jj) = 0.25_wp * r1_e1u(ji,jj) * &
- & ( zmf(ji ,jj) * ( e1v(ji ,jj) * v_ice(ji ,jj) + e1v(ji ,jj-1) * v_ice(ji ,jj-1) ) &
- & + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
- ! Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
- zTauE = zfU(ji,jj) + zTauU_ia(ji,jj) + zCorx(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
- ! ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
- u_ice(ji,jj) = ( ( zmU_t(ji,jj) * u_ice(ji,jj) + zTauE + zTauO * u_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
- & ) / MAX( zepsi, zmU_t(ji,jj) + zTauO ) * zswitchU(ji,jj) & ! m/dt + tau_io(only ice part)
- & + u_oce(ji,jj) * ( 1._wp - zswitchU(ji,jj) ) & ! v_ice = v_oce if mass < zmmin
- & ) * zmaskU(ji,jj)
- END DO
- END DO
- CALL lbc_lnk( u_ice, 'U', -1. )
-
- #if defined key_agrif && defined key_lim2
- CALL agrif_rhg_lim2( jter, nn_nevp, 'U' )
- #endif
- #if defined key_bdy
- CALL bdy_ice_lim_dyn( 'U' )
- #endif
- ELSE ! odd iterations
- DO jj = k_j1+1, k_jpj-1
- DO ji = fs_2, fs_jpim1
-
- ! tau_io/(u_oce - u_ice)
- zTauO = zaU(ji,jj) * rhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) ) &
- & + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
- ! Ocean-to-Ice stress
- ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
-
- ! Coriolis at U-points (energy conserving formulation)
- zCorx(ji,jj) = 0.25_wp * r1_e1u(ji,jj) * &
- & ( zmf(ji ,jj) * ( e1v(ji ,jj) * v_ice(ji ,jj) + e1v(ji ,jj-1) * v_ice(ji ,jj-1) ) &
- & + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
- ! Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
- zTauE = zfU(ji,jj) + zTauU_ia(ji,jj) + zCorx(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
- ! ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
- u_ice(ji,jj) = ( ( zmU_t(ji,jj) * u_ice(ji,jj) + zTauE + zTauO * u_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
- & ) / MAX( zepsi, zmU_t(ji,jj) + zTauO ) * zswitchU(ji,jj) & ! m/dt + tau_io(only ice part)
- & + u_oce(ji,jj) * ( 1._wp - zswitchU(ji,jj) ) & ! v_ice = v_oce if mass < zmmin
- & ) * zmaskU(ji,jj)
- END DO
- END DO
- CALL lbc_lnk( u_ice, 'U', -1. )
-
- #if defined key_agrif && defined key_lim2
- CALL agrif_rhg_lim2( jter, nn_nevp, 'U' )
- #endif
- #if defined key_bdy
- CALL bdy_ice_lim_dyn( 'U' )
- #endif
- DO jj = k_j1+1, k_jpj-1
- DO ji = fs_2, fs_jpim1
- ! tau_io/(v_oce - v_ice)
- zTauO = zaV(ji,jj) * rhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) ) &
- & + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
- ! Ocean-to-Ice stress
- ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
- ! Coriolis at V-points (energy conserving formulation)
- zCory(ji,jj) = - 0.25_wp * r1_e2v(ji,jj) * &
- & ( zmf(ji,jj ) * ( e2u(ji,jj ) * u_ice(ji,jj ) + e2u(ji-1,jj ) * u_ice(ji-1,jj ) ) &
- & + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
- ! Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
- zTauE = zfV(ji,jj) + zTauV_ia(ji,jj) + zCory(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
-
- ! ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
- v_ice(ji,jj) = ( ( zmV_t(ji,jj) * v_ice(ji,jj) + zTauE + zTauO * v_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
- & ) / MAX( zepsi, zmV_t(ji,jj) + zTauO ) * zswitchV(ji,jj) & ! m/dt + tau_io(only ice part)
- & + v_oce(ji,jj) * ( 1._wp - zswitchV(ji,jj) ) & ! v_ice = v_oce if mass < zmmin
- & ) * zmaskV(ji,jj)
- END DO
- END DO
- CALL lbc_lnk( v_ice, 'V', -1. )
-
- #if defined key_agrif && defined key_lim2
- CALL agrif_rhg_lim2( jter, nn_nevp, 'V' )
- #endif
- #if defined key_bdy
- CALL bdy_ice_lim_dyn( 'V' )
- #endif
- ENDIF
- IF(ln_ctl) THEN ! Convergence test
- DO jj = k_j1+1, k_jpj-1
- zresr(:,jj) = MAX( ABS( u_ice(:,jj) - zu_ice(:,jj) ), ABS( v_ice(:,jj) - zv_ice(:,jj) ) )
- END DO
- zresm = MAXVAL( zresr( 1:jpi, k_j1+1:k_jpj-1 ) )
- IF( lk_mpp ) CALL mpp_max( zresm ) ! max over the global domain
- ENDIF
- !
- ! ! ==================== !
- END DO ! end loop over jter !
- ! ! ==================== !
- !
- !------------------------------------------------------------------------------!
- ! 4) Recompute delta, shear and div (inputs for mechanical redistribution)
- !------------------------------------------------------------------------------!
- DO jj = k_j1, k_jpj-1
- DO ji = 1, jpim1
- ! shear at F points
- zds(ji,jj) = ( ( u_ice(ji,jj+1) * r1_e1u(ji,jj+1) - u_ice(ji,jj) * r1_e1u(ji,jj) ) * e1f(ji,jj) * e1f(ji,jj) &
- & + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj) &
- & ) * r1_e12f(ji,jj) * zfmask(ji,jj)
- END DO
- END DO
- CALL lbc_lnk( zds, 'F', 1. )
-
- DO jj = k_j1+1, k_jpj-1
- DO ji = 2, jpim1 ! no vector loop
-
- ! tension**2 at T points
- zdt = ( ( u_ice(ji,jj) * r1_e2u(ji,jj) - u_ice(ji-1,jj) * r1_e2u(ji-1,jj) ) * e2t(ji,jj) * e2t(ji,jj) &
- & - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj) &
- & ) * r1_e12t(ji,jj)
- zdt2 = zdt * zdt
-
- ! shear**2 at T points (doc eq. A16)
- zds2 = ( zds(ji,jj ) * zds(ji,jj ) * e12f(ji,jj ) + zds(ji-1,jj ) * zds(ji-1,jj ) * e12f(ji-1,jj ) &
- & + zds(ji,jj-1) * zds(ji,jj-1) * e12f(ji,jj-1) + zds(ji-1,jj-1) * zds(ji-1,jj-1) * e12f(ji-1,jj-1) &
- & ) * 0.25_wp * r1_e12t(ji,jj)
-
- ! shear at T points
- shear_i(ji,jj) = SQRT( zdt2 + zds2 )
- ! divergence at T points
- divu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj) &
- & + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1) &
- & ) * r1_e12t(ji,jj)
-
- ! delta at T points
- zdelta = SQRT( divu_i(ji,jj) * divu_i(ji,jj) + ( zdt2 + zds2 ) * usecc2 )
- rswitch = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zdelta ) ) ! 0 if delta=0
- delta_i(ji,jj) = zdelta + rn_creepl * rswitch
- END DO
- END DO
- CALL lbc_lnk_multi( shear_i, 'T', 1., divu_i, 'T', 1., delta_i, 'T', 1. )
-
- ! --- Store the stress tensor for the next time step --- !
- stress1_i (:,:) = zs1 (:,:)
- stress2_i (:,:) = zs2 (:,:)
- stress12_i(:,:) = zs12(:,:)
- !------------------------------------------------------------------------------!
- ! 5) SIMIP diagnostics
- !------------------------------------------------------------------------------!
- DO jj = 1, jpj
- DO ji = 1, jpi
- zswi(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) ) ! 1 if ice, 0 if no ice
- END DO
- END DO
- zmiss(:,:) = zmiss_val * ( 1. - zswi(:,:) )
-
- DO jj = k_j1+1, k_jpj-1
- DO ji = 2, jpim1
- ! Stress tensor invariants (normal and shear stress N/m)
- zdiag_sig1(ji,jj) = ( zs1(ji,jj) + zs2(ji,jj) ) * zswi(ji,jj) ! normal stress
- zdiag_sig2(ji,jj) = SQRT( ( zs1(ji,jj) - zs2(ji,jj) )**2 + 4*zs12(ji,jj)**2 ) * zswi(ji,jj) ! shear stress
- ! Stress terms of the momentum equation (N/m2)
- zdiag_dssh_dx(ji,jj) = zspgU(ji,jj) * zswi(ji,jj) ! sea surface slope stress term
- zdiag_dssh_dy(ji,jj) = zspgV(ji,jj) * zswi(ji,jj)
- zdiag_corstrx(ji,jj) = zCorx(ji,jj) * zswi(ji,jj) ! Coriolis stress term
- zdiag_corstry(ji,jj) = zCory(ji,jj) * zswi(ji,jj)
- zdiag_intstrx(ji,jj) = zfU(ji,jj) * zswi(ji,jj) ! internal stress term
- zdiag_intstry(ji,jj) = zfV(ji,jj) * zswi(ji,jj)
-
- zdiag_utau_oi(ji,jj) = ztaux_oi(ji,jj) * zswi(ji,jj) ! oceanic stress
- zdiag_vtau_oi(ji,jj) = ztauy_oi(ji,jj) * zswi(ji,jj)
- ! 2D ice mass, snow mass, area transport arrays (X, Y)
- zfac_x = 0.5 * u_ice(ji,jj) * e2u(ji,jj) * zswi(ji,jj)
- zfac_y = 0.5 * v_ice(ji,jj) * e1v(ji,jj) * zswi(ji,jj)
- zdiag_xmtrp_ice(ji,jj) = rhoic * zfac_x * ( vt_i(ji+1,jj) + vt_i(ji,jj) ) ! ice mass transport, X-component (kg/s)
- zdiag_ymtrp_ice(ji,jj) = rhoic * zfac_y * ( vt_i(ji,jj+1) + vt_i(ji,jj) ) ! '' Y- ''
- zdiag_xmtrp_snw(ji,jj) = rhosn * zfac_x * ( vt_s(ji+1,jj) + vt_s(ji,jj) ) ! snow mass transport, X-component
- zdiag_ymtrp_snw(ji,jj) = rhosn * zfac_y * ( vt_s(ji,jj+1) + vt_s(ji,jj) ) ! '' Y- ''
- zdiag_xatrp(ji,jj) = zfac_x * ( at_i(ji+1,jj) + at_i(ji,jj) ) ! area transport, X-component (m2/s)
- zdiag_yatrp(ji,jj) = zfac_y * ( at_i(ji,jj+1) + at_i(ji,jj) ) ! '' Y- ''
- END DO
- END DO
- CALL lbc_lnk_multi( zdiag_sig1 , 'T', 1., zdiag_sig2 , 'T', 1., &
- & zdiag_dssh_dx, 'U', -1., zdiag_dssh_dy, 'V', -1., &
- & zdiag_corstrx, 'U', -1., zdiag_corstry, 'V', -1., &
- & zdiag_intstrx, 'U', -1., zdiag_intstry, 'V', -1. )
- CALL lbc_lnk_multi( zdiag_utau_oi, 'U', -1., zdiag_vtau_oi, 'V', -1. )
- CALL lbc_lnk_multi( zdiag_xmtrp_ice, 'U', -1., zdiag_xmtrp_snw, 'U', -1., &
- & zdiag_xatrp , 'U', -1., zdiag_ymtrp_ice, 'V', -1., &
- & zdiag_ymtrp_snw, 'V', -1., zdiag_yatrp , 'V', -1. )
- IF ( iom_use( "xmtrpice" ) ) CALL iom_put( "xmtrpice" , zdiag_xmtrp_ice(:,:) ) ! X-component of sea-ice mass transport (kg/s)
- IF ( iom_use( "ymtrpice" ) ) CALL iom_put( "ymtrpice" , zdiag_ymtrp_ice(:,:) ) ! Y-component of sea-ice mass transport
- IF ( iom_use( "xmtrpsnw" ) ) CALL iom_put( "xmtrpsnw" , zdiag_xmtrp_snw(:,:) ) ! X-component of snow mass transport (kg/s)
- IF ( iom_use( "ymtrpsnw" ) ) CALL iom_put( "ymtrpsnw" , zdiag_ymtrp_snw(:,:) ) ! Y-component of snow mass transport
- IF ( iom_use( "xatrp" ) ) CALL iom_put( "xatrp" , zdiag_xatrp(:,:) ) ! X-component of ice area transport
- IF ( iom_use( "yatrp" ) ) CALL iom_put( "yatrp" , zdiag_yatrp(:,:) ) ! Y-component of ice area transport
- IF ( iom_use( "utau_ice" ) ) CALL iom_put( "utau_ice" , utau_ice(:,:) * zswi(:,:) + zmiss(:,:) ) ! Wind stress term in force balance (x)
- IF ( iom_use( "vtau_ice" ) ) CALL iom_put( "vtau_ice" , vtau_ice(:,:) * zswi(:,:) + zmiss(:,:) ) ! Wind stress term in force balance (y)
- IF ( iom_use( "utau_oi" ) ) CALL iom_put( "utau_oi" , zdiag_utau_oi(:,:) * zswi(:,:) + zmiss(:,:) ) ! Ocean stress term in force balance (x)
- IF ( iom_use( "vtau_oi" ) ) CALL iom_put( "vtau_oi" , zdiag_vtau_oi(:,:) * zswi(:,:) + zmiss(:,:) ) ! Ocean stress term in force balance (y)
- IF ( iom_use( "dssh_dx" ) ) CALL iom_put( "dssh_dx" , zdiag_dssh_dx(:,:) * zswi(:,:) + zmiss(:,:) ) ! Sea-surface tilt term in force balance (x)
- IF ( iom_use( "dssh_dy" ) ) CALL iom_put( "dssh_dy" , zdiag_dssh_dy(:,:) * zswi(:,:) + zmiss(:,:) ) ! Sea-surface tilt term in force balance (y)
- IF ( iom_use( "corstrx" ) ) CALL iom_put( "corstrx" , zdiag_corstrx(:,:) * zswi(:,:) + zmiss(:,:) ) ! Coriolis force term in force balance (x)
- IF ( iom_use( "corstry" ) ) CALL iom_put( "corstry" , zdiag_corstry(:,:) * zswi(:,:) + zmiss(:,:) ) ! Coriolis force term in force balance (y)
- IF ( iom_use( "intstrx" ) ) CALL iom_put( "intstrx" , zdiag_intstrx(:,:) * zswi(:,:) + zmiss(:,:) ) ! Internal force term in force balance (x)
- IF ( iom_use( "intstry" ) ) CALL iom_put( "intstry" , zdiag_intstry(:,:) * zswi(:,:) + zmiss(:,:) ) ! Internal force term in force balance (y)
- IF ( iom_use( "normstr" ) ) CALL iom_put( "normstr" , zdiag_sig1(:,:) * zswi(:,:) + zmiss(:,:) ) ! Normal stress
- IF ( iom_use( "sheastr" ) ) CALL iom_put( "sheastr" , zdiag_sig2(:,:) * zswi(:,:) + zmiss(:,:) ) ! Shear stress
- !
- !------------------------------------------------------------------------------!
- ! 6) Control prints of residual and charge ellipse
- !------------------------------------------------------------------------------!
- !
- ! print the residual for convergence
- IF(ln_ctl) THEN
- WRITE(charout,FMT="('lim_rhg : res =',D23.16, ' iter =',I4)") zresm, jter
- CALL prt_ctl_info(charout)
- CALL prt_ctl(tab2d_1=u_ice, clinfo1=' lim_rhg : u_ice :', tab2d_2=v_ice, clinfo2=' v_ice :')
- ENDIF
- ! print charge ellipse
- ! This can be desactivated once the user is sure that the stress state
- ! lie on the charge ellipse. See Bouillon et al. 08 for more details
- IF(ln_ctl) THEN
- CALL prt_ctl_info('lim_rhg : numit :',ivar1=numit)
- CALL prt_ctl_info('lim_rhg : nwrite :',ivar1=nwrite)
- CALL prt_ctl_info('lim_rhg : MOD :',ivar1=MOD(numit,nwrite))
- IF( MOD(numit,nwrite) .EQ. 0 ) THEN
- WRITE(charout,FMT="('lim_rhg :', I4, I6, I1, I1, A10)") 1000, numit, 0, 0, ' ch. ell. '
- CALL prt_ctl_info(charout)
- DO jj = k_j1+1, k_jpj-1
- DO ji = 2, jpim1
- IF (zpresh(ji,jj) > 1.0) THEN
- zsig1 = ( zs1(ji,jj) + (zs2(ji,jj)**2 + 4*zs12(ji,jj)**2 )**0.5 ) / ( 2*zpresh(ji,jj) )
- zsig2 = ( zs1(ji,jj) - (zs2(ji,jj)**2 + 4*zs12(ji,jj)**2 )**0.5 ) / ( 2*zpresh(ji,jj) )
- WRITE(charout,FMT="('lim_rhg :', I4, I4, D23.16, D23.16, D23.16, D23.16, A10)")
- CALL prt_ctl_info(charout)
- ENDIF
- END DO
- END DO
- WRITE(charout,FMT="('lim_rhg :', I4, I6, I1, I1, A10)") 2000, numit, 0, 0, ' ch. ell. '
- CALL prt_ctl_info(charout)
- ENDIF
- ENDIF
- !
- CALL wrk_dealloc( jpi,jpj, zpresh, z1_e1t0, z1_e2t0, zp_delt )
- CALL wrk_dealloc( jpi,jpj, zaU, zaV, zmU_t, zmV_t, zmf, zTauU_ia, ztauV_ia )
- CALL wrk_dealloc( jpi,jpj, zspgU, zspgV, v_oceU, u_oceV, v_iceU, u_iceV, zfU, zfV )
- CALL wrk_dealloc( jpi,jpj, zds, zs1, zs2, zs12, zu_ice, zv_ice, zresr, zpice )
- CALL wrk_dealloc( jpi,jpj, zswitchU, zswitchV, zmaskU, zmaskV, zfmask, zwf )
- CALL wrk_dealloc( jpi,jpj, zCorx, zCory )
- CALL wrk_dealloc( jpi,jpj, ztaux_oi, ztauy_oi )
-
- CALL wrk_dealloc( jpi,jpj, zdiag_xmtrp_ice, zdiag_ymtrp_ice )
- CALL wrk_dealloc( jpi,jpj, zdiag_xmtrp_snw, zdiag_ymtrp_snw )
- CALL wrk_dealloc( jpi,jpj, zdiag_xatrp , zdiag_yatrp )
- CALL wrk_dealloc( jpi,jpj, zdiag_utau_oi , zdiag_vtau_oi )
- CALL wrk_dealloc( jpi,jpj, zdiag_dssh_dx , zdiag_dssh_dy )
- CALL wrk_dealloc( jpi,jpj, zdiag_corstrx , zdiag_corstry )
- CALL wrk_dealloc( jpi,jpj, zdiag_intstrx , zdiag_intstry )
- CALL wrk_dealloc( jpi,jpj, zdiag_sig1 , zdiag_sig2 )
- CALL wrk_dealloc( jpi,jpj, zswi , zmiss )
- END SUBROUTINE lim_rhg
- #else
- !!----------------------------------------------------------------------
- !! Default option Dummy module NO LIM sea-ice model
- !!----------------------------------------------------------------------
- CONTAINS
- SUBROUTINE lim_rhg( k1 , k2 ) ! Dummy routine
- WRITE(*,*) 'lim_rhg: You should not have seen this print! error?', k1, k2
- END SUBROUTINE lim_rhg
- #endif
- !!==============================================================================
- END MODULE limrhg
|