1 | MODULE icestp |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE icestp *** |
---|
4 | !! sea ice : Master routine for all the sea ice model |
---|
5 | !!===================================================================== |
---|
6 | !! |
---|
7 | !! The sea ice model SI3 (Sea Ice modelling Integrated Initiative), |
---|
8 | !! aka Sea Ice cube for its nickname |
---|
9 | !! |
---|
10 | !! is originally based on LIM3, developed in Louvain-la-Neuve by: |
---|
11 | !! * Martin Vancoppenolle (UCL-ASTR, Belgium) |
---|
12 | !! * Sylvain Bouillon (UCL-ASTR, Belgium) |
---|
13 | !! * Miguel Angel Morales Maqueda (NOC-L, UK) |
---|
14 | !! thanks to valuable earlier work by |
---|
15 | !! * Thierry Fichefet |
---|
16 | !! * Hugues Goosse |
---|
17 | !! thanks also to the following persons who contributed |
---|
18 | !! * Gurvan Madec, Claude Talandier, Christian Ethe (LOCEAN, France) |
---|
19 | !! * Xavier Fettweis (UCL-ASTR), Ralph Timmermann (AWI, Germany) |
---|
20 | !! * Bill Lipscomb (LANL), Cecilia Bitz (UWa) and Elisabeth Hunke (LANL), USA. |
---|
21 | !! |
---|
22 | !! SI3 has been made possible by a handful of persons who met as working group |
---|
23 | !! (from France, Belgium, UK and Italy) |
---|
24 | !! * Clement Rousset, Martin Vancoppenolle & Gurvan Madec (LOCEAN, France) |
---|
25 | !! * Matthieu Chevalier & David Salas (Meteo France, France) |
---|
26 | !! * Gilles Garric (Mercator Ocean, France) |
---|
27 | !! * Thierry Fichefet & Francois Massonnet (UCL, Belgium) |
---|
28 | !! * Ed Blockley & Jeff Ridley (Met Office, UK) |
---|
29 | !! * Danny Feltham & David Schroeder (CPOM, UK) |
---|
30 | !! * Yevgeny Aksenov (NOC, UK) |
---|
31 | !! * Paul Holland (BAS, UK) |
---|
32 | !! * Dorotea Iovino (CMCC, Italy) |
---|
33 | !!====================================================================== |
---|
34 | !! History : 4.0 ! 2018 (C. Rousset) Original code SI3 |
---|
35 | !!---------------------------------------------------------------------- |
---|
36 | #if defined key_si3 |
---|
37 | !!---------------------------------------------------------------------- |
---|
38 | !! 'key_si3' SI3 sea-ice model |
---|
39 | !!---------------------------------------------------------------------- |
---|
40 | !! ice_stp : sea-ice model time-stepping and update ocean SBC over ice-covered area |
---|
41 | !! ice_init : initialize sea-ice |
---|
42 | !!---------------------------------------------------------------------- |
---|
43 | USE oce ! ocean dynamics and tracers |
---|
44 | USE dom_oce ! ocean space and time domain |
---|
45 | USE c1d ! 1D vertical configuration |
---|
46 | USE ice ! sea-ice: variables |
---|
47 | USE ice1D ! sea-ice: thermodynamical 1D variables |
---|
48 | ! |
---|
49 | USE phycst ! Define parameters for the routines |
---|
50 | USE eosbn2 ! equation of state |
---|
51 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
52 | USE sbc_ice ! Surface boundary condition: ice fields |
---|
53 | ! |
---|
54 | USE icesbc ! sea-ice: Surface boundary conditions |
---|
55 | USE icedyn ! sea-ice: dynamics |
---|
56 | USE icethd ! sea-ice: thermodynamics |
---|
57 | USE icecor ! sea-ice: corrections |
---|
58 | USE iceupdate ! sea-ice: sea surface boundary condition update |
---|
59 | USE icedia ! sea-ice: budget diagnostics |
---|
60 | USE icewri ! sea-ice: outputs |
---|
61 | USE icerst ! sea-ice: restarts |
---|
62 | USE icevar ! sea-ice: operations |
---|
63 | USE icectl ! sea-ice: control |
---|
64 | USE iceistate ! sea-ice: initial state |
---|
65 | USE iceitd ! sea-ice: remapping thickness distribution |
---|
66 | USE icealb ! sea-ice: albedo |
---|
67 | ! |
---|
68 | USE bdy_oce , ONLY : ln_bdy ! flag for bdy |
---|
69 | USE bdyice ! unstructured open boundary data for sea-ice |
---|
70 | # if defined key_agrif |
---|
71 | USE agrif_ice |
---|
72 | USE agrif_ice_interp |
---|
73 | # endif |
---|
74 | ! |
---|
75 | USE in_out_manager ! I/O manager |
---|
76 | USE iom ! I/O manager library |
---|
77 | USE lib_mpp ! MPP library |
---|
78 | USE lib_fortran ! fortran utilities (glob_sum + no signed zero) |
---|
79 | USE timing ! Timing |
---|
80 | USE prtctl ! Print control |
---|
81 | |
---|
82 | IMPLICIT NONE |
---|
83 | PRIVATE |
---|
84 | |
---|
85 | PUBLIC ice_stp ! called by sbcmod.F90 |
---|
86 | PUBLIC ice_init ! called by sbcmod.F90 |
---|
87 | |
---|
88 | !! * Substitutions |
---|
89 | # include "vectopt_loop_substitute.h90" |
---|
90 | !!---------------------------------------------------------------------- |
---|
91 | !! NEMO/ICE 4.0 , NEMO Consortium (2018) |
---|
92 | !! $Id$ |
---|
93 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
94 | !!---------------------------------------------------------------------- |
---|
95 | CONTAINS |
---|
96 | |
---|
97 | SUBROUTINE ice_stp( kt, ksbc ) |
---|
98 | !!--------------------------------------------------------------------- |
---|
99 | !! *** ROUTINE ice_stp *** |
---|
100 | !! |
---|
101 | !! ** Purpose : sea-ice model time-stepping and update ocean surface |
---|
102 | !! boundary condition over ice-covered area |
---|
103 | !! |
---|
104 | !! ** Method : ice model time stepping |
---|
105 | !! - call the ice dynamics routine |
---|
106 | !! - call the ice advection/diffusion routine |
---|
107 | !! - call the ice thermodynamics routine |
---|
108 | !! - call the routine that computes mass and |
---|
109 | !! heat fluxes at the ice/ocean interface |
---|
110 | !! - save the outputs |
---|
111 | !! - save the outputs for restart when necessary |
---|
112 | !! |
---|
113 | !! ** Action : - time evolution of the LIM sea-ice model |
---|
114 | !! - update all sbc variables below sea-ice: |
---|
115 | !! utau, vtau, taum, wndm, qns , qsr, emp , sfx |
---|
116 | !!--------------------------------------------------------------------- |
---|
117 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
118 | INTEGER, INTENT(in) :: ksbc ! flux formulation (user defined, bulk, or Pure Coupled) |
---|
119 | ! |
---|
120 | INTEGER :: jl ! dummy loop index |
---|
121 | !!---------------------------------------------------------------------- |
---|
122 | ! |
---|
123 | IF( ln_timing ) CALL timing_start('ice_stp') |
---|
124 | ! |
---|
125 | ! !-----------------------! |
---|
126 | IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN ! --- Ice time step --- ! |
---|
127 | ! !-----------------------! |
---|
128 | ! |
---|
129 | kt_ice = kt ! -- Ice model time step |
---|
130 | ! |
---|
131 | u_oce(:,:) = ssu_m(:,:) ! -- mean surface ocean current |
---|
132 | v_oce(:,:) = ssv_m(:,:) |
---|
133 | ! |
---|
134 | CALL eos_fzp( sss_m(:,:) , t_bo(:,:) ) ! -- freezing temperature [Kelvin] (set to rt0 over land) |
---|
135 | t_bo(:,:) = ( t_bo(:,:) + rt0 ) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) ) |
---|
136 | ! |
---|
137 | ! !== AGRIF Parent to Child ==! |
---|
138 | #if defined key_agrif |
---|
139 | ! ! nbstep_ice ranges from 1 to the nb of child ocean steps inside one parent ice step |
---|
140 | IF( .NOT. Agrif_Root() ) nbstep_ice = MOD( nbstep_ice, Agrif_irhot() * Agrif_Parent(nn_fsbc) / nn_fsbc ) + 1 |
---|
141 | ! ! these calls must remain here for restartability purposes |
---|
142 | CALL agrif_interp_ice( 'T' ) |
---|
143 | CALL agrif_interp_ice( 'U' ) |
---|
144 | CALL agrif_interp_ice( 'V' ) |
---|
145 | #endif |
---|
146 | CALL store_fields ! Store now ice values |
---|
147 | ! |
---|
148 | !------------------------------------------------! |
---|
149 | ! --- Dynamical coupling with the atmosphere --- ! |
---|
150 | !------------------------------------------------! |
---|
151 | ! It provides the following fields used in sea ice model: |
---|
152 | ! utau_ice, vtau_ice = surface ice stress [N/m2] |
---|
153 | !------------------------------------------------! |
---|
154 | CALL ice_sbc_tau( kt, ksbc, utau_ice, vtau_ice ) |
---|
155 | !-------------------------------------! |
---|
156 | ! --- ice dynamics and advection --- ! |
---|
157 | !-------------------------------------! |
---|
158 | CALL diag_set0 ! set diag of mass, heat and salt fluxes to 0 |
---|
159 | CALL ice_rst_opn( kt ) ! Open Ice restart file (if necessary) |
---|
160 | ! |
---|
161 | IF( ln_icedyn .AND. .NOT.lk_c1d ) & |
---|
162 | & CALL ice_dyn( kt ) ! -- Ice dynamics |
---|
163 | ! |
---|
164 | ! !== lateral boundary conditions ==! |
---|
165 | IF( ln_icethd .AND. ln_bdy ) CALL bdy_ice( kt ) ! -- bdy ice thermo |
---|
166 | ! |
---|
167 | ! !== previous lead fraction and ice volume for flux calculations |
---|
168 | CALL ice_var_glo2eqv ! h_i and h_s for ice albedo calculation |
---|
169 | CALL ice_var_agg(1) ! at_i for coupling |
---|
170 | CALL store_fields ! Store now ice values |
---|
171 | ! |
---|
172 | !------------------------------------------------------! |
---|
173 | ! --- Thermodynamical coupling with the atmosphere --- ! |
---|
174 | !------------------------------------------------------! |
---|
175 | ! It provides the following fields used in sea ice model: |
---|
176 | ! emp_oce , emp_ice = E-P over ocean and sea ice [Kg/m2/s] |
---|
177 | ! sprecip = solid precipitation [Kg/m2/s] |
---|
178 | ! evap_ice = sublimation [Kg/m2/s] |
---|
179 | ! qsr_tot , qns_tot = solar & non solar heat flux (total) [W/m2] |
---|
180 | ! qsr_ice , qns_ice = solar & non solar heat flux over ice [W/m2] |
---|
181 | ! dqns_ice = non solar heat sensistivity [W/m2] |
---|
182 | ! qemp_oce, qemp_ice, = sensible heat (associated with evap & precip) [W/m2] |
---|
183 | ! qprec_ice, qevap_ice |
---|
184 | !------------------------------------------------------! |
---|
185 | CALL ice_sbc_flx( kt, ksbc ) |
---|
186 | !----------------------------! |
---|
187 | ! --- ice thermodynamics --- ! |
---|
188 | !----------------------------! |
---|
189 | IF( ln_icethd ) CALL ice_thd( kt ) ! -- Ice thermodynamics |
---|
190 | ! |
---|
191 | CALL ice_cor( kt , 2 ) ! -- Corrections |
---|
192 | ! |
---|
193 | CALL ice_var_glo2eqv ! necessary calls (at least for coupling) |
---|
194 | CALL ice_var_agg( 2 ) ! necessary calls (at least for coupling) |
---|
195 | ! |
---|
196 | CALL ice_update_flx( kt ) ! -- Update ocean surface mass, heat and salt fluxes |
---|
197 | ! |
---|
198 | IF( ln_icediahsb ) CALL ice_dia( kt ) ! -- Diagnostics outputs |
---|
199 | ! |
---|
200 | CALL ice_wri( kt ) ! -- Ice outputs |
---|
201 | ! |
---|
202 | IF( lrst_ice ) CALL ice_rst_write( kt ) ! -- Ice restart file |
---|
203 | ! |
---|
204 | IF( ln_icectl ) CALL ice_ctl( kt ) ! -- alerts in case of model crash |
---|
205 | ! |
---|
206 | ENDIF ! End sea-ice time step only |
---|
207 | |
---|
208 | !-------------------------! |
---|
209 | ! --- Ocean time step --- ! |
---|
210 | !-------------------------! |
---|
211 | IF( ln_icedyn ) CALL ice_update_tau( kt, ub(:,:,1), vb(:,:,1) ) ! -- update surface ocean stresses |
---|
212 | !!gm remark, the ocean-ice stress is not saved in ice diag call above ..... find a solution!!! |
---|
213 | ! |
---|
214 | IF( ln_timing ) CALL timing_stop('ice_stp') |
---|
215 | ! |
---|
216 | END SUBROUTINE ice_stp |
---|
217 | |
---|
218 | |
---|
219 | SUBROUTINE ice_init |
---|
220 | !!---------------------------------------------------------------------- |
---|
221 | !! *** ROUTINE ice_init *** |
---|
222 | !! |
---|
223 | !! ** purpose : Initialize sea-ice parameters |
---|
224 | !!---------------------------------------------------------------------- |
---|
225 | INTEGER :: ji, jj, ierr |
---|
226 | !!---------------------------------------------------------------------- |
---|
227 | IF(lwp) WRITE(numout,*) |
---|
228 | IF(lwp) WRITE(numout,*) 'Sea Ice Model: SI3 (Sea Ice modelling Integrated Initiative)' |
---|
229 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~' |
---|
230 | IF(lwp) WRITE(numout,*) |
---|
231 | IF(lwp) WRITE(numout,*) 'ice_init: Arrays allocation & Initialization of all routines & init state' |
---|
232 | IF(lwp) WRITE(numout,*) '~~~~~~~~' |
---|
233 | ! |
---|
234 | ! ! Open the reference and configuration namelist files and namelist output file |
---|
235 | CALL ctl_opn( numnam_ice_ref, 'namelist_ice_ref', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) |
---|
236 | CALL ctl_opn( numnam_ice_cfg, 'namelist_ice_cfg', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) |
---|
237 | IF(lwm) CALL ctl_opn( numoni, 'output.namelist.ice', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, 1 ) |
---|
238 | ! |
---|
239 | CALL par_init ! set some ice run parameters |
---|
240 | ! |
---|
241 | ! ! Allocate the ice arrays (sbc_ice already allocated in sbc_init) |
---|
242 | ierr = ice_alloc () ! ice variables |
---|
243 | ierr = ierr + sbc_ice_alloc () ! surface boundary conditions |
---|
244 | ierr = ierr + ice1D_alloc () ! thermodynamics |
---|
245 | ! |
---|
246 | CALL mpp_sum( 'icestp', ierr ) |
---|
247 | IF( ierr /= 0 ) CALL ctl_stop('STOP', 'ice_init : unable to allocate ice arrays') |
---|
248 | ! |
---|
249 | CALL ice_itd_init ! ice thickness distribution initialization |
---|
250 | ! |
---|
251 | CALL ice_thd_init ! set ice thermodynics parameters (clem: important to call it first for melt ponds) |
---|
252 | ! |
---|
253 | ! ! Initial sea-ice state |
---|
254 | IF( .NOT. ln_rstart ) THEN ! start from rest: sea-ice deduced from sst |
---|
255 | CALL ice_istate_init |
---|
256 | CALL ice_istate( nit000 ) |
---|
257 | ELSE ! start from a restart file |
---|
258 | CALL ice_rst_read |
---|
259 | ENDIF |
---|
260 | CALL ice_var_glo2eqv |
---|
261 | CALL ice_var_agg(1) |
---|
262 | ! |
---|
263 | CALL ice_sbc_init ! set ice-ocean and ice-atm. coupling parameters |
---|
264 | ! |
---|
265 | CALL ice_dyn_init ! set ice dynamics parameters |
---|
266 | ! |
---|
267 | CALL ice_update_init ! ice surface boundary condition |
---|
268 | ! |
---|
269 | CALL ice_alb_init ! ice surface albedo |
---|
270 | ! |
---|
271 | CALL ice_dia_init ! initialization for diags |
---|
272 | ! |
---|
273 | fr_i (:,:) = at_i(:,:) ! initialisation of sea-ice fraction |
---|
274 | tn_ice(:,:,:) = t_su(:,:,:) ! initialisation of surface temp for coupled simu |
---|
275 | ! |
---|
276 | ! ! set max concentration in both hemispheres |
---|
277 | WHERE( gphit(:,:) > 0._wp ) ; rn_amax_2d(:,:) = rn_amax_n ! NH |
---|
278 | ELSEWHERE ; rn_amax_2d(:,:) = rn_amax_s ! SH |
---|
279 | END WHERE |
---|
280 | |
---|
281 | IF( ln_rstart ) CALL iom_close( numrir ) ! close input ice restart file |
---|
282 | ! |
---|
283 | END SUBROUTINE ice_init |
---|
284 | |
---|
285 | |
---|
286 | SUBROUTINE par_init |
---|
287 | !!------------------------------------------------------------------- |
---|
288 | !! *** ROUTINE par_init *** |
---|
289 | !! |
---|
290 | !! ** Purpose : Definition generic parameters for ice model |
---|
291 | !! |
---|
292 | !! ** Method : Read namelist and check the parameter |
---|
293 | !! values called at the first timestep (nit000) |
---|
294 | !! |
---|
295 | !! ** input : Namelist nampar |
---|
296 | !!------------------------------------------------------------------- |
---|
297 | INTEGER :: ios ! Local integer |
---|
298 | !! |
---|
299 | NAMELIST/nampar/ jpl, nlay_i, nlay_s, ln_virtual_itd, ln_icedyn, ln_icethd, rn_amax_n, rn_amax_s, & |
---|
300 | & cn_icerst_in, cn_icerst_indir, cn_icerst_out, cn_icerst_outdir |
---|
301 | !!------------------------------------------------------------------- |
---|
302 | ! |
---|
303 | REWIND( numnam_ice_ref ) ! Namelist nampar in reference namelist : Parameters for ice |
---|
304 | READ ( numnam_ice_ref, nampar, IOSTAT = ios, ERR = 901) |
---|
305 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampar in reference namelist' ) |
---|
306 | REWIND( numnam_ice_cfg ) ! Namelist nampar in configuration namelist : Parameters for ice |
---|
307 | READ ( numnam_ice_cfg, nampar, IOSTAT = ios, ERR = 902 ) |
---|
308 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'nampar in configuration namelist' ) |
---|
309 | IF(lwm) WRITE( numoni, nampar ) |
---|
310 | ! |
---|
311 | IF(lwp) THEN ! control print |
---|
312 | WRITE(numout,*) |
---|
313 | WRITE(numout,*) ' par_init: ice parameters shared among all the routines' |
---|
314 | WRITE(numout,*) ' ~~~~~~~~' |
---|
315 | WRITE(numout,*) ' Namelist nampar: ' |
---|
316 | WRITE(numout,*) ' number of ice categories jpl = ', jpl |
---|
317 | WRITE(numout,*) ' number of ice layers nlay_i = ', nlay_i |
---|
318 | WRITE(numout,*) ' number of snow layers nlay_s = ', nlay_s |
---|
319 | WRITE(numout,*) ' virtual ITD param for jpl=1 (T) or not (F) ln_virtual_itd = ', ln_virtual_itd |
---|
320 | WRITE(numout,*) ' Ice dynamics (T) or not (F) ln_icedyn = ', ln_icedyn |
---|
321 | WRITE(numout,*) ' Ice thermodynamics (T) or not (F) ln_icethd = ', ln_icethd |
---|
322 | WRITE(numout,*) ' maximum ice concentration for NH = ', rn_amax_n |
---|
323 | WRITE(numout,*) ' maximum ice concentration for SH = ', rn_amax_s |
---|
324 | ENDIF |
---|
325 | ! !--- change max ice concentration for roundoff errors |
---|
326 | rn_amax_n = MIN( rn_amax_n, 1._wp - epsi10 ) |
---|
327 | rn_amax_s = MIN( rn_amax_s, 1._wp - epsi10 ) |
---|
328 | ! !--- check consistency |
---|
329 | IF ( jpl > 1 .AND. ln_virtual_itd ) THEN |
---|
330 | ln_virtual_itd = .FALSE. |
---|
331 | IF(lwp) WRITE(numout,*) |
---|
332 | IF(lwp) WRITE(numout,*) ' ln_virtual_itd forced to false as jpl>1, no need with multiple categories to emulate them' |
---|
333 | ENDIF |
---|
334 | ! |
---|
335 | IF( ln_cpl .AND. nn_cats_cpl /= 1 .AND. nn_cats_cpl /= jpl ) THEN |
---|
336 | CALL ctl_stop( 'STOP', 'par_init: in coupled mode, nn_cats_cpl should be either 1 or jpl' ) |
---|
337 | ENDIF |
---|
338 | ! |
---|
339 | IF( ln_bdy .AND. ln_icediachk ) CALL ctl_warn('par_init: online conservation check does not work with BDY') |
---|
340 | ! |
---|
341 | rdt_ice = REAL(nn_fsbc) * rdt !--- sea-ice timestep and its inverse |
---|
342 | r1_rdtice = 1._wp / rdt_ice |
---|
343 | IF(lwp) WRITE(numout,*) |
---|
344 | IF(lwp) WRITE(numout,*) ' ice timestep rdt_ice = nn_fsbc*rdt = ', rdt_ice |
---|
345 | ! |
---|
346 | r1_nlay_i = 1._wp / REAL( nlay_i, wp ) !--- inverse of nlay_i and nlay_s |
---|
347 | r1_nlay_s = 1._wp / REAL( nlay_s, wp ) |
---|
348 | ! |
---|
349 | END SUBROUTINE par_init |
---|
350 | |
---|
351 | |
---|
352 | SUBROUTINE store_fields |
---|
353 | !!---------------------------------------------------------------------- |
---|
354 | !! *** ROUTINE store_fields *** |
---|
355 | !! |
---|
356 | !! ** purpose : store ice variables at "before" time step |
---|
357 | !!---------------------------------------------------------------------- |
---|
358 | INTEGER :: ji, jj, jl ! dummy loop index |
---|
359 | !!---------------------------------------------------------------------- |
---|
360 | ! |
---|
361 | a_i_b (:,:,:) = a_i (:,:,:) ! ice area |
---|
362 | v_i_b (:,:,:) = v_i (:,:,:) ! ice volume |
---|
363 | v_s_b (:,:,:) = v_s (:,:,:) ! snow volume |
---|
364 | sv_i_b(:,:,:) = sv_i(:,:,:) ! salt content |
---|
365 | oa_i_b(:,:,:) = oa_i(:,:,:) ! areal age content |
---|
366 | e_s_b (:,:,:,:) = e_s (:,:,:,:) ! snow thermal energy |
---|
367 | e_i_b (:,:,:,:) = e_i (:,:,:,:) ! ice thermal energy |
---|
368 | WHERE( a_i_b(:,:,:) >= epsi20 ) |
---|
369 | h_i_b(:,:,:) = v_i_b(:,:,:) / a_i_b(:,:,:) ! ice thickness |
---|
370 | h_s_b(:,:,:) = v_s_b(:,:,:) / a_i_b(:,:,:) ! snw thickness |
---|
371 | ELSEWHERE |
---|
372 | h_i_b(:,:,:) = 0._wp |
---|
373 | h_s_b(:,:,:) = 0._wp |
---|
374 | END WHERE |
---|
375 | |
---|
376 | WHERE( a_ip(:,:,:) >= epsi20 ) |
---|
377 | h_ip_b(:,:,:) = v_ip(:,:,:) / a_ip(:,:,:) ! ice pond thickness |
---|
378 | ELSEWHERE |
---|
379 | h_ip_b(:,:,:) = 0._wp |
---|
380 | END WHERE |
---|
381 | ! |
---|
382 | ! ice velocities & total concentration |
---|
383 | at_i_b(:,:) = SUM( a_i_b(:,:,:), dim=3 ) |
---|
384 | u_ice_b(:,:) = u_ice(:,:) |
---|
385 | v_ice_b(:,:) = v_ice(:,:) |
---|
386 | ! |
---|
387 | END SUBROUTINE store_fields |
---|
388 | |
---|
389 | |
---|
390 | SUBROUTINE diag_set0 |
---|
391 | !!---------------------------------------------------------------------- |
---|
392 | !! *** ROUTINE diag_set0 *** |
---|
393 | !! |
---|
394 | !! ** purpose : set ice-ocean and ice-atm. fluxes to zeros at the beggining |
---|
395 | !! of the time step |
---|
396 | !!---------------------------------------------------------------------- |
---|
397 | INTEGER :: ji, jj ! dummy loop index |
---|
398 | !!---------------------------------------------------------------------- |
---|
399 | sfx (:,:) = 0._wp ; |
---|
400 | sfx_bri(:,:) = 0._wp ; sfx_lam(:,:) = 0._wp |
---|
401 | sfx_sni(:,:) = 0._wp ; sfx_opw(:,:) = 0._wp |
---|
402 | sfx_bog(:,:) = 0._wp ; sfx_dyn(:,:) = 0._wp |
---|
403 | sfx_bom(:,:) = 0._wp ; sfx_sum(:,:) = 0._wp |
---|
404 | sfx_res(:,:) = 0._wp ; sfx_sub(:,:) = 0._wp |
---|
405 | ! |
---|
406 | wfx_snw(:,:) = 0._wp ; wfx_ice(:,:) = 0._wp |
---|
407 | wfx_sni(:,:) = 0._wp ; wfx_opw(:,:) = 0._wp |
---|
408 | wfx_bog(:,:) = 0._wp ; wfx_dyn(:,:) = 0._wp |
---|
409 | wfx_bom(:,:) = 0._wp ; wfx_sum(:,:) = 0._wp |
---|
410 | wfx_res(:,:) = 0._wp ; wfx_sub(:,:) = 0._wp |
---|
411 | wfx_spr(:,:) = 0._wp ; wfx_lam(:,:) = 0._wp |
---|
412 | wfx_snw_dyn(:,:) = 0._wp ; wfx_snw_sum(:,:) = 0._wp |
---|
413 | wfx_snw_sub(:,:) = 0._wp ; wfx_ice_sub(:,:) = 0._wp |
---|
414 | wfx_snw_sni(:,:) = 0._wp |
---|
415 | wfx_pnd(:,:) = 0._wp |
---|
416 | |
---|
417 | hfx_thd(:,:) = 0._wp ; |
---|
418 | hfx_snw(:,:) = 0._wp ; hfx_opw(:,:) = 0._wp |
---|
419 | hfx_bog(:,:) = 0._wp ; hfx_dyn(:,:) = 0._wp |
---|
420 | hfx_bom(:,:) = 0._wp ; hfx_sum(:,:) = 0._wp |
---|
421 | hfx_res(:,:) = 0._wp ; hfx_sub(:,:) = 0._wp |
---|
422 | hfx_spr(:,:) = 0._wp ; hfx_dif(:,:) = 0._wp |
---|
423 | hfx_err_rem(:,:) = 0._wp |
---|
424 | hfx_err_dif(:,:) = 0._wp |
---|
425 | wfx_err_sub(:,:) = 0._wp |
---|
426 | ! |
---|
427 | diag_heat(:,:) = 0._wp ; diag_sice(:,:) = 0._wp |
---|
428 | diag_vice(:,:) = 0._wp ; diag_vsnw(:,:) = 0._wp |
---|
429 | |
---|
430 | ! SIMIP diagnostics |
---|
431 | qcn_ice_bot(:,:,:) = 0._wp ; qcn_ice_top(:,:,:) = 0._wp ! conductive fluxes |
---|
432 | t_si (:,:,:) = rt0 ! temp at the ice-snow interface |
---|
433 | |
---|
434 | tau_icebfr (:,:) = 0._wp ! landfast ice param only (clem: important to keep the init here) |
---|
435 | cnd_ice (:,:,:) = 0._wp ! initialisation: effective conductivity at the top of ice/snow (ln_cndflx=T) |
---|
436 | qcn_ice (:,:,:) = 0._wp ! initialisation: conductive flux (ln_cndflx=T & ln_cndemule=T) |
---|
437 | qtr_ice_bot(:,:,:) = 0._wp ! initialization: part of solar radiation transmitted through the ice needed at least for outputs |
---|
438 | qsb_ice_bot(:,:) = 0._wp ! (needed if ln_icethd=F) |
---|
439 | ! |
---|
440 | ! for control checks (ln_icediachk) |
---|
441 | diag_trp_vi(:,:) = 0._wp ; diag_trp_vs(:,:) = 0._wp |
---|
442 | diag_trp_ei(:,:) = 0._wp ; diag_trp_es(:,:) = 0._wp |
---|
443 | diag_trp_sv(:,:) = 0._wp |
---|
444 | |
---|
445 | END SUBROUTINE diag_set0 |
---|
446 | |
---|
447 | #else |
---|
448 | !!---------------------------------------------------------------------- |
---|
449 | !! Default option Dummy module NO SI3 sea-ice model |
---|
450 | !!---------------------------------------------------------------------- |
---|
451 | CONTAINS |
---|
452 | SUBROUTINE ice_stp ( kt, ksbc ) ! Dummy routine |
---|
453 | INTEGER, INTENT(in) :: kt, ksbc |
---|
454 | WRITE(*,*) 'ice_stp: You should not have seen this print! error?', kt |
---|
455 | END SUBROUTINE ice_stp |
---|
456 | SUBROUTINE ice_init ! Dummy routine |
---|
457 | END SUBROUTINE ice_init |
---|
458 | #endif |
---|
459 | |
---|
460 | !!====================================================================== |
---|
461 | END MODULE icestp |
---|