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

Contents of /MITgcm/pkg/seaice/seaice_init.F

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


Revision 1.28 - (show annotations) (download)
Mon Dec 27 20:34:11 2004 UTC (19 years, 6 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint57t_post, checkpoint57m_post, checkpoint57v_post, checkpoint57f_post, checkpoint57s_post, checkpoint57j_post, checkpoint57f_pre, checkpoint57g_post, checkpoint57h_pre, checkpoint57y_post, checkpoint57x_post, checkpoint57g_pre, checkpoint57e_post, checkpoint57h_post, checkpoint57y_pre, checkpoint57c_pre, checkpoint57o_post, checkpoint57r_post, checkpoint57k_post, checkpoint57d_post, checkpoint57i_post, checkpoint58, checkpoint57h_done, checkpoint57n_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint58a_post, checkpoint57q_post, checkpoint57z_post, eckpoint57e_pre, checkpoint57c_post, checkpoint57l_post
Changes since 1.27: +1 -4 lines
o added seaice_summary.F and removed obsolete ALLOW_SEAICE's from pkg/seaice

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init.F,v 1.27 2004/12/22 00:49:36 dimitri Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5
6 CStartOfInterface
7 SUBROUTINE SEAICE_INIT( myThid )
8 C /==========================================================\
9 C | SUBROUTINE SEAICE_INIT |
10 C | o Initialization of sea ice model. |
11 C |==========================================================|
12 C \==========================================================/
13 IMPLICIT NONE
14
15 C === Global variables ===
16 #include "SIZE.h"
17 #include "EEPARAMS.h"
18 #include "PARAMS.h"
19 #include "GRID.h"
20 #include "SEAICE.h"
21 #include "SEAICE_GRID.h"
22 #include "SEAICE_DIAGS.h"
23 #include "SEAICE_PARAMS.h"
24 #ifdef ALLOW_EXCH2
25 #include "W2_EXCH2_TOPOLOGY.h"
26 #include "W2_EXCH2_PARAMS.h"
27 #endif /* ALLOW_EXCH2 */
28
29 C === Routine arguments ===
30 C myThid - Thread no. that called this routine.
31 INTEGER myThid
32 CEndOfInterface
33
34 C === Local variables ===
35 C i,j,k,bi,bj - Loop counters
36
37 INTEGER i, j, k, bi, bj
38 _RS mask_uice
39 INTEGER myIter, myTile
40
41 #ifdef ALLOW_TIMEAVE
42 C Initialize averages to zero
43 DO bj = myByLo(myThid), myByHi(myThid)
44 DO bi = myBxLo(myThid), myBxHi(myThid)
45 CALL TIMEAVE_RESET(FUtave ,1,bi,bj,myThid)
46 CALL TIMEAVE_RESET(FVtave ,1,bi,bj,myThid)
47 CALL TIMEAVE_RESET(EmPmRtave,1,bi,bj,myThid)
48 CALL TIMEAVE_RESET(QNETtave ,1,bi,bj,myThid)
49 CALL TIMEAVE_RESET(QSWtave ,1,bi,bj,myThid)
50 CALL TIMEAVE_RESET(UICEtave ,1,bi,bj,myThid)
51 CALL TIMEAVE_RESET(VICEtave ,1,bi,bj,myThid)
52 CALL TIMEAVE_RESET(HEFFtave ,1,bi,bj,myThid)
53 CALL TIMEAVE_RESET(AREAtave ,1,bi,bj,myThid)
54 DO k=1,Nr
55 SEAICE_TimeAve(k,bi,bj)=ZERO
56 ENDDO
57 ENDDO
58 ENDDO
59 #endif /* ALLOW_TIMEAVE */
60
61 cph(
62 cph make sure TAF sees proper initialisation
63 cph to avoid partial recomputation issues
64 DO bj=myByLo(myThid),myByHi(myThid)
65 DO bi=myBxLo(myThid),myBxHi(myThid)
66 c
67 DO K=1,3
68 DO J=1-OLy,sNy+OLy
69 DO I=1-OLx,sNx+OLx
70 HEFF(I,J,k,bi,bj)=ZERO
71 AREA(I,J,k,bi,bj)=ZERO
72 UICE(I,J,k,bi,bj)=ZERO
73 VICE(I,J,k,bi,bj)=ZERO
74 ENDDO
75 ENDDO
76 ENDDO
77 c
78 DO J=1-OLy,sNy+OLy
79 DO I=1-OLx,sNx+OLx
80 HSNOW(I,J,bi,bj)=ZERO
81 ZETA(I,J,bi,bj)=ZERO
82 ENDDO
83 ENDDO
84 c
85 ENDDO
86 ENDDO
87 cph)
88
89 C-- Initialize grid info
90 DO bj=myByLo(myThid),myByHi(myThid)
91 DO bi=myBxLo(myThid),myBxHi(myThid)
92 DO J=1-OLy,sNy+OLy
93 DO I=1-OLx,sNx+OLx
94 CSTICE(i,j,bi,bj) =cos(yC(I,J,bi,bj)*deg2rad)
95 CSUICE(i,j,bi,bj) =cos(yG(I,J,bi,bj)*deg2rad)
96 CML(
97 C inverses of CSTICE and CSUICE. Lets hope we are never
98 C at the poles
99 IF ( CSTICE(I,J,bi,bj) .ne. 0. _d 0 ) THEN
100 RECIP_CSTICE(I,J,bi,bj) = 1./CSTICE(I,J,bi,bj)
101 ELSE
102 RECIP_CSTICE(I,J,bi,bj) =0. _d 0
103 ENDIF
104 IF ( CSUICE(I,J,bi,bj) .ne. 0. _d 0 ) THEN
105 RECIP_CSUICE(I,J,bi,bj) = 1./CSUICE(I,J,bi,bj)
106 ELSE
107 RECIP_CSUICE(I,J,bi,bj) =0. _d 0
108 ENDIF
109 CML)
110 SINEICE(i,j,bi,bj)=sin(yC(I,J,bi,bj)*deg2rad)
111 TNGTICE(i,j,bi,bj)=SINEICE(i,j,bi,bj)*RECIP_CSTICE(i,j,bi,bj)
112 SINEICE(i,j,bi,bj)=sin(yG(I,J,bi,bj)*deg2rad)
113 TNGICE(i,j,bi,bj) =SINEICE(i,j,bi,bj)*RECIP_CSUICE(i,j,bi,bj)
114 DXTICE(i,j,bi,bj)=dxF(i,j,bi,bj)*RECIP_CSTICE(i,j,bi,bj)
115 DXUICE(i,j,bi,bj)=dxV(i,j,bi,bj)*RECIP_CSUICE(i,j,bi,bj)
116 DYTICE(i,j,bi,bj)=dyF(i,j,bi,bj)
117 DYUICE(i,j,bi,bj)=dyU(i,j,bi,bj)
118 ENDDO
119 ENDDO
120 DO j=1-OLy,sNy+OLy
121 DO i=1-OLx,sNx+OLx
122 HEFFM(i,j,bi,bj)=ONE
123 IF (_hFacC(i,j,1,bi,bj).eq.0.) HEFFM(i,j,bi,bj)=ZERO
124 ENDDO
125 ENDDO
126 DO J=2-OLy,sNy+OLy
127 DO I=2-OLx,sNx+OLx
128 UVM(i,j,bi,bj)=ZERO
129 mask_uice=HEFFM(I,J, bi,bj)+HEFFM(I-1,J-1,bi,bj)
130 & +HEFFM(I,J-1,bi,bj)+HEFFM(I-1,J, bi,bj)
131 IF(mask_uice.GT.3.5) UVM(I,J,bi,bj)=ONE
132 ENDDO
133 ENDDO
134
135 #ifdef ALLOW_EXCH2
136 C-- Special stuff for cubed sphere: assume grid is rectangular and
137 C set UV mask to zero except for Arctic and Antarctic cube faces.
138 IF (useCubedSphereExchange) THEN
139 DO J=1-OLy,sNy+OLy
140 DO I=1-OLx,sNx+OLx
141 CSTICE(i,j,bi,bj) = ONE
142 CSUICE(i,j,bi,bj) = ONE
143 CML(
144 RECIP_CSTICE(I,J,bi,bj) = ONE
145 RECIP_CSUICE(I,J,bi,bj) = ONE
146 CML)
147 TNGTICE(i,j,bi,bj)= ZERO
148 TNGICE(i,j,bi,bj) = ZERO
149 DXTICE(i,j,bi,bj) = dxF(i,j,bi,bj)
150 DXUICE(i,j,bi,bj) = dxV(i,j,bi,bj)
151 ENDDO
152 ENDDO
153 myTile = W2_myTileList(bi)
154 IF ( exch2_myFace(myTile) .EQ. 1 .OR.
155 & exch2_myFace(myTile) .EQ. 2 .OR.
156 & exch2_myFace(myTile) .EQ. 4 .OR.
157 & exch2_myFace(myTile) .EQ. 5 ) THEN
158 DO J=1-OLy,sNy+OLy
159 DO I=1-OLx,sNx+OLx
160 UVM(i,j,bi,bj)=ZERO
161 ENDDO
162 ENDDO
163 ELSEIF ( exch2_isWedge(myTile) .EQ. 1 ) THEN
164 I=1
165 DO J=1-OLy,sNy+OLy
166 UVM(i,j,bi,bj)=ZERO
167 ENDDO
168 ELSEIF ( exch2_isSedge(myTile) .EQ. 1 ) THEN
169 J=1
170 DO I=1-OLx,sNx+OLx
171 UVM(i,j,bi,bj)=ZERO
172 ENDDO
173 ENDIF
174 C-- Make sure that DXTICE and DYTICE do not contain any zeros
175 DO J=1,sNy
176 DO I=1-OLx,0
177 IF(DXTICE(i,j,bi,bj).EQ.0)
178 & DXTICE(i,j,bi,bj)=DXTICE(1,j,bi,bj)
179 IF(DYTICE(i,j,bi,bj).EQ.0)
180 & DYTICE(i,j,bi,bj)=DYTICE(1,j,bi,bj)
181 ENDDO
182 DO I=sNx+1,sNx+OLx
183 IF(DXTICE(i,j,bi,bj).EQ.0)
184 & DXTICE(i,j,bi,bj)=DXTICE(sNx,j,bi,bj)
185 IF(DYTICE(i,j,bi,bj).EQ.0)
186 & DYTICE(i,j,bi,bj)=DYTICE(sNx,j,bi,bj)
187 ENDDO
188 ENDDO
189 DO J=1-OLy,0
190 DO I=1-OLx,sNx+OLx
191 IF(DXTICE(i,j,bi,bj).EQ.0)
192 & DXTICE(i,j,bi,bj)=DXTICE(i,1,bi,bj)
193 IF(DYTICE(i,j,bi,bj).EQ.0)
194 & DYTICE(i,j,bi,bj)=DYTICE(i,1,bi,bj)
195 ENDDO
196 ENDDO
197 DO J=sNy+1,sNy+OLy
198 DO I=1-OLx,sNx+OLx
199 IF(DXTICE(i,j,bi,bj).EQ.0)
200 & DXTICE(i,j,bi,bj)=DXTICE(i,sNy,bi,bj)
201 IF(DYTICE(i,j,bi,bj).EQ.0)
202 & DYTICE(i,j,bi,bj)=DYTICE(i,sNy,bi,bj)
203 ENDDO
204 ENDDO
205 ENDIF
206 #endif /* ALLOW_EXCH2 */
207
208 DO j=1-OLy,sNy+OLy
209 DO i=1-OLx,sNx+OLx
210 TICE(I,J,bi,bj)=273.0 _d 0
211 #ifdef SEAICE_MULTILEVEL
212 DO k=1,MULTDIM
213 TICES(I,J,k,bi,bj)=273.0 _d 0
214 ENDDO
215 #endif
216 UICEC(I,J,bi,bj)=ZERO
217 VICEC(I,J,bi,bj)=ZERO
218 AMASS(I,J,bi,bj)=1000.0 _d 0
219 ENDDO
220 ENDDO
221
222 C-- Choose a proxy level for geostrophic velocity,
223 DO j=1-OLy,sNy+OLy
224 DO i=1-OLx,sNx+OLx
225 #ifdef SEAICE_TEST_ICE_STRESS_1
226 KGEO(I,J,bi,bj) = 1
227 #else /* SEAICE_TEST_ICE_STRESS_1 */
228 IF (klowc(i,j,bi,bj) .LT. 2) THEN
229 KGEO(I,J,bi,bj) = 1
230 ELSE
231 KGEO(I,J,bi,bj) = 2
232 DO WHILE ( abs(rC(KGEO(I,J,bi,bj))) .LT. 50.0 .AND.
233 & KGEO(I,J,bi,bj) .LT. (klowc(i,j,bi,bj)-1) )
234 KGEO(I,J,bi,bj) = KGEO(I,J,bi,bj) + 1
235 ENDDO
236 ENDIF
237 #endif /* SEAICE_TEST_ICE_STRESS_1 */
238 ENDDO
239 ENDDO
240
241 ENDDO
242 ENDDO
243
244 C-- Update overlap regions
245 _EXCH_XY_R8(UVM, myThid)
246
247 C-- Now lets look at all these beasts
248 IF ( debugLevel .GE. debLevB ) THEN
249 myIter=0
250 CALL PLOT_FIELD_XYRL( CSTICE , 'Current CSTICE ' ,
251 & myIter, myThid )
252 CALL PLOT_FIELD_XYRL( CSUICE , 'Current CSUICE ' ,
253 & myIter, myThid )
254 CALL PLOT_FIELD_XYRL( TNGTICE , 'Current TNGTICE ' ,
255 & myIter, myThid )
256 CALL PLOT_FIELD_XYRL( TNGICE , 'Current TNGICE ' ,
257 & myIter, myThid )
258 CALL PLOT_FIELD_XYRL( SINEICE , 'Current SINEICE ' ,
259 & myIter, myThid )
260 CALL PLOT_FIELD_XYRL( DXTICE , 'Current DXTICE ' ,
261 & myIter, myThid )
262 CALL PLOT_FIELD_XYRL( DXUICE , 'Current DXUICE ' ,
263 & myIter, myThid )
264 CALL PLOT_FIELD_XYRL( DYTICE , 'Current DYTICE ' ,
265 & myIter, myThid )
266 CALL PLOT_FIELD_XYRL( DYUICE , 'Current DYUICE ' ,
267 & myIter, myThid )
268 CALL PLOT_FIELD_XYRL( HEFFM , 'Current HEFFM ' ,
269 & myIter, myThid )
270 CALL PLOT_FIELD_XYRL( UVM , 'Current UVM ' ,
271 & myIter, myThid )
272 ENDIF
273
274 C-- Set model variables to initial/restart conditions
275 IF ( nIter0 .NE. 0 ) THEN
276
277 CALL SEAICE_READ_PICKUP ( myThid )
278
279 ELSE
280
281 DO bj=myByLo(myThid),myByHi(myThid)
282 DO bi=myBxLo(myThid),myBxHi(myThid)
283 DO j=1-OLy,sNy+OLy
284 DO i=1-OLx,sNx+OLx
285 HSNOW(I,J,bi,bj)=0.2*HEFFM(i,j,bi,bj)
286 YNEG(I,J,bi,bj)=ZERO
287 TMIX(I,J,bi,bj)=TICE(I,J,bi,bj)
288 DO k=1,3
289 HEFF(I,J,k,bi,bj)=SEAICE_initialHEFF*HEFFM(i,j,bi,bj)
290 UICE(I,J,k,bi,bj)=ZERO
291 VICE(I,J,k,bi,bj)=ZERO
292 ENDDO
293 ENDDO
294 ENDDO
295 ENDDO
296 ENDDO
297
298 C-- Read initial sea-ice thickness from file if available.
299 IF ( HeffFile .NE. ' ' ) THEN
300 _BEGIN_MASTER( myThid )
301 CALL READ_FLD_XY_RL( HeffFile, ' ', ZETA, 0, myThid )
302 _END_MASTER(myThid)
303 _EXCH_XY_R8(ZETA,myThid)
304 DO bj=myByLo(myThid),myByHi(myThid)
305 DO bi=myBxLo(myThid),myBxHi(myThid)
306 DO j=1-OLy,sNy+OLy
307 DO i=1-OLx,sNx+OLx
308 DO k=1,3
309 HEFF(I,J,k,bi,bj) = MAX(ZETA(i,j,bi,bj),ZERO)
310 ENDDO
311 ENDDO
312 ENDDO
313 ENDDO
314 ENDDO
315 ENDIF
316
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 DO k=1,3
322 IF(HEFF(I,J,k,bi,bj).GT.ZERO) AREA(I,J,k,bi,bj)=ONE
323 ENDDO
324 ENDDO
325 ENDDO
326 ENDDO
327 ENDDO
328
329 ENDIF
330
331 C--- Complete initialization
332 DO bj=myByLo(myThid),myByHi(myThid)
333 DO bi=myBxLo(myThid),myBxHi(myThid)
334 DO j=1-OLy,sNy+OLy
335 DO i=1-OLx,sNx+OLx
336 ZETA(I,J,bi,bj)=HEFF(I,J,1,bi,bj)*(1.0 _d 11)
337 ETA(I,J,bi,bj)=ZETA(I,J,bi,bj)/4.0 _d 0
338 ENDDO
339 ENDDO
340 ENDDO
341 ENDDO
342
343 RETURN
344 END

  ViewVC Help
Powered by ViewVC 1.1.22