/[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.82 - (hide annotations) (download)
Wed Apr 23 12:41:09 2014 UTC (10 years, 1 month ago) by mlosch
Branch: MAIN
Changes since 1.81: +2 -7 lines
call seaice_itd_redist once in the initialisation phase to have the initial
conditions in the correct category

1 mlosch 1.82 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_varia.F,v 1.81 2014/04/09 16:24:51 mlosch Exp $
2 heimbach 1.1 C $Name: $
3    
4     #include "SEAICE_OPTIONS.h"
5 jmc 1.63 #ifdef ALLOW_OBCS
6     # include "OBCS_OPTIONS.h"
7     #endif
8 heimbach 1.1
9     CStartOfInterface
10     SUBROUTINE SEAICE_INIT_VARIA( myThid )
11 jmc 1.47 C *==========================================================*
12 heimbach 1.1 C | SUBROUTINE SEAICE_INIT_VARIA |
13     C | o Initialization of sea ice model. |
14 jmc 1.47 C *==========================================================*
15     C *==========================================================*
16 heimbach 1.1 IMPLICIT NONE
17 jmc 1.23
18 heimbach 1.1 C === Global variables ===
19     #include "SIZE.h"
20     #include "EEPARAMS.h"
21     #include "PARAMS.h"
22     #include "GRID.h"
23 dimitri 1.21 #include "DYNVARS.h"
24 jmc 1.50 #include "FFIELDS.h"
25 heimbach 1.54 #include "SEAICE_SIZE.h"
26 jmc 1.50 #include "SEAICE_PARAMS.h"
27 heimbach 1.1 #include "SEAICE.h"
28 heimbach 1.54 #include "SEAICE_TRACER.h"
29 jmc 1.50 #include "SEAICE_TAVE.h"
30 jmc 1.63 #ifdef OBCS_UVICE_OLD
31 jmc 1.57 # include "OBCS_GRID.h"
32 dimitri 1.24 #endif
33 heimbach 1.1
34     C === Routine arguments ===
35 jmc 1.62 C myThid :: Thread no. that called this routine.
36 heimbach 1.1 INTEGER myThid
37     CEndOfInterface
38 jmc 1.23
39 heimbach 1.1 C === Local variables ===
40 jmc 1.62 C i,j,k,bi,bj :: Loop counters
41 heimbach 1.1
42 jmc 1.55 INTEGER i, j, bi, bj
43 mlosch 1.8 _RL PSTAR
44 jmc 1.59 INTEGER kSurface
45 jmc 1.55 #ifdef SEAICE_CGRID
46 heimbach 1.1 _RS mask_uice
47 jmc 1.55 #endif
48 jmc 1.50 INTEGER k
49 gforget 1.60 #ifdef ALLOW_SITRACER
50     INTEGER iTr, jTh
51     #endif
52 jmc 1.63 #ifdef OBCS_UVICE_OLD
53 dimitri 1.53 INTEGER I_obc, J_obc
54     #endif /* ALLOW_OBCS */
55 heimbach 1.1
56 jmc 1.23 IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
57 mlosch 1.44 kSurface = Nr
58 jmc 1.23 ELSE
59 mlosch 1.44 kSurface = 1
60 jmc 1.23 ENDIF
61 dimitri 1.11
62 jmc 1.23 C-- Initialise all variables in common blocks:
63 heimbach 1.1 DO bj=myByLo(myThid),myByHi(myThid)
64     DO bi=myBxLo(myThid),myBxHi(myThid)
65 mlosch 1.42 DO j=1-OLy,sNy+OLy
66     DO i=1-OLx,sNx+OLx
67 mlosch 1.43 HEFF(i,j,bi,bj)=0. _d 0
68     AREA(i,j,bi,bj)=0. _d 0
69 heimbach 1.75 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 mlosch 1.43 UICE(i,j,bi,bj)=0. _d 0
78     VICE(i,j,bi,bj)=0. _d 0
79 heimbach 1.65 #ifdef SEAICE_ALLOW_FREEDRIFT
80 heimbach 1.64 uice_fd(i,j,bi,bj)=0. _d 0
81     vice_fd(i,j,bi,bj)=0. _d 0
82 heimbach 1.65 #endif
83 mlosch 1.43 C
84 mlosch 1.42 uIceNm1(i,j,bi,bj)=0. _d 0
85     vIceNm1(i,j,bi,bj)=0. _d 0
86 jmc 1.23 ETA (i,j,bi,bj) = 0. _d 0
87 mlosch 1.77 etaZ(i,j,bi,bj) = 0. _d 0
88 jmc 1.23 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 jmc 1.68 uIceC(i,j,bi,bj) = 0. _d 0
92     vIceC(i,j,bi,bj) = 0. _d 0
93 jmc 1.23 #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 heimbach 1.39 # ifdef SEAICE_ALLOW_EVP
100 jmc 1.23 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 dimitri 1.25 GWATX(i,j,bi,bj) = 0. _d 0
110     GWATY(i,j,bi,bj) = 0. _d 0
111 jmc 1.23 #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 heimbach 1.75 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 ifenty 1.56 #ifdef SEAICE_VARIABLE_SALINITY
127 jmc 1.23 HSALT(i,j,bi,bj) = 0. _d 0
128     #endif
129 gforget 1.60 #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 gforget 1.61 DO jTh = 1, 3
140     SItrAREA (i,j,bi,bj,jTh) = 0. _d 0
141     ENDDO
142 gforget 1.60 #endif
143 jmc 1.23 TICE(i,j,bi,bj) = 0. _d 0
144 gforget 1.70 DO k=1,MULTDIM
145     TICES(i,j,k,bi,bj)=0. _d 0
146     ENDDO
147 jmc 1.23 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 heimbach 1.76 icevolMeanCell(i,j,bi,bj) = 0. _d 0
153 dimitri 1.18 #endif
154 heimbach 1.32 saltWtrIce(i,j,bi,bj) = 0. _d 0
155     frWtrIce(i,j,bi,bj) = 0. _d 0
156 heimbach 1.72 #if (defined (ALLOW_MEAN_SFLUX_COST_CONTRIBUTION) || defined (ALLOW_SSH_GLOBMEAN_COST_CONTRIBUTION))
157 heimbach 1.49 frWtrAtm(i,j,bi,bj) = 0. _d 0
158 heimbach 1.72 AREAforAtmFW(i,j,bi,bj)=0. _d 0
159 heimbach 1.49 #endif
160 heimbach 1.1 ENDDO
161     ENDDO
162     ENDDO
163     ENDDO
164    
165 jmc 1.50 #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 mlosch 1.40 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 mlosch 1.66 #ifdef SEAICE_CGRID
187 heimbach 1.1 DO bj=myByLo(myThid),myByHi(myThid)
188     DO bi=myBxLo(myThid),myBxHi(myThid)
189 jmc 1.23 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 heimbach 1.45 IF(mask_uice.GT.1.5 _d 0) seaiceMaskU(i,j,bi,bj)=1.0 _d 0
195 jmc 1.23 mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i ,j-1,bi,bj)
196 heimbach 1.45 IF(mask_uice.GT.1.5 _d 0) seaiceMaskV(i,j,bi,bj)=1.0 _d 0
197 dimitri 1.22 ENDDO
198     ENDDO
199 mlosch 1.40 ENDDO
200     ENDDO
201 mlosch 1.66 #endif /* SEAICE_CGRID */
202 dimitri 1.11
203 mlosch 1.40 DO bj=myByLo(myThid),myByHi(myThid)
204     DO bi=myBxLo(myThid),myBxHi(myThid)
205 jmc 1.63 #ifdef OBCS_UVICE_OLD
206     #ifdef SEAICE_CGRID
207 dimitri 1.16 IF (useOBCS) THEN
208 dimitri 1.22 C-- If OBCS is turned on, close southern and western boundaries
209 jmc 1.67 DO i=1-OLx,sNx+OLx
210 dimitri 1.11 C Southern boundary
211 dimitri 1.53 J_obc = OB_Js(i,bi,bj)
212 jmc 1.73 IF ( J_obc.NE.OB_indexNone ) THEN
213 dimitri 1.53 seaiceMaskU(i,J_obc,bi,bj)= 0.0 _d 0
214     seaiceMaskV(i,J_obc,bi,bj)= 0.0 _d 0
215 dimitri 1.16 ENDIF
216     ENDDO
217 jmc 1.67 DO j=1-OLy,sNy+OLy
218 dimitri 1.11 C Western boundary
219 jmc 1.73 I_obc = OB_Iw(j,bi,bj)
220     IF ( I_obc.NE.OB_indexNone ) THEN
221 dimitri 1.53 seaiceMaskU(I_obc,j,bi,bj)= 0.0 _d 0
222     seaiceMaskV(I_obc,j,bi,bj)= 0.0 _d 0
223 dimitri 1.16 ENDIF
224     ENDDO
225     ENDIF
226 dimitri 1.11 #endif /* SEAICE_CGRID */
227 jmc 1.63 #endif /* OBCS_UVICE_OLD */
228 heimbach 1.1
229     DO j=1-OLy,sNy+OLy
230     DO i=1-OLx,sNx+OLx
231 jmc 1.23 TICE(i,j,bi,bj)=273.0 _d 0
232 heimbach 1.1 DO k=1,MULTDIM
233 jmc 1.23 TICES(i,j,k,bi,bj)=273.0 _d 0
234 heimbach 1.1 ENDDO
235     #ifndef SEAICE_CGRID
236 jmc 1.23 AMASS (i,j,bi,bj)=1000.0 _d 0
237 heimbach 1.1 #else
238 jmc 1.23 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 heimbach 1.1 #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 jmc 1.67 _EXCH_XY_RS(UVM, myThid)
253 heimbach 1.1 #endif
254    
255     C-- Now lets look at all these beasts
256 jmc 1.59 IF ( debugLevel .GE. debLevC ) THEN
257 heimbach 1.1 CALL PLOT_FIELD_XYRL( HEFFM , 'Current HEFFM ' ,
258 jmc 1.59 & nIter0, myThid )
259 heimbach 1.1 #ifdef SEAICE_CGRID
260     CALL PLOT_FIELD_XYRL( seaiceMaskU, 'Current seaiceMaskU',
261 jmc 1.59 & nIter0, myThid )
262 heimbach 1.1 CALL PLOT_FIELD_XYRL( seaiceMaskV, 'Current seaiceMaskV',
263 jmc 1.59 & nIter0, myThid )
264 heimbach 1.1 #else
265 jmc 1.67 CALL PLOT_FIELD_XYRS( UVM , 'Current UVM ' ,
266 jmc 1.59 & nIter0, myThid )
267 heimbach 1.1 #endif
268     ENDIF
269    
270     C-- Set model variables to initial/restart conditions
271 mlosch 1.26 IF ( .NOT. ( startTime .EQ. baseTime .AND. nIter0 .EQ. 0
272     & .AND. pickupSuff .EQ. ' ') ) THEN
273 heimbach 1.1
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 mlosch 1.43 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 heimbach 1.1 ENDDO
286     ENDDO
287     ENDDO
288     ENDDO
289    
290 jmc 1.62 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 heimbach 1.1 C-- Read initial sea-ice thickness from file if available.
312     IF ( HeffFile .NE. ' ' ) THEN
313 mlosch 1.43 CALL READ_FLD_XY_RL( HeffFile, ' ', HEFF, 0, myThid )
314     _EXCH_XY_RL(HEFF,myThid)
315 heimbach 1.1 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 mlosch 1.43 HEFF(i,j,bi,bj) = MAX(HEFF(i,j,bi,bj),ZERO)
320 heimbach 1.1 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 mlosch 1.43 IF(HEFF(i,j,bi,bj).GT.ZERO) AREA(i,j,bi,bj)=ONE
331 heimbach 1.1 ENDDO
332     ENDDO
333     ENDDO
334     ENDDO
335 jmc 1.23
336 dimitri 1.27 C-- Read initial sea-ice area from file if available.
337 mlosch 1.7 IF ( AreaFile .NE. ' ' ) THEN
338 mlosch 1.43 CALL READ_FLD_XY_RL( AreaFile, ' ', AREA, 0, myThid )
339     _EXCH_XY_RL(AREA,myThid)
340 mlosch 1.7 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 mlosch 1.43 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 mlosch 1.7 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 jmc 1.47 HSNOW(i,j,bi,bj) = 0.2 _d 0 * AREA(i,j,bi,bj)
359 mlosch 1.7 ENDDO
360     ENDDO
361     ENDDO
362     ENDDO
363 dimitri 1.9
364     C-- Read initial snow thickness from file if available.
365     IF ( HsnowFile .NE. ' ' ) THEN
366 mlosch 1.43 CALL READ_FLD_XY_RL( HsnowFile, ' ', HSNOW, 0, myThid )
367     _EXCH_XY_RL(HSNOW,myThid)
368 dimitri 1.9 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 mlosch 1.43 HSNOW(i,j,bi,bj) = MAX(HSNOW(i,j,bi,bj),ZERO)
373 dimitri 1.9 ENDDO
374     ENDDO
375     ENDDO
376     ENDDO
377     ENDIF
378 dimitri 1.18
379 heimbach 1.75 #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 mlosch 1.81 opnWtrFrac(I,J,bi,bj)=1. _d 0 - AREA(I,J,bi,bj)
388 heimbach 1.75 ENDDO
389     ENDDO
390 mlosch 1.82 CALL SEAICE_ITD_REDIST(bi, bj, baseTime, nIter0, myThid)
391 heimbach 1.75 ENDDO
392     ENDDO
393     #endif
394    
395 ifenty 1.56 #ifdef SEAICE_VARIABLE_SALINITY
396 dimitri 1.20 DO bj=myByLo(myThid),myByHi(myThid)
397     DO bi=myBxLo(myThid),myBxHi(myThid)
398 jmc 1.23 DO j=1-OLy,sNy+OLy
399     DO i=1-OLx,sNx+OLx
400 mlosch 1.44 HSALT(i,j,bi,bj)=HEFF(i,j,bi,bj)*salt(i,j,kSurface,bi,bj)*
401 jmc 1.71 & SEAICE_rhoIce*SEAICE_saltFrac
402     cif & ICE2WATR*rhoConstFresh*SEAICE_saltFrac
403 ifenty 1.56
404 dimitri 1.20 ENDDO
405     ENDDO
406     ENDDO
407     ENDDO
408    
409 dimitri 1.18 C-- Read initial sea ice salinity from file if available.
410     IF ( HsaltFile .NE. ' ' ) THEN
411 mlosch 1.43 CALL READ_FLD_XY_RL( HsaltFile, ' ', HSALT, 0, myThid )
412     _EXCH_XY_RL(HSALT,myThid)
413 dimitri 1.18 ENDIF
414 ifenty 1.56 #endif /* SEAICE_VARIABLE_SALINITY */
415 jmc 1.23
416 gforget 1.69 #ifdef ALLOW_SITRACER
417 dimitri 1.33 C-- Read initial sea ice age from file if available.
418 gforget 1.69 DO iTr = 1, SItrMaxNum
419     IF ( SItrFile(iTr) .NE. ' ' ) THEN
420     CALL READ_FLD_XY_RL( siTrFile(iTr), ' ',
421     & SItracer(1-OLx,1-OLy,1,1,iTr), 0, myThid )
422     _EXCH_XY_RL(SItracer(1-OLx,1-OLy,1,1,iTr),myThid)
423 heimbach 1.54 ENDIF
424     ENDDO
425 gforget 1.69 #endif /* ALLOW_SITRACER */
426 dimitri 1.33
427 heimbach 1.1 ENDIF
428    
429 heimbach 1.72 #if (defined (ALLOW_MEAN_SFLUX_COST_CONTRIBUTION) || defined (ALLOW_SSH_GLOBMEAN_COST_CONTRIBUTION))
430     DO bj=myByLo(myThid),myByHi(myThid)
431     DO bi=myBxLo(myThid),myBxHi(myThid)
432     DO j=1-OLy,sNy+OLy
433     DO i=1-OLx,sNx+OLx
434     AREAforAtmFW(i,j,bi,bj) = AREA(i,j,bi,bj)
435     ENDDO
436     ENDDO
437     ENDDO
438     ENDDO
439     #endif
440    
441 jmc 1.58 #ifdef ALLOW_OBCS
442 jmc 1.52 C-- In case we use scheme with a large stencil that extends into overlap:
443 jmc 1.58 C no longer needed with the right masking in advection & diffusion S/R.
444     c IF ( useOBCS ) THEN
445     c DO bj=myByLo(myThid),myByHi(myThid)
446     c DO bi=myBxLo(myThid),myBxHi(myThid)
447 jmc 1.67 c CALL OBCS_COPY_TRACER( HEFF(1-OLx,1-OLy,bi,bj),
448 jmc 1.58 c I 1, bi, bj, myThid )
449 jmc 1.67 c CALL OBCS_COPY_TRACER( AREA(1-OLx,1-OLy,bi,bj),
450 jmc 1.58 c I 1, bi, bj, myThid )
451 jmc 1.67 c CALL OBCS_COPY_TRACER( HSNOW(1-OLx,1-OLy,bi,bj),
452 jmc 1.58 c I 1, bi, bj, myThid )
453 ifenty 1.56 #ifdef SEAICE_VARIABLE_SALINITY
454 jmc 1.67 c CALL OBCS_COPY_TRACER( HSALT(1-OLx,1-OLy,bi,bj),
455 jmc 1.58 c I 1, bi, bj, myThid )
456 jmc 1.52 #endif
457 jmc 1.58 c ENDDO
458     c ENDDO
459     c ENDIF
460 jmc 1.52 #endif /* ALLOW_OBCS */
461    
462 mlosch 1.78 #ifdef SEAICE_ALLOW_JFNK
463 jmc 1.80 C Computing this metric cannot be done in S/R SEAICE_INIT_FIXED
464     C where it belongs, because globalArea is only defined later after
465 mlosch 1.78 C S/R PACKAGES_INIT_FIXED, so we move this computation here.
466 jmc 1.80 CALL SEAICE_MAP_RS2VEC( nVec, rAw, rAs,
467 mlosch 1.78 & scalarProductMetric, .TRUE., myThid )
468     DO bj=myByLo(myThid),myByHi(myThid)
469     DO bi=myBxLo(myThid),myBxHi(myThid)
470     DO i=1,nVec
471     scalarProductMetric(i,1,bi,bj) =
472     & scalarProductMetric(i,1,bi,bj)/globalArea
473     ENDDO
474     ENDDO
475     ENDDO
476     #endif /* SEAICE_ALLOW_JFNK */
477    
478 heimbach 1.1 C--- Complete initialization
479 mlosch 1.8 PSTAR = SEAICE_strength
480 heimbach 1.1 DO bj=myByLo(myThid),myByHi(myThid)
481     DO bi=myBxLo(myThid),myBxHi(myThid)
482     DO j=1-OLy,sNy+OLy
483     DO i=1-OLx,sNx+OLx
484 mlosch 1.48 ZETA(i,j,bi,bj) = HEFF(i,j,bi,bj)*(1.0 _d 11)
485     ETA(i,j,bi,bj) = ZETA(i,j,bi,bj)/SEAICE_eccen**2
486     PRESS0(i,j,bi,bj) = PSTAR*HEFF(i,j,bi,bj)
487 mlosch 1.43 & *EXP(-20.0 _d 0*(ONE-AREA(i,j,bi,bj)))
488 mlosch 1.48 ZMAX(I,J,bi,bj) = SEAICE_zetaMaxFac*PRESS0(I,J,bi,bj)
489     ZMIN(i,j,bi,bj) = SEAICE_zetaMin
490     PRESS0(i,j,bi,bj) = PRESS0(i,j,bi,bj)*HEFFM(i,j,bi,bj)
491 heimbach 1.1 ENDDO
492     ENDDO
493     IF ( useRealFreshWaterFlux .AND. .NOT.useThSIce ) THEN
494     DO j=1-OLy,sNy+OLy
495     DO i=1-OLx,sNx+OLx
496 mlosch 1.43 sIceLoad(i,j,bi,bj) = HEFF(i,j,bi,bj)*SEAICE_rhoIce
497 mlosch 1.36 & + HSNOW(i,j,bi,bj)*SEAICE_rhoSnow
498 jmc 1.23
499 heimbach 1.1 ENDDO
500     ENDDO
501     ENDIF
502     ENDDO
503     ENDDO
504    
505     RETURN
506     END

  ViewVC Help
Powered by ViewVC 1.1.22