/[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.75 - (show annotations) (download)
Mon Oct 22 21:14:06 2012 UTC (11 years, 7 months ago) by heimbach
Branch: MAIN
Changes since 1.74: +36 -1 lines
Step 1 of merging ice-thickness distribution (ITD) code from
MITgcm_contrib/torge/itd/code/ to main repository
(author: Torge Martin)

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

  ViewVC Help
Powered by ViewVC 1.1.22