1 | MODULE limsbc |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE limsbc *** |
---|
4 | !! computation of the flux at the sea ice/ocean interface |
---|
5 | !!====================================================================== |
---|
6 | !! History : - ! 2006-07 (M. Vancoppelle) LIM3 original code |
---|
7 | !! 3.0 ! 2008-03 (C. Tallandier) surface module |
---|
8 | !! - ! 2008-04 (C. Tallandier) split in 2 + new ice-ocean coupling |
---|
9 | !! 3.3 ! 2010-05 (G. Madec) decrease ocean & ice reference salinities in the Baltic sea |
---|
10 | !! ! + simplification of the ice-ocean stress calculation |
---|
11 | !! 3.4 ! 2011-02 (G. Madec) dynamical allocation |
---|
12 | !! - ! 2012 (D. Iovino) salt flux change |
---|
13 | !! - ! 2012-05 (C. Rousset) add penetration solar flux |
---|
14 | !! 3.5 ! 2012-10 (A. Coward, G. Madec) salt fluxes ; ice+snow mass |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | #if defined key_lim3 |
---|
17 | !!---------------------------------------------------------------------- |
---|
18 | !! 'key_lim3' LIM 3.0 sea-ice model |
---|
19 | !!---------------------------------------------------------------------- |
---|
20 | !! lim_sbc_alloc : allocate the limsbc arrays |
---|
21 | !! lim_sbc_init : initialisation |
---|
22 | !! lim_sbc_flx : updates mass, heat and salt fluxes at the ocean surface |
---|
23 | !! lim_sbc_tau : update i- and j-stresses, and its modulus at the ocean surface |
---|
24 | !!---------------------------------------------------------------------- |
---|
25 | USE par_oce ! ocean parameters |
---|
26 | USE phycst ! physical constants |
---|
27 | USE dom_oce ! ocean domain |
---|
28 | USE ice ! LIM sea-ice variables |
---|
29 | USE sbc_ice ! Surface boundary condition: sea-ice fields |
---|
30 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
31 | USE sbccpl |
---|
32 | USE oce , ONLY : sshn, sshb, snwice_mass, snwice_mass_b, snwice_fmass |
---|
33 | USE albedo ! albedo parameters |
---|
34 | USE lbclnk ! ocean lateral boundary condition - MPP exchanges |
---|
35 | USE lib_mpp ! MPP library |
---|
36 | USE wrk_nemo ! work arrays |
---|
37 | USE in_out_manager ! I/O manager |
---|
38 | USE prtctl ! Print control |
---|
39 | USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) |
---|
40 | USE traqsr ! add penetration of solar flux in the calculation of heat budget |
---|
41 | USE iom |
---|
42 | USE domvvl ! Variable volume |
---|
43 | USE limctl |
---|
44 | USE limcons |
---|
45 | |
---|
46 | IMPLICIT NONE |
---|
47 | PRIVATE |
---|
48 | |
---|
49 | PUBLIC lim_sbc_init ! called by sbc_lim_init |
---|
50 | PUBLIC lim_sbc_flx ! called by sbc_ice_lim |
---|
51 | PUBLIC lim_sbc_tau ! called by sbc_ice_lim |
---|
52 | |
---|
53 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: utau_oce, vtau_oce ! air-ocean surface i- & j-stress [N/m2] |
---|
54 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: tmod_io ! modulus of the ice-ocean velocity [m/s] |
---|
55 | |
---|
56 | !! * Substitutions |
---|
57 | # include "vectopt_loop_substitute.h90" |
---|
58 | # include "domzgr_substitute.h90" |
---|
59 | !!---------------------------------------------------------------------- |
---|
60 | !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011) |
---|
61 | !! $Id$ |
---|
62 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
63 | !!---------------------------------------------------------------------- |
---|
64 | CONTAINS |
---|
65 | |
---|
66 | INTEGER FUNCTION lim_sbc_alloc() |
---|
67 | !!------------------------------------------------------------------- |
---|
68 | !! *** ROUTINE lim_sbc_alloc *** |
---|
69 | !!------------------------------------------------------------------- |
---|
70 | ALLOCATE( utau_oce(jpi,jpj) , vtau_oce(jpi,jpj) , tmod_io(jpi,jpj), STAT=lim_sbc_alloc ) |
---|
71 | ! |
---|
72 | IF( lk_mpp ) CALL mpp_sum( lim_sbc_alloc ) |
---|
73 | IF( lim_sbc_alloc /= 0 ) CALL ctl_warn('lim_sbc_alloc: failed to allocate arrays') |
---|
74 | END FUNCTION lim_sbc_alloc |
---|
75 | |
---|
76 | |
---|
77 | SUBROUTINE lim_sbc_flx( kt ) |
---|
78 | !!------------------------------------------------------------------- |
---|
79 | !! *** ROUTINE lim_sbc_flx *** |
---|
80 | !! |
---|
81 | !! ** Purpose : Update the surface ocean boundary condition for heat |
---|
82 | !! salt and mass over areas where sea-ice is non-zero |
---|
83 | !! |
---|
84 | !! ** Action : - computes the heat and freshwater/salt fluxes |
---|
85 | !! at the ice-ocean interface. |
---|
86 | !! - Update the ocean sbc |
---|
87 | !! |
---|
88 | !! ** Outputs : - qsr : sea heat flux: solar |
---|
89 | !! - qns : sea heat flux: non solar |
---|
90 | !! - emp : freshwater budget: volume flux |
---|
91 | !! - sfx : salt flux |
---|
92 | !! - fr_i : ice fraction |
---|
93 | !! - tn_ice : sea-ice surface temperature |
---|
94 | !! - alb_ice : sea-ice albedo (recomputed only for coupled mode) |
---|
95 | !! |
---|
96 | !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90. |
---|
97 | !! Tartinville et al. 2001 Ocean Modelling, 3, 95-108. |
---|
98 | !! These refs are now obsolete since everything has been revised |
---|
99 | !! The ref should be Rousset et al., 2015 |
---|
100 | !!--------------------------------------------------------------------- |
---|
101 | INTEGER, INTENT(in) :: kt ! number of iteration |
---|
102 | INTEGER :: ji, jj, jl, jk ! dummy loop indices |
---|
103 | REAL(wp) :: zqmass ! Heat flux associated with mass exchange ice->ocean (W.m-2) |
---|
104 | REAL(wp) :: zqsr ! New solar flux received by the ocean |
---|
105 | ! |
---|
106 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_cs, zalb_os ! 3D workspace |
---|
107 | REAL(wp), POINTER, DIMENSION(:,:) :: zalb ! 2D workspace |
---|
108 | !!--------------------------------------------------------------------- |
---|
109 | |
---|
110 | ! --- case we bypass ice thermodynamics --- ! |
---|
111 | IF( .NOT. ln_limthd ) THEN ! we suppose ice is impermeable => ocean is isolated from atmosphere |
---|
112 | hfx_in (:,:) = pfrld(:,:) * ( qns_oce(:,:) + qsr_oce(:,:) ) + qemp_oce(:,:) |
---|
113 | hfx_out (:,:) = pfrld(:,:) * qns_oce(:,:) + qemp_oce(:,:) |
---|
114 | ftr_ice (:,:,:) = 0._wp |
---|
115 | emp_ice (:,:) = 0._wp |
---|
116 | qemp_ice (:,:) = 0._wp |
---|
117 | qevap_ice(:,:,:) = 0._wp |
---|
118 | ENDIF |
---|
119 | |
---|
120 | ! albedo output |
---|
121 | CALL wrk_alloc( jpi,jpj, zalb ) |
---|
122 | |
---|
123 | zalb(:,:) = 0._wp |
---|
124 | WHERE ( at_i_b <= epsi06 ) ; zalb(:,:) = 0.066_wp |
---|
125 | ELSEWHERE ; zalb(:,:) = SUM( alb_ice * a_i_b, dim=3 ) / at_i_b |
---|
126 | END WHERE |
---|
127 | IF( iom_use('alb_ice' ) ) CALL iom_put( "alb_ice" , zalb(:,:) ) ! ice albedo output |
---|
128 | |
---|
129 | zalb(:,:) = SUM( alb_ice * a_i_b, dim=3 ) + 0.066_wp * ( 1._wp - at_i_b ) |
---|
130 | IF( iom_use('albedo' ) ) CALL iom_put( "albedo" , zalb(:,:) ) ! ice albedo output |
---|
131 | |
---|
132 | CALL wrk_dealloc( jpi,jpj, zalb ) |
---|
133 | ! |
---|
134 | |
---|
135 | DO jj = 1, jpj |
---|
136 | DO ji = 1, jpi |
---|
137 | |
---|
138 | !------------------------------------------! |
---|
139 | ! heat flux at the ocean surface ! |
---|
140 | !------------------------------------------! |
---|
141 | ! Solar heat flux reaching the ocean = zqsr (W.m-2) |
---|
142 | !--------------------------------------------------- |
---|
143 | zqsr = qsr_tot(ji,jj) |
---|
144 | DO jl = 1, jpl |
---|
145 | zqsr = zqsr - a_i_b(ji,jj,jl) * ( qsr_ice(ji,jj,jl) - ftr_ice(ji,jj,jl) ) |
---|
146 | END DO |
---|
147 | |
---|
148 | ! Total heat flux reaching the ocean = hfx_out (W.m-2) |
---|
149 | !--------------------------------------------------- |
---|
150 | zqmass = hfx_thd(ji,jj) + hfx_dyn(ji,jj) + hfx_res(ji,jj) ! heat flux from snow is 0 (T=0 degC) |
---|
151 | hfx_out(ji,jj) = hfx_out(ji,jj) + zqmass + zqsr |
---|
152 | |
---|
153 | ! Add the residual from heat diffusion equation and sublimation (W.m-2) |
---|
154 | !---------------------------------------------------------------------- |
---|
155 | hfx_out(ji,jj) = hfx_out(ji,jj) + hfx_err_dif(ji,jj) + & |
---|
156 | & ( hfx_sub(ji,jj) - SUM( qevap_ice(ji,jj,:) * a_i_b(ji,jj,:) ) ) |
---|
157 | |
---|
158 | ! New qsr and qns used to compute the oceanic heat flux at the next time step |
---|
159 | !---------------------------------------------------------------------------- |
---|
160 | qsr(ji,jj) = zqsr |
---|
161 | qns(ji,jj) = hfx_out(ji,jj) - zqsr |
---|
162 | |
---|
163 | !------------------------------------------! |
---|
164 | ! mass flux at the ocean surface ! |
---|
165 | !------------------------------------------! |
---|
166 | ! case of realistic freshwater flux (Tartinville et al., 2001) (presently ACTIVATED) |
---|
167 | ! ------------------------------------------------------------------------------------- |
---|
168 | ! The idea of this approach is that the system that we consider is the ICE-OCEAN system |
---|
169 | ! Thus FW flux = External ( E-P+snow melt) |
---|
170 | ! Salt flux = Exchanges in the ice-ocean system then converted into FW |
---|
171 | ! Associated to Ice formation AND Ice melting |
---|
172 | ! Even if i see Ice melting as a FW and SALT flux |
---|
173 | ! |
---|
174 | ! mass flux from ice/ocean |
---|
175 | wfx_ice(ji,jj) = wfx_bog(ji,jj) + wfx_bom(ji,jj) + wfx_sum(ji,jj) + wfx_sni(ji,jj) & |
---|
176 | + wfx_opw(ji,jj) + wfx_dyn(ji,jj) + wfx_res(ji,jj) + wfx_lam(ji,jj) |
---|
177 | |
---|
178 | ! mass flux at the ocean/ice interface |
---|
179 | fmmflx(ji,jj) = - ( wfx_ice(ji,jj) + wfx_snw(ji,jj) + wfx_err_sub(ji,jj) ) ! F/M mass flux save at least for biogeochemical model |
---|
180 | emp(ji,jj) = emp_oce(ji,jj) - wfx_ice(ji,jj) - wfx_snw(ji,jj) - wfx_err_sub(ji,jj) ! mass flux + F/M mass flux (always ice/ocean mass exchange) |
---|
181 | END DO |
---|
182 | END DO |
---|
183 | |
---|
184 | !------------------------------------------! |
---|
185 | ! salt flux at the ocean surface ! |
---|
186 | !------------------------------------------! |
---|
187 | sfx(:,:) = sfx_bog(:,:) + sfx_bom(:,:) + sfx_sum(:,:) + sfx_sni(:,:) + sfx_opw(:,:) & |
---|
188 | & + sfx_res(:,:) + sfx_dyn(:,:) + sfx_bri(:,:) + sfx_sub(:,:) + sfx_lam(:,:) |
---|
189 | |
---|
190 | !-------------------------------------------------------------! |
---|
191 | ! mass of snow and ice per unit area for embedded sea-ice ! |
---|
192 | !-------------------------------------------------------------! |
---|
193 | IF( nn_ice_embd /= 0 ) THEN |
---|
194 | ! save mass from the previous ice time step |
---|
195 | snwice_mass_b(:,:) = snwice_mass(:,:) |
---|
196 | ! new mass per unit area |
---|
197 | snwice_mass (:,:) = tmask(:,:,1) * ( rhosn * vt_s(:,:) + rhoic * vt_i(:,:) ) |
---|
198 | ! time evolution of snow+ice mass |
---|
199 | snwice_fmass (:,:) = ( snwice_mass(:,:) - snwice_mass_b(:,:) ) * r1_rdtice |
---|
200 | ENDIF |
---|
201 | |
---|
202 | !-----------------------------------------------! |
---|
203 | ! Storing the transmitted variables ! |
---|
204 | !-----------------------------------------------! |
---|
205 | fr_i (:,:) = at_i(:,:) ! Sea-ice fraction |
---|
206 | tn_ice(:,:,:) = t_su(:,:,:) ! Ice surface temperature |
---|
207 | |
---|
208 | !------------------------------------------------------------------------! |
---|
209 | ! Snow/ice albedo (only if sent to coupler, useless in forced mode) ! |
---|
210 | !------------------------------------------------------------------------! |
---|
211 | CALL wrk_alloc( jpi, jpj, jpl, zalb_cs, zalb_os ) |
---|
212 | CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os ) ! cloud-sky and overcast-sky ice albedos |
---|
213 | alb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) |
---|
214 | CALL wrk_dealloc( jpi, jpj, jpl, zalb_cs, zalb_os ) |
---|
215 | |
---|
216 | ! conservation test |
---|
217 | IF( ln_limdiahsb ) CALL lim_cons_final( 'limsbc' ) |
---|
218 | |
---|
219 | ! control prints |
---|
220 | IF( ln_icectl ) CALL lim_prt( kt, iiceprt, jiceprt, 3, ' - Final state lim_sbc - ' ) |
---|
221 | |
---|
222 | IF(ln_ctl) THEN |
---|
223 | CALL prt_ctl( tab2d_1=qsr , clinfo1=' lim_sbc: qsr : ', tab2d_2=qns , clinfo2=' qns : ' ) |
---|
224 | CALL prt_ctl( tab2d_1=emp , clinfo1=' lim_sbc: emp : ', tab2d_2=sfx , clinfo2=' sfx : ' ) |
---|
225 | CALL prt_ctl( tab2d_1=fr_i , clinfo1=' lim_sbc: fr_i : ' ) |
---|
226 | CALL prt_ctl( tab3d_1=tn_ice, clinfo1=' lim_sbc: tn_ice : ', kdim=jpl ) |
---|
227 | ENDIF |
---|
228 | |
---|
229 | END SUBROUTINE lim_sbc_flx |
---|
230 | |
---|
231 | |
---|
232 | SUBROUTINE lim_sbc_tau( kt , pu_oce, pv_oce ) |
---|
233 | !!------------------------------------------------------------------- |
---|
234 | !! *** ROUTINE lim_sbc_tau *** |
---|
235 | !! |
---|
236 | !! ** Purpose : Update the ocean surface stresses due to the ice |
---|
237 | !! |
---|
238 | !! ** Action : * at each ice time step (every nn_fsbc time step): |
---|
239 | !! - compute the modulus of ice-ocean relative velocity |
---|
240 | !! (*rho*Cd) at T-point (C-grid) or I-point (B-grid) |
---|
241 | !! tmod_io = rhoco * | U_ice-U_oce | |
---|
242 | !! - update the modulus of stress at ocean surface |
---|
243 | !! taum = frld * taum + (1-frld) * tmod_io * | U_ice-U_oce | |
---|
244 | !! * at each ocean time step (every kt): |
---|
245 | !! compute linearized ice-ocean stresses as |
---|
246 | !! Utau = tmod_io * | U_ice - pU_oce | |
---|
247 | !! using instantaneous current ocean velocity (usually before) |
---|
248 | !! |
---|
249 | !! NB: - ice-ocean rotation angle no more allowed |
---|
250 | !! - here we make an approximation: taum is only computed every ice time step |
---|
251 | !! This avoids mutiple average to pass from T -> U,V grids and next from U,V grids |
---|
252 | !! to T grid. taum is used in TKE and GLS, which should not be too sensitive to this approximaton... |
---|
253 | !! |
---|
254 | !! ** Outputs : - utau, vtau : surface ocean i- and j-stress (u- & v-pts) updated with ice-ocean fluxes |
---|
255 | !! - taum : modulus of the surface ocean stress (T-point) updated with ice-ocean fluxes |
---|
256 | !!--------------------------------------------------------------------- |
---|
257 | INTEGER , INTENT(in) :: kt ! ocean time-step index |
---|
258 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pu_oce, pv_oce ! surface ocean currents |
---|
259 | !! |
---|
260 | INTEGER :: ji, jj ! dummy loop indices |
---|
261 | REAL(wp) :: zat_u, zutau_ice, zu_t, zmodt ! local scalar |
---|
262 | REAL(wp) :: zat_v, zvtau_ice, zv_t, zrhoco ! - - |
---|
263 | !!--------------------------------------------------------------------- |
---|
264 | zrhoco = rau0 * rn_cio |
---|
265 | ! |
---|
266 | IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN !== Ice time-step only ==! (i.e. surface module time-step) |
---|
267 | DO jj = 2, jpjm1 !* update the modulus of stress at ocean surface (T-point) |
---|
268 | DO ji = fs_2, fs_jpim1 |
---|
269 | ! ! 2*(U_ice-U_oce) at T-point |
---|
270 | zu_t = u_ice(ji,jj) + u_ice(ji-1,jj) - u_oce(ji,jj) - u_oce(ji-1,jj) |
---|
271 | zv_t = v_ice(ji,jj) + v_ice(ji,jj-1) - v_oce(ji,jj) - v_oce(ji,jj-1) |
---|
272 | ! ! |U_ice-U_oce|^2 |
---|
273 | zmodt = 0.25_wp * ( zu_t * zu_t + zv_t * zv_t ) |
---|
274 | ! ! update the ocean stress modulus |
---|
275 | taum(ji,jj) = ( 1._wp - at_i(ji,jj) ) * taum(ji,jj) + at_i(ji,jj) * zrhoco * zmodt |
---|
276 | tmod_io(ji,jj) = zrhoco * SQRT( zmodt ) ! rhoco * |U_ice-U_oce| at T-point |
---|
277 | END DO |
---|
278 | END DO |
---|
279 | CALL lbc_lnk_multi( taum, 'T', 1., tmod_io, 'T', 1. ) |
---|
280 | ! |
---|
281 | utau_oce(:,:) = utau(:,:) !* save the air-ocean stresses at ice time-step |
---|
282 | vtau_oce(:,:) = vtau(:,:) |
---|
283 | ! |
---|
284 | ENDIF |
---|
285 | ! |
---|
286 | ! !== every ocean time-step ==! |
---|
287 | ! |
---|
288 | DO jj = 2, jpjm1 !* update the stress WITHOUT a ice-ocean rotation angle |
---|
289 | DO ji = fs_2, fs_jpim1 ! Vect. Opt. |
---|
290 | zat_u = ( at_i(ji,jj) + at_i(ji+1,jj) ) * 0.5_wp ! ice area at u and V-points |
---|
291 | zat_v = ( at_i(ji,jj) + at_i(ji,jj+1) ) * 0.5_wp |
---|
292 | ! ! linearized quadratic drag formulation |
---|
293 | zutau_ice = 0.5_wp * ( tmod_io(ji,jj) + tmod_io(ji+1,jj) ) * ( u_ice(ji,jj) - pu_oce(ji,jj) ) |
---|
294 | zvtau_ice = 0.5_wp * ( tmod_io(ji,jj) + tmod_io(ji,jj+1) ) * ( v_ice(ji,jj) - pv_oce(ji,jj) ) |
---|
295 | ! ! stresses at the ocean surface |
---|
296 | utau(ji,jj) = ( 1._wp - zat_u ) * utau_oce(ji,jj) + zat_u * zutau_ice |
---|
297 | vtau(ji,jj) = ( 1._wp - zat_v ) * vtau_oce(ji,jj) + zat_v * zvtau_ice |
---|
298 | END DO |
---|
299 | END DO |
---|
300 | CALL lbc_lnk_multi( utau, 'U', -1., vtau, 'V', -1. ) ! lateral boundary condition |
---|
301 | ! |
---|
302 | IF(ln_ctl) CALL prt_ctl( tab2d_1=utau, clinfo1=' lim_sbc: utau : ', mask1=umask, & |
---|
303 | & tab2d_2=vtau, clinfo2=' vtau : ' , mask2=vmask ) |
---|
304 | ! |
---|
305 | END SUBROUTINE lim_sbc_tau |
---|
306 | |
---|
307 | |
---|
308 | SUBROUTINE lim_sbc_init |
---|
309 | !!------------------------------------------------------------------- |
---|
310 | !! *** ROUTINE lim_sbc_init *** |
---|
311 | !! |
---|
312 | !! ** Purpose : Preparation of the file ice_evolu for the output of |
---|
313 | !! the temporal evolution of key variables |
---|
314 | !! |
---|
315 | !! ** input : Namelist namicedia |
---|
316 | !!------------------------------------------------------------------- |
---|
317 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
318 | REAL(wp) :: zcoefu, zcoefv, zcoeff ! local scalar |
---|
319 | IF(lwp) WRITE(numout,*) |
---|
320 | IF(lwp) WRITE(numout,*) 'lim_sbc_init : LIM-3 sea-ice - surface boundary condition' |
---|
321 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~ ' |
---|
322 | |
---|
323 | ! ! allocate lim_sbc array |
---|
324 | IF( lim_sbc_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'lim_sbc_init : unable to allocate standard arrays' ) |
---|
325 | ! |
---|
326 | IF( .NOT. ln_rstart ) THEN |
---|
327 | ! ! embedded sea ice |
---|
328 | IF( nn_ice_embd /= 0 ) THEN ! mass exchanges between ice and ocean (case 1 or 2) set the snow+ice mass |
---|
329 | snwice_mass (:,:) = tmask(:,:,1) * ( rhosn * vt_s(:,:) + rhoic * vt_i(:,:) ) |
---|
330 | snwice_mass_b(:,:) = snwice_mass(:,:) |
---|
331 | ELSE |
---|
332 | snwice_mass (:,:) = 0.0_wp ! no mass exchanges |
---|
333 | snwice_mass_b(:,:) = 0.0_wp ! no mass exchanges |
---|
334 | ENDIF |
---|
335 | IF( nn_ice_embd == 2 ) THEN ! full embedment (case 2) deplete the initial ssh below sea-ice area |
---|
336 | sshn(:,:) = sshn(:,:) - snwice_mass(:,:) * r1_rau0 |
---|
337 | sshb(:,:) = sshb(:,:) - snwice_mass(:,:) * r1_rau0 |
---|
338 | #if defined key_vvl |
---|
339 | ! key_vvl necessary? clem: yes for compilation purpose |
---|
340 | DO jk = 1,jpkm1 ! adjust initial vertical scale factors |
---|
341 | fse3t_n(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshn(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) ) |
---|
342 | fse3t_b(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshb(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) ) |
---|
343 | ENDDO |
---|
344 | fse3t_a(:,:,:) = fse3t_b(:,:,:) |
---|
345 | ! Reconstruction of all vertical scale factors at now and before time |
---|
346 | ! steps |
---|
347 | ! ============================================================================= |
---|
348 | ! Horizontal scale factor interpolations |
---|
349 | ! -------------------------------------- |
---|
350 | CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3u_b(:,:,:), 'U' ) |
---|
351 | CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3v_b(:,:,:), 'V' ) |
---|
352 | CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3u_n(:,:,:), 'U' ) |
---|
353 | CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3v_n(:,:,:), 'V' ) |
---|
354 | CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3f_n(:,:,:), 'F' ) |
---|
355 | ! Vertical scale factor interpolations |
---|
356 | ! ------------------------------------ |
---|
357 | CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n (:,:,:), 'W' ) |
---|
358 | CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3uw_n(:,:,:), 'UW' ) |
---|
359 | CALL dom_vvl_interpol( fse3v_n(:,:,:), fse3vw_n(:,:,:), 'VW' ) |
---|
360 | CALL dom_vvl_interpol( fse3u_b(:,:,:), fse3uw_b(:,:,:), 'UW' ) |
---|
361 | CALL dom_vvl_interpol( fse3v_b(:,:,:), fse3vw_b(:,:,:), 'VW' ) |
---|
362 | ! t- and w- points depth |
---|
363 | ! ---------------------- |
---|
364 | fsdept_n(:,:,1) = 0.5_wp * fse3w_n(:,:,1) |
---|
365 | fsdepw_n(:,:,1) = 0.0_wp |
---|
366 | fsde3w_n(:,:,1) = fsdept_n(:,:,1) - sshn(:,:) |
---|
367 | DO jk = 2, jpk |
---|
368 | fsdept_n(:,:,jk) = fsdept_n(:,:,jk-1) + fse3w_n(:,:,jk) |
---|
369 | fsdepw_n(:,:,jk) = fsdepw_n(:,:,jk-1) + fse3t_n(:,:,jk-1) |
---|
370 | fsde3w_n(:,:,jk) = fsdept_n(:,:,jk ) - sshn (:,:) |
---|
371 | END DO |
---|
372 | #endif |
---|
373 | ENDIF |
---|
374 | ENDIF ! .NOT. ln_rstart |
---|
375 | ! |
---|
376 | |
---|
377 | END SUBROUTINE lim_sbc_init |
---|
378 | |
---|
379 | #else |
---|
380 | !!---------------------------------------------------------------------- |
---|
381 | !! Default option : Dummy module NO LIM 3.0 sea-ice model |
---|
382 | !!---------------------------------------------------------------------- |
---|
383 | CONTAINS |
---|
384 | SUBROUTINE lim_sbc ! Dummy routine |
---|
385 | END SUBROUTINE lim_sbc |
---|
386 | #endif |
---|
387 | |
---|
388 | !!====================================================================== |
---|
389 | END MODULE limsbc |
---|