1 | MODULE icedyn |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE icedyn *** |
---|
4 | !! Sea-Ice dynamics : master routine for sea ice dynamics |
---|
5 | !!====================================================================== |
---|
6 | !! history : 4.0 ! 2018 (C. Rousset) original code SI3 [aka Sea Ice cube] |
---|
7 | !!---------------------------------------------------------------------- |
---|
8 | #if defined key_si3 |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | !! 'key_si3' SI3 sea-ice model |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | !! ice_dyn : dynamics of sea ice |
---|
13 | !! ice_dyn_init : initialization and namelist read |
---|
14 | !!---------------------------------------------------------------------- |
---|
15 | USE phycst ! physical constants |
---|
16 | USE dom_oce ! ocean space and time domain |
---|
17 | USE ice ! sea-ice: variables |
---|
18 | USE icedyn_rhg ! sea-ice: rheology |
---|
19 | USE icedyn_adv ! sea-ice: advection |
---|
20 | USE icedyn_rdgrft ! sea-ice: ridging/rafting |
---|
21 | USE icecor ! sea-ice: corrections |
---|
22 | USE icevar ! sea-ice: operations |
---|
23 | USE icectl ! sea-ice: control prints |
---|
24 | ! |
---|
25 | USE in_out_manager ! I/O manager |
---|
26 | USE iom ! I/O manager library |
---|
27 | USE lib_mpp ! MPP library |
---|
28 | USE lib_fortran ! fortran utilities (glob_sum + no signed zero) |
---|
29 | USE lbclnk ! lateral boundary conditions (or mpp links) |
---|
30 | USE timing ! Timing |
---|
31 | |
---|
32 | IMPLICIT NONE |
---|
33 | PRIVATE |
---|
34 | |
---|
35 | PUBLIC ice_dyn ! called by icestp.F90 |
---|
36 | PUBLIC ice_dyn_init ! called by icestp.F90 |
---|
37 | |
---|
38 | INTEGER :: nice_dyn ! choice of the type of dynamics |
---|
39 | ! ! associated indices: |
---|
40 | INTEGER, PARAMETER :: np_dynALL = 1 ! full ice dynamics (rheology + advection + ridging/rafting + correction) |
---|
41 | INTEGER, PARAMETER :: np_dynRHGADV = 2 ! pure dynamics (rheology + advection) |
---|
42 | INTEGER, PARAMETER :: np_dynADV1D = 3 ! only advection 1D - test case from Schar & Smolarkiewicz 1996 |
---|
43 | INTEGER, PARAMETER :: np_dynADV2D = 4 ! only advection 2D w prescribed vel.(rn_uvice + advection) |
---|
44 | ! |
---|
45 | ! ** namelist (namdyn) ** |
---|
46 | LOGICAL :: ln_dynALL ! full ice dynamics (rheology + advection + ridging/rafting + correction) |
---|
47 | LOGICAL :: ln_dynRHGADV ! no ridge/raft & no corrections (rheology + advection) |
---|
48 | LOGICAL :: ln_dynADV1D ! only advection in 1D w ice convergence (test case from Schar & Smolarkiewicz 1996) |
---|
49 | LOGICAL :: ln_dynADV2D ! only advection in 2D w prescribed vel. (rn_uvice + advection) |
---|
50 | REAL(wp) :: rn_uice ! prescribed u-vel (case np_dynADV1D & np_dynADV2D) |
---|
51 | REAL(wp) :: rn_vice ! prescribed v-vel (case np_dynADV1D & np_dynADV2D) |
---|
52 | |
---|
53 | !! * Substitutions |
---|
54 | # include "vectopt_loop_substitute.h90" |
---|
55 | !!---------------------------------------------------------------------- |
---|
56 | !! NEMO/ICE 4.0 , NEMO Consortium (2018) |
---|
57 | !! $Id: icedyn.F90 8378 2017-07-26 13:55:59Z clem $ |
---|
58 | !! Software governed by the CeCILL licence (./LICENSE) |
---|
59 | !!---------------------------------------------------------------------- |
---|
60 | CONTAINS |
---|
61 | |
---|
62 | SUBROUTINE ice_dyn( kt ) |
---|
63 | !!------------------------------------------------------------------- |
---|
64 | !! *** ROUTINE ice_dyn *** |
---|
65 | !! |
---|
66 | !! ** Purpose : this routine manages sea ice dynamics |
---|
67 | !! |
---|
68 | !! ** Action : - calculation of friction in case of landfast ice |
---|
69 | !! - call ice_dyn_rhg = rheology |
---|
70 | !! - call ice_dyn_adv = advection |
---|
71 | !! - call ice_dyn_rdgrft = ridging/rafting |
---|
72 | !! - call ice_cor = corrections if fields are out of bounds |
---|
73 | !!-------------------------------------------------------------------- |
---|
74 | INTEGER, INTENT(in) :: kt ! ice time step |
---|
75 | !! |
---|
76 | INTEGER :: ji, jj, jl ! dummy loop indices |
---|
77 | REAL(wp) :: zcoefu, zcoefv |
---|
78 | REAL(wp), DIMENSION(jpi,jpj,jpl) :: zhi_max, zhs_max |
---|
79 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zdivu_i |
---|
80 | !!-------------------------------------------------------------------- |
---|
81 | ! |
---|
82 | ! controls |
---|
83 | IF( ln_timing ) CALL timing_start('icedyn') ! timing |
---|
84 | IF( ln_icediachk ) CALL ice_cons_hsm(0, 'icedyn', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation |
---|
85 | ! |
---|
86 | IF( kt == nit000 .AND. lwp ) THEN |
---|
87 | WRITE(numout,*) |
---|
88 | WRITE(numout,*)'ice_dyn: sea-ice dynamics' |
---|
89 | WRITE(numout,*)'~~~~~~~' |
---|
90 | ENDIF |
---|
91 | ! |
---|
92 | IF( ln_landfast_home ) THEN !-- Landfast ice parameterization |
---|
93 | tau_icebfr(:,:) = 0._wp |
---|
94 | DO jl = 1, jpl |
---|
95 | WHERE( h_i_b(:,:,jl) > ht_n(:,:) * rn_depfra ) tau_icebfr(:,:) = tau_icebfr(:,:) + a_i(:,:,jl) * rn_icebfr |
---|
96 | END DO |
---|
97 | ENDIF |
---|
98 | ! |
---|
99 | ! !-- Record max of the surrounding 9-pts ice thick. (for CALL Hbig) |
---|
100 | DO jl = 1, jpl |
---|
101 | DO jj = 2, jpjm1 |
---|
102 | DO ji = fs_2, fs_jpim1 |
---|
103 | zhi_max(ji,jj,jl) = MAX( epsi20, h_i_b(ji,jj,jl), h_i_b(ji+1,jj ,jl), h_i_b(ji ,jj+1,jl), & |
---|
104 | & h_i_b(ji-1,jj ,jl), h_i_b(ji ,jj-1,jl), & |
---|
105 | & h_i_b(ji+1,jj+1,jl), h_i_b(ji-1,jj-1,jl), & |
---|
106 | & h_i_b(ji+1,jj-1,jl), h_i_b(ji-1,jj+1,jl) ) |
---|
107 | zhs_max(ji,jj,jl) = MAX( epsi20, h_s_b(ji,jj,jl), h_s_b(ji+1,jj ,jl), h_s_b(ji ,jj+1,jl), & |
---|
108 | & h_s_b(ji-1,jj ,jl), h_s_b(ji ,jj-1,jl), & |
---|
109 | & h_s_b(ji+1,jj+1,jl), h_s_b(ji-1,jj-1,jl), & |
---|
110 | & h_s_b(ji+1,jj-1,jl), h_s_b(ji-1,jj+1,jl) ) |
---|
111 | END DO |
---|
112 | END DO |
---|
113 | END DO |
---|
114 | CALL lbc_lnk_multi( zhi_max(:,:,:), 'T', 1., zhs_max(:,:,:), 'T', 1. ) |
---|
115 | ! |
---|
116 | ! |
---|
117 | SELECT CASE( nice_dyn ) !-- Set which dynamics is running |
---|
118 | |
---|
119 | CASE ( np_dynALL ) !== all dynamical processes ==! |
---|
120 | CALL ice_dyn_rhg ( kt ) ! -- rheology |
---|
121 | CALL ice_dyn_adv ( kt ) ; CALL Hbig( zhi_max, zhs_max ) ! -- advection of ice + correction on ice thickness |
---|
122 | CALL ice_dyn_rdgrft( kt ) ! -- ridging/rafting |
---|
123 | CALL ice_cor ( kt , 1 ) ! -- Corrections |
---|
124 | |
---|
125 | CASE ( np_dynRHGADV ) !== no ridge/raft & no corrections ==! |
---|
126 | CALL ice_dyn_rhg ( kt ) ! -- rheology |
---|
127 | CALL ice_dyn_adv ( kt ) ; CALL Hbig( zhi_max, zhs_max ) ! -- advection of ice + correction on ice thickness |
---|
128 | CALL Hpiling ! -- simple pile-up (replaces ridging/rafting) |
---|
129 | |
---|
130 | CASE ( np_dynADV1D ) !== pure advection ==! (1D) |
---|
131 | ALLOCATE( zdivu_i(jpi,jpj) ) |
---|
132 | ! --- monotonicity test from Schar & Smolarkiewicz 1996 --- ! |
---|
133 | ! CFL = 0.5 at a distance from the bound of 1/6 of the basin length |
---|
134 | ! Then for dx = 2m and dt = 1s => rn_uice = u (1/6th) = 1m/s |
---|
135 | DO jj = 1, jpj |
---|
136 | DO ji = 1, jpi |
---|
137 | zcoefu = ( REAL(jpiglo+1)*0.5 - REAL(ji+nimpp-1) ) / ( REAL(jpiglo+1)*0.5 - 1. ) |
---|
138 | zcoefv = ( REAL(jpjglo+1)*0.5 - REAL(jj+njmpp-1) ) / ( REAL(jpjglo+1)*0.5 - 1. ) |
---|
139 | u_ice(ji,jj) = rn_uice * 1.5 * SIGN( 1., zcoefu ) * ABS( zcoefu ) * umask(ji,jj,1) |
---|
140 | v_ice(ji,jj) = rn_vice * 1.5 * SIGN( 1., zcoefv ) * ABS( zcoefv ) * vmask(ji,jj,1) |
---|
141 | END DO |
---|
142 | END DO |
---|
143 | ! --- |
---|
144 | CALL ice_dyn_adv ( kt ) ! -- advection of ice + correction on ice thickness |
---|
145 | |
---|
146 | ! diagnostics: divergence at T points |
---|
147 | DO jj = 2, jpjm1 |
---|
148 | DO ji = 2, jpim1 |
---|
149 | zdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj) & |
---|
150 | & + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1) ) * r1_e1e2t(ji,jj) |
---|
151 | END DO |
---|
152 | END DO |
---|
153 | CALL lbc_lnk( zdivu_i, 'T', 1. ) |
---|
154 | IF( iom_use('icediv') ) CALL iom_put( "icediv" , zdivu_i(:,:) ) |
---|
155 | |
---|
156 | DEALLOCATE( zdivu_i ) |
---|
157 | |
---|
158 | CASE ( np_dynADV2D ) !== pure advection ==! (2D w prescribed velocities) |
---|
159 | ALLOCATE( zdivu_i(jpi,jpj) ) |
---|
160 | u_ice(:,:) = rn_uice * umask(:,:,1) |
---|
161 | v_ice(:,:) = rn_vice * vmask(:,:,1) |
---|
162 | !CALL RANDOM_NUMBER(u_ice(:,:)) ; u_ice(:,:) = u_ice(:,:) * 0.1 + rn_uice * 0.9 * umask(:,:,1) |
---|
163 | !CALL RANDOM_NUMBER(v_ice(:,:)) ; v_ice(:,:) = v_ice(:,:) * 0.1 + rn_vice * 0.9 * vmask(:,:,1) |
---|
164 | ! --- |
---|
165 | CALL ice_dyn_adv ( kt ) ; CALL Hbig( zhi_max, zhs_max ) ! -- advection of ice + correction on ice thickness |
---|
166 | |
---|
167 | ! diagnostics: divergence at T points |
---|
168 | DO jj = 2, jpjm1 |
---|
169 | DO ji = 2, jpim1 |
---|
170 | zdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj) & |
---|
171 | & + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1) ) * r1_e1e2t(ji,jj) |
---|
172 | END DO |
---|
173 | END DO |
---|
174 | CALL lbc_lnk( zdivu_i, 'T', 1. ) |
---|
175 | IF( iom_use('icediv') ) CALL iom_put( "icediv" , zdivu_i(:,:) ) |
---|
176 | |
---|
177 | DEALLOCATE( zdivu_i ) |
---|
178 | |
---|
179 | END SELECT |
---|
180 | ! |
---|
181 | ! controls |
---|
182 | IF( ln_icediachk ) CALL ice_cons_hsm(1, 'icedyn', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation |
---|
183 | IF( ln_timing ) CALL timing_stop ('icedyn') ! timing |
---|
184 | ! |
---|
185 | END SUBROUTINE ice_dyn |
---|
186 | |
---|
187 | |
---|
188 | SUBROUTINE Hbig( phi_max, phs_max ) |
---|
189 | !!------------------------------------------------------------------- |
---|
190 | !! *** ROUTINE Hbig *** |
---|
191 | !! |
---|
192 | !! ** Purpose : Thickness correction in case advection scheme creates |
---|
193 | !! abnormally tick ice or snow |
---|
194 | !! |
---|
195 | !! ** Method : 1- check whether ice thickness is larger than the surrounding 9-points |
---|
196 | !! (before advection) and reduce it by adapting ice concentration |
---|
197 | !! 2- check whether snow thickness is larger than the surrounding 9-points |
---|
198 | !! (before advection) and reduce it by sending the excess in the ocean |
---|
199 | !! 3- check whether snow load deplets the snow-ice interface below sea level$ |
---|
200 | !! and reduce it by sending the excess in the ocean |
---|
201 | !! 4- correct pond fraction to avoid a_ip > a_i |
---|
202 | !! |
---|
203 | !! ** input : Max thickness of the surrounding 9-points |
---|
204 | !!------------------------------------------------------------------- |
---|
205 | REAL(wp), DIMENSION(:,:,:), INTENT(in) :: phi_max, phs_max ! max ice thick from surrounding 9-pts |
---|
206 | ! |
---|
207 | INTEGER :: ji, jj, jl ! dummy loop indices |
---|
208 | REAL(wp) :: zhi, zhs, zvs_excess, zfra |
---|
209 | !!------------------------------------------------------------------- |
---|
210 | ! |
---|
211 | CALL ice_var_zapsmall !-- zap small areas |
---|
212 | ! |
---|
213 | DO jl = 1, jpl |
---|
214 | DO jj = 1, jpj |
---|
215 | DO ji = 1, jpi |
---|
216 | IF ( v_i(ji,jj,jl) > 0._wp ) THEN |
---|
217 | ! |
---|
218 | ! ! -- check h_i -- ! |
---|
219 | ! if h_i is larger than the surrounding 9 pts => reduce h_i and increase a_i |
---|
220 | zhi = v_i (ji,jj,jl) / a_i(ji,jj,jl) |
---|
221 | !!clem zdv = v_i(ji,jj,jl) - v_i_b(ji,jj,jl) |
---|
222 | !!clem IF ( ( zdv > 0.0 .AND. zh > phmax(ji,jj,jl) .AND. at_i_b(ji,jj) < 0.80 ) .OR. & |
---|
223 | !!clem & ( zdv <= 0.0 .AND. zh > phmax(ji,jj,jl) ) ) THEN |
---|
224 | IF( zhi > phi_max(ji,jj,jl) .AND. a_i(ji,jj,jl) < 0.15 ) THEN |
---|
225 | a_i(ji,jj,jl) = v_i(ji,jj,jl) / MIN( phi_max(ji,jj,jl), hi_max(jpl) ) !-- bound h_i to hi_max (99 m) |
---|
226 | ENDIF |
---|
227 | ! |
---|
228 | ! ! -- check h_s -- ! |
---|
229 | ! if h_s is larger than the surrounding 9 pts => put the snow excess in the ocean |
---|
230 | zhs = v_s (ji,jj,jl) / a_i(ji,jj,jl) |
---|
231 | IF( v_s(ji,jj,jl) > 0._wp .AND. zhs > phs_max(ji,jj,jl) .AND. a_i(ji,jj,jl) < 0.15 ) THEN |
---|
232 | zfra = a_i(ji,jj,jl) * phs_max(ji,jj,jl) / MAX( v_s(ji,jj,jl), epsi20 ) |
---|
233 | ! |
---|
234 | wfx_res(ji,jj) = wfx_res(ji,jj) + ( v_s(ji,jj,jl) - a_i(ji,jj,jl) * phs_max(ji,jj,jl) ) * rhos * r1_rdtice |
---|
235 | hfx_res(ji,jj) = hfx_res(ji,jj) - e_s(ji,jj,1,jl) * ( 1._wp - zfra ) * r1_rdtice ! W.m-2 <0 |
---|
236 | ! |
---|
237 | e_s(ji,jj,1,jl) = e_s(ji,jj,1,jl) * zfra |
---|
238 | v_s(ji,jj,jl) = a_i(ji,jj,jl) * phs_max(ji,jj,jl) |
---|
239 | ENDIF |
---|
240 | ! |
---|
241 | ! ! -- check snow load -- ! |
---|
242 | ! if snow load makes snow-ice interface to deplet below the ocean surface => put the snow excess in the ocean |
---|
243 | ! this correction is crucial because of the call to routine icecor afterwards which imposes a mini of ice thick. (rn_himin) |
---|
244 | ! this imposed mini can artificially make the snow thickness very high (if concentration decreases drastically) |
---|
245 | zvs_excess = MAX( 0._wp, v_s(ji,jj,jl) - v_i(ji,jj,jl) * (rau0-rhoi) * r1_rhos ) |
---|
246 | IF( zvs_excess > 0._wp ) THEN |
---|
247 | zfra = zvs_excess / MAX( v_s(ji,jj,jl), epsi20 ) |
---|
248 | wfx_res(ji,jj) = wfx_res(ji,jj) + zvs_excess * rhos * r1_rdtice |
---|
249 | hfx_res(ji,jj) = hfx_res(ji,jj) - e_s(ji,jj,1,jl) * ( 1._wp - zfra ) * r1_rdtice ! W.m-2 <0 |
---|
250 | ! |
---|
251 | e_s(ji,jj,1,jl) = e_s(ji,jj,1,jl) * zfra |
---|
252 | v_s(ji,jj,jl) = v_s(ji,jj,jl) - zvs_excess |
---|
253 | ENDIF |
---|
254 | |
---|
255 | ENDIF |
---|
256 | END DO |
---|
257 | END DO |
---|
258 | END DO |
---|
259 | ! !-- correct pond fraction to avoid a_ip > a_i |
---|
260 | WHERE( a_ip(:,:,:) > a_i(:,:,:) ) a_ip(:,:,:) = a_i(:,:,:) |
---|
261 | ! |
---|
262 | END SUBROUTINE Hbig |
---|
263 | |
---|
264 | |
---|
265 | SUBROUTINE Hpiling |
---|
266 | !!------------------------------------------------------------------- |
---|
267 | !! *** ROUTINE Hpiling *** |
---|
268 | !! |
---|
269 | !! ** Purpose : Simple conservative piling comparable with 1-cat models |
---|
270 | !! |
---|
271 | !! ** Method : pile-up ice when no ridging/rafting |
---|
272 | !! |
---|
273 | !! ** input : a_i |
---|
274 | !!------------------------------------------------------------------- |
---|
275 | INTEGER :: jl ! dummy loop indices |
---|
276 | !!------------------------------------------------------------------- |
---|
277 | ! |
---|
278 | CALL ice_var_zapsmall !-- zap small areas |
---|
279 | ! |
---|
280 | at_i(:,:) = SUM( a_i(:,:,:), dim=3 ) |
---|
281 | DO jl = 1, jpl |
---|
282 | WHERE( at_i(:,:) > epsi20 ) |
---|
283 | a_i(:,:,jl) = a_i(:,:,jl) * ( 1._wp + MIN( rn_amax_2d(:,:) - at_i(:,:) , 0._wp ) / at_i(:,:) ) |
---|
284 | END WHERE |
---|
285 | END DO |
---|
286 | ! |
---|
287 | END SUBROUTINE Hpiling |
---|
288 | |
---|
289 | |
---|
290 | SUBROUTINE ice_dyn_init |
---|
291 | !!------------------------------------------------------------------- |
---|
292 | !! *** ROUTINE ice_dyn_init *** |
---|
293 | !! |
---|
294 | !! ** Purpose : Physical constants and parameters linked to the ice |
---|
295 | !! dynamics |
---|
296 | !! |
---|
297 | !! ** Method : Read the namdyn namelist and check the ice-dynamic |
---|
298 | !! parameter values called at the first timestep (nit000) |
---|
299 | !! |
---|
300 | !! ** input : Namelist namdyn |
---|
301 | !!------------------------------------------------------------------- |
---|
302 | INTEGER :: ios, ioptio ! Local integer output status for namelist read |
---|
303 | !! |
---|
304 | NAMELIST/namdyn/ ln_dynALL, ln_dynRHGADV, ln_dynADV1D, ln_dynADV2D, rn_uice, rn_vice, & |
---|
305 | & rn_ishlat , & |
---|
306 | & ln_landfast_L16, ln_landfast_home, rn_depfra, rn_icebfr, rn_lfrelax, rn_tensile |
---|
307 | !!------------------------------------------------------------------- |
---|
308 | ! |
---|
309 | REWIND( numnam_ice_ref ) ! Namelist namdyn in reference namelist : Ice dynamics |
---|
310 | READ ( numnam_ice_ref, namdyn, IOSTAT = ios, ERR = 901) |
---|
311 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdyn in reference namelist', lwp ) |
---|
312 | REWIND( numnam_ice_cfg ) ! Namelist namdyn in configuration namelist : Ice dynamics |
---|
313 | READ ( numnam_ice_cfg, namdyn, IOSTAT = ios, ERR = 902 ) |
---|
314 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namdyn in configuration namelist', lwp ) |
---|
315 | IF(lwm) WRITE( numoni, namdyn ) |
---|
316 | ! |
---|
317 | IF(lwp) THEN ! control print |
---|
318 | WRITE(numout,*) |
---|
319 | WRITE(numout,*) 'ice_dyn_init: ice parameters for ice dynamics ' |
---|
320 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
321 | WRITE(numout,*) ' Namelist namdyn:' |
---|
322 | WRITE(numout,*) ' Full ice dynamics (rhg + adv + ridge/raft + corr) ln_dynALL = ', ln_dynALL |
---|
323 | WRITE(numout,*) ' No ridge/raft & No cor (rhg + adv) ln_dynRHGADV = ', ln_dynRHGADV |
---|
324 | WRITE(numout,*) ' Advection 1D only (Schar & Smolarkiewicz 1996) ln_dynADV1D = ', ln_dynADV1D |
---|
325 | WRITE(numout,*) ' Advection 2D only (rn_uvice + adv) ln_dynADV2D = ', ln_dynADV2D |
---|
326 | WRITE(numout,*) ' with prescribed velocity given by (u,v)_ice = (rn_uice,rn_vice) = (', rn_uice,',', rn_vice,')' |
---|
327 | WRITE(numout,*) ' lateral boundary condition for sea ice dynamics rn_ishlat = ', rn_ishlat |
---|
328 | WRITE(numout,*) ' Landfast: param from Lemieux 2016 ln_landfast_L16 = ', ln_landfast_L16 |
---|
329 | WRITE(numout,*) ' Landfast: param from home made ln_landfast_home= ', ln_landfast_home |
---|
330 | WRITE(numout,*) ' fraction of ocean depth that ice must reach rn_depfra = ', rn_depfra |
---|
331 | WRITE(numout,*) ' maximum bottom stress per unit area of contact rn_icebfr = ', rn_icebfr |
---|
332 | WRITE(numout,*) ' relax time scale (s-1) to reach static friction rn_lfrelax = ', rn_lfrelax |
---|
333 | WRITE(numout,*) ' isotropic tensile strength rn_tensile = ', rn_tensile |
---|
334 | WRITE(numout,*) |
---|
335 | ENDIF |
---|
336 | ! !== set the choice of ice dynamics ==! |
---|
337 | ioptio = 0 |
---|
338 | ! !--- full dynamics (rheology + advection + ridging/rafting + correction) |
---|
339 | IF( ln_dynALL ) THEN ; ioptio = ioptio + 1 ; nice_dyn = np_dynALL ; ENDIF |
---|
340 | ! !--- dynamics without ridging/rafting and corr (rheology + advection) |
---|
341 | IF( ln_dynRHGADV ) THEN ; ioptio = ioptio + 1 ; nice_dyn = np_dynRHGADV ; ENDIF |
---|
342 | ! !--- advection 1D only - test case from Schar & Smolarkiewicz 1996 |
---|
343 | IF( ln_dynADV1D ) THEN ; ioptio = ioptio + 1 ; nice_dyn = np_dynADV1D ; ENDIF |
---|
344 | ! !--- advection 2D only with prescribed ice velocities (from namelist) |
---|
345 | IF( ln_dynADV2D ) THEN ; ioptio = ioptio + 1 ; nice_dyn = np_dynADV2D ; ENDIF |
---|
346 | ! |
---|
347 | IF( ioptio /= 1 ) CALL ctl_stop( 'ice_dyn_init: one and only one ice dynamics option has to be defined ' ) |
---|
348 | ! |
---|
349 | ! !--- Lateral boundary conditions |
---|
350 | IF ( rn_ishlat == 0. ) THEN ; IF(lwp) WRITE(numout,*) ' ===>>> ice lateral free-slip' |
---|
351 | ELSEIF ( rn_ishlat == 2. ) THEN ; IF(lwp) WRITE(numout,*) ' ===>>> ice lateral no-slip' |
---|
352 | ELSEIF ( 0. < rn_ishlat .AND. rn_ishlat < 2. ) THEN ; IF(lwp) WRITE(numout,*) ' ===>>> ice lateral partial-slip' |
---|
353 | ELSEIF ( 2. < rn_ishlat ) THEN ; IF(lwp) WRITE(numout,*) ' ===>>> ice lateral strong-slip' |
---|
354 | ENDIF |
---|
355 | ! !--- Landfast ice |
---|
356 | IF( .NOT.ln_landfast_L16 .AND. .NOT.ln_landfast_home ) tau_icebfr(:,:) = 0._wp |
---|
357 | ! |
---|
358 | IF ( ln_landfast_L16 .AND. ln_landfast_home ) THEN |
---|
359 | CALL ctl_stop( 'ice_dyn_init: choose one and only one landfast parameterization (ln_landfast_L16 or ln_landfast_home)' ) |
---|
360 | ENDIF |
---|
361 | ! |
---|
362 | CALL ice_dyn_rdgrft_init ! set ice ridging/rafting parameters |
---|
363 | CALL ice_dyn_rhg_init ! set ice rheology parameters |
---|
364 | CALL ice_dyn_adv_init ! set ice advection parameters |
---|
365 | ! |
---|
366 | END SUBROUTINE ice_dyn_init |
---|
367 | |
---|
368 | #else |
---|
369 | !!---------------------------------------------------------------------- |
---|
370 | !! Default option Empty module NO SI3 sea-ice model |
---|
371 | !!---------------------------------------------------------------------- |
---|
372 | #endif |
---|
373 | |
---|
374 | !!====================================================================== |
---|
375 | END MODULE icedyn |
---|