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

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

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


Revision 1.53 - (hide annotations) (download)
Mon Nov 8 21:50:56 2010 UTC (13 years, 6 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint62o, checkpoint62n, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62t
Changes since 1.52: +14 -9 lines
undoing 1.50 to 1.51 modification

1 dimitri 1.53 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_varia.F,v 1.52 2010/11/08 17:38:04 jmc Exp $
2 heimbach 1.1 C $Name: $
3    
4     #include "SEAICE_OPTIONS.h"
5    
6     CStartOfInterface
7     SUBROUTINE SEAICE_INIT_VARIA( myThid )
8 jmc 1.47 C *==========================================================*
9 heimbach 1.1 C | SUBROUTINE SEAICE_INIT_VARIA |
10     C | o Initialization of sea ice model. |
11 jmc 1.47 C *==========================================================*
12     C *==========================================================*
13 heimbach 1.1 IMPLICIT NONE
14 jmc 1.23
15 heimbach 1.1 C === Global variables ===
16     #include "SIZE.h"
17     #include "EEPARAMS.h"
18     #include "PARAMS.h"
19     #include "GRID.h"
20 dimitri 1.21 #include "DYNVARS.h"
21 jmc 1.50 #include "FFIELDS.h"
22     #include "SEAICE_PARAMS.h"
23 heimbach 1.1 #include "SEAICE.h"
24 jmc 1.50 #include "SEAICE_TAVE.h"
25 heimbach 1.1 #ifdef ALLOW_EXCH2
26 jmc 1.38 # include "W2_EXCH2_SIZE.h"
27 dimitri 1.24 # include "W2_EXCH2_TOPOLOGY.h"
28     #endif
29 dimitri 1.11 #ifdef ALLOW_OBCS
30 dimitri 1.24 # include "OBCS_OPTIONS.h"
31     # include "OBCS.h"
32     #endif
33 heimbach 1.1
34     C === Routine arguments ===
35     C myThid - Thread no. that called this routine.
36     INTEGER myThid
37     CEndOfInterface
38 jmc 1.23
39 heimbach 1.1 C === Local variables ===
40     C i,j,k,bi,bj - Loop counters
41    
42 jmc 1.50 INTEGER i, j, bi, bj
43 mlosch 1.8 _RL PSTAR
44 heimbach 1.1 _RS mask_uice
45 jmc 1.50 INTEGER myIter, kSurface
46     #ifdef SEAICE_MULTICATEGORY
47     INTEGER k
48     #endif
49 dimitri 1.53 #ifdef ALLOW_OBCS
50     INTEGER I_obc, J_obc
51     #endif /* ALLOW_OBCS */
52 jmc 1.50 #ifdef ALLOW_EXCH2
53 dimitri 1.51 # ifndef SEAICE_CGRID
54 jmc 1.50 INTEGER myTile
55 dimitri 1.51 # endif
56 jmc 1.50 #endif
57 heimbach 1.1
58 jmc 1.23 IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
59 mlosch 1.44 kSurface = Nr
60 jmc 1.23 ELSE
61 mlosch 1.44 kSurface = 1
62 jmc 1.23 ENDIF
63 dimitri 1.11
64 jmc 1.23 C-- Initialise all variables in common blocks:
65 heimbach 1.1 DO bj=myByLo(myThid),myByHi(myThid)
66     DO bi=myBxLo(myThid),myBxHi(myThid)
67 mlosch 1.42 DO j=1-OLy,sNy+OLy
68     DO i=1-OLx,sNx+OLx
69 mlosch 1.43 HEFF(i,j,bi,bj)=0. _d 0
70     AREA(i,j,bi,bj)=0. _d 0
71     UICE(i,j,bi,bj)=0. _d 0
72     VICE(i,j,bi,bj)=0. _d 0
73     C
74 mlosch 1.42 uIceNm1(i,j,bi,bj)=0. _d 0
75     vIceNm1(i,j,bi,bj)=0. _d 0
76     areaNm1(i,j,bi,bj)=0. _d 0
77     hEffNm1(i,j,bi,bj)=0. _d 0
78     ENDDO
79     ENDDO
80 jmc 1.23 #ifdef SEAICE_MULTICATEGORY
81     DO k=1,MULTDIM
82     DO j=1-OLy,sNy+OLy
83     DO i=1-OLx,sNx+OLx
84     TICES(i,j,k,bi,bj)=0. _d 0
85 heimbach 1.1 ENDDO
86     ENDDO
87     ENDDO
88 jmc 1.23 #endif
89     DO j=1-OLy,sNy+OLy
90     DO i=1-OLx,sNx+OLx
91     ETA (i,j,bi,bj) = 0. _d 0
92     ZETA(i,j,bi,bj) = 0. _d 0
93     DRAGS(i,j,bi,bj) = 0. _d 0
94     DRAGA(i,j,bi,bj) = 0. _d 0
95     FORCEX(i,j,bi,bj) = 0. _d 0
96     FORCEY(i,j,bi,bj) = 0. _d 0
97     UICEC(i,j,bi,bj) = 0. _d 0
98     VICEC(i,j,bi,bj) = 0. _d 0
99     #ifdef SEAICE_CGRID
100     seaiceMassC(i,j,bi,bj)=0. _d 0
101     seaiceMassU(i,j,bi,bj)=0. _d 0
102     seaiceMassV(i,j,bi,bj)=0. _d 0
103     stressDivergenceX(i,j,bi,bj) = 0. _d 0
104     stressDivergenceY(i,j,bi,bj) = 0. _d 0
105 heimbach 1.39 # ifdef SEAICE_ALLOW_EVP
106 jmc 1.23 seaice_sigma1 (i,j,bi,bj) = 0. _d 0
107     seaice_sigma2 (i,j,bi,bj) = 0. _d 0
108     seaice_sigma12(i,j,bi,bj) = 0. _d 0
109     # endif /* SEAICE_ALLOW_EVP */
110     #else /* SEAICE_CGRID */
111     AMASS(i,j,bi,bj) = 0. _d 0
112     DAIRN(i,j,bi,bj) = 0. _d 0
113     WINDX(i,j,bi,bj) = 0. _d 0
114     WINDY(i,j,bi,bj) = 0. _d 0
115 dimitri 1.25 GWATX(i,j,bi,bj) = 0. _d 0
116     GWATY(i,j,bi,bj) = 0. _d 0
117 jmc 1.23 #endif /* SEAICE_CGRID */
118     DWATN(i,j,bi,bj) = 0. _d 0
119     PRESS0(i,j,bi,bj) = 0. _d 0
120     FORCEX0(i,j,bi,bj)= 0. _d 0
121     FORCEY0(i,j,bi,bj)= 0. _d 0
122     ZMAX(i,j,bi,bj) = 0. _d 0
123     ZMIN(i,j,bi,bj) = 0. _d 0
124     HSNOW(i,j,bi,bj) = 0. _d 0
125 dimitri 1.18 #ifdef SEAICE_SALINITY
126 jmc 1.23 HSALT(i,j,bi,bj) = 0. _d 0
127     #endif
128 dimitri 1.33 #ifdef SEAICE_AGE
129     ICEAGE(i,j,bi,bj) = 0. _d 0
130     #endif
131 jmc 1.23 YNEG (i,j,bi,bj) = 0. _d 0
132     RIVER(i,j,bi,bj) = 0. _d 0
133     TMIX(i,j,bi,bj) = 0. _d 0
134     TICE(i,j,bi,bj) = 0. _d 0
135     TAUX(i,j,bi,bj) = 0. _d 0
136     TAUY(i,j,bi,bj) = 0. _d 0
137     #ifdef ALLOW_SEAICE_COST_EXPORT
138     uHeffExportCell(i,j,bi,bj) = 0. _d 0
139     vHeffExportCell(i,j,bi,bj) = 0. _d 0
140 dimitri 1.18 #endif
141 heimbach 1.32 saltWtrIce(i,j,bi,bj) = 0. _d 0
142     frWtrIce(i,j,bi,bj) = 0. _d 0
143 heimbach 1.49 #ifdef ALLOW_MEAN_SFLUX_COST_CONTRIBUTION
144     frWtrAtm(i,j,bi,bj) = 0. _d 0
145     #endif
146 heimbach 1.1 ENDDO
147     ENDDO
148     ENDDO
149     ENDDO
150    
151 jmc 1.50 #ifdef ALLOW_TIMEAVE
152     C Initialize averages to zero
153     DO bj = myByLo(myThid), myByHi(myThid)
154     DO bi = myBxLo(myThid), myBxHi(myThid)
155     CALL TIMEAVE_RESET( FUtave , 1, bi, bj, myThid )
156     CALL TIMEAVE_RESET( FVtave , 1, bi, bj, myThid )
157     CALL TIMEAVE_RESET( EmPmRtave, 1, bi, bj, myThid )
158     CALL TIMEAVE_RESET( QNETtave , 1, bi, bj, myThid )
159     CALL TIMEAVE_RESET( QSWtave , 1, bi, bj, myThid )
160     CALL TIMEAVE_RESET( UICEtave , 1, bi, bj, myThid )
161     CALL TIMEAVE_RESET( VICEtave , 1, bi, bj, myThid )
162     CALL TIMEAVE_RESET( HEFFtave , 1, bi, bj, myThid )
163     CALL TIMEAVE_RESET( AREAtave , 1, bi, bj, myThid )
164     SEAICE_timeAve(bi,bj) = ZERO
165     ENDDO
166     ENDDO
167     #endif /* ALLOW_TIMEAVE */
168    
169 mlosch 1.40 C-- Initialize (variable) grid info. As long as we allow masking of
170     C-- velocities outside of ice covered areas (in seaice_dynsolver)
171     C-- we need to re-initialize seaiceMaskU/V here for TAF/TAMC
172 heimbach 1.1 DO bj=myByLo(myThid),myByHi(myThid)
173     DO bi=myBxLo(myThid),myBxHi(myThid)
174 jmc 1.23 DO j=1-OLy+1,sNy+OLy
175     DO i=1-OLx+1,sNx+OLx
176 heimbach 1.49 #ifdef SEAICE_CGRID
177 jmc 1.23 seaiceMaskU(i,j,bi,bj)= 0.0 _d 0
178     seaiceMaskV(i,j,bi,bj)= 0.0 _d 0
179     mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i-1,j ,bi,bj)
180 heimbach 1.45 IF(mask_uice.GT.1.5 _d 0) seaiceMaskU(i,j,bi,bj)=1.0 _d 0
181 jmc 1.23 mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i ,j-1,bi,bj)
182 heimbach 1.45 IF(mask_uice.GT.1.5 _d 0) seaiceMaskV(i,j,bi,bj)=1.0 _d 0
183 heimbach 1.49 #else
184     UVM(i,j,bi,bj)= 0.0 _d 0
185     #endif /* SEAICE_CGRID */
186 dimitri 1.22 ENDDO
187     ENDDO
188 mlosch 1.40 ENDDO
189     ENDDO
190 dimitri 1.11
191 mlosch 1.40 DO bj=myByLo(myThid),myByHi(myThid)
192     DO bi=myBxLo(myThid),myBxHi(myThid)
193 dimitri 1.11 #ifdef ALLOW_OBCS
194 dimitri 1.16 IF (useOBCS) THEN
195 dimitri 1.22 C-- If OBCS is turned on, close southern and western boundaries
196 dimitri 1.11 #ifdef ALLOW_OBCS_SOUTH
197 jmc 1.23 DO i=1-Olx,sNx+Olx
198 dimitri 1.11 C Southern boundary
199 dimitri 1.53 J_obc = OB_Js(i,bi,bj)
200     IF (J_obc.NE.0) THEN
201 dimitri 1.22 #ifdef SEAICE_CGRID
202 dimitri 1.53 seaiceMaskU(i,J_obc,bi,bj)= 0.0 _d 0
203     seaiceMaskV(i,J_obc,bi,bj)= 0.0 _d 0
204 dimitri 1.22 #else
205 dimitri 1.53 UVM(i,J_obc,bi,bj)=0. _d 0
206 dimitri 1.22 #endif /* SEAICE_CGRID */
207 dimitri 1.16 ENDIF
208     ENDDO
209 dimitri 1.22 #endif /* ALLOW_OBCS_SOUTH */
210 dimitri 1.11 #ifdef ALLOW_OBCS_WEST
211 jmc 1.23 DO j=1-Oly,sNy+Oly
212 dimitri 1.11 C Western boundary
213 dimitri 1.53 I_obc=OB_Iw(j,bi,bj)
214     IF (I_obc.NE.0) THEN
215 dimitri 1.22 #ifdef SEAICE_CGRID
216 dimitri 1.53 seaiceMaskU(I_obc,j,bi,bj)= 0.0 _d 0
217     seaiceMaskV(I_obc,j,bi,bj)= 0.0 _d 0
218 dimitri 1.22 #else
219 dimitri 1.53 UVM(I_obc,j,bi,bj)=0. _d 0
220 dimitri 1.22 #endif /* SEAICE_CGRID */
221 dimitri 1.16 ENDIF
222     ENDDO
223 dimitri 1.22 #endif /* ALLOW_OBCS_WEST */
224 dimitri 1.16 ENDIF
225 dimitri 1.11 #endif /* ALLOW_OBCS */
226    
227 heimbach 1.1 #ifdef ALLOW_EXCH2
228 dimitri 1.11 #ifndef SEAICE_CGRID
229 heimbach 1.1 C-- Special stuff for cubed sphere: assume grid is rectangular and
230     C set UV mask to zero except for Arctic and Antarctic cube faces.
231     IF (useCubedSphereExchange) THEN
232 jmc 1.46 myTile = W2_myTileList(bi,bj)
233 heimbach 1.1 IF ( exch2_myFace(myTile) .EQ. 1 .OR.
234     & exch2_myFace(myTile) .EQ. 2 .OR.
235     & exch2_myFace(myTile) .EQ. 4 .OR.
236     & exch2_myFace(myTile) .EQ. 5 ) THEN
237 jmc 1.23 DO j=1-OLy,sNy+OLy
238     DO i=1-OLx,sNx+OLx
239 heimbach 1.10 UVM(i,j,bi,bj)=0. _d 0
240 heimbach 1.1 ENDDO
241     ENDDO
242     ELSEIF ( exch2_isWedge(myTile) .EQ. 1 ) THEN
243 jmc 1.23 i=1
244     DO j=1-OLy,sNy+OLy
245 heimbach 1.10 UVM(i,j,bi,bj)=0. _d 0
246 heimbach 1.1 ENDDO
247     ELSEIF ( exch2_isSedge(myTile) .EQ. 1 ) THEN
248 jmc 1.23 j=1
249     DO i=1-OLx,sNx+OLx
250 heimbach 1.10 UVM(i,j,bi,bj)=0. _d 0
251 heimbach 1.1 ENDDO
252     ENDIF
253     ENDIF
254 dimitri 1.11 #endif /* SEAICE_CGRID */
255 heimbach 1.1 #endif /* ALLOW_EXCH2 */
256    
257     DO j=1-OLy,sNy+OLy
258     DO i=1-OLx,sNx+OLx
259 jmc 1.23 TICE(i,j,bi,bj)=273.0 _d 0
260 mlosch 1.2 #ifdef SEAICE_MULTICATEGORY
261 heimbach 1.1 DO k=1,MULTDIM
262 jmc 1.23 TICES(i,j,k,bi,bj)=273.0 _d 0
263 heimbach 1.1 ENDDO
264 mlosch 1.2 #endif /* SEAICE_MULTICATEGORY */
265 heimbach 1.1 #ifndef SEAICE_CGRID
266 jmc 1.23 AMASS (i,j,bi,bj)=1000.0 _d 0
267 heimbach 1.1 #else
268 jmc 1.23 seaiceMassC(i,j,bi,bj)=1000.0 _d 0
269     seaiceMassU(i,j,bi,bj)=1000.0 _d 0
270     seaiceMassV(i,j,bi,bj)=1000.0 _d 0
271 heimbach 1.1 #endif
272     ENDDO
273     ENDDO
274    
275     ENDDO
276     ENDDO
277    
278     C-- Update overlap regions
279     #ifdef SEAICE_CGRID
280     CALL EXCH_UV_XY_RL(seaiceMaskU,seaiceMaskV,.FALSE.,myThid)
281     #else
282 jmc 1.37 _EXCH_XY_RL(UVM, myThid)
283 heimbach 1.1 #endif
284    
285     C-- Now lets look at all these beasts
286     IF ( debugLevel .GE. debLevB ) THEN
287     myIter=0
288     CALL PLOT_FIELD_XYRL( HEFFM , 'Current HEFFM ' ,
289     & myIter, myThid )
290     #ifdef SEAICE_CGRID
291     CALL PLOT_FIELD_XYRL( seaiceMaskU, 'Current seaiceMaskU',
292     & myIter, myThid )
293     CALL PLOT_FIELD_XYRL( seaiceMaskV, 'Current seaiceMaskV',
294     & myIter, myThid )
295     #else
296     CALL PLOT_FIELD_XYRL( UVM , 'Current UVM ' ,
297     & myIter, myThid )
298     #endif
299     ENDIF
300    
301     C-- Set model variables to initial/restart conditions
302 mlosch 1.26 IF ( .NOT. ( startTime .EQ. baseTime .AND. nIter0 .EQ. 0
303     & .AND. pickupSuff .EQ. ' ') ) THEN
304 heimbach 1.1
305     CALL SEAICE_READ_PICKUP ( myThid )
306    
307     ELSE
308    
309     DO bj=myByLo(myThid),myByHi(myThid)
310     DO bi=myBxLo(myThid),myBxHi(myThid)
311     DO j=1-OLy,sNy+OLy
312     DO i=1-OLx,sNx+OLx
313 jmc 1.23 TMIX(i,j,bi,bj)=TICE(i,j,bi,bj)
314 mlosch 1.43 HEFF(i,j,bi,bj)=SEAICE_initialHEFF*HEFFM(i,j,bi,bj)
315     UICE(i,j,bi,bj)=ZERO
316     VICE(i,j,bi,bj)=ZERO
317 heimbach 1.1 ENDDO
318     ENDDO
319     ENDDO
320     ENDDO
321    
322     C-- Read initial sea-ice thickness from file if available.
323     IF ( HeffFile .NE. ' ' ) THEN
324 mlosch 1.43 CALL READ_FLD_XY_RL( HeffFile, ' ', HEFF, 0, myThid )
325     _EXCH_XY_RL(HEFF,myThid)
326 heimbach 1.1 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 mlosch 1.43 HEFF(i,j,bi,bj) = MAX(HEFF(i,j,bi,bj),ZERO)
331 heimbach 1.1 ENDDO
332     ENDDO
333     ENDDO
334     ENDDO
335     ENDIF
336    
337     DO bj=myByLo(myThid),myByHi(myThid)
338     DO bi=myBxLo(myThid),myBxHi(myThid)
339     DO j=1-OLy,sNy+OLy
340     DO i=1-OLx,sNx+OLx
341 mlosch 1.43 IF(HEFF(i,j,bi,bj).GT.ZERO) AREA(i,j,bi,bj)=ONE
342 heimbach 1.1 ENDDO
343     ENDDO
344     ENDDO
345     ENDDO
346 jmc 1.23
347 dimitri 1.27 C-- Read initial sea-ice area from file if available.
348 mlosch 1.7 IF ( AreaFile .NE. ' ' ) THEN
349 mlosch 1.43 CALL READ_FLD_XY_RL( AreaFile, ' ', AREA, 0, myThid )
350     _EXCH_XY_RL(AREA,myThid)
351 mlosch 1.7 DO bj=myByLo(myThid),myByHi(myThid)
352     DO bi=myBxLo(myThid),myBxHi(myThid)
353     DO j=1-OLy,sNy+OLy
354     DO i=1-OLx,sNx+OLx
355 mlosch 1.43 AREA(i,j,bi,bj) = MAX(AREA(i,j,bi,bj),ZERO)
356     AREA(i,j,bi,bj) = MIN(AREA(i,j,bi,bj),ONE)
357     IF ( AREA(i,j,bi,bj) .LE. ZERO ) HEFF(i,j,bi,bj) = ZERO
358     IF ( HEFF(i,j,bi,bj) .LE. ZERO ) AREA(i,j,bi,bj) = ZERO
359 mlosch 1.7 ENDDO
360     ENDDO
361     ENDDO
362     ENDDO
363     ENDIF
364    
365     DO bj=myByLo(myThid),myByHi(myThid)
366     DO bi=myBxLo(myThid),myBxHi(myThid)
367     DO j=1-OLy,sNy+OLy
368     DO i=1-OLx,sNx+OLx
369 jmc 1.47 HSNOW(i,j,bi,bj) = 0.2 _d 0 * AREA(i,j,bi,bj)
370 mlosch 1.7 ENDDO
371     ENDDO
372     ENDDO
373     ENDDO
374 dimitri 1.9
375     C-- Read initial snow thickness from file if available.
376     IF ( HsnowFile .NE. ' ' ) THEN
377 mlosch 1.43 CALL READ_FLD_XY_RL( HsnowFile, ' ', HSNOW, 0, myThid )
378     _EXCH_XY_RL(HSNOW,myThid)
379 dimitri 1.9 DO bj=myByLo(myThid),myByHi(myThid)
380     DO bi=myBxLo(myThid),myBxHi(myThid)
381     DO j=1-OLy,sNy+OLy
382     DO i=1-OLx,sNx+OLx
383 mlosch 1.43 HSNOW(i,j,bi,bj) = MAX(HSNOW(i,j,bi,bj),ZERO)
384 dimitri 1.9 ENDDO
385     ENDDO
386     ENDDO
387     ENDDO
388     ENDIF
389 dimitri 1.18
390     #ifdef SEAICE_SALINITY
391 dimitri 1.20 DO bj=myByLo(myThid),myByHi(myThid)
392     DO bi=myBxLo(myThid),myBxHi(myThid)
393 jmc 1.23 DO j=1-OLy,sNy+OLy
394     DO i=1-OLx,sNx+OLx
395 mlosch 1.44 HSALT(i,j,bi,bj)=HEFF(i,j,bi,bj)*salt(i,j,kSurface,bi,bj)*
396 dimitri 1.20 & ICE2WATR*rhoConstFresh*SEAICE_salinity
397     ENDDO
398     ENDDO
399     ENDDO
400     ENDDO
401    
402 dimitri 1.18 C-- Read initial sea ice salinity from file if available.
403     IF ( HsaltFile .NE. ' ' ) THEN
404 mlosch 1.43 CALL READ_FLD_XY_RL( HsaltFile, ' ', HSALT, 0, myThid )
405     _EXCH_XY_RL(HSALT,myThid)
406 dimitri 1.18 ENDIF
407     #endif /* SEAICE_SALINITY */
408 jmc 1.23
409 dimitri 1.33 #ifdef SEAICE_AGE
410     C-- Read initial sea ice age from file if available.
411     IF ( IceAgeFile .NE. ' ' ) THEN
412 mlosch 1.43 CALL READ_FLD_XY_RL( IceAgeFile, ' ', ICEAGE, 0, myThid )
413     _EXCH_XY_RL(ICEAGE,myThid)
414 dimitri 1.33 ENDIF
415     #endif /* SEAICE_AGE */
416    
417 heimbach 1.1 ENDIF
418    
419 jmc 1.52 C-- In case we use scheme with a large stencil that extends into overlap:
420     #ifdef ALLOW_OBCS
421     IF ( useOBCS ) THEN
422     DO bj=myByLo(myThid),myByHi(myThid)
423     DO bi=myBxLo(myThid),myBxHi(myThid)
424     CALL OBCS_COPY_TRACER( HEFF(1-Olx,1-Oly,bi,bj),
425     I 1, bi, bj, myThid )
426     CALL OBCS_COPY_TRACER( AREA(1-Olx,1-Oly,bi,bj),
427     I 1, bi, bj, myThid )
428     CALL OBCS_COPY_TRACER( HSNOW(1-Olx,1-Oly,bi,bj),
429     I 1, bi, bj, myThid )
430     #ifdef SEAICE_SALINITY
431     CALL OBCS_COPY_TRACER( HSALT(1-Olx,1-Oly,bi,bj),
432     I 1, bi, bj, myThid )
433     #endif
434     #ifdef SEAICE_AGE
435     CALL OBCS_COPY_TRACER( ICEAGE(1-Olx,1-Oly,bi,bj),
436     I 1, bi, bj, myThid )
437     #endif
438     ENDDO
439     ENDDO
440     ENDIF
441     #endif /* ALLOW_OBCS */
442    
443 heimbach 1.1 C--- Complete initialization
444 mlosch 1.8 PSTAR = SEAICE_strength
445 heimbach 1.1 DO bj=myByLo(myThid),myByHi(myThid)
446     DO bi=myBxLo(myThid),myBxHi(myThid)
447     DO j=1-OLy,sNy+OLy
448     DO i=1-OLx,sNx+OLx
449 mlosch 1.48 ZETA(i,j,bi,bj) = HEFF(i,j,bi,bj)*(1.0 _d 11)
450     ETA(i,j,bi,bj) = ZETA(i,j,bi,bj)/SEAICE_eccen**2
451     PRESS0(i,j,bi,bj) = PSTAR*HEFF(i,j,bi,bj)
452 mlosch 1.43 & *EXP(-20.0 _d 0*(ONE-AREA(i,j,bi,bj)))
453 mlosch 1.48 ZMAX(I,J,bi,bj) = SEAICE_zetaMaxFac*PRESS0(I,J,bi,bj)
454     ZMIN(i,j,bi,bj) = SEAICE_zetaMin
455     PRESS0(i,j,bi,bj) = PRESS0(i,j,bi,bj)*HEFFM(i,j,bi,bj)
456 heimbach 1.1 ENDDO
457     ENDDO
458     IF ( useRealFreshWaterFlux .AND. .NOT.useThSIce ) THEN
459     DO j=1-OLy,sNy+OLy
460     DO i=1-OLx,sNx+OLx
461 mlosch 1.43 sIceLoad(i,j,bi,bj) = HEFF(i,j,bi,bj)*SEAICE_rhoIce
462 mlosch 1.36 & + HSNOW(i,j,bi,bj)*SEAICE_rhoSnow
463 jmc 1.23
464 heimbach 1.1 ENDDO
465     ENDDO
466     ENDIF
467     ENDDO
468     ENDDO
469    
470     RETURN
471     END

  ViewVC Help
Powered by ViewVC 1.1.22