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