[888] | 1 | MODULE limsbc_2 |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE limsbc_2 *** |
---|
[2528] | 4 | !! LIM-2 : updates the fluxes at the ocean surface with ice-ocean fluxes |
---|
[888] | 5 | !!====================================================================== |
---|
[2528] | 6 | !! History : LIM ! 2000-01 (H. Goosse) Original code |
---|
| 7 | !! 1.0 ! 2002-07 (C. Ethe, G. Madec) re-writing F90 |
---|
| 8 | !! 3.0 ! 2006-07 (G. Madec) surface module |
---|
| 9 | !! 3.3 ! 2009-05 (G. Garric, C. Bricaud) addition of the lim2_evp case |
---|
| 10 | !! - ! 2010-11 (G. Madec) ice-ocean stress computed at each ocean time-step |
---|
[3625] | 11 | !! 3.3.1 ! 2011-01 (A. R. Porter, STFC Daresbury) dynamical allocation |
---|
| 12 | !! 3.5 ! 2012-11 ((G. Madec, Y. Aksenov, A. Coward) salt and heat fluxes associated with e-p |
---|
[888] | 13 | !!---------------------------------------------------------------------- |
---|
| 14 | #if defined key_lim2 |
---|
| 15 | !!---------------------------------------------------------------------- |
---|
| 16 | !! 'key_lim2' LIM 2.0 sea-ice model |
---|
| 17 | !!---------------------------------------------------------------------- |
---|
[2715] | 18 | !! lim_sbc_alloc_2 : allocate the limsbc arrays |
---|
| 19 | !! lim_sbc_init : initialisation |
---|
| 20 | !! lim_sbc_flx_2 : update mass, heat and salt fluxes at the ocean surface |
---|
| 21 | !! lim_sbc_tau_2 : update i- and j-stresses, and its modulus at the ocean surface |
---|
[888] | 22 | !!---------------------------------------------------------------------- |
---|
| 23 | USE par_oce ! ocean parameters |
---|
[2528] | 24 | USE phycst ! physical constants |
---|
[888] | 25 | USE dom_oce ! ocean domain |
---|
[4292] | 26 | USE domvvl ! ocean vertical scale factors |
---|
[2528] | 27 | USE dom_ice_2 ! LIM-2: ice domain |
---|
| 28 | USE ice_2 ! LIM-2: ice variables |
---|
| 29 | USE sbc_ice ! surface boundary condition: ice |
---|
| 30 | USE sbc_oce ! surface boundary condition: ocean |
---|
[6004] | 31 | USE sbccpl ! surface boundary condition: coupled interface |
---|
[3625] | 32 | USE oce , ONLY : sshn, sshb, snwice_mass, snwice_mass_b, snwice_fmass |
---|
[2566] | 33 | USE albedo ! albedo parameters |
---|
[6004] | 34 | ! |
---|
[2528] | 35 | USE lbclnk ! ocean lateral boundary condition - MPP exchanges |
---|
[2715] | 36 | USE lib_mpp ! MPP library |
---|
[3294] | 37 | USE wrk_nemo ! work arrays |
---|
[888] | 38 | USE in_out_manager ! I/O manager |
---|
[2528] | 39 | USE iom ! I/O library |
---|
[888] | 40 | USE prtctl ! Print control |
---|
[3625] | 41 | USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) |
---|
[888] | 42 | |
---|
| 43 | IMPLICIT NONE |
---|
| 44 | PRIVATE |
---|
| 45 | |
---|
[6004] | 46 | PUBLIC lim_sbc_init_2 ! called by ice_init_2 |
---|
| 47 | PUBLIC lim_sbc_flx_2 ! called by sbc_ice_lim_2 |
---|
| 48 | PUBLIC lim_sbc_tau_2 ! called by sbc_ice_lim_2 |
---|
[888] | 49 | |
---|
[2528] | 50 | REAL(wp) :: r1_rdtice ! = 1. / rdt_ice |
---|
| 51 | REAL(wp) :: epsi16 = 1.e-16_wp ! constant values |
---|
| 52 | REAL(wp) :: rzero = 0._wp ! - - |
---|
| 53 | REAL(wp) :: rone = 1._wp ! - - |
---|
| 54 | ! |
---|
[6004] | 55 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: soce_0, sice_0 ! fix SSS and ice salinity used in levitating case 0 |
---|
[2715] | 56 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: utau_oce, vtau_oce ! air-ocean surface i- & j-stress [N/m2] |
---|
| 57 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: tmod_io ! modulus of the ice-ocean relative velocity [m/s] |
---|
[2528] | 58 | |
---|
[888] | 59 | !! * Substitutions |
---|
| 60 | # include "vectopt_loop_substitute.h90" |
---|
| 61 | !!---------------------------------------------------------------------- |
---|
[2715] | 62 | !! NEMO/LIM2 4.0 , UCL - NEMO Consortium (2011) |
---|
[1156] | 63 | !! $Id$ |
---|
[2528] | 64 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[888] | 65 | !!---------------------------------------------------------------------- |
---|
| 66 | CONTAINS |
---|
| 67 | |
---|
[2715] | 68 | INTEGER FUNCTION lim_sbc_alloc_2() |
---|
| 69 | !!------------------------------------------------------------------- |
---|
| 70 | !! *** ROUTINE lim_sbc_alloc_2 *** |
---|
| 71 | !!------------------------------------------------------------------- |
---|
| 72 | ALLOCATE( soce_0(jpi,jpj) , utau_oce(jpi,jpj) , & |
---|
| 73 | & sice_0(jpi,jpj) , vtau_oce(jpi,jpj) , tmod_io(jpi,jpj), STAT=lim_sbc_alloc_2) |
---|
| 74 | ! |
---|
| 75 | IF( lk_mpp ) CALL mpp_sum( lim_sbc_alloc_2 ) |
---|
| 76 | IF( lim_sbc_alloc_2 /= 0 ) CALL ctl_warn('lim_sbc_alloc_2: failed to allocate arrays.') |
---|
| 77 | ! |
---|
| 78 | END FUNCTION lim_sbc_alloc_2 |
---|
| 79 | |
---|
| 80 | |
---|
[2528] | 81 | SUBROUTINE lim_sbc_flx_2( kt ) |
---|
[888] | 82 | !!------------------------------------------------------------------- |
---|
| 83 | !! *** ROUTINE lim_sbc_2 *** |
---|
| 84 | !! |
---|
[2566] | 85 | !! ** Purpose : Update surface ocean boundary condition over areas |
---|
| 86 | !! that are at least partially covered by sea-ice |
---|
[888] | 87 | !! |
---|
| 88 | !! ** Action : - comput. of the momentum, heat and freshwater/salt |
---|
[2566] | 89 | !! fluxes at the ice-ocean interface. |
---|
| 90 | !! - Update the fluxes provided to the ocean |
---|
[888] | 91 | !! |
---|
[3625] | 92 | !! ** Outputs : - qsr : sea heat flux : solar |
---|
| 93 | !! - qns : sea heat flux : non solar (including heat content of the mass flux) |
---|
| 94 | !! - emp : freshwater budget: mass flux |
---|
| 95 | !! - sfx : freshwater budget: salt flux due to Freezing/Melting |
---|
[1037] | 96 | !! - fr_i : ice fraction |
---|
| 97 | !! - tn_ice : sea-ice surface temperature |
---|
[5407] | 98 | !! - alb_ice : sea-ice albedo (ln_cpl=T) |
---|
[888] | 99 | !! |
---|
| 100 | !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90. |
---|
| 101 | !! Tartinville et al. 2001 Ocean Modelling, 3, 95-108. |
---|
| 102 | !!--------------------------------------------------------------------- |
---|
[2528] | 103 | INTEGER, INTENT(in) :: kt ! number of iteration |
---|
[6004] | 104 | ! |
---|
[2528] | 105 | INTEGER :: ji, jj ! dummy loop indices |
---|
| 106 | INTEGER :: ii0, ii1, ij0, ij1 ! local integers |
---|
| 107 | INTEGER :: ifvt, i1mfr, idfr, iflt ! - - |
---|
| 108 | INTEGER :: ial, iadv, ifral, ifrdv ! - - |
---|
[3625] | 109 | REAL(wp) :: zqsr, zqns, zfmm ! local scalars |
---|
| 110 | REAL(wp) :: zinda, zfsalt, zemp ! - - |
---|
| 111 | REAL(wp) :: zemp_snw, zqhc, zcd ! - - |
---|
| 112 | REAL(wp) :: zswitch ! - - |
---|
[3294] | 113 | REAL(wp), POINTER, DIMENSION(:,:) :: zqnsoce ! 2D workspace |
---|
[2715] | 114 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb, zalbp ! 2D/3D workspace |
---|
[888] | 115 | !!--------------------------------------------------------------------- |
---|
[6004] | 116 | ! |
---|
[3294] | 117 | CALL wrk_alloc( jpi, jpj, zqnsoce ) |
---|
| 118 | CALL wrk_alloc( jpi, jpj, 1, zalb, zalbp ) |
---|
[6004] | 119 | ! |
---|
| 120 | SELECT CASE( nn_ice_embd ) ! levitating or embedded sea-ice option |
---|
| 121 | CASE( 0 ) ; zswitch = 1 ! (0) old levitating sea-ice : salt exchange only |
---|
| 122 | CASE( 1, 2 ) ; zswitch = 0 ! (1) levitating sea-ice: salt and volume exchange but no pressure effect |
---|
| 123 | ! ! (2) embedded sea-ice : salt and volume fluxes and pressure |
---|
| 124 | END SELECT |
---|
[888] | 125 | |
---|
| 126 | !------------------------------------------! |
---|
| 127 | ! heat flux at the ocean surface ! |
---|
| 128 | !------------------------------------------! |
---|
| 129 | |
---|
[1482] | 130 | zqnsoce(:,:) = qns(:,:) |
---|
[888] | 131 | DO jj = 1, jpj |
---|
| 132 | DO ji = 1, jpi |
---|
| 133 | zinda = 1.0 - MAX( rzero , SIGN( rone, - ( 1.0 - pfrld(ji,jj) ) ) ) |
---|
| 134 | ifvt = zinda * MAX( rzero , SIGN( rone, - phicif(ji,jj) ) ) |
---|
| 135 | i1mfr = 1.0 - MAX( rzero , SIGN( rone, - ( 1.0 - frld(ji,jj) ) ) ) |
---|
| 136 | idfr = 1.0 - MAX( rzero , SIGN( rone, frld(ji,jj) - pfrld(ji,jj) ) ) |
---|
| 137 | iflt = zinda * (1 - i1mfr) * (1 - ifvt ) |
---|
| 138 | ial = ifvt * i1mfr + ( 1 - ifvt ) * idfr |
---|
| 139 | iadv = ( 1 - i1mfr ) * zinda |
---|
| 140 | ifral = ( 1 - i1mfr * ( 1 - ial ) ) |
---|
| 141 | ifrdv = ( 1 - ifral * ( 1 - ial ) ) * iadv |
---|
[1218] | 142 | |
---|
[3625] | 143 | !!$ attempt to explain the tricky flags set above.... |
---|
| 144 | !!$ zinda = 1.0 - AINT( pfrld(ji,jj) ) ! = 0. if ice-free ocean else 1. (after ice adv, but before ice thermo) |
---|
| 145 | !!$ i1mfr = 1.0 - AINT( frld(ji,jj) ) ! = 0. if ice-free ocean else 1. (after ice thermo) |
---|
[1218] | 146 | !!$ |
---|
[3625] | 147 | !!$ IF( phicif(ji,jj) <= 0. ) THEN ; ifvt = zinda ! = zinda if previous thermodynamic step overmelted the ice??? |
---|
| 148 | !!$ ELSE ; ifvt = 0. ! |
---|
[1218] | 149 | !!$ ENDIF |
---|
| 150 | !!$ |
---|
[3625] | 151 | !!$ IF( frld(ji,jj) >= pfrld(ji,jj) ) THEN ; idfr = 0. ! = 0. if lead fraction increases due to ice thermodynamics |
---|
[1218] | 152 | !!$ ELSE ; idfr = 1. |
---|
| 153 | !!$ ENDIF |
---|
| 154 | !!$ |
---|
[3625] | 155 | !!$ iflt = zinda * (1 - i1mfr) * (1 - ifvt ) ! = 1. if ice (not only snow) at previous time and ice-free ocean currently |
---|
[1218] | 156 | !!$ |
---|
| 157 | !!$ ial = ifvt * i1mfr + ( 1 - ifvt ) * idfr |
---|
[3625] | 158 | !!$ = i1mfr if ifvt = 1 i.e. |
---|
| 159 | !!$ = idfr if ifvt = 0 |
---|
[1218] | 160 | !!$! snow no ice ice ice or nothing lead fraction increases |
---|
| 161 | !!$! at previous now at previous |
---|
[3625] | 162 | !!$! -> ice area increases ??? -> ice area decreases ??? |
---|
[1218] | 163 | !!$ |
---|
[2715] | 164 | !!$ iadv = ( 1 - i1mfr ) * zinda |
---|
[1218] | 165 | !!$! pure ocean ice at |
---|
| 166 | !!$! at current previous |
---|
| 167 | !!$! -> = 1. if ice disapear between previous and current |
---|
| 168 | !!$ |
---|
| 169 | !!$ ifral = ( 1 - i1mfr * ( 1 - ial ) ) |
---|
| 170 | !!$! ice at ??? |
---|
| 171 | !!$! current |
---|
| 172 | !!$! -> ??? |
---|
| 173 | !!$ |
---|
[2715] | 174 | !!$ ifrdv = ( 1 - ifral * ( 1 - ial ) ) * iadv |
---|
| 175 | !!$! ice disapear |
---|
[1218] | 176 | !!$ |
---|
[2715] | 177 | !!$ |
---|
[1218] | 178 | |
---|
[888] | 179 | ! computation the solar flux at ocean surface |
---|
[5407] | 180 | IF( ln_cpl ) THEN |
---|
[4990] | 181 | zqsr = qsr_tot(ji,jj) + ( fstric(ji,jj) - qsr_ice(ji,jj,1) ) * ( 1.0 - pfrld(ji,jj) ) |
---|
| 182 | ELSE |
---|
| 183 | zqsr = pfrld(ji,jj) * qsr(ji,jj) + ( 1. - pfrld(ji,jj) ) * fstric(ji,jj) |
---|
| 184 | ENDIF |
---|
[888] | 185 | ! computation the non solar heat flux at ocean surface |
---|
[3625] | 186 | zqns = - ( 1. - thcm(ji,jj) ) * zqsr & ! part of the solar energy used in leads |
---|
| 187 | & + iflt * ( fscmbq(ji,jj) + ffltbif(ji,jj) ) & |
---|
| 188 | & + ifral * ( ial * qcmif(ji,jj) + (1 - ial) * qldif(ji,jj) ) * r1_rdtice & |
---|
| 189 | & + ifrdv * ( qfvbq(ji,jj) + qdtcn(ji,jj) ) * r1_rdtice |
---|
[888] | 190 | |
---|
[3625] | 191 | fsbbq(ji,jj) = ( 1.0 - ( ifvt + iflt ) ) * fscmbq(ji,jj) ! store residual heat flux (to put into the ocean at the next time-step) |
---|
| 192 | zqhc = ( rdq_snw(ji,jj) & |
---|
| 193 | & + rdq_ice(ji,jj) * ( 1.- zswitch) ) * r1_rdtice ! heat flux due to snow ( & ice heat content, |
---|
| 194 | ! ! if ice/ocean mass exchange active) |
---|
[888] | 195 | qsr (ji,jj) = zqsr ! solar heat flux |
---|
[3625] | 196 | qns (ji,jj) = zqns - fdtcn(ji,jj) + zqhc ! non solar heat flux |
---|
[2528] | 197 | ! |
---|
[3625] | 198 | ! !------------------------------------------! |
---|
| 199 | ! ! mass and salt flux at the ocean surface ! |
---|
| 200 | ! !------------------------------------------! |
---|
| 201 | ! |
---|
| 202 | ! mass flux at the ocean-atmosphere interface (open ocean fraction = leads area) |
---|
| 203 | ! ! coupled mode: |
---|
[5407] | 204 | IF( ln_cpl ) THEN |
---|
[4990] | 205 | zemp = + emp_tot(ji,jj) & ! net mass flux over the grid cell (ice+ocean area) |
---|
| 206 | & - emp_ice(ji,jj) * ( 1. - pfrld(ji,jj) ) ! minus the mass flux intercepted by sea-ice |
---|
| 207 | ELSE |
---|
| 208 | ! ! forced mode: |
---|
| 209 | zemp = + emp(ji,jj) * frld(ji,jj) & ! mass flux over open ocean fraction |
---|
| 210 | & - tprecip(ji,jj) * ( 1. - frld(ji,jj) ) & ! liquid precip. over ice reaches directly the ocean |
---|
| 211 | & + sprecip(ji,jj) * ( 1. - pfrld(ji,jj) ) ! snow is intercepted by sea-ice (previous frld) |
---|
| 212 | ENDIF |
---|
[2528] | 213 | ! |
---|
[3625] | 214 | ! mass flux at the ocean/ice interface (sea ice fraction) |
---|
| 215 | zemp_snw = rdm_snw(ji,jj) * r1_rdtice ! snow melting = pure water that enters the ocean |
---|
| 216 | zfmm = rdm_ice(ji,jj) * r1_rdtice ! Freezing minus Melting (F-M) |
---|
| 217 | |
---|
[4148] | 218 | fmmflx(ji,jj) = zfmm ! F/M mass flux save at least for biogeochemical model |
---|
| 219 | |
---|
[3625] | 220 | ! salt flux at the ice/ocean interface (sea ice fraction) [PSU*kg/m2/s] |
---|
| 221 | zfsalt = - sice_0(ji,jj) * zfmm ! F-M salt exchange |
---|
| 222 | zcd = soce_0(ji,jj) * zfmm ! concentration/dilution term due to F-M |
---|
[2528] | 223 | ! |
---|
[3625] | 224 | ! salt flux only : add concentration dilution term in salt flux and no F-M term in volume flux |
---|
| 225 | ! salt and mass fluxes : non concentration dilution term in salt flux and add F-M term in volume flux |
---|
| 226 | sfx (ji,jj) = zfsalt + zswitch * zcd ! salt flux (+ C/D if no ice/ocean mass exchange) |
---|
| 227 | emp (ji,jj) = zemp + zemp_snw + ( 1.- zswitch) * zfmm ! mass flux (+ F/M mass flux if ice/ocean mass exchange) |
---|
[2528] | 228 | ! |
---|
[888] | 229 | END DO |
---|
| 230 | END DO |
---|
[3625] | 231 | ! !------------------------------------------! |
---|
| 232 | ! ! mass of snow and ice per unit area ! |
---|
| 233 | ! !------------------------------------------! |
---|
| 234 | IF( nn_ice_embd /= 0 ) THEN ! embedded sea-ice (mass required) |
---|
| 235 | snwice_mass_b(:,:) = snwice_mass(:,:) ! save mass from the previous ice time step |
---|
| 236 | ! ! new mass per unit area |
---|
| 237 | snwice_mass (:,:) = tms(:,:) * ( rhosn * hsnif(:,:) + rhoic * hicif(:,:) ) * ( 1.0 - frld(:,:) ) |
---|
| 238 | ! ! time evolution of snow+ice mass |
---|
| 239 | snwice_fmass (:,:) = ( snwice_mass(:,:) - snwice_mass_b(:,:) ) / rdt_ice |
---|
| 240 | ENDIF |
---|
[888] | 241 | |
---|
[4990] | 242 | IF( iom_use('hflx_ice_cea' ) ) CALL iom_put( 'hflx_ice_cea', - fdtcn(:,:) ) |
---|
| 243 | IF( iom_use('qns_io_cea' ) ) CALL iom_put( 'qns_io_cea', qns(:,:) - zqnsoce(:,:) * pfrld(:,:) ) |
---|
| 244 | IF( iom_use('qsr_io_cea' ) ) CALL iom_put( 'qsr_io_cea', fstric(:,:) * (1.e0 - pfrld(:,:)) ) |
---|
[3625] | 245 | |
---|
[4990] | 246 | IF( iom_use('isnwmlt_cea' ) ) CALL iom_put( 'isnwmlt_cea' , rdm_snw(:,:) * r1_rdtice ) |
---|
| 247 | IF( iom_use('fsal_virt_cea') ) CALL iom_put( 'fsal_virt_cea', soce_0(:,:) * rdm_ice(:,:) * r1_rdtice ) |
---|
| 248 | IF( iom_use('fsal_real_cea') ) CALL iom_put( 'fsal_real_cea', - sice_0(:,:) * rdm_ice(:,:) * r1_rdtice ) |
---|
[1756] | 249 | |
---|
[888] | 250 | !-----------------------------------------------! |
---|
[1482] | 251 | ! Coupling variables ! |
---|
[888] | 252 | !-----------------------------------------------! |
---|
| 253 | |
---|
[5407] | 254 | IF( ln_cpl) THEN |
---|
[4990] | 255 | tn_ice(:,:,1) = sist(:,:) ! sea-ice surface temperature |
---|
| 256 | ht_i(:,:,1) = hicif(:,:) |
---|
| 257 | ht_s(:,:,1) = hsnif(:,:) |
---|
| 258 | a_i(:,:,1) = fr_i(:,:) |
---|
| 259 | ! ! Computation of snow/ice and ocean albedo |
---|
| 260 | CALL albedo_ice( tn_ice, ht_i, ht_s, zalbp, zalb ) |
---|
| 261 | alb_ice(:,:,1) = 0.5 * ( zalbp(:,:,1) + zalb (:,:,1) ) ! Ice albedo (mean clear and overcast skys) |
---|
| 262 | IF( iom_use('icealb_cea' ) ) CALL iom_put( 'icealb_cea', alb_ice(:,:,1) * fr_i(:,:) ) ! ice albedo |
---|
| 263 | ENDIF |
---|
[888] | 264 | |
---|
[2528] | 265 | IF(ln_ctl) THEN ! control print |
---|
[888] | 266 | CALL prt_ctl(tab2d_1=qsr , clinfo1=' lim_sbc: qsr : ', tab2d_2=qns , clinfo2=' qns : ') |
---|
[3625] | 267 | CALL prt_ctl(tab2d_1=emp , clinfo1=' lim_sbc: emp : ', tab2d_2=sfx , clinfo2=' sfx : ') |
---|
[1463] | 268 | CALL prt_ctl(tab2d_1=fr_i , clinfo1=' lim_sbc: fr_i : ', tab2d_2=tn_ice(:,:,1), clinfo2=' tn_ice : ') |
---|
[888] | 269 | ENDIF |
---|
[2528] | 270 | ! |
---|
[3294] | 271 | CALL wrk_dealloc( jpi, jpj, zqnsoce ) |
---|
| 272 | CALL wrk_dealloc( jpi, jpj, 1, zalb, zalbp ) |
---|
[2715] | 273 | ! |
---|
[2528] | 274 | END SUBROUTINE lim_sbc_flx_2 |
---|
[888] | 275 | |
---|
[2528] | 276 | |
---|
| 277 | SUBROUTINE lim_sbc_tau_2( kt , pu_oce, pv_oce ) |
---|
| 278 | !!------------------------------------------------------------------- |
---|
[2566] | 279 | !! *** ROUTINE lim_sbc_tau_2 *** |
---|
[2528] | 280 | !! |
---|
| 281 | !! ** Purpose : Update the ocean surface stresses due to the ice |
---|
| 282 | !! |
---|
| 283 | !! ** Action : * at each ice time step (every nn_fsbc time step): |
---|
| 284 | !! - compute the modulus of ice-ocean relative velocity |
---|
| 285 | !! at T-point (C-grid) or I-point (B-grid) |
---|
| 286 | !! tmod_io = rhoco * | U_ice-U_oce | |
---|
| 287 | !! - update the modulus of stress at ocean surface |
---|
| 288 | !! taum = frld * taum + (1-frld) * tmod_io * | U_ice-U_oce | |
---|
| 289 | !! * at each ocean time step (each kt): |
---|
| 290 | !! compute linearized ice-ocean stresses as |
---|
| 291 | !! Utau = tmod_io * | U_ice - pU_oce | |
---|
| 292 | !! using instantaneous current ocean velocity (usually before) |
---|
| 293 | !! |
---|
| 294 | !! NB: - the averaging operator used depends on the ice dynamics grid (cp_ice_msh='I' or 'C') |
---|
| 295 | !! - ice-ocean rotation angle only allowed in cp_ice_msh='I' case |
---|
| 296 | !! - here we make an approximation: taum is only computed every ice time step |
---|
| 297 | !! This avoids mutiple average to pass from T -> U,V grids and next from U,V grids |
---|
[2566] | 298 | !! to T grid. taum is used in TKE and GLS, which should not be too sensitive to this approximation... |
---|
[2528] | 299 | !! |
---|
[2566] | 300 | !! ** Outputs : - utau, vtau : surface ocean i- and j-stress (u- & v-pts) updated with ice-ocean fluxes |
---|
| 301 | !! - taum : modulus of the surface ocean stress (T-point) updated with ice-ocean fluxes |
---|
[2528] | 302 | !!--------------------------------------------------------------------- |
---|
| 303 | INTEGER , INTENT(in) :: kt ! ocean time-step index |
---|
| 304 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pu_oce, pv_oce ! surface ocean currents |
---|
[6004] | 305 | ! |
---|
[2528] | 306 | INTEGER :: ji, jj ! dummy loop indices |
---|
| 307 | REAL(wp) :: zfrldu, zat_u, zu_i, zutau_ice, zu_t, zmodt ! local scalar |
---|
| 308 | REAL(wp) :: zfrldv, zat_v, zv_i, zvtau_ice, zv_t, zmodi ! - - |
---|
[2566] | 309 | REAL(wp) :: zsang, zumt ! - - |
---|
[3294] | 310 | REAL(wp), POINTER, DIMENSION(:,:) :: ztio_u, ztio_v ! ocean stress below sea-ice |
---|
[2528] | 311 | !!--------------------------------------------------------------------- |
---|
| 312 | ! |
---|
[3294] | 313 | CALL wrk_alloc( jpi, jpj, ztio_u, ztio_v ) |
---|
[2528] | 314 | ! |
---|
| 315 | SELECT CASE( cp_ice_msh ) |
---|
| 316 | ! !-----------------------! |
---|
| 317 | CASE( 'I' ) ! B-grid ice dynamics ! I-point (i.e. F-point with sea-ice indexation) |
---|
| 318 | ! !--=--------------------! |
---|
| 319 | ! |
---|
| 320 | IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN !== Ice time-step only ==! (i.e. surface module time-step) |
---|
[5836] | 321 | ! |
---|
[2528] | 322 | DO jj = 1, jpj !* modulus of ice-ocean relative velocity at I-point |
---|
| 323 | DO ji = 1, jpi |
---|
| 324 | zu_i = u_ice(ji,jj) - u_oce(ji,jj) ! ice-ocean relative velocity at I-point |
---|
| 325 | zv_i = v_ice(ji,jj) - v_oce(ji,jj) |
---|
| 326 | tmod_io(ji,jj) = SQRT( zu_i * zu_i + zv_i * zv_i ) ! modulus of this velocity (at I-point) |
---|
| 327 | END DO |
---|
| 328 | END DO |
---|
| 329 | DO jj = 1, jpjm1 !* update the modulus of stress at ocean surface (T-point) |
---|
| 330 | DO ji = 1, jpim1 ! NO vector opt. |
---|
| 331 | ! ! modulus of U_ice-U_oce at T-point |
---|
[2566] | 332 | zumt = 0.25_wp * ( tmod_io(ji+1,jj) + tmod_io(ji ,jj ) & |
---|
[2528] | 333 | & + tmod_io(ji,jj+1) + tmod_io(ji+1,jj+1) ) |
---|
| 334 | ! ! update the modulus of stress at ocean surface |
---|
| 335 | taum(ji,jj) = frld(ji,jj) * taum(ji,jj) + ( 1._wp - frld(ji,jj) ) * rhoco * zumt * zumt |
---|
| 336 | END DO |
---|
| 337 | END DO |
---|
| 338 | CALL lbc_lnk( taum, 'T', 1. ) |
---|
| 339 | ! |
---|
| 340 | utau_oce(:,:) = utau(:,:) !* save the air-ocean stresses at ice time-step |
---|
| 341 | vtau_oce(:,:) = vtau(:,:) |
---|
| 342 | ! |
---|
| 343 | ENDIF |
---|
| 344 | ! |
---|
| 345 | ! !== at each ocean time-step ==! |
---|
| 346 | ! |
---|
| 347 | ! !* ice/ocean stress WITH a ice-ocean rotation angle at I-point |
---|
| 348 | DO jj = 2, jpj |
---|
| 349 | zsang = SIGN( 1._wp, gphif(1,jj) ) * sangvg ! change the cosine angle sign in the SH |
---|
| 350 | DO ji = 2, jpi ! NO vect. opt. possible |
---|
| 351 | ! ... ice-ocean relative velocity at I-point using instantaneous surface ocean current at u- & v-pts |
---|
| 352 | zu_i = u_ice(ji,jj) - 0.5_wp * ( pu_oce(ji-1,jj ) + pu_oce(ji-1,jj-1) ) * tmu(ji,jj) |
---|
| 353 | zv_i = v_ice(ji,jj) - 0.5_wp * ( pv_oce(ji ,jj-1) + pv_oce(ji-1,jj-1) ) * tmu(ji,jj) |
---|
| 354 | ! ... components of stress with a ice-ocean rotation angle |
---|
| 355 | zmodi = rhoco * tmod_io(ji,jj) |
---|
| 356 | ztio_u(ji,jj) = zmodi * ( cangvg * zu_i - zsang * zv_i ) |
---|
| 357 | ztio_v(ji,jj) = zmodi * ( cangvg * zv_i + zsang * zu_i ) |
---|
| 358 | END DO |
---|
| 359 | END DO |
---|
| 360 | ! !* surface ocean stresses at u- and v-points |
---|
| 361 | DO jj = 2, jpjm1 |
---|
| 362 | DO ji = 2, jpim1 ! NO vector opt. |
---|
| 363 | ! ! ice-ocean stress at U and V-points (from I-point values) |
---|
| 364 | zutau_ice = 0.5_wp * ( ztio_u(ji+1,jj) + ztio_u(ji+1,jj+1) ) |
---|
| 365 | zvtau_ice = 0.5_wp * ( ztio_v(ji,jj+1) + ztio_v(ji+1,jj+1) ) |
---|
| 366 | ! ! open-ocean (lead) fraction at U- & V-points (from T-point values) |
---|
[2566] | 367 | zfrldu = 0.5_wp * ( frld(ji,jj) + frld(ji+1,jj) ) |
---|
| 368 | zfrldv = 0.5_wp * ( frld(ji,jj) + frld(ji,jj+1) ) |
---|
[2528] | 369 | ! ! update the surface ocean stress (ice-cover wheighted) |
---|
| 370 | utau(ji,jj) = zfrldu * utau_oce(ji,jj) + ( 1._wp - zfrldu ) * zutau_ice |
---|
| 371 | vtau(ji,jj) = zfrldv * vtau_oce(ji,jj) + ( 1._wp - zfrldv ) * zvtau_ice |
---|
| 372 | END DO |
---|
| 373 | END DO |
---|
| 374 | CALL lbc_lnk( utau, 'U', -1. ) ; CALL lbc_lnk( vtau, 'V', -1. ) ! lateral boundary condition |
---|
| 375 | ! |
---|
| 376 | ! |
---|
| 377 | ! !-----------------------! |
---|
| 378 | CASE( 'C' ) ! C-grid ice dynamics ! U & V-points (same as in the ocean) |
---|
| 379 | ! !--=--------------------! |
---|
| 380 | ! |
---|
| 381 | IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN !== Ice time-step only ==! (i.e. surface module time-step) |
---|
[5836] | 382 | ! |
---|
[2528] | 383 | DO jj = 2, jpjm1 !* modulus of the ice-ocean velocity at T-point |
---|
| 384 | DO ji = fs_2, fs_jpim1 |
---|
| 385 | zu_t = u_ice(ji,jj) + u_ice(ji-1,jj) - u_oce(ji,jj) - u_oce(ji-1,jj) ! 2*(U_ice-U_oce) at T-point |
---|
| 386 | zv_t = v_ice(ji,jj) + v_ice(ji,jj-1) - v_oce(ji,jj) - v_oce(ji,jj-1) |
---|
| 387 | zmodt = 0.25_wp * ( zu_t * zu_t + zv_t * zv_t ) ! |U_ice-U_oce|^2 |
---|
| 388 | ! ! update the modulus of stress at ocean surface |
---|
| 389 | taum (ji,jj) = frld(ji,jj) * taum(ji,jj) + ( 1._wp - frld(ji,jj) ) * rhoco * zmodt |
---|
| 390 | tmod_io(ji,jj) = SQRT( zmodt ) * rhoco ! rhoco*|Uice-Uoce| |
---|
| 391 | END DO |
---|
| 392 | END DO |
---|
| 393 | CALL lbc_lnk( taum, 'T', 1. ) ; CALL lbc_lnk( tmod_io, 'T', 1. ) |
---|
| 394 | ! |
---|
| 395 | utau_oce(:,:) = utau(:,:) !* save the air-ocean stresses at ice time-step |
---|
| 396 | vtau_oce(:,:) = vtau(:,:) |
---|
| 397 | ! |
---|
| 398 | ENDIF |
---|
| 399 | ! |
---|
| 400 | ! !== at each ocean time-step ==! |
---|
| 401 | ! |
---|
| 402 | DO jj = 2, jpjm1 !* ice stress over ocean WITHOUT a ice-ocean rotation angle |
---|
| 403 | DO ji = fs_2, fs_jpim1 |
---|
| 404 | ! ! ocean area at u- & v-points |
---|
[2566] | 405 | zfrldu = 0.5_wp * ( frld(ji,jj) + frld(ji+1,jj) ) |
---|
| 406 | zfrldv = 0.5_wp * ( frld(ji,jj) + frld(ji,jj+1) ) |
---|
[2528] | 407 | ! ! quadratic drag formulation without rotation |
---|
| 408 | ! ! using instantaneous surface ocean current |
---|
| 409 | zutau_ice = 0.5 * ( tmod_io(ji,jj) + tmod_io(ji+1,jj) ) * ( u_ice(ji,jj) - pu_oce(ji,jj) ) |
---|
| 410 | zvtau_ice = 0.5 * ( tmod_io(ji,jj) + tmod_io(ji,jj+1) ) * ( v_ice(ji,jj) - pv_oce(ji,jj) ) |
---|
| 411 | ! ! update the surface ocean stress (ice-cover wheighted) |
---|
| 412 | utau(ji,jj) = zfrldu * utau_oce(ji,jj) + ( 1._wp - zfrldu ) * zutau_ice |
---|
| 413 | vtau(ji,jj) = zfrldv * vtau_oce(ji,jj) + ( 1._wp - zfrldv ) * zvtau_ice |
---|
| 414 | END DO |
---|
| 415 | END DO |
---|
| 416 | CALL lbc_lnk( utau, 'U', -1. ) ; CALL lbc_lnk( vtau, 'V', -1. ) ! lateral boundary condition |
---|
| 417 | ! |
---|
| 418 | END SELECT |
---|
| 419 | |
---|
| 420 | IF(ln_ctl) CALL prt_ctl( tab2d_1=utau, clinfo1=' lim_sbc: utau : ', mask1=umask, & |
---|
| 421 | & tab2d_2=vtau, clinfo2=' vtau : ' , mask2=vmask ) |
---|
| 422 | ! |
---|
[3294] | 423 | CALL wrk_dealloc( jpi, jpj, ztio_u, ztio_v ) |
---|
[2715] | 424 | ! |
---|
[2528] | 425 | END SUBROUTINE lim_sbc_tau_2 |
---|
| 426 | |
---|
[2715] | 427 | |
---|
| 428 | SUBROUTINE lim_sbc_init_2 |
---|
| 429 | !!------------------------------------------------------------------- |
---|
| 430 | !! *** ROUTINE lim_sbc_init *** |
---|
| 431 | !! |
---|
| 432 | !! ** Purpose : Preparation of the file ice_evolu for the output of |
---|
| 433 | !! the temporal evolution of key variables |
---|
| 434 | !! |
---|
| 435 | !! ** input : Namelist namicedia |
---|
| 436 | !!------------------------------------------------------------------- |
---|
[6004] | 437 | INTEGER :: jk ! local integer |
---|
| 438 | !!------------------------------------------------------------------- |
---|
[2715] | 439 | ! |
---|
| 440 | IF(lwp) WRITE(numout,*) |
---|
| 441 | IF(lwp) WRITE(numout,*) 'lim_sbc_init_2 : LIM-2 sea-ice - surface boundary condition' |
---|
| 442 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~ ' |
---|
| 443 | |
---|
| 444 | ! ! allocate lim_sbc arrays |
---|
| 445 | IF( lim_sbc_alloc_2() /= 0 ) CALL ctl_stop( 'STOP', 'lim_sbc_flx_2 : unable to allocate arrays' ) |
---|
| 446 | ! |
---|
| 447 | r1_rdtice = 1._wp / rdt_ice |
---|
| 448 | ! |
---|
| 449 | soce_0(:,:) = soce ! constant SSS and ice salinity used in levitating sea-ice case |
---|
| 450 | sice_0(:,:) = sice |
---|
| 451 | ! |
---|
| 452 | IF( cp_cfg == "orca" ) THEN ! decrease ocean & ice reference salinities in the Baltic sea |
---|
| 453 | WHERE( 14._wp <= glamt(:,:) .AND. glamt(:,:) <= 32._wp .AND. & |
---|
| 454 | & 54._wp <= gphit(:,:) .AND. gphit(:,:) <= 66._wp ) |
---|
| 455 | soce_0(:,:) = 4._wp |
---|
| 456 | sice_0(:,:) = 2._wp |
---|
| 457 | END WHERE |
---|
| 458 | ENDIF |
---|
[3625] | 459 | ! ! embedded sea ice |
---|
| 460 | IF( nn_ice_embd /= 0 ) THEN ! mass exchanges between ice and ocean (case 1 or 2) set the snow+ice mass |
---|
| 461 | snwice_mass (:,:) = tms(:,:) * ( rhosn * hsnif(:,:) + rhoic * hicif(:,:) ) * ( 1.0 - frld(:,:) ) |
---|
| 462 | snwice_mass_b(:,:) = snwice_mass(:,:) |
---|
| 463 | ELSE |
---|
| 464 | snwice_mass (:,:) = 0.e0 ! no mass exchanges |
---|
| 465 | snwice_mass_b(:,:) = 0.e0 ! no mass exchanges |
---|
[4292] | 466 | snwice_fmass (:,:) = 0.e0 ! no mass exchanges |
---|
[3625] | 467 | ENDIF |
---|
| 468 | IF( nn_ice_embd == 2 .AND. & ! full embedment (case 2) & no restart : |
---|
| 469 | & .NOT.ln_rstart ) THEN ! deplete the initial ssh below sea-ice area |
---|
| 470 | sshn(:,:) = sshn(:,:) - snwice_mass(:,:) * r1_rau0 |
---|
| 471 | sshb(:,:) = sshb(:,:) - snwice_mass(:,:) * r1_rau0 |
---|
[5866] | 472 | !!gm I really don't like this staff here... Find a way to put that elsewhere or differently |
---|
[5862] | 473 | !!gm |
---|
[6004] | 474 | IF( .NOT.ln_linssh ) THEN |
---|
[5866] | 475 | |
---|
| 476 | do jk = 1,jpkm1 ! adjust initial vertical scale factors |
---|
| 477 | e3t_n(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshn(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) ) |
---|
| 478 | e3t_b(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshb(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) ) |
---|
| 479 | end do |
---|
| 480 | e3t_a(:,:,:) = e3t_b(:,:,:) |
---|
| 481 | ! Reconstruction of all vertical scale factors at now and before time steps |
---|
| 482 | ! ! Horizontal scale factor interpolations |
---|
| 483 | CALL dom_vvl_interpol( e3t_b(:,:,:), e3u_b(:,:,:), 'U' ) |
---|
| 484 | CALL dom_vvl_interpol( e3t_b(:,:,:), e3v_b(:,:,:), 'V' ) |
---|
| 485 | CALL dom_vvl_interpol( e3t_n(:,:,:), e3u_n(:,:,:), 'U' ) |
---|
| 486 | CALL dom_vvl_interpol( e3t_n(:,:,:), e3v_n(:,:,:), 'V' ) |
---|
| 487 | CALL dom_vvl_interpol( e3u_n(:,:,:), e3f_n(:,:,:), 'F' ) |
---|
| 488 | ! ! Vertical scale factor interpolations |
---|
| 489 | CALL dom_vvl_interpol( e3t_n(:,:,:), e3w_n (:,:,:), 'W' ) |
---|
| 490 | CALL dom_vvl_interpol( e3u_n(:,:,:), e3uw_n(:,:,:), 'UW' ) |
---|
| 491 | CALL dom_vvl_interpol( e3v_n(:,:,:), e3vw_n(:,:,:), 'VW' ) |
---|
| 492 | CALL dom_vvl_interpol( e3u_b(:,:,:), e3uw_b(:,:,:), 'UW' ) |
---|
| 493 | CALL dom_vvl_interpol( e3v_b(:,:,:), e3vw_b(:,:,:), 'VW' ) |
---|
| 494 | ! ! t- and w- points depth |
---|
| 495 | gdept_n(:,:,1) = 0.5_wp * e3w_n(:,:,1) |
---|
| 496 | gdepw_n(:,:,1) = 0.0_wp |
---|
| 497 | gde3w_n(:,:,1) = gdept_n(:,:,1) - sshn(:,:) |
---|
| 498 | DO jk = 2, jpk |
---|
| 499 | gdept_n(:,:,jk) = gdept_n(:,:,jk-1) + e3w_n(:,:,jk) |
---|
| 500 | gdepw_n(:,:,jk) = gdepw_n(:,:,jk-1) + e3t_n(:,:,jk-1) |
---|
| 501 | gde3w_n(:,:,jk) = gdept_n(:,:,jk ) - sshn (:,:) |
---|
| 502 | END DO |
---|
| 503 | ENDIF |
---|
| 504 | !!gm end |
---|
[3625] | 505 | ENDIF |
---|
[2715] | 506 | ! |
---|
| 507 | END SUBROUTINE lim_sbc_init_2 |
---|
| 508 | |
---|
[888] | 509 | #else |
---|
| 510 | !!---------------------------------------------------------------------- |
---|
[2528] | 511 | !! Default option Empty module NO LIM 2.0 sea-ice model |
---|
[888] | 512 | !!---------------------------------------------------------------------- |
---|
| 513 | #endif |
---|
| 514 | |
---|
| 515 | !!====================================================================== |
---|
| 516 | END MODULE limsbc_2 |
---|