/[MITgcm]/MITgcm/pkg/seaice/seaice_init_varia.F
ViewVC logotype

Contents of /MITgcm/pkg/seaice/seaice_init_varia.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.94 - (show annotations) (download)
Tue Apr 4 23:31:27 2017 UTC (7 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.93: +2 -2 lines
- add specific run-time param to select level of printed plot-field-maps,
  set by default to debugLevel.

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_varia.F,v 1.93 2016/11/30 00:17:16 jmc Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5 #ifdef ALLOW_OBCS
6 # include "OBCS_OPTIONS.h"
7 #endif
8
9 CStartOfInterface
10 SUBROUTINE SEAICE_INIT_VARIA( myThid )
11 C *==========================================================*
12 C | SUBROUTINE SEAICE_INIT_VARIA |
13 C | o Initialization of sea ice model. |
14 C *==========================================================*
15 C *==========================================================*
16 IMPLICIT NONE
17
18 C === Global variables ===
19 #include "SIZE.h"
20 #include "EEPARAMS.h"
21 #include "PARAMS.h"
22 #include "GRID.h"
23 #include "DYNVARS.h"
24 #include "FFIELDS.h"
25 #include "SEAICE_SIZE.h"
26 #include "SEAICE_PARAMS.h"
27 #include "SEAICE.h"
28 #include "SEAICE_TRACER.h"
29 #include "SEAICE_TAVE.h"
30 #ifdef OBCS_UVICE_OLD
31 # include "OBCS_GRID.h"
32 #endif
33
34 C === Routine arguments ===
35 C myThid :: Thread no. that called this routine.
36 INTEGER myThid
37 CEndOfInterface
38
39 C === Local variables ===
40 C i,j,k,bi,bj :: Loop counters
41
42 INTEGER i, j, bi, bj
43 INTEGER kSurface
44 _RS mask_uice
45 INTEGER k
46 #ifdef ALLOW_SITRACER
47 INTEGER iTr, jTh
48 #endif
49 #ifdef OBCS_UVICE_OLD
50 INTEGER I_obc, J_obc
51 #endif /* ALLOW_OBCS */
52 #ifdef SEAICE_CGRID
53 _RL recip_tensilDepth
54 #endif
55
56 IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
57 kSurface = Nr
58 ELSE
59 kSurface = 1
60 ENDIF
61
62 C-- Initialize grid info
63 DO bj=myByLo(myThid),myByHi(myThid)
64 DO bi=myBxLo(myThid),myBxHi(myThid)
65 DO j=1-OLy,sNy+OLy
66 DO i=1-OLx,sNx+OLx
67 HEFFM(i,j,bi,bj) = 0. _d 0
68 ENDDO
69 ENDDO
70 DO j=1-OLy,sNy+OLy
71 DO i=1-OLx,sNx+OLx
72 HEFFM(i,j,bi,bj)= 1. _d 0
73 IF (_hFacC(i,j,kSurface,bi,bj).eq.0.)
74 & HEFFM(i,j,bi,bj)= 0. _d 0
75 ENDDO
76 ENDDO
77 #ifndef SEAICE_CGRID
78 DO j=1-OLy+1,sNy+OLy
79 DO i=1-OLx+1,sNx+OLx
80 UVM(i,j,bi,bj)=0. _d 0
81 mask_uice=HEFFM(i,j, bi,bj)+HEFFM(i-1,j-1,bi,bj)
82 & +HEFFM(i,j-1,bi,bj)+HEFFM(i-1,j, bi,bj)
83 IF(mask_uice.GT.3.5 _d 0) UVM(i,j,bi,bj)=1. _d 0
84 ENDDO
85 ENDDO
86 #endif /* SEAICE_CGRID */
87 ENDDO
88 ENDDO
89
90 C coefficients for metric terms
91 DO bj=myByLo(myThid),myByHi(myThid)
92 DO bi=myBxLo(myThid),myBxHi(myThid)
93 #ifdef SEAICE_CGRID
94 DO j=1-OLy,sNy+OLy
95 DO i=1-OLx,sNx+OLx
96 k1AtC(I,J,bi,bj) = 0.0 _d 0
97 k1AtZ(I,J,bi,bj) = 0.0 _d 0
98 k2AtC(I,J,bi,bj) = 0.0 _d 0
99 k2AtZ(I,J,bi,bj) = 0.0 _d 0
100 ENDDO
101 ENDDO
102 IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN
103 C This is the only case where tan(phi) is not zero. In this case
104 C C and U points, and Z and V points have the same phi, so that we
105 C only need a copy here. Do not use tan(YC) and tan(YG), because
106 C these
107 C can be the geographical coordinates and not the correct grid
108 C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0)
109 DO j=1-OLy,sNy+OLy
110 DO i=1-OLx,sNx+OLx
111 k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere
112 k2AtZ(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere
113 ENDDO
114 ENDDO
115 ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN
116 C compute metric term coefficients from finite difference
117 C approximation
118 DO j=1-OLy,sNy+OLy
119 DO i=1-OLx,sNx+OLx-1
120 k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj)
121 & * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) )
122 & * _recip_dxF(I,J,bi,bj)
123 ENDDO
124 ENDDO
125 DO j=1-OLy,sNy+OLy
126 DO i=1-OLx+1,sNx+OLx
127 k1AtZ(I,J,bi,bj) = _recip_dyU(I,J,bi,bj)
128 & * ( _dyC(I,J,bi,bj) - _dyC(I-1,J,bi,bj) )
129 & * _recip_dxV(I,J,bi,bj)
130 ENDDO
131 ENDDO
132 DO j=1-OLy,sNy+OLy-1
133 DO i=1-OLx,sNx+OLx
134 k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj)
135 & * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) )
136 & * _recip_dyF(I,J,bi,bj)
137 ENDDO
138 ENDDO
139 DO j=1-OLy+1,sNy+OLy
140 DO i=1-OLx,sNx+OLx
141 k2AtZ(I,J,bi,bj) = _recip_dxV(I,J,bi,bj)
142 & * ( _dxC(I,J,bi,bj) - _dxC(I,J-1,bi,bj) )
143 & * _recip_dyU(I,J,bi,bj)
144 ENDDO
145 ENDDO
146 ENDIF
147 #else /* not SEAICE_CGRID */
148 DO j=1-OLy,sNy+OLy
149 DO i=1-OLx,sNx+OLx
150 k1AtC(I,J,bi,bj) = 0.0 _d 0
151 k1AtU(I,J,bi,bj) = 0.0 _d 0
152 k1AtV(I,J,bi,bj) = 0.0 _d 0
153 k2AtC(I,J,bi,bj) = 0.0 _d 0
154 k2AtU(I,J,bi,bj) = 0.0 _d 0
155 k2AtV(I,J,bi,bj) = 0.0 _d 0
156 ENDDO
157 ENDDO
158 IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN
159 C This is the only case where tan(phi) is not zero. In this case
160 C C and U points, and Z and V points have the same phi, so that we
161 C only need a copy here. Do not use tan(YC) and tan(YG), because
162 C these
163 C can be the geographical coordinates and not the correct grid
164 C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0)
165 DO j=1-OLy,sNy+OLy
166 DO i=1-OLx,sNx+OLx
167 k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere
168 k2AtU(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere
169 k2AtV(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere
170 ENDDO
171 ENDDO
172 ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN
173 C compute metric term coefficients from finite difference
174 C approximation
175 DO j=1-OLy,sNy+OLy
176 DO i=1-OLx,sNx+OLx-1
177 k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj)
178 & * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) )
179 & * _recip_dxF(I,J,bi,bj)
180 ENDDO
181 ENDDO
182 DO j=1-OLy,sNy+OLy
183 DO i=1-OLx+1,sNx+OLx
184 k1AtU(I,J,bi,bj) = _recip_dyG(I,J,bi,bj)
185 & * ( _dyF(I,J,bi,bj) - _dyF(I-1,J,bi,bj) )
186 & * _recip_dxC(I,J,bi,bj)
187 ENDDO
188 ENDDO
189 DO j=1-OLy,sNy+OLy
190 DO i=1-OLx,sNx+OLx-1
191 k1AtV(I,J,bi,bj) = _recip_dyC(I,J,bi,bj)
192 & * ( _dyU(I+1,J,bi,bj) - _dyU(I,J,bi,bj) )
193 & * _recip_dxG(I,J,bi,bj)
194 ENDDO
195 ENDDO
196 DO j=1-OLy,sNy+OLy-1
197 DO i=1-OLx,sNx+OLx
198 k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj)
199 & * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) )
200 & * _recip_dyF(I,J,bi,bj)
201 ENDDO
202 ENDDO
203 DO j=1-OLy,sNy+OLy-1
204 DO i=1-OLx,sNx+OLx
205 k2AtU(I,J,bi,bj) = _recip_dxC(I,J,bi,bj)
206 & * ( _dxV(I,J+1,bi,bj) - _dxV(I,J,bi,bj) )
207 & * _recip_dyG(I,J,bi,bj)
208 ENDDO
209 ENDDO
210 DO j=1-OLy+1,sNy+OLy
211 DO i=1-OLx,sNx+OLx
212 k2AtV(I,J,bi,bj) = _recip_dxG(I,J,bi,bj)
213 & * ( _dxF(I,J,bi,bj) - _dxF(I,J-1,bi,bj) )
214 & * _recip_dyC(I,J,bi,bj)
215 ENDDO
216 ENDDO
217 ENDIF
218 #endif /* not SEAICE_CGRID */
219 ENDDO
220 ENDDO
221
222 #ifndef SEAICE_CGRID
223 C-- Choose a proxy level for geostrophic velocity,
224 DO bj=myByLo(myThid),myByHi(myThid)
225 DO bi=myBxLo(myThid),myBxHi(myThid)
226 DO j=1-OLy,sNy+OLy
227 DO i=1-OLx,sNx+OLx
228 KGEO(i,j,bi,bj) = 0
229 ENDDO
230 ENDDO
231 DO j=1-OLy,sNy+OLy
232 DO i=1-OLx,sNx+OLx
233 #ifdef SEAICE_BICE_STRESS
234 KGEO(i,j,bi,bj) = 1
235 #else /* SEAICE_BICE_STRESS */
236 IF (klowc(i,j,bi,bj) .LT. 2) THEN
237 KGEO(i,j,bi,bj) = 1
238 ELSE
239 KGEO(i,j,bi,bj) = 2
240 DO WHILE ( abs(rC(KGEO(i,j,bi,bj))) .LT. 50.0 _d 0 .AND.
241 & KGEO(i,j,bi,bj) .LT. (klowc(i,j,bi,bj)-1) )
242 KGEO(i,j,bi,bj) = KGEO(i,j,bi,bj) + 1
243 ENDDO
244 ENDIF
245 #endif /* SEAICE_BICE_STRESS */
246 ENDDO
247 ENDDO
248 ENDDO
249 ENDDO
250 #endif /* SEAICE_CGRID */
251
252 C-- Initialise all variables in common blocks:
253 DO bj=myByLo(myThid),myByHi(myThid)
254 DO bi=myBxLo(myThid),myBxHi(myThid)
255 DO j=1-OLy,sNy+OLy
256 DO i=1-OLx,sNx+OLx
257 HEFF(i,j,bi,bj)=0. _d 0
258 AREA(i,j,bi,bj)=0. _d 0
259 CToM<<<
260 #ifdef SEAICE_ITD
261 DO k=1,nITD
262 AREAITD(i,j,k,bi,bj) =0. _d 0
263 HEFFITD(i,j,k,bi,bj) =0. _d 0
264 ENDDO
265 #endif
266 C>>>ToM
267 UICE(i,j,bi,bj)=0. _d 0
268 VICE(i,j,bi,bj)=0. _d 0
269 #ifdef SEAICE_ALLOW_FREEDRIFT
270 uice_fd(i,j,bi,bj)=0. _d 0
271 vice_fd(i,j,bi,bj)=0. _d 0
272 #endif
273 C
274 uIceNm1(i,j,bi,bj)=0. _d 0
275 vIceNm1(i,j,bi,bj)=0. _d 0
276 ETA (i,j,bi,bj) = 0. _d 0
277 etaZ(i,j,bi,bj) = 0. _d 0
278 ZETA(i,j,bi,bj) = 0. _d 0
279 FORCEX(i,j,bi,bj) = 0. _d 0
280 FORCEY(i,j,bi,bj) = 0. _d 0
281 #ifdef SEAICE_CGRID
282 seaiceMassC(i,j,bi,bj)=0. _d 0
283 seaiceMassU(i,j,bi,bj)=0. _d 0
284 seaiceMassV(i,j,bi,bj)=0. _d 0
285 stressDivergenceX(i,j,bi,bj) = 0. _d 0
286 stressDivergenceY(i,j,bi,bj) = 0. _d 0
287 # ifdef SEAICE_ALLOW_EVP
288 seaice_sigma1 (i,j,bi,bj) = 0. _d 0
289 seaice_sigma2 (i,j,bi,bj) = 0. _d 0
290 seaice_sigma12(i,j,bi,bj) = 0. _d 0
291 # endif /* SEAICE_ALLOW_EVP */
292 #else /* SEAICE_CGRID */
293 uIceC(i,j,bi,bj) = 0. _d 0
294 vIceC(i,j,bi,bj) = 0. _d 0
295 AMASS(i,j,bi,bj) = 0. _d 0
296 DAIRN(i,j,bi,bj) = 0. _d 0
297 WINDX(i,j,bi,bj) = 0. _d 0
298 WINDY(i,j,bi,bj) = 0. _d 0
299 GWATX(i,j,bi,bj) = 0. _d 0
300 GWATY(i,j,bi,bj) = 0. _d 0
301 #endif /* SEAICE_CGRID */
302 DWATN(i,j,bi,bj) = 0. _d 0
303 #ifdef SEAICE_ALLOW_BOTTOMDRAG
304 CbotC(i,j,bi,bj) = 0. _d 0
305 #endif /* SEAICE_ALLOW_BOTTOMDRAG */
306 PRESS0(i,j,bi,bj) = 0. _d 0
307 FORCEX0(i,j,bi,bj)= 0. _d 0
308 FORCEY0(i,j,bi,bj)= 0. _d 0
309 ZMAX(i,j,bi,bj) = 0. _d 0
310 ZMIN(i,j,bi,bj) = 0. _d 0
311 HSNOW(i,j,bi,bj) = 0. _d 0
312 tensileStrFac(i,j,bi,bj) = 0. _d 0
313 CToM<<<
314 #ifdef SEAICE_ITD
315 DO k=1,nITD
316 HSNOWITD(i,j,k,bi,bj)=0. _d 0
317 ENDDO
318 #endif
319 C>>>ToM
320 #ifdef SEAICE_VARIABLE_SALINITY
321 HSALT(i,j,bi,bj) = 0. _d 0
322 #endif
323 #ifdef ALLOW_SITRACER
324 DO iTr = 1, SItrMaxNum
325 SItracer(i,j,bi,bj,iTr) = 0. _d 0
326 SItrBucket(i,j,bi,bj,iTr) = 0. _d 0
327 c "ice concentration" tracer that should remain .EQ.1.
328 if (SItrName(iTr).EQ.'one') SItracer(i,j,bi,bj,iTr)=1. _d 0
329 ENDDO
330 DO jTh = 1, 5
331 SItrHEFF (i,j,bi,bj,jTh) = 0. _d 0
332 ENDDO
333 DO jTh = 1, 3
334 SItrAREA (i,j,bi,bj,jTh) = 0. _d 0
335 ENDDO
336 #endif
337 DO k=1,nITD
338 TICES(i,j,k,bi,bj)=0. _d 0
339 ENDDO
340 TAUX(i,j,bi,bj) = 0. _d 0
341 TAUY(i,j,bi,bj) = 0. _d 0
342 #ifdef ALLOW_SEAICE_COST_EXPORT
343 uHeffExportCell(i,j,bi,bj) = 0. _d 0
344 vHeffExportCell(i,j,bi,bj) = 0. _d 0
345 icevolMeanCell(i,j,bi,bj) = 0. _d 0
346 #endif
347 saltWtrIce(i,j,bi,bj) = 0. _d 0
348 frWtrIce(i,j,bi,bj) = 0. _d 0
349 #if (defined (ALLOW_MEAN_SFLUX_COST_CONTRIBUTION) || defined (ALLOW_SSH_GLOBMEAN_COST_CONTRIBUTION))
350 frWtrAtm(i,j,bi,bj) = 0. _d 0
351 AREAforAtmFW(i,j,bi,bj)=0. _d 0
352 #endif
353 ENDDO
354 ENDDO
355 ENDDO
356 ENDDO
357
358 #ifdef ALLOW_TIMEAVE
359 C Initialize averages to zero
360 DO bj = myByLo(myThid), myByHi(myThid)
361 DO bi = myBxLo(myThid), myBxHi(myThid)
362 CALL TIMEAVE_RESET( FUtave , 1, bi, bj, myThid )
363 CALL TIMEAVE_RESET( FVtave , 1, bi, bj, myThid )
364 CALL TIMEAVE_RESET( EmPmRtave, 1, bi, bj, myThid )
365 CALL TIMEAVE_RESET( QNETtave , 1, bi, bj, myThid )
366 CALL TIMEAVE_RESET( QSWtave , 1, bi, bj, myThid )
367 CALL TIMEAVE_RESET( UICEtave , 1, bi, bj, myThid )
368 CALL TIMEAVE_RESET( VICEtave , 1, bi, bj, myThid )
369 CALL TIMEAVE_RESET( HEFFtave , 1, bi, bj, myThid )
370 CALL TIMEAVE_RESET( AREAtave , 1, bi, bj, myThid )
371 SEAICE_timeAve(bi,bj) = ZERO
372 ENDDO
373 ENDDO
374 #endif /* ALLOW_TIMEAVE */
375
376 C-- Initialize (variable) grid info. As long as we allow masking of
377 C-- velocities outside of ice covered areas (in seaice_dynsolver)
378 C-- we need to re-initialize seaiceMaskU/V here for TAF/TAMC
379 #ifdef SEAICE_CGRID
380 DO bj=myByLo(myThid),myByHi(myThid)
381 DO bi=myBxLo(myThid),myBxHi(myThid)
382 DO j=1-OLy+1,sNy+OLy
383 DO i=1-OLx+1,sNx+OLx
384 seaiceMaskU(i,j,bi,bj)= 0.0 _d 0
385 seaiceMaskV(i,j,bi,bj)= 0.0 _d 0
386 mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i-1,j ,bi,bj)
387 IF(mask_uice.GT.1.5 _d 0) seaiceMaskU(i,j,bi,bj)=1.0 _d 0
388 mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i ,j-1,bi,bj)
389 IF(mask_uice.GT.1.5 _d 0) seaiceMaskV(i,j,bi,bj)=1.0 _d 0
390 ENDDO
391 ENDDO
392 ENDDO
393 ENDDO
394 #endif /* SEAICE_CGRID */
395
396 DO bj=myByLo(myThid),myByHi(myThid)
397 DO bi=myBxLo(myThid),myBxHi(myThid)
398 #ifdef OBCS_UVICE_OLD
399 #ifdef SEAICE_CGRID
400 IF (useOBCS) THEN
401 C-- If OBCS is turned on, close southern and western boundaries
402 DO i=1-OLx,sNx+OLx
403 C Southern boundary
404 J_obc = OB_Js(i,bi,bj)
405 IF ( J_obc.NE.OB_indexNone ) THEN
406 seaiceMaskU(i,J_obc,bi,bj)= 0.0 _d 0
407 seaiceMaskV(i,J_obc,bi,bj)= 0.0 _d 0
408 ENDIF
409 ENDDO
410 DO j=1-OLy,sNy+OLy
411 C Western boundary
412 I_obc = OB_Iw(j,bi,bj)
413 IF ( I_obc.NE.OB_indexNone ) THEN
414 seaiceMaskU(I_obc,j,bi,bj)= 0.0 _d 0
415 seaiceMaskV(I_obc,j,bi,bj)= 0.0 _d 0
416 ENDIF
417 ENDDO
418 ENDIF
419 #endif /* SEAICE_CGRID */
420 #endif /* OBCS_UVICE_OLD */
421
422 DO j=1-OLy,sNy+OLy
423 DO i=1-OLx,sNx+OLx
424 DO k=1,nITD
425 TICES(i,j,k,bi,bj)=273.0 _d 0
426 ENDDO
427 #ifndef SEAICE_CGRID
428 AMASS (i,j,bi,bj)=1000.0 _d 0
429 #else
430 seaiceMassC(i,j,bi,bj)=1000.0 _d 0
431 seaiceMassU(i,j,bi,bj)=1000.0 _d 0
432 seaiceMassV(i,j,bi,bj)=1000.0 _d 0
433 #endif
434 ENDDO
435 ENDDO
436
437 ENDDO
438 ENDDO
439
440 C-- Update overlap regions
441 #ifdef SEAICE_CGRID
442 CALL EXCH_UV_XY_RL(seaiceMaskU,seaiceMaskV,.FALSE.,myThid)
443 #else
444 _EXCH_XY_RS(UVM, myThid)
445 #endif
446
447 C-- Now lets look at all these beasts
448 IF ( plotLevel.GE.debLevC ) THEN
449 CALL PLOT_FIELD_XYRL( HEFFM , 'Current HEFFM ' ,
450 & nIter0, myThid )
451 #ifdef SEAICE_CGRID
452 CALL PLOT_FIELD_XYRL( seaiceMaskU, 'Current seaiceMaskU',
453 & nIter0, myThid )
454 CALL PLOT_FIELD_XYRL( seaiceMaskV, 'Current seaiceMaskV',
455 & nIter0, myThid )
456 #else
457 CALL PLOT_FIELD_XYRS( UVM , 'Current UVM ' ,
458 & nIter0, myThid )
459 #endif
460 ENDIF
461
462 C-- Set model variables to initial/restart conditions
463 IF ( .NOT. ( startTime .EQ. baseTime .AND. nIter0 .EQ. 0
464 & .AND. pickupSuff .EQ. ' ') ) THEN
465
466 CALL SEAICE_READ_PICKUP ( myThid )
467
468 ELSE
469
470 DO bj=myByLo(myThid),myByHi(myThid)
471 DO bi=myBxLo(myThid),myBxHi(myThid)
472 DO j=1-OLy,sNy+OLy
473 DO i=1-OLx,sNx+OLx
474 HEFF(i,j,bi,bj)=SEAICE_initialHEFF*HEFFM(i,j,bi,bj)
475 UICE(i,j,bi,bj)=ZERO
476 VICE(i,j,bi,bj)=ZERO
477 ENDDO
478 ENDDO
479 ENDDO
480 ENDDO
481
482 C-- Read initial sea-ice velocity from file (if available)
483 IF ( uIceFile .NE. ' ' )
484 & CALL READ_FLD_XY_RL( uIceFile, ' ', uIce, 0, myThid )
485 IF ( vIceFile .NE. ' ' )
486 & CALL READ_FLD_XY_RL( vIceFile, ' ', vIce, 0, myThid )
487 IF ( uIceFile .NE. ' ' .OR. vIceFile .NE. ' ' ) THEN
488 #ifdef SEAICE_CGRID
489 DO bj=myByLo(myThid),myByHi(myThid)
490 DO bi=myBxLo(myThid),myBxHi(myThid)
491 DO j=1-OLy,sNy+OLy
492 DO i=1-OLx,sNx+OLx
493 uIce(i,j,bi,bj) = uIce(i,j,bi,bj)*seaiceMaskU(i,j,bi,bj)
494 vIce(i,j,bi,bj) = vIce(i,j,bi,bj)*seaiceMaskV(i,j,bi,bj)
495 ENDDO
496 ENDDO
497 ENDDO
498 ENDDO
499 #endif /* SEAICE_CGRID */
500 CALL EXCH_UV_XY_RL( uIce, vIce, .TRUE., myThid )
501 ENDIF
502
503 C-- Read initial sea-ice thickness from file if available.
504 IF ( HeffFile .NE. ' ' ) THEN
505 CALL READ_FLD_XY_RL( HeffFile, ' ', HEFF, 0, myThid )
506 _EXCH_XY_RL(HEFF,myThid)
507 DO bj=myByLo(myThid),myByHi(myThid)
508 DO bi=myBxLo(myThid),myBxHi(myThid)
509 DO j=1-OLy,sNy+OLy
510 DO i=1-OLx,sNx+OLx
511 HEFF(i,j,bi,bj) = MAX(HEFF(i,j,bi,bj),ZERO)
512 ENDDO
513 ENDDO
514 ENDDO
515 ENDDO
516 ENDIF
517
518 DO bj=myByLo(myThid),myByHi(myThid)
519 DO bi=myBxLo(myThid),myBxHi(myThid)
520 DO j=1-OLy,sNy+OLy
521 DO i=1-OLx,sNx+OLx
522 IF(HEFF(i,j,bi,bj).GT.ZERO) AREA(i,j,bi,bj)=ONE
523 ENDDO
524 ENDDO
525 ENDDO
526 ENDDO
527
528 C-- Read initial sea-ice area from file if available.
529 IF ( AreaFile .NE. ' ' ) THEN
530 CALL READ_FLD_XY_RL( AreaFile, ' ', AREA, 0, myThid )
531 _EXCH_XY_RL(AREA,myThid)
532 DO bj=myByLo(myThid),myByHi(myThid)
533 DO bi=myBxLo(myThid),myBxHi(myThid)
534 DO j=1-OLy,sNy+OLy
535 DO i=1-OLx,sNx+OLx
536 AREA(i,j,bi,bj) = MAX(AREA(i,j,bi,bj),ZERO)
537 AREA(i,j,bi,bj) = MIN(AREA(i,j,bi,bj),ONE)
538 IF ( AREA(i,j,bi,bj) .LE. ZERO ) HEFF(i,j,bi,bj) = ZERO
539 IF ( HEFF(i,j,bi,bj) .LE. ZERO ) AREA(i,j,bi,bj) = ZERO
540 ENDDO
541 ENDDO
542 ENDDO
543 ENDDO
544 ENDIF
545
546 DO bj=myByLo(myThid),myByHi(myThid)
547 DO bi=myBxLo(myThid),myBxHi(myThid)
548 DO j=1-OLy,sNy+OLy
549 DO i=1-OLx,sNx+OLx
550 HSNOW(i,j,bi,bj) = 0.2 _d 0 * AREA(i,j,bi,bj)
551 ENDDO
552 ENDDO
553 ENDDO
554 ENDDO
555
556 C-- Read initial snow thickness from file if available.
557 IF ( HsnowFile .NE. ' ' ) THEN
558 CALL READ_FLD_XY_RL( HsnowFile, ' ', HSNOW, 0, myThid )
559 _EXCH_XY_RL(HSNOW,myThid)
560 DO bj=myByLo(myThid),myByHi(myThid)
561 DO bi=myBxLo(myThid),myBxHi(myThid)
562 DO j=1-OLy,sNy+OLy
563 DO i=1-OLx,sNx+OLx
564 HSNOW(i,j,bi,bj) = MAX(HSNOW(i,j,bi,bj),ZERO)
565 ENDDO
566 ENDDO
567 ENDDO
568 ENDDO
569 ENDIF
570
571 #ifdef SEAICE_ITD
572 DO bj=myByLo(myThid),myByHi(myThid)
573 DO bi=myBxLo(myThid),myBxHi(myThid)
574 DO j=1-OLy,sNy+OLy
575 DO i=1-OLx,sNx+OLx
576 AREAITD(I,J,1,bi,bj) = AREA(I,J,bi,bj)
577 HEFFITD(I,J,1,bi,bj) = HEFF(I,J,bi,bj)
578 HSNOWITD(I,J,1,bi,bj) = HSNOW(I,J,bi,bj)
579 opnWtrFrac(I,J,bi,bj) = 1. _d 0 - AREA(I,J,bi,bj)
580 fw2ObyRidge(I,J,bi,bj) = 0. _d 0
581 ENDDO
582 ENDDO
583 CALL SEAICE_ITD_REDIST(bi, bj, baseTime, nIter0, myThid)
584 ENDDO
585 ENDDO
586 #endif
587
588 #ifdef SEAICE_VARIABLE_SALINITY
589 DO bj=myByLo(myThid),myByHi(myThid)
590 DO bi=myBxLo(myThid),myBxHi(myThid)
591 DO j=1-OLy,sNy+OLy
592 DO i=1-OLx,sNx+OLx
593 HSALT(i,j,bi,bj)=HEFF(i,j,bi,bj)*salt(i,j,kSurface,bi,bj)*
594 & SEAICE_rhoIce*SEAICE_saltFrac
595 cif & ICE2WATR*rhoConstFresh*SEAICE_saltFrac
596
597 ENDDO
598 ENDDO
599 ENDDO
600 ENDDO
601
602 C-- Read initial sea ice salinity from file if available.
603 IF ( HsaltFile .NE. ' ' ) THEN
604 CALL READ_FLD_XY_RL( HsaltFile, ' ', HSALT, 0, myThid )
605 _EXCH_XY_RL(HSALT,myThid)
606 ENDIF
607 #endif /* SEAICE_VARIABLE_SALINITY */
608
609 #ifdef ALLOW_SITRACER
610 C-- Read initial sea ice age from file if available.
611 DO iTr = 1, SItrMaxNum
612 IF ( SItrFile(iTr) .NE. ' ' ) THEN
613 CALL READ_FLD_XY_RL( siTrFile(iTr), ' ',
614 & SItracer(1-OLx,1-OLy,1,1,iTr), 0, myThid )
615 _EXCH_XY_RL(SItracer(1-OLx,1-OLy,1,1,iTr),myThid)
616 ENDIF
617 ENDDO
618 #endif /* ALLOW_SITRACER */
619
620 ENDIF
621
622 #if (defined (ALLOW_MEAN_SFLUX_COST_CONTRIBUTION) || defined (ALLOW_SSH_GLOBMEAN_COST_CONTRIBUTION))
623 DO bj=myByLo(myThid),myByHi(myThid)
624 DO bi=myBxLo(myThid),myBxHi(myThid)
625 DO j=1-OLy,sNy+OLy
626 DO i=1-OLx,sNx+OLx
627 AREAforAtmFW(i,j,bi,bj) = AREA(i,j,bi,bj)
628 ENDDO
629 ENDDO
630 ENDDO
631 ENDDO
632 #endif
633
634 #ifdef ALLOW_OBCS
635 C-- In case we use scheme with a large stencil that extends into overlap:
636 C no longer needed with the right masking in advection & diffusion S/R.
637 c IF ( useOBCS ) THEN
638 c DO bj=myByLo(myThid),myByHi(myThid)
639 c DO bi=myBxLo(myThid),myBxHi(myThid)
640 c CALL OBCS_COPY_TRACER( HEFF(1-OLx,1-OLy,bi,bj),
641 c I 1, bi, bj, myThid )
642 c CALL OBCS_COPY_TRACER( AREA(1-OLx,1-OLy,bi,bj),
643 c I 1, bi, bj, myThid )
644 c CALL OBCS_COPY_TRACER( HSNOW(1-OLx,1-OLy,bi,bj),
645 c I 1, bi, bj, myThid )
646 #ifdef SEAICE_VARIABLE_SALINITY
647 c CALL OBCS_COPY_TRACER( HSALT(1-OLx,1-OLy,bi,bj),
648 c I 1, bi, bj, myThid )
649 #endif
650 c ENDDO
651 c ENDDO
652 c ENDIF
653 #endif /* ALLOW_OBCS */
654
655 #ifdef SEAICE_ALLOW_JFNK
656 C Computing this metric cannot be done in S/R SEAICE_INIT_FIXED
657 C where it belongs, because globalArea is only defined later after
658 C S/R PACKAGES_INIT_FIXED, so we move this computation here.
659 CALL SEAICE_MAP_RS2VEC( nVec, rAw, rAs,
660 & scalarProductMetric, .TRUE., myThid )
661 DO bj=myByLo(myThid),myByHi(myThid)
662 DO bi=myBxLo(myThid),myBxHi(myThid)
663 DO i=1,nVec
664 scalarProductMetric(i,1,bi,bj) =
665 & scalarProductMetric(i,1,bi,bj)/globalArea
666 ENDDO
667 ENDDO
668 ENDDO
669 #endif /* SEAICE_ALLOW_JFNK */
670
671 C--- Complete initialization
672 DO bj=myByLo(myThid),myByHi(myThid)
673 DO bi=myBxLo(myThid),myBxHi(myThid)
674 DO j=1-OLy,sNy+OLy
675 DO i=1-OLx,sNx+OLx
676 ZETA(i,j,bi,bj) = HEFF(i,j,bi,bj)*(1.0 _d 11)
677 ETA(i,j,bi,bj) = ZETA(i,j,bi,bj)/SEAICE_eccen**2
678 PRESS0(i,j,bi,bj) = SEAICE_strength*HEFF(i,j,bi,bj)
679 & *EXP(-SEAICE_cStar*(ONE-AREA(i,j,bi,bj)))
680 ZMAX(I,J,bi,bj) = SEAICE_zetaMaxFac*PRESS0(I,J,bi,bj)
681 ZMIN(i,j,bi,bj) = SEAICE_zetaMin
682 PRESS0(i,j,bi,bj) = PRESS0(i,j,bi,bj)*HEFFM(i,j,bi,bj)
683 ENDDO
684 ENDDO
685 IF ( useRealFreshWaterFlux .AND. .NOT.useThSIce ) THEN
686 DO j=1-OLy,sNy+OLy
687 DO i=1-OLx,sNx+OLx
688 sIceLoad(i,j,bi,bj) = HEFF(i,j,bi,bj)*SEAICE_rhoIce
689 & + HSNOW(i,j,bi,bj)*SEAICE_rhoSnow
690
691 ENDDO
692 ENDDO
693 ENDIF
694 ENDDO
695 ENDDO
696 #ifdef SEAICE_CGRID
697 C compute tensile strength factor k: tensileStrength = k*PRESS
698 C can be done in initialisation phase as long as it depends only
699 C on depth
700 IF ( SEAICE_tensilFac .NE. 0. _d 0 ) THEN
701 recip_tensilDepth = 0. _d 0
702 IF ( SEAICE_tensilDepth .GT. 0. _d 0 )
703 & recip_tensilDepth = 1. _d 0 / SEAICE_tensilDepth
704 DO bj=myByLo(myThid),myByHi(myThid)
705 DO bi=myBxLo(myThid),myBxHi(myThid)
706 DO j=1-OLy,sNy+OLy
707 DO i=1-OLx,sNx+OLx
708 tensileStrFac(i,j,bi,bj) = SEAICE_tensilFac*HEFFM(i,j,bi,bj)
709 & *exp(-ABS(R_low(i,j,bi,bj))*recip_tensilDepth)
710 ENDDO
711 ENDDO
712 ENDDO
713 ENDDO
714 ENDIF
715 #endif /* SEAICE_CGRID */
716
717 RETURN
718 END

  ViewVC Help
Powered by ViewVC 1.1.22