/[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.79 - (show annotations) (download)
Thu Dec 27 23:05:47 2012 UTC (11 years, 5 months ago) by gforget
Branch: MAIN
Changes since 1.78: +1 -5 lines
- remove legacy branch code.
- retire SEAICE_GROWTH_LEGACY and
  SEAICE_CAP_HEFF accordingly.

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_varia.F,v 1.78 2012/12/17 15:06:02 mlosch 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 #ifdef SEAICE_CGRID
46 _RS mask_uice
47 #endif
48 INTEGER k
49 #ifdef ALLOW_SITRACER
50 INTEGER iTr, jTh
51 #endif
52 #ifdef OBCS_UVICE_OLD
53 INTEGER I_obc, J_obc
54 #endif /* ALLOW_OBCS */
55
56 IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
57 kSurface = Nr
58 ELSE
59 kSurface = 1
60 ENDIF
61
62 C-- Initialise all variables in common blocks:
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 HEFF(i,j,bi,bj)=0. _d 0
68 AREA(i,j,bi,bj)=0. _d 0
69 CToM<<<
70 #ifdef SEAICE_ITD
71 DO k=1,nITD
72 AREAITD(i,j,k,bi,bj) =0. _d 0
73 HEFFITD(i,j,k,bi,bj) =0. _d 0
74 ENDDO
75 #endif
76 C>>>ToM
77 UICE(i,j,bi,bj)=0. _d 0
78 VICE(i,j,bi,bj)=0. _d 0
79 #ifdef SEAICE_ALLOW_FREEDRIFT
80 uice_fd(i,j,bi,bj)=0. _d 0
81 vice_fd(i,j,bi,bj)=0. _d 0
82 #endif
83 C
84 uIceNm1(i,j,bi,bj)=0. _d 0
85 vIceNm1(i,j,bi,bj)=0. _d 0
86 ETA (i,j,bi,bj) = 0. _d 0
87 etaZ(i,j,bi,bj) = 0. _d 0
88 ZETA(i,j,bi,bj) = 0. _d 0
89 FORCEX(i,j,bi,bj) = 0. _d 0
90 FORCEY(i,j,bi,bj) = 0. _d 0
91 uIceC(i,j,bi,bj) = 0. _d 0
92 vIceC(i,j,bi,bj) = 0. _d 0
93 #ifdef SEAICE_CGRID
94 seaiceMassC(i,j,bi,bj)=0. _d 0
95 seaiceMassU(i,j,bi,bj)=0. _d 0
96 seaiceMassV(i,j,bi,bj)=0. _d 0
97 stressDivergenceX(i,j,bi,bj) = 0. _d 0
98 stressDivergenceY(i,j,bi,bj) = 0. _d 0
99 # ifdef SEAICE_ALLOW_EVP
100 seaice_sigma1 (i,j,bi,bj) = 0. _d 0
101 seaice_sigma2 (i,j,bi,bj) = 0. _d 0
102 seaice_sigma12(i,j,bi,bj) = 0. _d 0
103 # endif /* SEAICE_ALLOW_EVP */
104 #else /* SEAICE_CGRID */
105 AMASS(i,j,bi,bj) = 0. _d 0
106 DAIRN(i,j,bi,bj) = 0. _d 0
107 WINDX(i,j,bi,bj) = 0. _d 0
108 WINDY(i,j,bi,bj) = 0. _d 0
109 GWATX(i,j,bi,bj) = 0. _d 0
110 GWATY(i,j,bi,bj) = 0. _d 0
111 #endif /* SEAICE_CGRID */
112 DWATN(i,j,bi,bj) = 0. _d 0
113 PRESS0(i,j,bi,bj) = 0. _d 0
114 FORCEX0(i,j,bi,bj)= 0. _d 0
115 FORCEY0(i,j,bi,bj)= 0. _d 0
116 ZMAX(i,j,bi,bj) = 0. _d 0
117 ZMIN(i,j,bi,bj) = 0. _d 0
118 HSNOW(i,j,bi,bj) = 0. _d 0
119 CToM<<<
120 #ifdef SEAICE_ITD
121 DO k=1,nITD
122 HSNOWITD(i,j,k,bi,bj)=0. _d 0
123 ENDDO
124 #endif
125 C>>>ToM
126 #ifdef SEAICE_VARIABLE_SALINITY
127 HSALT(i,j,bi,bj) = 0. _d 0
128 #endif
129 #ifdef ALLOW_SITRACER
130 DO iTr = 1, SItrMaxNum
131 SItracer(i,j,bi,bj,iTr) = 0. _d 0
132 SItrBucket(i,j,bi,bj,iTr) = 0. _d 0
133 c "ice concentration" tracer that should remain .EQ.1.
134 if (SItrName(iTr).EQ.'one') SItracer(i,j,bi,bj,iTr)=1. _d 0
135 ENDDO
136 DO jTh = 1, 5
137 SItrHEFF (i,j,bi,bj,jTh) = 0. _d 0
138 ENDDO
139 DO jTh = 1, 3
140 SItrAREA (i,j,bi,bj,jTh) = 0. _d 0
141 ENDDO
142 #endif
143 TICE(i,j,bi,bj) = 0. _d 0
144 DO k=1,MULTDIM
145 TICES(i,j,k,bi,bj)=0. _d 0
146 ENDDO
147 TAUX(i,j,bi,bj) = 0. _d 0
148 TAUY(i,j,bi,bj) = 0. _d 0
149 #ifdef ALLOW_SEAICE_COST_EXPORT
150 uHeffExportCell(i,j,bi,bj) = 0. _d 0
151 vHeffExportCell(i,j,bi,bj) = 0. _d 0
152 icevolMeanCell(i,j,bi,bj) = 0. _d 0
153 #endif
154 saltWtrIce(i,j,bi,bj) = 0. _d 0
155 frWtrIce(i,j,bi,bj) = 0. _d 0
156 #if (defined (ALLOW_MEAN_SFLUX_COST_CONTRIBUTION) || defined (ALLOW_SSH_GLOBMEAN_COST_CONTRIBUTION))
157 frWtrAtm(i,j,bi,bj) = 0. _d 0
158 AREAforAtmFW(i,j,bi,bj)=0. _d 0
159 #endif
160 ENDDO
161 ENDDO
162 ENDDO
163 ENDDO
164
165 #ifdef ALLOW_TIMEAVE
166 C Initialize averages to zero
167 DO bj = myByLo(myThid), myByHi(myThid)
168 DO bi = myBxLo(myThid), myBxHi(myThid)
169 CALL TIMEAVE_RESET( FUtave , 1, bi, bj, myThid )
170 CALL TIMEAVE_RESET( FVtave , 1, bi, bj, myThid )
171 CALL TIMEAVE_RESET( EmPmRtave, 1, bi, bj, myThid )
172 CALL TIMEAVE_RESET( QNETtave , 1, bi, bj, myThid )
173 CALL TIMEAVE_RESET( QSWtave , 1, bi, bj, myThid )
174 CALL TIMEAVE_RESET( UICEtave , 1, bi, bj, myThid )
175 CALL TIMEAVE_RESET( VICEtave , 1, bi, bj, myThid )
176 CALL TIMEAVE_RESET( HEFFtave , 1, bi, bj, myThid )
177 CALL TIMEAVE_RESET( AREAtave , 1, bi, bj, myThid )
178 SEAICE_timeAve(bi,bj) = ZERO
179 ENDDO
180 ENDDO
181 #endif /* ALLOW_TIMEAVE */
182
183 C-- Initialize (variable) grid info. As long as we allow masking of
184 C-- velocities outside of ice covered areas (in seaice_dynsolver)
185 C-- we need to re-initialize seaiceMaskU/V here for TAF/TAMC
186 #ifdef SEAICE_CGRID
187 DO bj=myByLo(myThid),myByHi(myThid)
188 DO bi=myBxLo(myThid),myBxHi(myThid)
189 DO j=1-OLy+1,sNy+OLy
190 DO i=1-OLx+1,sNx+OLx
191 seaiceMaskU(i,j,bi,bj)= 0.0 _d 0
192 seaiceMaskV(i,j,bi,bj)= 0.0 _d 0
193 mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i-1,j ,bi,bj)
194 IF(mask_uice.GT.1.5 _d 0) seaiceMaskU(i,j,bi,bj)=1.0 _d 0
195 mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i ,j-1,bi,bj)
196 IF(mask_uice.GT.1.5 _d 0) seaiceMaskV(i,j,bi,bj)=1.0 _d 0
197 ENDDO
198 ENDDO
199 ENDDO
200 ENDDO
201 #endif /* SEAICE_CGRID */
202
203 DO bj=myByLo(myThid),myByHi(myThid)
204 DO bi=myBxLo(myThid),myBxHi(myThid)
205 #ifdef OBCS_UVICE_OLD
206 #ifdef SEAICE_CGRID
207 IF (useOBCS) THEN
208 C-- If OBCS is turned on, close southern and western boundaries
209 DO i=1-OLx,sNx+OLx
210 C Southern boundary
211 J_obc = OB_Js(i,bi,bj)
212 IF ( J_obc.NE.OB_indexNone ) THEN
213 seaiceMaskU(i,J_obc,bi,bj)= 0.0 _d 0
214 seaiceMaskV(i,J_obc,bi,bj)= 0.0 _d 0
215 ENDIF
216 ENDDO
217 DO j=1-OLy,sNy+OLy
218 C Western boundary
219 I_obc = OB_Iw(j,bi,bj)
220 IF ( I_obc.NE.OB_indexNone ) THEN
221 seaiceMaskU(I_obc,j,bi,bj)= 0.0 _d 0
222 seaiceMaskV(I_obc,j,bi,bj)= 0.0 _d 0
223 ENDIF
224 ENDDO
225 ENDIF
226 #endif /* SEAICE_CGRID */
227 #endif /* OBCS_UVICE_OLD */
228
229 DO j=1-OLy,sNy+OLy
230 DO i=1-OLx,sNx+OLx
231 TICE(i,j,bi,bj)=273.0 _d 0
232 DO k=1,MULTDIM
233 TICES(i,j,k,bi,bj)=273.0 _d 0
234 ENDDO
235 #ifndef SEAICE_CGRID
236 AMASS (i,j,bi,bj)=1000.0 _d 0
237 #else
238 seaiceMassC(i,j,bi,bj)=1000.0 _d 0
239 seaiceMassU(i,j,bi,bj)=1000.0 _d 0
240 seaiceMassV(i,j,bi,bj)=1000.0 _d 0
241 #endif
242 ENDDO
243 ENDDO
244
245 ENDDO
246 ENDDO
247
248 C-- Update overlap regions
249 #ifdef SEAICE_CGRID
250 CALL EXCH_UV_XY_RL(seaiceMaskU,seaiceMaskV,.FALSE.,myThid)
251 #else
252 _EXCH_XY_RS(UVM, myThid)
253 #endif
254
255 C-- Now lets look at all these beasts
256 IF ( debugLevel .GE. debLevC ) THEN
257 CALL PLOT_FIELD_XYRL( HEFFM , 'Current HEFFM ' ,
258 & nIter0, myThid )
259 #ifdef SEAICE_CGRID
260 CALL PLOT_FIELD_XYRL( seaiceMaskU, 'Current seaiceMaskU',
261 & nIter0, myThid )
262 CALL PLOT_FIELD_XYRL( seaiceMaskV, 'Current seaiceMaskV',
263 & nIter0, myThid )
264 #else
265 CALL PLOT_FIELD_XYRS( UVM , 'Current UVM ' ,
266 & nIter0, myThid )
267 #endif
268 ENDIF
269
270 C-- Set model variables to initial/restart conditions
271 IF ( .NOT. ( startTime .EQ. baseTime .AND. nIter0 .EQ. 0
272 & .AND. pickupSuff .EQ. ' ') ) THEN
273
274 CALL SEAICE_READ_PICKUP ( myThid )
275
276 ELSE
277
278 DO bj=myByLo(myThid),myByHi(myThid)
279 DO bi=myBxLo(myThid),myBxHi(myThid)
280 DO j=1-OLy,sNy+OLy
281 DO i=1-OLx,sNx+OLx
282 HEFF(i,j,bi,bj)=SEAICE_initialHEFF*HEFFM(i,j,bi,bj)
283 UICE(i,j,bi,bj)=ZERO
284 VICE(i,j,bi,bj)=ZERO
285 ENDDO
286 ENDDO
287 ENDDO
288 ENDDO
289
290 C-- Read initial sea-ice velocity from file (if available)
291 IF ( uIceFile .NE. ' ' )
292 & CALL READ_FLD_XY_RL( uIceFile, ' ', uIce, 0, myThid )
293 IF ( vIceFile .NE. ' ' )
294 & CALL READ_FLD_XY_RL( vIceFile, ' ', vIce, 0, myThid )
295 IF ( uIceFile .NE. ' ' .OR. vIceFile .NE. ' ' ) THEN
296 #ifdef SEAICE_CGRID
297 DO bj=myByLo(myThid),myByHi(myThid)
298 DO bi=myBxLo(myThid),myBxHi(myThid)
299 DO j=1-OLy,sNy+OLy
300 DO i=1-OLx,sNx+OLx
301 uIce(i,j,bi,bj) = uIce(i,j,bi,bj)*seaiceMaskU(i,j,bi,bj)
302 vIce(i,j,bi,bj) = vIce(i,j,bi,bj)*seaiceMaskV(i,j,bi,bj)
303 ENDDO
304 ENDDO
305 ENDDO
306 ENDDO
307 #endif /* SEAICE_CGRID */
308 CALL EXCH_UV_XY_RL( uIce, vIce, .TRUE., myThid )
309 ENDIF
310
311 C-- Read initial sea-ice thickness from file if available.
312 IF ( HeffFile .NE. ' ' ) THEN
313 CALL READ_FLD_XY_RL( HeffFile, ' ', HEFF, 0, myThid )
314 _EXCH_XY_RL(HEFF,myThid)
315 DO bj=myByLo(myThid),myByHi(myThid)
316 DO bi=myBxLo(myThid),myBxHi(myThid)
317 DO j=1-OLy,sNy+OLy
318 DO i=1-OLx,sNx+OLx
319 HEFF(i,j,bi,bj) = MAX(HEFF(i,j,bi,bj),ZERO)
320 ENDDO
321 ENDDO
322 ENDDO
323 ENDDO
324 ENDIF
325
326 DO bj=myByLo(myThid),myByHi(myThid)
327 DO bi=myBxLo(myThid),myBxHi(myThid)
328 DO j=1-OLy,sNy+OLy
329 DO i=1-OLx,sNx+OLx
330 IF(HEFF(i,j,bi,bj).GT.ZERO) AREA(i,j,bi,bj)=ONE
331 ENDDO
332 ENDDO
333 ENDDO
334 ENDDO
335
336 C-- Read initial sea-ice area from file if available.
337 IF ( AreaFile .NE. ' ' ) THEN
338 CALL READ_FLD_XY_RL( AreaFile, ' ', AREA, 0, myThid )
339 _EXCH_XY_RL(AREA,myThid)
340 DO bj=myByLo(myThid),myByHi(myThid)
341 DO bi=myBxLo(myThid),myBxHi(myThid)
342 DO j=1-OLy,sNy+OLy
343 DO i=1-OLx,sNx+OLx
344 AREA(i,j,bi,bj) = MAX(AREA(i,j,bi,bj),ZERO)
345 AREA(i,j,bi,bj) = MIN(AREA(i,j,bi,bj),ONE)
346 IF ( AREA(i,j,bi,bj) .LE. ZERO ) HEFF(i,j,bi,bj) = ZERO
347 IF ( HEFF(i,j,bi,bj) .LE. ZERO ) AREA(i,j,bi,bj) = ZERO
348 ENDDO
349 ENDDO
350 ENDDO
351 ENDDO
352 ENDIF
353
354 DO bj=myByLo(myThid),myByHi(myThid)
355 DO bi=myBxLo(myThid),myBxHi(myThid)
356 DO j=1-OLy,sNy+OLy
357 DO i=1-OLx,sNx+OLx
358 HSNOW(i,j,bi,bj) = 0.2 _d 0 * AREA(i,j,bi,bj)
359 ENDDO
360 ENDDO
361 ENDDO
362 ENDDO
363
364 C-- Read initial snow thickness from file if available.
365 IF ( HsnowFile .NE. ' ' ) THEN
366 CALL READ_FLD_XY_RL( HsnowFile, ' ', HSNOW, 0, myThid )
367 _EXCH_XY_RL(HSNOW,myThid)
368 DO bj=myByLo(myThid),myByHi(myThid)
369 DO bi=myBxLo(myThid),myBxHi(myThid)
370 DO j=1-OLy,sNy+OLy
371 DO i=1-OLx,sNx+OLx
372 HSNOW(i,j,bi,bj) = MAX(HSNOW(i,j,bi,bj),ZERO)
373 ENDDO
374 ENDDO
375 ENDDO
376 ENDDO
377 ENDIF
378
379 #ifdef SEAICE_ITD
380 DO bj=myByLo(myThid),myByHi(myThid)
381 DO bi=myBxLo(myThid),myBxHi(myThid)
382 DO j=1-OLy,sNy+OLy
383 DO i=1-OLx,sNx+OLx
384 AREAITD(I,J,1,bi,bj) =AREA(I,J,bi,bj)
385 HEFFITD(I,J,1,bi,bj) =HEFF(I,J,bi,bj)
386 HSNOWITD(I,J,1,bi,bj)=HSNOW(I,J,bi,bj)
387 ENDDO
388 ENDDO
389 ENDDO
390 ENDDO
391 C DO bj=myByLo(myThid),myByHi(myThid)
392 C DO bi=myBxLo(myThid),myBxHi(myThid)
393 C CALL SEAICE_ITD_REDIST(bi, bj, myTime, myIter, myThid)
394 C ENDDO
395 C ENDDO
396 C it is sufficient to have first call of SEAICE_ITD_REDIST before SEAICE_GROWTH
397 #endif
398
399 #ifdef SEAICE_VARIABLE_SALINITY
400 DO bj=myByLo(myThid),myByHi(myThid)
401 DO bi=myBxLo(myThid),myBxHi(myThid)
402 DO j=1-OLy,sNy+OLy
403 DO i=1-OLx,sNx+OLx
404 HSALT(i,j,bi,bj)=HEFF(i,j,bi,bj)*salt(i,j,kSurface,bi,bj)*
405 & SEAICE_rhoIce*SEAICE_saltFrac
406 cif & ICE2WATR*rhoConstFresh*SEAICE_saltFrac
407
408 ENDDO
409 ENDDO
410 ENDDO
411 ENDDO
412
413 C-- Read initial sea ice salinity from file if available.
414 IF ( HsaltFile .NE. ' ' ) THEN
415 CALL READ_FLD_XY_RL( HsaltFile, ' ', HSALT, 0, myThid )
416 _EXCH_XY_RL(HSALT,myThid)
417 ENDIF
418 #endif /* SEAICE_VARIABLE_SALINITY */
419
420 #ifdef ALLOW_SITRACER
421 C-- Read initial sea ice age from file if available.
422 DO iTr = 1, SItrMaxNum
423 IF ( SItrFile(iTr) .NE. ' ' ) THEN
424 CALL READ_FLD_XY_RL( siTrFile(iTr), ' ',
425 & SItracer(1-OLx,1-OLy,1,1,iTr), 0, myThid )
426 _EXCH_XY_RL(SItracer(1-OLx,1-OLy,1,1,iTr),myThid)
427 ENDIF
428 ENDDO
429 #endif /* ALLOW_SITRACER */
430
431 ENDIF
432
433 #if (defined (ALLOW_MEAN_SFLUX_COST_CONTRIBUTION) || defined (ALLOW_SSH_GLOBMEAN_COST_CONTRIBUTION))
434 DO bj=myByLo(myThid),myByHi(myThid)
435 DO bi=myBxLo(myThid),myBxHi(myThid)
436 DO j=1-OLy,sNy+OLy
437 DO i=1-OLx,sNx+OLx
438 AREAforAtmFW(i,j,bi,bj) = AREA(i,j,bi,bj)
439 ENDDO
440 ENDDO
441 ENDDO
442 ENDDO
443 #endif
444
445 #ifdef ALLOW_OBCS
446 C-- In case we use scheme with a large stencil that extends into overlap:
447 C no longer needed with the right masking in advection & diffusion S/R.
448 c IF ( useOBCS ) THEN
449 c DO bj=myByLo(myThid),myByHi(myThid)
450 c DO bi=myBxLo(myThid),myBxHi(myThid)
451 c CALL OBCS_COPY_TRACER( HEFF(1-OLx,1-OLy,bi,bj),
452 c I 1, bi, bj, myThid )
453 c CALL OBCS_COPY_TRACER( AREA(1-OLx,1-OLy,bi,bj),
454 c I 1, bi, bj, myThid )
455 c CALL OBCS_COPY_TRACER( HSNOW(1-OLx,1-OLy,bi,bj),
456 c I 1, bi, bj, myThid )
457 #ifdef SEAICE_VARIABLE_SALINITY
458 c CALL OBCS_COPY_TRACER( HSALT(1-OLx,1-OLy,bi,bj),
459 c I 1, bi, bj, myThid )
460 #endif
461 c ENDDO
462 c ENDDO
463 c ENDIF
464 #endif /* ALLOW_OBCS */
465
466 #ifdef SEAICE_ALLOW_JFNK
467 C Computing this metric cannot be done in S/R SEAICE_INIT_FIXED
468 C where it belongs, because globalArea is only defined later after
469 C S/R PACKAGES_INIT_FIXED, so we move this computation here.
470 CALL SEAICE_MAP2VEC( nVec, rAw, rAs,
471 & scalarProductMetric, .TRUE., myThid )
472 DO bj=myByLo(myThid),myByHi(myThid)
473 DO bi=myBxLo(myThid),myBxHi(myThid)
474 DO i=1,nVec
475 scalarProductMetric(i,1,bi,bj) =
476 & scalarProductMetric(i,1,bi,bj)/globalArea
477 ENDDO
478 ENDDO
479 ENDDO
480 #endif /* SEAICE_ALLOW_JFNK */
481
482 C--- Complete initialization
483 PSTAR = SEAICE_strength
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 ZETA(i,j,bi,bj) = HEFF(i,j,bi,bj)*(1.0 _d 11)
489 ETA(i,j,bi,bj) = ZETA(i,j,bi,bj)/SEAICE_eccen**2
490 PRESS0(i,j,bi,bj) = PSTAR*HEFF(i,j,bi,bj)
491 & *EXP(-20.0 _d 0*(ONE-AREA(i,j,bi,bj)))
492 ZMAX(I,J,bi,bj) = SEAICE_zetaMaxFac*PRESS0(I,J,bi,bj)
493 ZMIN(i,j,bi,bj) = SEAICE_zetaMin
494 PRESS0(i,j,bi,bj) = PRESS0(i,j,bi,bj)*HEFFM(i,j,bi,bj)
495 ENDDO
496 ENDDO
497 IF ( useRealFreshWaterFlux .AND. .NOT.useThSIce ) THEN
498 DO j=1-OLy,sNy+OLy
499 DO i=1-OLx,sNx+OLx
500 sIceLoad(i,j,bi,bj) = HEFF(i,j,bi,bj)*SEAICE_rhoIce
501 & + HSNOW(i,j,bi,bj)*SEAICE_rhoSnow
502
503 ENDDO
504 ENDDO
505 ENDIF
506 ENDDO
507 ENDDO
508
509 RETURN
510 END

  ViewVC Help
Powered by ViewVC 1.1.22