1 | MODULE limsbc_2 |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE limsbc_2 *** |
---|
4 | !! LIM-2 : updates the fluxes at the ocean surface with ice-ocean fluxes |
---|
5 | !!====================================================================== |
---|
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 |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | #if defined key_lim2 |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | !! 'key_lim2' LIM 2.0 sea-ice model |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | !! lim_sbc_flx_2 : update mass, heat and salt fluxes at the ocean surface |
---|
17 | !! lim_sbc_tau_2 : update i- and j-stresses, and its modulus at the ocean surface |
---|
18 | !!---------------------------------------------------------------------- |
---|
19 | USE par_oce ! ocean parameters |
---|
20 | USE phycst ! physical constants |
---|
21 | USE dom_oce ! ocean domain |
---|
22 | USE dom_ice_2 ! LIM-2: ice domain |
---|
23 | USE ice_2 ! LIM-2: ice variables |
---|
24 | USE sbc_ice ! surface boundary condition: ice |
---|
25 | USE sbc_oce ! surface boundary condition: ocean |
---|
26 | |
---|
27 | USE albedo ! albedo parameters |
---|
28 | USE lbclnk ! ocean lateral boundary condition - MPP exchanges |
---|
29 | USE in_out_manager ! I/O manager |
---|
30 | USE diaar5, ONLY : lk_diaar5 |
---|
31 | USE iom ! I/O library |
---|
32 | USE prtctl ! Print control |
---|
33 | USE cpl_oasis3, ONLY : lk_cpl |
---|
34 | |
---|
35 | IMPLICIT NONE |
---|
36 | PRIVATE |
---|
37 | |
---|
38 | PUBLIC lim_sbc_flx_2 ! called by sbc_ice_lim_2 |
---|
39 | PUBLIC lim_sbc_tau_2 ! called by sbc_ice_lim_2 |
---|
40 | |
---|
41 | REAL(wp) :: r1_rdtice ! = 1. / rdt_ice |
---|
42 | REAL(wp) :: epsi16 = 1.e-16_wp ! constant values |
---|
43 | REAL(wp) :: rzero = 0._wp ! - - |
---|
44 | REAL(wp) :: rone = 1._wp ! - - |
---|
45 | ! |
---|
46 | REAL(wp), DIMENSION(jpi,jpj) :: soce_0, sice_0 ! constant SSS and ice salinity used in levitating sea-ice case |
---|
47 | |
---|
48 | REAL(wp), DIMENSION(jpi,jpj) :: utau_oce, vtau_oce ! air-ocean surface i- & j-stress [N/m2] |
---|
49 | REAL(wp), DIMENSION(jpi,jpj) :: tmod_io ! modulus of the ice-ocean relative velocity [m/s] |
---|
50 | |
---|
51 | !! * Substitutions |
---|
52 | # include "vectopt_loop_substitute.h90" |
---|
53 | !!---------------------------------------------------------------------- |
---|
54 | !! NEMO/LIM2 3.3 , UCL - NEMO Consortium (2010) |
---|
55 | !! $Id$ |
---|
56 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
57 | !!---------------------------------------------------------------------- |
---|
58 | CONTAINS |
---|
59 | |
---|
60 | SUBROUTINE lim_sbc_flx_2( kt ) |
---|
61 | !!------------------------------------------------------------------- |
---|
62 | !! *** ROUTINE lim_sbc_2 *** |
---|
63 | !! |
---|
64 | !! ** Purpose : Update surface ocean boundary condition over areas |
---|
65 | !! that are at least partially covered by sea-ice |
---|
66 | !! |
---|
67 | !! ** Action : - comput. of the momentum, heat and freshwater/salt |
---|
68 | !! fluxes at the ice-ocean interface. |
---|
69 | !! - Update the fluxes provided to the ocean |
---|
70 | !! |
---|
71 | !! ** Outputs : - qsr : sea heat flux: solar |
---|
72 | !! - qns : sea heat flux: non solar |
---|
73 | !! - emp : freshwater budget: volume flux |
---|
74 | !! - emps : freshwater budget: concentration/dillution |
---|
75 | !! - utau : sea surface i-stress (ocean referential) |
---|
76 | !! - vtau : sea surface j-stress (ocean referential) |
---|
77 | !! - fr_i : ice fraction |
---|
78 | !! - tn_ice : sea-ice surface temperature |
---|
79 | !! - alb_ice : sea-ice alberdo (lk_cpl=T) |
---|
80 | !! |
---|
81 | !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90. |
---|
82 | !! Tartinville et al. 2001 Ocean Modelling, 3, 95-108. |
---|
83 | !!--------------------------------------------------------------------- |
---|
84 | INTEGER, INTENT(in) :: kt ! number of iteration |
---|
85 | !! |
---|
86 | INTEGER :: ji, jj ! dummy loop indices |
---|
87 | INTEGER :: ii0, ii1, ij0, ij1 ! local integers |
---|
88 | INTEGER :: ifvt, i1mfr, idfr, iflt ! - - |
---|
89 | INTEGER :: ial, iadv, ifral, ifrdv ! - - |
---|
90 | REAL(wp) :: zqsr, zqns, zfm ! local scalars |
---|
91 | REAL(wp) :: zinda, zfons, zemp ! - - |
---|
92 | REAL(wp), DIMENSION(jpi,jpj) :: zqnsoce ! 2D workspace |
---|
93 | REAL(wp), DIMENSION(jpi,jpj,1) :: zalb, zalbp ! 2D/3D workspace |
---|
94 | !!--------------------------------------------------------------------- |
---|
95 | |
---|
96 | IF( kt == nit000 ) THEN |
---|
97 | IF(lwp) WRITE(numout,*) |
---|
98 | IF(lwp) WRITE(numout,*) 'lim_sbc_flx_2 : LIM-2 sea-ice - surface boundary condition - Mass, heat & salt fluxes' |
---|
99 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~ ' |
---|
100 | ! |
---|
101 | r1_rdtice = 1._wp / rdt_ice |
---|
102 | ! |
---|
103 | soce_0(:,:) = soce ! constant SSS and ice salinity used in levitating sea-ice case |
---|
104 | sice_0(:,:) = sice |
---|
105 | ! |
---|
106 | IF( cp_cfg == "orca" ) THEN ! decrease ocean & ice reference salinities in the Baltic sea |
---|
107 | WHERE( 14._wp <= glamt(:,:) .AND. glamt(:,:) <= 32._wp .AND. & |
---|
108 | & 54._wp <= gphit(:,:) .AND. gphit(:,:) <= 66._wp ) |
---|
109 | soce_0(:,:) = 4._wp |
---|
110 | sice_0(:,:) = 2._wp |
---|
111 | END WHERE |
---|
112 | ENDIF |
---|
113 | ! |
---|
114 | ENDIF |
---|
115 | |
---|
116 | !------------------------------------------! |
---|
117 | ! heat flux at the ocean surface ! |
---|
118 | !------------------------------------------! |
---|
119 | |
---|
120 | zqnsoce(:,:) = qns(:,:) |
---|
121 | DO jj = 1, jpj |
---|
122 | DO ji = 1, jpi |
---|
123 | zinda = 1.0 - MAX( rzero , SIGN( rone, - ( 1.0 - pfrld(ji,jj) ) ) ) |
---|
124 | ifvt = zinda * MAX( rzero , SIGN( rone, - phicif(ji,jj) ) ) |
---|
125 | i1mfr = 1.0 - MAX( rzero , SIGN( rone, - ( 1.0 - frld(ji,jj) ) ) ) |
---|
126 | idfr = 1.0 - MAX( rzero , SIGN( rone, frld(ji,jj) - pfrld(ji,jj) ) ) |
---|
127 | iflt = zinda * (1 - i1mfr) * (1 - ifvt ) |
---|
128 | ial = ifvt * i1mfr + ( 1 - ifvt ) * idfr |
---|
129 | iadv = ( 1 - i1mfr ) * zinda |
---|
130 | ifral = ( 1 - i1mfr * ( 1 - ial ) ) |
---|
131 | ifrdv = ( 1 - ifral * ( 1 - ial ) ) * iadv |
---|
132 | |
---|
133 | !!$ zinda = 1.0 - AINT( pfrld(ji,jj) ) ! = 0. if pure ocean else 1. (at previous time) |
---|
134 | !!$ |
---|
135 | !!$ i1mfr = 1.0 - AINT( frld(ji,jj) ) ! = 0. if pure ocean else 1. (at current time) |
---|
136 | !!$ |
---|
137 | !!$ IF( phicif(ji,jj) <= 0. ) THEN ; ifvt = zinda ! = 1. if (snow and no ice at previous time) else 0. ??? |
---|
138 | !!$ ELSE ; ifvt = 0. |
---|
139 | !!$ ENDIF |
---|
140 | !!$ |
---|
141 | !!$ IF( frld(ji,jj) >= pfrld(ji,jj) ) THEN ; idfr = 0. ! = 0. if lead fraction increases from previous to current |
---|
142 | !!$ ELSE ; idfr = 1. |
---|
143 | !!$ ENDIF |
---|
144 | !!$ |
---|
145 | !!$ iflt = zinda * (1 - i1mfr) * (1 - ifvt ) ! = 1. if ice (not only snow) at previous and pure ocean at current |
---|
146 | !!$ |
---|
147 | !!$ ial = ifvt * i1mfr + ( 1 - ifvt ) * idfr |
---|
148 | !!$! snow no ice ice ice or nothing lead fraction increases |
---|
149 | !!$! at previous now at previous |
---|
150 | !!$! -> ice aera increases ??? -> ice aera decreases ??? |
---|
151 | !!$ |
---|
152 | !!$ iadv = ( 1 - i1mfr ) * zinda |
---|
153 | !!$! pure ocean ice at |
---|
154 | !!$! at current previous |
---|
155 | !!$! -> = 1. if ice disapear between previous and current |
---|
156 | !!$ |
---|
157 | !!$ ifral = ( 1 - i1mfr * ( 1 - ial ) ) |
---|
158 | !!$! ice at ??? |
---|
159 | !!$! current |
---|
160 | !!$! -> ??? |
---|
161 | !!$ |
---|
162 | !!$ ifrdv = ( 1 - ifral * ( 1 - ial ) ) * iadv |
---|
163 | !!$! ice disapear |
---|
164 | !!$ |
---|
165 | !!$ |
---|
166 | |
---|
167 | ! computation the solar flux at ocean surface |
---|
168 | #if defined key_coupled |
---|
169 | zqsr = qsr_tot(ji,jj) + ( fstric(ji,jj) - qsr_ice(ji,jj,1) ) * ( 1.0 - pfrld(ji,jj) ) |
---|
170 | #else |
---|
171 | zqsr = pfrld(ji,jj) * qsr(ji,jj) + ( 1. - pfrld(ji,jj) ) * fstric(ji,jj) |
---|
172 | #endif |
---|
173 | ! computation the non solar heat flux at ocean surface |
---|
174 | zqns = - ( 1. - thcm(ji,jj) ) * zqsr & ! part of the solar energy used in leads |
---|
175 | & + iflt * ( fscmbq(ji,jj) + ffltbif(ji,jj) ) & |
---|
176 | & + ifral * ( ial * qcmif(ji,jj) + (1 - ial) * qldif(ji,jj) ) * r1_rdtice & |
---|
177 | & + ifrdv * ( qfvbq(ji,jj) + qdtcn(ji,jj) ) * r1_rdtice |
---|
178 | |
---|
179 | fsbbq(ji,jj) = ( 1.0 - ( ifvt + iflt ) ) * fscmbq(ji,jj) ! ??? |
---|
180 | ! |
---|
181 | qsr (ji,jj) = zqsr ! solar heat flux |
---|
182 | qns (ji,jj) = zqns - fdtcn(ji,jj) ! non solar heat flux |
---|
183 | END DO |
---|
184 | END DO |
---|
185 | |
---|
186 | CALL iom_put( 'hflx_ice_cea', - fdtcn(:,:) ) |
---|
187 | CALL iom_put( 'qns_io_cea', qns(:,:) - zqnsoce(:,:) * pfrld(:,:) ) |
---|
188 | CALL iom_put( 'qsr_io_cea', fstric(:,:) * (1.e0 - pfrld(:,:)) ) |
---|
189 | |
---|
190 | !------------------------------------------! |
---|
191 | ! mass flux at the ocean surface ! |
---|
192 | !------------------------------------------! |
---|
193 | DO jj = 1, jpj |
---|
194 | DO ji = 1, jpi |
---|
195 | ! |
---|
196 | #if defined key_coupled |
---|
197 | ! freshwater exchanges at the ice-atmosphere / ocean interface (coupled mode) |
---|
198 | zemp = emp_tot(ji,jj) - emp_ice(ji,jj) * ( 1. - pfrld(ji,jj) ) & ! |
---|
199 | & + rdmsnif(ji,jj) * r1_rdtice ! freshwaterflux due to snow melting |
---|
200 | #else |
---|
201 | ! computing freshwater exchanges at the ice/ocean interface |
---|
202 | zemp = + emp(ji,jj) * frld(ji,jj) & ! e-p budget over open ocean fraction |
---|
203 | & - tprecip(ji,jj) * ( 1. - frld(ji,jj) ) & ! liquid precipitation reaches directly the ocean |
---|
204 | & + sprecip(ji,jj) * ( 1. - pfrld(ji,jj) ) & ! change in ice cover within the time step |
---|
205 | & + rdmsnif(ji,jj) * r1_rdtice ! freshwater flux due to snow melting |
---|
206 | #endif |
---|
207 | ! |
---|
208 | ! computing salt exchanges at the ice/ocean interface |
---|
209 | zfons = ( soce_0(ji,jj) - sice_0(ji,jj) ) * ( rdmicif(ji,jj) * r1_rdtice ) |
---|
210 | ! |
---|
211 | ! converting the salt flux from ice to a freshwater flux from ocean |
---|
212 | zfm = zfons / ( sss_m(ji,jj) + epsi16 ) |
---|
213 | ! |
---|
214 | emps(ji,jj) = zemp + zfm ! surface ocean concentration/dilution effect (use on SSS evolution) |
---|
215 | emp (ji,jj) = zemp ! surface ocean volume flux (use on sea-surface height evolution) |
---|
216 | ! |
---|
217 | END DO |
---|
218 | END DO |
---|
219 | |
---|
220 | IF( lk_diaar5 ) THEN ! AR5 diagnostics |
---|
221 | CALL iom_put( 'isnwmlt_cea' , rdmsnif(:,:) * r1_rdtice ) |
---|
222 | CALL iom_put( 'fsal_virt_cea', soce_0(:,:) * rdmicif(:,:) * r1_rdtice ) |
---|
223 | CALL iom_put( 'fsal_real_cea', - sice_0(:,:) * rdmicif(:,:) * r1_rdtice ) |
---|
224 | ENDIF |
---|
225 | |
---|
226 | !-----------------------------------------------! |
---|
227 | ! Coupling variables ! |
---|
228 | !-----------------------------------------------! |
---|
229 | |
---|
230 | IF( lk_cpl ) THEN ! coupled case |
---|
231 | ! Ice surface temperature |
---|
232 | tn_ice(:,:,1) = sist(:,:) ! sea-ice surface temperature |
---|
233 | ! Computation of snow/ice and ocean albedo |
---|
234 | CALL albedo_ice( tn_ice, reshape( hicif, (/jpi,jpj,1/) ), reshape( hsnif, (/jpi,jpj,1/) ), zalbp, zalb ) |
---|
235 | alb_ice(:,:,1) = 0.5 * ( zalbp(:,:,1) + zalb (:,:,1) ) ! Ice albedo (mean clear and overcast skys) |
---|
236 | CALL iom_put( "icealb_cea", alb_ice(:,:,1) * fr_i(:,:) ) ! ice albedo |
---|
237 | ENDIF |
---|
238 | |
---|
239 | IF(ln_ctl) THEN ! control print |
---|
240 | CALL prt_ctl(tab2d_1=qsr , clinfo1=' lim_sbc: qsr : ', tab2d_2=qns , clinfo2=' qns : ') |
---|
241 | CALL prt_ctl(tab2d_1=emp , clinfo1=' lim_sbc: emp : ', tab2d_2=emps , clinfo2=' emps : ') |
---|
242 | CALL prt_ctl(tab2d_1=utau , clinfo1=' lim_sbc: utau : ', mask1=umask, & |
---|
243 | & tab2d_2=vtau , clinfo2=' vtau : ' , mask2=vmask ) |
---|
244 | CALL prt_ctl(tab2d_1=fr_i , clinfo1=' lim_sbc: fr_i : ', tab2d_2=tn_ice(:,:,1), clinfo2=' tn_ice : ') |
---|
245 | ENDIF |
---|
246 | ! |
---|
247 | END SUBROUTINE lim_sbc_flx_2 |
---|
248 | |
---|
249 | |
---|
250 | SUBROUTINE lim_sbc_tau_2( kt , pu_oce, pv_oce ) |
---|
251 | !!------------------------------------------------------------------- |
---|
252 | !! *** ROUTINE lim_sbc_tau_2 *** |
---|
253 | !! |
---|
254 | !! ** Purpose : Update the ocean surface stresses due to the ice |
---|
255 | !! |
---|
256 | !! ** Action : * at each ice time step (every nn_fsbc time step): |
---|
257 | !! - compute the modulus of ice-ocean relative velocity |
---|
258 | !! at T-point (C-grid) or I-point (B-grid) |
---|
259 | !! tmod_io = rhoco * | U_ice-U_oce | |
---|
260 | !! - update the modulus of stress at ocean surface |
---|
261 | !! taum = frld * taum + (1-frld) * tmod_io * | U_ice-U_oce | |
---|
262 | !! * at each ocean time step (each kt): |
---|
263 | !! compute linearized ice-ocean stresses as |
---|
264 | !! Utau = tmod_io * | U_ice - pU_oce | |
---|
265 | !! using instantaneous current ocean velocity (usually before) |
---|
266 | !! |
---|
267 | !! NB: - the averaging operator used depends on the ice dynamics grid (cp_ice_msh='I' or 'C') |
---|
268 | !! - ice-ocean rotation angle only allowed in cp_ice_msh='I' case |
---|
269 | !! - here we make an approximation: taum is only computed every ice time step |
---|
270 | !! This avoids mutiple average to pass from T -> U,V grids and next from U,V grids |
---|
271 | !! to T grid. taum is used in TKE and GLS, which should not be too sensitive to this approximation... |
---|
272 | !! |
---|
273 | !! ** Outputs : - utau, vtau : surface ocean i- and j-stress (u- & v-pts) updated with ice-ocean fluxes |
---|
274 | !! - taum : modulus of the surface ocean stress (T-point) updated with ice-ocean fluxes |
---|
275 | !!--------------------------------------------------------------------- |
---|
276 | INTEGER , INTENT(in) :: kt ! ocean time-step index |
---|
277 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pu_oce, pv_oce ! surface ocean currents |
---|
278 | !! |
---|
279 | INTEGER :: ji, jj ! dummy loop indices |
---|
280 | REAL(wp) :: zfrldu, zat_u, zu_i, zutau_ice, zu_t, zmodt ! local scalar |
---|
281 | REAL(wp) :: zfrldv, zat_v, zv_i, zvtau_ice, zv_t, zmodi ! - - |
---|
282 | REAL(wp) :: zsang, zumt ! - - |
---|
283 | REAL(wp), DIMENSION(jpi,jpj) :: ztio_u, ztio_v ! ocean stress below sea-ice |
---|
284 | !!--------------------------------------------------------------------- |
---|
285 | ! |
---|
286 | IF( kt == nit000 .AND. lwp ) THEN ! control print |
---|
287 | WRITE(numout,*) |
---|
288 | WRITE(numout,*) 'lim_sbc_tau_2 : LIM 2.0 sea-ice - surface ocean momentum fluxes' |
---|
289 | WRITE(numout,*) '~~~~~~~~~~~~~ ' |
---|
290 | IF( lk_lim2_vp ) THEN ; WRITE(numout,*) ' VP rheology - B-grid case' |
---|
291 | ELSE ; WRITE(numout,*) ' EVP rheology - C-grid case' |
---|
292 | ENDIF |
---|
293 | ENDIF |
---|
294 | ! |
---|
295 | SELECT CASE( cp_ice_msh ) |
---|
296 | ! !-----------------------! |
---|
297 | CASE( 'I' ) ! B-grid ice dynamics ! I-point (i.e. F-point with sea-ice indexation) |
---|
298 | ! !--=--------------------! |
---|
299 | ! |
---|
300 | IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN !== Ice time-step only ==! (i.e. surface module time-step) |
---|
301 | !CDIR NOVERRCHK |
---|
302 | DO jj = 1, jpj !* modulus of ice-ocean relative velocity at I-point |
---|
303 | !CDIR NOVERRCHK |
---|
304 | DO ji = 1, jpi |
---|
305 | zu_i = u_ice(ji,jj) - u_oce(ji,jj) ! ice-ocean relative velocity at I-point |
---|
306 | zv_i = v_ice(ji,jj) - v_oce(ji,jj) |
---|
307 | tmod_io(ji,jj) = SQRT( zu_i * zu_i + zv_i * zv_i ) ! modulus of this velocity (at I-point) |
---|
308 | END DO |
---|
309 | END DO |
---|
310 | !CDIR NOVERRCHK |
---|
311 | DO jj = 1, jpjm1 !* update the modulus of stress at ocean surface (T-point) |
---|
312 | !CDIR NOVERRCHK |
---|
313 | DO ji = 1, jpim1 ! NO vector opt. |
---|
314 | ! ! modulus of U_ice-U_oce at T-point |
---|
315 | zumt = 0.25_wp * ( tmod_io(ji+1,jj) + tmod_io(ji ,jj ) & |
---|
316 | & + tmod_io(ji,jj+1) + tmod_io(ji+1,jj+1) ) |
---|
317 | ! ! update the modulus of stress at ocean surface |
---|
318 | taum(ji,jj) = frld(ji,jj) * taum(ji,jj) + ( 1._wp - frld(ji,jj) ) * rhoco * zumt * zumt |
---|
319 | END DO |
---|
320 | END DO |
---|
321 | CALL lbc_lnk( taum, 'T', 1. ) |
---|
322 | ! |
---|
323 | utau_oce(:,:) = utau(:,:) !* save the air-ocean stresses at ice time-step |
---|
324 | vtau_oce(:,:) = vtau(:,:) |
---|
325 | ! |
---|
326 | ENDIF |
---|
327 | ! |
---|
328 | ! !== at each ocean time-step ==! |
---|
329 | ! |
---|
330 | ! !* ice/ocean stress WITH a ice-ocean rotation angle at I-point |
---|
331 | DO jj = 2, jpj |
---|
332 | zsang = SIGN( 1._wp, gphif(1,jj) ) * sangvg ! change the cosine angle sign in the SH |
---|
333 | DO ji = 2, jpi ! NO vect. opt. possible |
---|
334 | ! ... ice-ocean relative velocity at I-point using instantaneous surface ocean current at u- & v-pts |
---|
335 | zu_i = u_ice(ji,jj) - 0.5_wp * ( pu_oce(ji-1,jj ) + pu_oce(ji-1,jj-1) ) * tmu(ji,jj) |
---|
336 | zv_i = v_ice(ji,jj) - 0.5_wp * ( pv_oce(ji ,jj-1) + pv_oce(ji-1,jj-1) ) * tmu(ji,jj) |
---|
337 | ! ... components of stress with a ice-ocean rotation angle |
---|
338 | zmodi = rhoco * tmod_io(ji,jj) |
---|
339 | ztio_u(ji,jj) = zmodi * ( cangvg * zu_i - zsang * zv_i ) |
---|
340 | ztio_v(ji,jj) = zmodi * ( cangvg * zv_i + zsang * zu_i ) |
---|
341 | END DO |
---|
342 | END DO |
---|
343 | ! !* surface ocean stresses at u- and v-points |
---|
344 | DO jj = 2, jpjm1 |
---|
345 | DO ji = 2, jpim1 ! NO vector opt. |
---|
346 | ! ! ice-ocean stress at U and V-points (from I-point values) |
---|
347 | zutau_ice = 0.5_wp * ( ztio_u(ji+1,jj) + ztio_u(ji+1,jj+1) ) |
---|
348 | zvtau_ice = 0.5_wp * ( ztio_v(ji,jj+1) + ztio_v(ji+1,jj+1) ) |
---|
349 | ! ! open-ocean (lead) fraction at U- & V-points (from T-point values) |
---|
350 | zfrldu = 0.5_wp * ( frld(ji,jj) + frld(ji+1,jj) ) |
---|
351 | zfrldv = 0.5_wp * ( frld(ji,jj) + frld(ji,jj+1) ) |
---|
352 | ! ! update the surface ocean stress (ice-cover wheighted) |
---|
353 | utau(ji,jj) = zfrldu * utau_oce(ji,jj) + ( 1._wp - zfrldu ) * zutau_ice |
---|
354 | vtau(ji,jj) = zfrldv * vtau_oce(ji,jj) + ( 1._wp - zfrldv ) * zvtau_ice |
---|
355 | END DO |
---|
356 | END DO |
---|
357 | CALL lbc_lnk( utau, 'U', -1. ) ; CALL lbc_lnk( vtau, 'V', -1. ) ! lateral boundary condition |
---|
358 | ! |
---|
359 | ! |
---|
360 | ! !-----------------------! |
---|
361 | CASE( 'C' ) ! C-grid ice dynamics ! U & V-points (same as in the ocean) |
---|
362 | ! !--=--------------------! |
---|
363 | ! |
---|
364 | IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN !== Ice time-step only ==! (i.e. surface module time-step) |
---|
365 | !CDIR NOVERRCHK |
---|
366 | DO jj = 2, jpjm1 !* modulus of the ice-ocean velocity at T-point |
---|
367 | !CDIR NOVERRCHK |
---|
368 | DO ji = fs_2, fs_jpim1 |
---|
369 | 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 |
---|
370 | zv_t = v_ice(ji,jj) + v_ice(ji,jj-1) - v_oce(ji,jj) - v_oce(ji,jj-1) |
---|
371 | zmodt = 0.25_wp * ( zu_t * zu_t + zv_t * zv_t ) ! |U_ice-U_oce|^2 |
---|
372 | ! ! update the modulus of stress at ocean surface |
---|
373 | taum (ji,jj) = frld(ji,jj) * taum(ji,jj) + ( 1._wp - frld(ji,jj) ) * rhoco * zmodt |
---|
374 | tmod_io(ji,jj) = SQRT( zmodt ) * rhoco ! rhoco*|Uice-Uoce| |
---|
375 | END DO |
---|
376 | END DO |
---|
377 | CALL lbc_lnk( taum, 'T', 1. ) ; CALL lbc_lnk( tmod_io, 'T', 1. ) |
---|
378 | ! |
---|
379 | utau_oce(:,:) = utau(:,:) !* save the air-ocean stresses at ice time-step |
---|
380 | vtau_oce(:,:) = vtau(:,:) |
---|
381 | ! |
---|
382 | ENDIF |
---|
383 | ! |
---|
384 | ! !== at each ocean time-step ==! |
---|
385 | ! |
---|
386 | DO jj = 2, jpjm1 !* ice stress over ocean WITHOUT a ice-ocean rotation angle |
---|
387 | DO ji = fs_2, fs_jpim1 |
---|
388 | ! ! ocean area at u- & v-points |
---|
389 | zfrldu = 0.5_wp * ( frld(ji,jj) + frld(ji+1,jj) ) |
---|
390 | zfrldv = 0.5_wp * ( frld(ji,jj) + frld(ji,jj+1) ) |
---|
391 | ! ! quadratic drag formulation without rotation |
---|
392 | ! ! using instantaneous surface ocean current |
---|
393 | zutau_ice = 0.5 * ( tmod_io(ji,jj) + tmod_io(ji+1,jj) ) * ( u_ice(ji,jj) - pu_oce(ji,jj) ) |
---|
394 | zvtau_ice = 0.5 * ( tmod_io(ji,jj) + tmod_io(ji,jj+1) ) * ( v_ice(ji,jj) - pv_oce(ji,jj) ) |
---|
395 | ! ! update the surface ocean stress (ice-cover wheighted) |
---|
396 | utau(ji,jj) = zfrldu * utau_oce(ji,jj) + ( 1._wp - zfrldu ) * zutau_ice |
---|
397 | vtau(ji,jj) = zfrldv * vtau_oce(ji,jj) + ( 1._wp - zfrldv ) * zvtau_ice |
---|
398 | END DO |
---|
399 | END DO |
---|
400 | CALL lbc_lnk( utau, 'U', -1. ) ; CALL lbc_lnk( vtau, 'V', -1. ) ! lateral boundary condition |
---|
401 | ! |
---|
402 | END SELECT |
---|
403 | |
---|
404 | IF(ln_ctl) CALL prt_ctl( tab2d_1=utau, clinfo1=' lim_sbc: utau : ', mask1=umask, & |
---|
405 | & tab2d_2=vtau, clinfo2=' vtau : ' , mask2=vmask ) |
---|
406 | ! |
---|
407 | END SUBROUTINE lim_sbc_tau_2 |
---|
408 | |
---|
409 | #else |
---|
410 | !!---------------------------------------------------------------------- |
---|
411 | !! Default option Empty module NO LIM 2.0 sea-ice model |
---|
412 | !!---------------------------------------------------------------------- |
---|
413 | #endif |
---|
414 | |
---|
415 | !!====================================================================== |
---|
416 | END MODULE limsbc_2 |
---|