1 | MODULE sbcice_lim |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE sbcice_lim *** |
---|
4 | !! Surface module : update the ocean surface boundary condition over ice |
---|
5 | !! & covered area using LIM sea-ice model |
---|
6 | !! Sea-Ice model : LIM-3 Sea ice model time-stepping |
---|
7 | !!===================================================================== |
---|
8 | !! History : 2.0 ! 2006-12 (M. Vancoppenolle) Original code |
---|
9 | !! 3.0 ! 2008-02 (C. Talandier) Surface module from icestp.F90 |
---|
10 | !! - ! 2008-04 (G. Madec) sltyle and lim_ctl routine |
---|
11 | !! 3.3 ! 2010-11 (G. Madec) ice-ocean stress always computed at each ocean time-step |
---|
12 | !! 3.4 ! 2011-01 (A Porter) dynamical allocation |
---|
13 | !! - ! 2012-10 (C. Rousset) add lim_diahsb |
---|
14 | !! 3.6 ! 2014-07 (M. Vancoppenolle, G. Madec, O. Marti) revise coupled interface |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | #if defined key_lim3 |
---|
17 | !!---------------------------------------------------------------------- |
---|
18 | !! 'key_lim3' : LIM 3.0 sea-ice model |
---|
19 | !!---------------------------------------------------------------------- |
---|
20 | !! sbc_ice_lim : sea-ice model time-stepping and update ocean sbc over ice-covered area |
---|
21 | !!---------------------------------------------------------------------- |
---|
22 | USE oce ! ocean dynamics and tracers |
---|
23 | USE dom_oce ! ocean space and time domain |
---|
24 | USE ice ! LIM-3: ice variables |
---|
25 | USE thd_ice ! LIM-3: thermodynamical variables |
---|
26 | |
---|
27 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
28 | USE sbc_ice ! Surface boundary condition: ice fields |
---|
29 | USE sbcblk_core ! Surface boundary condition: CORE bulk |
---|
30 | USE sbcblk_clio ! Surface boundary condition: CLIO bulk |
---|
31 | USE sbccpl ! Surface boundary condition: coupled interface |
---|
32 | USE albedo ! ocean & ice albedo |
---|
33 | |
---|
34 | USE phycst ! Define parameters for the routines |
---|
35 | USE eosbn2 ! equation of state |
---|
36 | USE limdyn ! Ice dynamics |
---|
37 | USE limtrp ! Ice transport |
---|
38 | USE limhdf ! Ice horizontal diffusion |
---|
39 | USE limthd ! Ice thermodynamics |
---|
40 | USE limitd_me ! Mechanics on ice thickness distribution |
---|
41 | USE limsbc ! sea surface boundary condition |
---|
42 | USE limdiahsb ! Ice budget diagnostics |
---|
43 | USE limwri ! Ice outputs |
---|
44 | USE limrst ! Ice restarts |
---|
45 | USE limupdate1 ! update of global variables |
---|
46 | USE limupdate2 ! update of global variables |
---|
47 | USE limvar ! Ice variables switch |
---|
48 | |
---|
49 | USE limistate ! LIM initial state |
---|
50 | USE limthd_sal ! LIM ice thermodynamics: salinity |
---|
51 | |
---|
52 | USE c1d ! 1D vertical configuration |
---|
53 | USE lbclnk ! lateral boundary condition - MPP link |
---|
54 | USE lib_mpp ! MPP library |
---|
55 | USE wrk_nemo ! work arrays |
---|
56 | USE timing ! Timing |
---|
57 | USE iom ! I/O manager library |
---|
58 | USE in_out_manager ! I/O manager |
---|
59 | USE prtctl ! Print control |
---|
60 | USE lib_fortran ! |
---|
61 | USE limctl |
---|
62 | |
---|
63 | #if defined key_bdy |
---|
64 | USE bdyice_lim ! unstructured open boundary data (bdy_ice_lim routine) |
---|
65 | #endif |
---|
66 | # if defined key_agrif |
---|
67 | USE agrif_ice |
---|
68 | USE agrif_lim3_update |
---|
69 | USE agrif_lim3_interp |
---|
70 | # endif |
---|
71 | |
---|
72 | IMPLICIT NONE |
---|
73 | PRIVATE |
---|
74 | |
---|
75 | PUBLIC sbc_ice_lim ! routine called by sbcmod.F90 |
---|
76 | PUBLIC sbc_lim_init ! routine called by sbcmod.F90 |
---|
77 | |
---|
78 | !! * Substitutions |
---|
79 | # include "domzgr_substitute.h90" |
---|
80 | # include "vectopt_loop_substitute.h90" |
---|
81 | !!---------------------------------------------------------------------- |
---|
82 | !! NEMO/OPA 4.0 , UCL NEMO Consortium (2011) |
---|
83 | !! $Id$ |
---|
84 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
85 | !!---------------------------------------------------------------------- |
---|
86 | CONTAINS |
---|
87 | |
---|
88 | !!====================================================================== |
---|
89 | |
---|
90 | SUBROUTINE sbc_ice_lim( kt, kblk ) |
---|
91 | !!--------------------------------------------------------------------- |
---|
92 | !! *** ROUTINE sbc_ice_lim *** |
---|
93 | !! |
---|
94 | !! ** Purpose : update the ocean surface boundary condition via the |
---|
95 | !! Louvain la Neuve Sea Ice Model time stepping |
---|
96 | !! |
---|
97 | !! ** Method : ice model time stepping |
---|
98 | !! - call the ice dynamics routine |
---|
99 | !! - call the ice advection/diffusion routine |
---|
100 | !! - call the ice thermodynamics routine |
---|
101 | !! - call the routine that computes mass and |
---|
102 | !! heat fluxes at the ice/ocean interface |
---|
103 | !! - save the outputs |
---|
104 | !! - save the outputs for restart when necessary |
---|
105 | !! |
---|
106 | !! ** Action : - time evolution of the LIM sea-ice model |
---|
107 | !! - update all sbc variables below sea-ice: |
---|
108 | !! utau, vtau, taum, wndm, qns , qsr, emp , sfx |
---|
109 | !!--------------------------------------------------------------------- |
---|
110 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
111 | INTEGER, INTENT(in) :: kblk ! type of bulk (=3 CLIO, =4 CORE, =5 COUPLED) |
---|
112 | !! |
---|
113 | INTEGER :: jl ! dummy loop index |
---|
114 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_os, zalb_cs ! ice albedo under overcast/clear sky |
---|
115 | REAL(wp), POINTER, DIMENSION(:,: ) :: zutau_ice, zvtau_ice |
---|
116 | !!---------------------------------------------------------------------- |
---|
117 | |
---|
118 | IF( nn_timing == 1 ) CALL timing_start('sbc_ice_lim') |
---|
119 | |
---|
120 | !-----------------------! |
---|
121 | ! --- Ice time step --- ! |
---|
122 | !-----------------------! |
---|
123 | IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN |
---|
124 | |
---|
125 | # if defined key_agrif |
---|
126 | IF( .NOT. Agrif_Root() ) THEN |
---|
127 | lim_nbstep = MOD( lim_nbstep, Agrif_irhot() * Agrif_Parent(nn_fsbc) / nn_fsbc ) + 1 |
---|
128 | ENDIF |
---|
129 | # endif |
---|
130 | |
---|
131 | ! mean surface ocean current at ice velocity point (C-grid dynamics : U- & V-points as the ocean) |
---|
132 | u_oce(:,:) = ssu_m(:,:) * umask(:,:,1) |
---|
133 | v_oce(:,:) = ssv_m(:,:) * vmask(:,:,1) |
---|
134 | |
---|
135 | ! masked sea surface freezing temperature [Kelvin] (set to rt0 over land) |
---|
136 | CALL eos_fzp( sss_m(:,:) , t_bo(:,:) ) |
---|
137 | t_bo(:,:) = ( t_bo(:,:) + rt0 ) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) ) |
---|
138 | |
---|
139 | ! Mask sea ice surface temperature (set to rt0 over land) |
---|
140 | DO jl = 1, jpl |
---|
141 | t_su(:,:,jl) = t_su(:,:,jl) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) ) |
---|
142 | END DO |
---|
143 | ! |
---|
144 | !------------------------------------------------! |
---|
145 | ! --- Dynamical coupling with the atmosphere --- ! |
---|
146 | !------------------------------------------------! |
---|
147 | ! It provides the following fields: |
---|
148 | ! utau_ice, vtau_ice : surface ice stress (U- & V-points) [N/m2] |
---|
149 | !----------------------------------------------------------------- |
---|
150 | SELECT CASE( kblk ) |
---|
151 | CASE( jp_clio ) ; CALL blk_ice_clio_tau ! CLIO bulk formulation |
---|
152 | CASE( jp_core ) ; CALL blk_ice_core_tau ! CORE bulk formulation |
---|
153 | CASE( jp_purecpl ) ; CALL sbc_cpl_ice_tau( utau_ice , vtau_ice ) ! Coupled formulation |
---|
154 | END SELECT |
---|
155 | |
---|
156 | IF( ln_mixcpl) THEN ! Case of a mixed Bulk/Coupled formulation |
---|
157 | CALL wrk_alloc( jpi,jpj , zutau_ice, zvtau_ice) |
---|
158 | CALL sbc_cpl_ice_tau( zutau_ice , zvtau_ice ) |
---|
159 | utau_ice(:,:) = utau_ice(:,:) * xcplmask(:,:,0) + zutau_ice(:,:) * ( 1. - xcplmask(:,:,0) ) |
---|
160 | vtau_ice(:,:) = vtau_ice(:,:) * xcplmask(:,:,0) + zvtau_ice(:,:) * ( 1. - xcplmask(:,:,0) ) |
---|
161 | CALL wrk_dealloc( jpi,jpj , zutau_ice, zvtau_ice) |
---|
162 | ENDIF |
---|
163 | |
---|
164 | !-------------------------------------------------------! |
---|
165 | ! --- ice dynamics and transport (except in 1D case) ---! |
---|
166 | !-------------------------------------------------------! |
---|
167 | numit = numit + nn_fsbc ! Ice model time step |
---|
168 | ! |
---|
169 | CALL sbc_lim_bef ! Store previous ice values |
---|
170 | CALL sbc_lim_diag0 ! set diag of mass, heat and salt fluxes to 0 |
---|
171 | CALL lim_rst_opn( kt ) ! Open Ice restart file |
---|
172 | ! |
---|
173 | #if defined key_agrif |
---|
174 | IF( .NOT. Agrif_Root() ) CALL agrif_interp_lim3('T') |
---|
175 | #endif |
---|
176 | ! --- zap this if no ice dynamics --- ! |
---|
177 | IF( .NOT. lk_c1d .AND. ln_limdyn ) THEN |
---|
178 | ! |
---|
179 | IF( nn_limdyn /= 0 ) THEN ! -- Ice dynamics |
---|
180 | CALL lim_dyn( kt ) ! rheology |
---|
181 | ELSE |
---|
182 | u_ice(:,:) = rn_uice * umask(:,:,1) ! or prescribed velocity |
---|
183 | v_ice(:,:) = rn_vice * vmask(:,:,1) |
---|
184 | ENDIF |
---|
185 | CALL lim_trp( kt ) ! -- Ice transport (Advection/diffusion) |
---|
186 | IF( nn_limdyn == 2 .AND. nn_monocat /= 2 ) & ! -- Mechanical redistribution (ridging/rafting) |
---|
187 | & CALL lim_itd_me |
---|
188 | IF( nn_limdyn == 2 ) CALL lim_update1( kt ) ! -- Corrections |
---|
189 | ! |
---|
190 | ENDIF |
---|
191 | ! --- |
---|
192 | #if defined key_agrif |
---|
193 | IF( .NOT. Agrif_Root() ) CALL agrif_interp_lim3('T') |
---|
194 | #endif |
---|
195 | #if defined key_bdy |
---|
196 | IF( ln_limthd ) CALL bdy_ice_lim( kt ) ! -- bdy ice thermo |
---|
197 | IF( ln_icectl ) CALL lim_prt( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' ) |
---|
198 | #endif |
---|
199 | |
---|
200 | ! previous lead fraction and ice volume for flux calculations |
---|
201 | CALL sbc_lim_bef |
---|
202 | CALL lim_var_glo2eqv ! ht_i and ht_s for ice albedo calculation |
---|
203 | CALL lim_var_agg(1) ! at_i for coupling (via pfrld) |
---|
204 | ! |
---|
205 | pfrld(:,:) = 1._wp - at_i(:,:) |
---|
206 | phicif(:,:) = vt_i(:,:) |
---|
207 | |
---|
208 | !------------------------------------------------------! |
---|
209 | ! --- Thermodynamical coupling with the atmosphere --- ! |
---|
210 | !------------------------------------------------------! |
---|
211 | ! It provides the following fields: |
---|
212 | ! qsr_ice , qns_ice : solar & non solar heat flux over ice (T-point) [W/m2] |
---|
213 | ! qla_ice : latent heat flux over ice (T-point) [W/m2] |
---|
214 | ! dqns_ice, dqla_ice : non solar & latent heat sensistivity (T-point) [W/m2] |
---|
215 | ! tprecip , sprecip : total & solid precipitation (T-point) [Kg/m2/s] |
---|
216 | ! fr1_i0 , fr2_i0 : 1sr & 2nd fraction of qsr penetration in ice [%] |
---|
217 | !---------------------------------------------------------------------------------------- |
---|
218 | CALL wrk_alloc( jpi,jpj,jpl, zalb_os, zalb_cs ) |
---|
219 | |
---|
220 | CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os ) ! cloud-sky and overcast-sky ice albedos |
---|
221 | SELECT CASE( kblk ) |
---|
222 | |
---|
223 | CASE( jp_clio ) ! CLIO bulk formulation |
---|
224 | ! In CLIO the cloud fraction is read in the climatology and the all-sky albedo |
---|
225 | ! (alb_ice) is computed within the bulk routine |
---|
226 | CALL blk_ice_clio_flx( t_su, zalb_cs, zalb_os, alb_ice ) |
---|
227 | IF( ln_mixcpl ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su ) |
---|
228 | IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) |
---|
229 | |
---|
230 | CASE( jp_core ) ! CORE bulk formulation |
---|
231 | ! albedo depends on cloud fraction because of non-linear spectral effects |
---|
232 | alb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) |
---|
233 | CALL blk_ice_core_flx( t_su, alb_ice ) |
---|
234 | IF( ln_mixcpl ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su ) |
---|
235 | IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) |
---|
236 | |
---|
237 | CASE ( jp_purecpl ) ! Coupled formulation |
---|
238 | ! albedo depends on cloud fraction because of non-linear spectral effects |
---|
239 | alb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) |
---|
240 | CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su ) |
---|
241 | IF( nn_limflx == 2 ) CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) |
---|
242 | |
---|
243 | END SELECT |
---|
244 | |
---|
245 | CALL wrk_dealloc( jpi,jpj,jpl, zalb_os, zalb_cs ) |
---|
246 | |
---|
247 | !----------------------------! |
---|
248 | ! --- ice thermodynamics --- ! |
---|
249 | !----------------------------! |
---|
250 | ! --- zap this if no ice thermo --- ! |
---|
251 | IF( ln_limthd ) CALL lim_thd( kt ) ! -- Ice thermodynamics |
---|
252 | IF( ln_limthd ) CALL lim_update2( kt ) ! -- Corrections |
---|
253 | ! --- |
---|
254 | # if defined key_agrif |
---|
255 | IF( .NOT. Agrif_Root() ) CALL agrif_update_lim3( kt ) |
---|
256 | # endif |
---|
257 | CALL lim_var_glo2eqv ! necessary calls (at least for coupling) |
---|
258 | CALL lim_var_agg( 2 ) ! necessary calls (at least for coupling) |
---|
259 | ! |
---|
260 | CALL lim_sbc_flx( kt ) ! -- Update surface ocean mass, heat and salt fluxes |
---|
261 | ! |
---|
262 | IF(ln_limdiaout) CALL lim_diahsb( kt ) ! -- Diagnostics and outputs |
---|
263 | ! |
---|
264 | CALL lim_wri( 1 ) ! -- Ice outputs |
---|
265 | ! |
---|
266 | IF( kt == nit000 .AND. ln_rstart ) & |
---|
267 | & CALL iom_close( numrir ) ! close input ice restart file |
---|
268 | ! |
---|
269 | IF( lrst_ice ) CALL lim_rst_write( kt ) ! -- Ice restart file |
---|
270 | ! |
---|
271 | IF( ln_icectl ) CALL lim_ctl( kt ) ! alerts in case of model crash |
---|
272 | ! |
---|
273 | ENDIF ! End sea-ice time step only |
---|
274 | |
---|
275 | !-------------------------! |
---|
276 | ! --- Ocean time step --- ! |
---|
277 | !-------------------------! |
---|
278 | ! Update surface ocean stresses (only in ice-dynamic case) otherwise the atm.-ocean stresses are used everywhere |
---|
279 | ! using before instantaneous surf. currents |
---|
280 | IF( ln_limdyn ) CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) ) |
---|
281 | !!gm remark, the ocean-ice stress is not saved in ice diag call above ..... find a solution!!! |
---|
282 | ! |
---|
283 | IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim') |
---|
284 | ! |
---|
285 | END SUBROUTINE sbc_ice_lim |
---|
286 | |
---|
287 | |
---|
288 | SUBROUTINE sbc_lim_init |
---|
289 | !!---------------------------------------------------------------------- |
---|
290 | !! *** ROUTINE sbc_lim_init *** |
---|
291 | !! |
---|
292 | !! ** purpose : Allocate all the dynamic arrays of the LIM-3 modules |
---|
293 | !!---------------------------------------------------------------------- |
---|
294 | INTEGER :: ierr |
---|
295 | INTEGER :: ji, jj |
---|
296 | !!---------------------------------------------------------------------- |
---|
297 | IF(lwp) WRITE(numout,*) |
---|
298 | IF(lwp) WRITE(numout,*) 'sbc_lim_init : update ocean surface boudary condition' |
---|
299 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ via Louvain la Neuve Ice Model (LIM-3) time stepping' |
---|
300 | ! |
---|
301 | ! Open the reference and configuration namelist files and namelist output file |
---|
302 | CALL ctl_opn( numnam_ice_ref, 'namelist_ice_ref', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) |
---|
303 | CALL ctl_opn( numnam_ice_cfg, 'namelist_ice_cfg', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) |
---|
304 | IF(lwm) CALL ctl_opn( numoni, 'output.namelist.ice', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, 1 ) |
---|
305 | |
---|
306 | CALL ice_run ! set some ice run parameters |
---|
307 | ! |
---|
308 | ! ! Allocate the ice arrays |
---|
309 | ierr = ice_alloc () ! ice variables |
---|
310 | ierr = ierr + sbc_ice_alloc () ! surface forcing |
---|
311 | ierr = ierr + thd_ice_alloc () ! thermodynamics |
---|
312 | IF( ln_limdyn ) ierr = ierr + lim_itd_me_alloc () ! ice thickness distribution - mechanics |
---|
313 | ! |
---|
314 | IF( lk_mpp ) CALL mpp_sum( ierr ) |
---|
315 | IF( ierr /= 0 ) CALL ctl_stop('STOP', 'sbc_lim_init : unable to allocate ice arrays') |
---|
316 | ! |
---|
317 | ! ! adequation jpk versus ice/snow layers/categories |
---|
318 | IF( jpl > jpk .OR. (nlay_i+1) > jpk .OR. nlay_s > jpk ) & |
---|
319 | & CALL ctl_stop( 'STOP', & |
---|
320 | & 'sbc_lim_init: the 3rd dimension of workspace arrays is too small.', & |
---|
321 | & 'use more ocean levels or less ice/snow layers/categories.' ) |
---|
322 | ! |
---|
323 | CALL lim_dyn_init ! set ice dynamics parameters |
---|
324 | ! |
---|
325 | CALL lim_itd_init ! ice thickness distribution initialization |
---|
326 | ! |
---|
327 | CALL lim_hdf_init ! set ice horizontal diffusion computation parameters |
---|
328 | ! |
---|
329 | CALL lim_thd_init ! set ice thermodynics parameters |
---|
330 | ! |
---|
331 | CALL lim_thd_sal_init ! set ice salinity parameters |
---|
332 | ! |
---|
333 | IF( ln_limdyn ) CALL lim_itd_me_init ! ice thickness distribution initialization for mecanical deformation |
---|
334 | ! ! Initial sea-ice state |
---|
335 | IF( .NOT. ln_rstart ) THEN ! start from rest: sea-ice deduced from sst |
---|
336 | numit = 0 |
---|
337 | numit = nit000 - 1 |
---|
338 | CALL lim_istate |
---|
339 | ELSE ! start from a restart file |
---|
340 | CALL lim_rst_read |
---|
341 | numit = nit000 - 1 |
---|
342 | ENDIF |
---|
343 | CALL lim_var_agg(2) |
---|
344 | CALL lim_var_glo2eqv |
---|
345 | ! |
---|
346 | CALL lim_sbc_init ! ice surface boundary condition |
---|
347 | ! |
---|
348 | IF( ln_limdiaout) CALL lim_diahsb_init ! initialization for diags |
---|
349 | ! |
---|
350 | fr_i(:,:) = at_i(:,:) ! initialisation of sea-ice fraction |
---|
351 | tn_ice(:,:,:) = t_su(:,:,:) ! initialisation of surface temp for coupled simu |
---|
352 | ! |
---|
353 | DO jj = 1, jpj |
---|
354 | DO ji = 1, jpi |
---|
355 | IF( gphit(ji,jj) > 0._wp ) THEN ; rn_amax_2d(ji,jj) = rn_amax_n ! NH |
---|
356 | ELSE ; rn_amax_2d(ji,jj) = rn_amax_s ! SH |
---|
357 | ENDIF |
---|
358 | ENDDO |
---|
359 | ENDDO |
---|
360 | ! |
---|
361 | nstart = numit + nn_fsbc |
---|
362 | nitrun = nitend - nit000 + 1 |
---|
363 | nlast = numit + nitrun |
---|
364 | ! |
---|
365 | IF( nstock == 0 ) nstock = nlast + 1 |
---|
366 | ! |
---|
367 | # if defined key_agrif |
---|
368 | IF( .NOT. Agrif_Root() ) CALL Agrif_InitValues_cont_lim3 |
---|
369 | # endif |
---|
370 | ! |
---|
371 | END SUBROUTINE sbc_lim_init |
---|
372 | |
---|
373 | |
---|
374 | SUBROUTINE ice_run |
---|
375 | !!------------------------------------------------------------------- |
---|
376 | !! *** ROUTINE ice_run *** |
---|
377 | !! |
---|
378 | !! ** Purpose : Definition some run parameter for ice model |
---|
379 | !! |
---|
380 | !! ** Method : Read the namicerun namelist and check the parameter |
---|
381 | !! values called at the first timestep (nit000) |
---|
382 | !! |
---|
383 | !! ** input : Namelist namicerun |
---|
384 | !!------------------------------------------------------------------- |
---|
385 | INTEGER :: ios ! Local integer output status for namelist read |
---|
386 | NAMELIST/namicerun/ jpl, nlay_i, nlay_s, rn_amax_n, rn_amax_s, cn_icerst_in, cn_icerst_indir, & |
---|
387 | & cn_icerst_out, cn_icerst_outdir, ln_limthd, ln_limdyn, nn_limdyn, rn_uice, rn_vice |
---|
388 | NAMELIST/namicediag/ ln_limdiahsb, ln_limdiaout, ln_icectl, iiceprt, jiceprt |
---|
389 | !!------------------------------------------------------------------- |
---|
390 | ! |
---|
391 | REWIND( numnam_ice_ref ) ! Namelist namicerun in reference namelist : Parameters for ice |
---|
392 | READ ( numnam_ice_ref, namicerun, IOSTAT = ios, ERR = 901) |
---|
393 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in reference namelist', lwp ) |
---|
394 | |
---|
395 | REWIND( numnam_ice_cfg ) ! Namelist namicerun in configuration namelist : Parameters for ice |
---|
396 | READ ( numnam_ice_cfg, namicerun, IOSTAT = ios, ERR = 902 ) |
---|
397 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in configuration namelist', lwp ) |
---|
398 | IF(lwm) WRITE ( numoni, namicerun ) |
---|
399 | ! |
---|
400 | REWIND( numnam_ice_ref ) ! Namelist namicediag in reference namelist : Parameters for ice |
---|
401 | READ ( numnam_ice_ref, namicediag, IOSTAT = ios, ERR = 903) |
---|
402 | 903 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicediag in reference namelist', lwp ) |
---|
403 | |
---|
404 | REWIND( numnam_ice_cfg ) ! Namelist namicediag in configuration namelist : Parameters for ice |
---|
405 | READ ( numnam_ice_cfg, namicediag, IOSTAT = ios, ERR = 904 ) |
---|
406 | 904 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicediag in configuration namelist', lwp ) |
---|
407 | IF(lwm) WRITE ( numoni, namicediag ) |
---|
408 | ! |
---|
409 | IF(lwp) THEN ! control print |
---|
410 | WRITE(numout,*) |
---|
411 | WRITE(numout,*) 'ice_run : ice share~d parameters for dynamics/advection/thermo of sea-ice' |
---|
412 | WRITE(numout,*) ' ~~~~~~' |
---|
413 | WRITE(numout,*) ' number of ice categories = ', jpl |
---|
414 | WRITE(numout,*) ' number of ice layers = ', nlay_i |
---|
415 | WRITE(numout,*) ' number of snow layers = ', nlay_s |
---|
416 | WRITE(numout,*) ' maximum ice concentration for NH = ', rn_amax_n |
---|
417 | WRITE(numout,*) ' maximum ice concentration for SH = ', rn_amax_s |
---|
418 | WRITE(numout,*) ' Ice thermodynamics (T) or not (F) ln_limthd = ', ln_limthd |
---|
419 | WRITE(numout,*) ' Ice dynamics (T) or not (F) ln_limdyn = ', ln_limdyn |
---|
420 | WRITE(numout,*) ' (ln_limdyn=T) Ice dynamics switch nn_limdyn = ', nn_limdyn |
---|
421 | WRITE(numout,*) ' 2: total' |
---|
422 | WRITE(numout,*) ' 1: advection only (no diffusion, no ridging/rafting)' |
---|
423 | WRITE(numout,*) ' 0: advection only (as 1 + prescribed velocity, bypass rheology)' |
---|
424 | WRITE(numout,*) ' (ln_limdyn=T) prescribed u-vel (case nn_limdyn=0) = ', rn_uice |
---|
425 | WRITE(numout,*) ' (ln_limdyn=T) prescribed v-vel (case nn_limdyn=0) = ', rn_vice |
---|
426 | WRITE(numout,*) |
---|
427 | WRITE(numout,*) '...and ice diagnostics' |
---|
428 | WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~' |
---|
429 | WRITE(numout,*) ' Diagnose heat/mass/salt budget or not ln_limdiahsb = ', ln_limdiahsb |
---|
430 | WRITE(numout,*) ' Output heat/mass/salt budget or not ln_limdiaout = ', ln_limdiaout |
---|
431 | WRITE(numout,*) ' control prints in ocean.out for (i,j)=(iiceprt,jiceprt) = ', ln_icectl |
---|
432 | WRITE(numout,*) ' i-index for control prints (ln_icectl=true) = ', iiceprt |
---|
433 | WRITE(numout,*) ' j-index for control prints (ln_icectl=true) = ', jiceprt |
---|
434 | ENDIF |
---|
435 | ! |
---|
436 | ! sea-ice timestep and inverse |
---|
437 | rdt_ice = nn_fsbc * rdttra(1) |
---|
438 | r1_rdtice = 1._wp / rdt_ice |
---|
439 | |
---|
440 | ! inverse of nlay_i and nlay_s |
---|
441 | r1_nlay_i = 1._wp / REAL( nlay_i, wp ) |
---|
442 | r1_nlay_s = 1._wp / REAL( nlay_s, wp ) |
---|
443 | ! |
---|
444 | #if defined key_bdy |
---|
445 | IF( lwp .AND. ln_limdiahsb ) CALL ctl_warn('online conservation check activated but it does not work with BDY') |
---|
446 | #endif |
---|
447 | ! |
---|
448 | END SUBROUTINE ice_run |
---|
449 | |
---|
450 | |
---|
451 | SUBROUTINE lim_itd_init |
---|
452 | !!------------------------------------------------------------------ |
---|
453 | !! *** ROUTINE lim_itd_init *** |
---|
454 | !! |
---|
455 | !! ** Purpose : Initializes the ice thickness distribution |
---|
456 | !! ** Method : ... |
---|
457 | !! ** input : Namelist namiceitd |
---|
458 | !!------------------------------------------------------------------- |
---|
459 | INTEGER :: ios ! Local integer output status for namelist read |
---|
460 | NAMELIST/namiceitd/ nn_catbnd, rn_himean |
---|
461 | ! |
---|
462 | INTEGER :: jl ! dummy loop index |
---|
463 | REAL(wp) :: zc1, zc2, zc3, zx1 ! local scalars |
---|
464 | REAL(wp) :: zhmax, znum, zden, zalpha ! |
---|
465 | !!------------------------------------------------------------------ |
---|
466 | ! |
---|
467 | REWIND( numnam_ice_ref ) ! Namelist namiceitd in reference namelist : Parameters for ice |
---|
468 | READ ( numnam_ice_ref, namiceitd, IOSTAT = ios, ERR = 905) |
---|
469 | 905 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in reference namelist', lwp ) |
---|
470 | |
---|
471 | REWIND( numnam_ice_cfg ) ! Namelist namiceitd in configuration namelist : Parameters for ice |
---|
472 | READ ( numnam_ice_cfg, namiceitd, IOSTAT = ios, ERR = 906 ) |
---|
473 | 906 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in configuration namelist', lwp ) |
---|
474 | IF(lwm) WRITE ( numoni, namiceitd ) |
---|
475 | ! |
---|
476 | IF(lwp) THEN ! control print |
---|
477 | WRITE(numout,*) |
---|
478 | WRITE(numout,*) 'lim_itd_init : Initialization of ice cat distribution ' |
---|
479 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
480 | WRITE(numout,*) ' shape of ice categories distribution nn_catbnd = ', nn_catbnd |
---|
481 | WRITE(numout,*) ' mean ice thickness in the domain (only active if nn_catbnd=2) rn_himean = ', rn_himean |
---|
482 | ENDIF |
---|
483 | |
---|
484 | !---------------------------------- |
---|
485 | !- Thickness categories boundaries |
---|
486 | !---------------------------------- |
---|
487 | hi_max(:) = 0._wp |
---|
488 | |
---|
489 | SELECT CASE ( nn_catbnd ) |
---|
490 | !---------------------- |
---|
491 | CASE (1) ! tanh function (CICE) |
---|
492 | !---------------------- |
---|
493 | zc1 = 3._wp / REAL( jpl, wp ) |
---|
494 | zc2 = 10._wp * zc1 |
---|
495 | zc3 = 3._wp |
---|
496 | |
---|
497 | DO jl = 1, jpl |
---|
498 | zx1 = REAL( jl-1, wp ) / REAL( jpl, wp ) |
---|
499 | hi_max(jl) = hi_max(jl-1) + zc1 + zc2 * (1._wp + TANH( zc3 * (zx1 - 1._wp ) ) ) |
---|
500 | END DO |
---|
501 | |
---|
502 | !---------------------- |
---|
503 | CASE (2) ! h^(-alpha) function |
---|
504 | !---------------------- |
---|
505 | zalpha = 0.05 ! exponent of the transform function |
---|
506 | |
---|
507 | zhmax = 3.*rn_himean |
---|
508 | |
---|
509 | DO jl = 1, jpl |
---|
510 | znum = jpl * ( zhmax+1 )**zalpha |
---|
511 | zden = ( jpl - jl ) * ( zhmax+1 )**zalpha + jl |
---|
512 | hi_max(jl) = ( znum / zden )**(1./zalpha) - 1 |
---|
513 | END DO |
---|
514 | |
---|
515 | END SELECT |
---|
516 | |
---|
517 | DO jl = 1, jpl |
---|
518 | hi_mean(jl) = ( hi_max(jl) + hi_max(jl-1) ) * 0.5_wp |
---|
519 | END DO |
---|
520 | |
---|
521 | ! Set hi_max(jpl) to a big value to ensure that all ice is thinner than hi_max(jpl) |
---|
522 | hi_max(jpl) = 99._wp |
---|
523 | |
---|
524 | IF(lwp) WRITE(numout,*) ' Thickness category boundaries ' |
---|
525 | IF(lwp) WRITE(numout,*) ' hi_max ', hi_max(0:jpl) |
---|
526 | ! |
---|
527 | END SUBROUTINE lim_itd_init |
---|
528 | |
---|
529 | |
---|
530 | SUBROUTINE ice_lim_flx( ptn_ice, palb_ice, pqns_ice, pqsr_ice, pdqn_ice, pevap_ice, pdevap_ice, k_limflx ) |
---|
531 | !!--------------------------------------------------------------------- |
---|
532 | !! *** ROUTINE ice_lim_flx *** |
---|
533 | !! |
---|
534 | !! ** Purpose : update the ice surface boundary condition by averaging and / or |
---|
535 | !! redistributing fluxes on ice categories |
---|
536 | !! |
---|
537 | !! ** Method : average then redistribute |
---|
538 | !! |
---|
539 | !! ** Action : |
---|
540 | !!--------------------------------------------------------------------- |
---|
541 | INTEGER , INTENT(in ) :: k_limflx ! =-1 do nothing; =0 average ; |
---|
542 | ! =1 average and redistribute ; =2 redistribute |
---|
543 | REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: ptn_ice ! ice surface temperature |
---|
544 | REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: palb_ice ! ice albedo |
---|
545 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqns_ice ! non solar flux |
---|
546 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqsr_ice ! net solar flux |
---|
547 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdqn_ice ! non solar flux sensitivity |
---|
548 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pevap_ice ! sublimation |
---|
549 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdevap_ice ! sublimation sensitivity |
---|
550 | ! |
---|
551 | INTEGER :: jl ! dummy loop index |
---|
552 | ! |
---|
553 | REAL(wp), POINTER, DIMENSION(:,:) :: zalb_m ! Mean albedo over all categories |
---|
554 | REAL(wp), POINTER, DIMENSION(:,:) :: ztem_m ! Mean temperature over all categories |
---|
555 | ! |
---|
556 | REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_m ! Mean solar heat flux over all categories |
---|
557 | REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_m ! Mean non solar heat flux over all categories |
---|
558 | REAL(wp), POINTER, DIMENSION(:,:) :: z_evap_m ! Mean sublimation over all categories |
---|
559 | REAL(wp), POINTER, DIMENSION(:,:) :: z_dqn_m ! Mean d(qns)/dT over all categories |
---|
560 | REAL(wp), POINTER, DIMENSION(:,:) :: z_devap_m ! Mean d(evap)/dT over all categories |
---|
561 | !!---------------------------------------------------------------------- |
---|
562 | |
---|
563 | IF( nn_timing == 1 ) CALL timing_start('ice_lim_flx') |
---|
564 | ! |
---|
565 | ! |
---|
566 | SELECT CASE( k_limflx ) !== averaged on all ice categories ==! |
---|
567 | CASE( 0 , 1 ) |
---|
568 | CALL wrk_alloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m) |
---|
569 | ! |
---|
570 | z_qns_m (:,:) = fice_ice_ave ( pqns_ice (:,:,:) ) |
---|
571 | z_qsr_m (:,:) = fice_ice_ave ( pqsr_ice (:,:,:) ) |
---|
572 | z_dqn_m (:,:) = fice_ice_ave ( pdqn_ice (:,:,:) ) |
---|
573 | z_evap_m (:,:) = fice_ice_ave ( pevap_ice (:,:,:) ) |
---|
574 | z_devap_m(:,:) = fice_ice_ave ( pdevap_ice (:,:,:) ) |
---|
575 | DO jl = 1, jpl |
---|
576 | pdqn_ice (:,:,jl) = z_dqn_m(:,:) |
---|
577 | pdevap_ice(:,:,jl) = z_devap_m(:,:) |
---|
578 | END DO |
---|
579 | ! |
---|
580 | DO jl = 1, jpl |
---|
581 | pqns_ice (:,:,jl) = z_qns_m(:,:) |
---|
582 | pqsr_ice (:,:,jl) = z_qsr_m(:,:) |
---|
583 | pevap_ice(:,:,jl) = z_evap_m(:,:) |
---|
584 | END DO |
---|
585 | ! |
---|
586 | CALL wrk_dealloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m) |
---|
587 | END SELECT |
---|
588 | |
---|
589 | SELECT CASE( k_limflx ) !== redistribution on all ice categories ==! |
---|
590 | CASE( 1 , 2 ) |
---|
591 | CALL wrk_alloc( jpi,jpj, zalb_m, ztem_m ) |
---|
592 | ! |
---|
593 | zalb_m(:,:) = fice_ice_ave ( palb_ice (:,:,:) ) |
---|
594 | ztem_m(:,:) = fice_ice_ave ( ptn_ice (:,:,:) ) |
---|
595 | DO jl = 1, jpl |
---|
596 | pqns_ice (:,:,jl) = pqns_ice (:,:,jl) + pdqn_ice (:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) ) |
---|
597 | pevap_ice(:,:,jl) = pevap_ice(:,:,jl) + pdevap_ice(:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) ) |
---|
598 | pqsr_ice (:,:,jl) = pqsr_ice (:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) ) |
---|
599 | END DO |
---|
600 | ! |
---|
601 | CALL wrk_dealloc( jpi,jpj, zalb_m, ztem_m ) |
---|
602 | END SELECT |
---|
603 | ! |
---|
604 | IF( nn_timing == 1 ) CALL timing_stop('ice_lim_flx') |
---|
605 | ! |
---|
606 | END SUBROUTINE ice_lim_flx |
---|
607 | |
---|
608 | SUBROUTINE sbc_lim_bef |
---|
609 | !!---------------------------------------------------------------------- |
---|
610 | !! *** ROUTINE sbc_lim_bef *** |
---|
611 | !! |
---|
612 | !! ** purpose : store ice variables at "before" time step |
---|
613 | !!---------------------------------------------------------------------- |
---|
614 | a_i_b (:,:,:) = a_i (:,:,:) ! ice area |
---|
615 | e_i_b (:,:,:,:) = e_i (:,:,:,:) ! ice thermal energy |
---|
616 | v_i_b (:,:,:) = v_i (:,:,:) ! ice volume |
---|
617 | v_s_b (:,:,:) = v_s (:,:,:) ! snow volume |
---|
618 | e_s_b (:,:,:,:) = e_s (:,:,:,:) ! snow thermal energy |
---|
619 | smv_i_b(:,:,:) = smv_i(:,:,:) ! salt content |
---|
620 | oa_i_b (:,:,:) = oa_i (:,:,:) ! areal age content |
---|
621 | u_ice_b(:,:) = u_ice(:,:) |
---|
622 | v_ice_b(:,:) = v_ice(:,:) |
---|
623 | at_i_b (:,:) = SUM( a_i_b(:,:,:), dim=3 ) |
---|
624 | |
---|
625 | END SUBROUTINE sbc_lim_bef |
---|
626 | |
---|
627 | SUBROUTINE sbc_lim_diag0 |
---|
628 | !!---------------------------------------------------------------------- |
---|
629 | !! *** ROUTINE sbc_lim_diag0 *** |
---|
630 | !! |
---|
631 | !! ** purpose : set ice-ocean and ice-atm. fluxes to zeros at the beggining |
---|
632 | !! of the time step |
---|
633 | !!---------------------------------------------------------------------- |
---|
634 | sfx (:,:) = 0._wp ; |
---|
635 | sfx_bri(:,:) = 0._wp ; sfx_lam(:,:) = 0._wp |
---|
636 | sfx_sni(:,:) = 0._wp ; sfx_opw(:,:) = 0._wp |
---|
637 | sfx_bog(:,:) = 0._wp ; sfx_dyn(:,:) = 0._wp |
---|
638 | sfx_bom(:,:) = 0._wp ; sfx_sum(:,:) = 0._wp |
---|
639 | sfx_res(:,:) = 0._wp ; sfx_sub(:,:) = 0._wp |
---|
640 | |
---|
641 | wfx_snw(:,:) = 0._wp ; wfx_ice(:,:) = 0._wp |
---|
642 | wfx_sni(:,:) = 0._wp ; wfx_opw(:,:) = 0._wp |
---|
643 | wfx_bog(:,:) = 0._wp ; wfx_dyn(:,:) = 0._wp |
---|
644 | wfx_bom(:,:) = 0._wp ; wfx_sum(:,:) = 0._wp |
---|
645 | wfx_res(:,:) = 0._wp ; wfx_sub(:,:) = 0._wp |
---|
646 | wfx_spr(:,:) = 0._wp ; wfx_lam(:,:) = 0._wp |
---|
647 | |
---|
648 | hfx_thd(:,:) = 0._wp ; |
---|
649 | hfx_snw(:,:) = 0._wp ; hfx_opw(:,:) = 0._wp |
---|
650 | hfx_bog(:,:) = 0._wp ; hfx_dyn(:,:) = 0._wp |
---|
651 | hfx_bom(:,:) = 0._wp ; hfx_sum(:,:) = 0._wp |
---|
652 | hfx_res(:,:) = 0._wp ; hfx_sub(:,:) = 0._wp |
---|
653 | hfx_spr(:,:) = 0._wp ; hfx_dif(:,:) = 0._wp |
---|
654 | hfx_err(:,:) = 0._wp ; hfx_err_rem(:,:) = 0._wp |
---|
655 | hfx_err_dif(:,:) = 0._wp |
---|
656 | wfx_err_sub(:,:) = 0._wp |
---|
657 | |
---|
658 | afx_tot(:,:) = 0._wp ; |
---|
659 | afx_dyn(:,:) = 0._wp ; afx_thd(:,:) = 0._wp |
---|
660 | |
---|
661 | diag_heat(:,:) = 0._wp ; diag_smvi(:,:) = 0._wp ; |
---|
662 | diag_vice(:,:) = 0._wp ; diag_vsnw(:,:) = 0._wp ; |
---|
663 | |
---|
664 | tau_icebfr(:,:) = 0._wp; ! landfast ice param only (clem: important to keep the init here) |
---|
665 | |
---|
666 | END SUBROUTINE sbc_lim_diag0 |
---|
667 | |
---|
668 | |
---|
669 | FUNCTION fice_cell_ave ( ptab ) |
---|
670 | !!-------------------------------------------------------------------------- |
---|
671 | !! * Compute average over categories, for grid cell (ice covered and free ocean) |
---|
672 | !!-------------------------------------------------------------------------- |
---|
673 | REAL (wp), DIMENSION (jpi,jpj) :: fice_cell_ave |
---|
674 | REAL (wp), DIMENSION (jpi,jpj,jpl), INTENT (in) :: ptab |
---|
675 | INTEGER :: jl ! Dummy loop index |
---|
676 | |
---|
677 | fice_cell_ave (:,:) = 0.0_wp |
---|
678 | DO jl = 1, jpl |
---|
679 | fice_cell_ave (:,:) = fice_cell_ave (:,:) + a_i (:,:,jl) * ptab (:,:,jl) |
---|
680 | END DO |
---|
681 | |
---|
682 | END FUNCTION fice_cell_ave |
---|
683 | |
---|
684 | |
---|
685 | FUNCTION fice_ice_ave ( ptab ) |
---|
686 | !!-------------------------------------------------------------------------- |
---|
687 | !! * Compute average over categories, for ice covered part of grid cell |
---|
688 | !!-------------------------------------------------------------------------- |
---|
689 | REAL (kind=wp), DIMENSION (jpi,jpj) :: fice_ice_ave |
---|
690 | REAL (kind=wp), DIMENSION (jpi,jpj,jpl), INTENT(in) :: ptab |
---|
691 | |
---|
692 | fice_ice_ave (:,:) = 0.0_wp |
---|
693 | WHERE ( at_i (:,:) > 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:) |
---|
694 | |
---|
695 | END FUNCTION fice_ice_ave |
---|
696 | |
---|
697 | |
---|
698 | #else |
---|
699 | !!---------------------------------------------------------------------- |
---|
700 | !! Default option Dummy module NO LIM 3.0 sea-ice model |
---|
701 | !!---------------------------------------------------------------------- |
---|
702 | CONTAINS |
---|
703 | SUBROUTINE sbc_ice_lim ( kt, kblk ) ! Dummy routine |
---|
704 | WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk |
---|
705 | END SUBROUTINE sbc_ice_lim |
---|
706 | SUBROUTINE sbc_lim_init ! Dummy routine |
---|
707 | END SUBROUTINE sbc_lim_init |
---|
708 | #endif |
---|
709 | |
---|
710 | !!====================================================================== |
---|
711 | END MODULE sbcice_lim |
---|