/[MITgcm]/MITgcm_contrib/arctic40km/code/seaice_init.F
ViewVC logotype

Contents of /MITgcm_contrib/arctic40km/code/seaice_init.F

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


Revision 1.2 - (show annotations) (download)
Sat May 6 13:46:39 2006 UTC (19 years, 2 months ago) by dimitri
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +1 -1 lines
FILE REMOVED
Removed arctic40km/code/SEAICE_OPTIONS.h and seaice_init.F.
These files are not needed since Martin's pkg/seaice changes
that make use of model's grid definition directly.  However,
this change introduces some numerical noise in the coupled
ice-ocean results, whose cause has not yet been identified.

1 C $Header: /u/gcmpack/MITgcm_contrib/arctic40km/code/seaice_init.F,v 1.1 2005/09/01 14:16:40 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 #ifdef ALLOW_SEAICE
35 C === Local variables ===
36 C i,j,k,bi,bj - Loop counters
37
38 INTEGER i, j, k, bi, bj
39 _RS mask_uice
40 INTEGER myIter, myTile
41 INTEGER iG, jG
42 _RL lat, dlat, dlon, xG0, yG0
43 _RL xGloc(1-Olx:sNx+Olx+1,1-Oly:sNy+Oly+1)
44 _RL yGloc(1-Olx:sNx+Olx+1,1-Oly:sNy+Oly+1)
45 INTEGER iGl,jGl
46 iGl(I,bi) = 1+mod(myXGlobalLo-1+(bi-1)*sNx+I+Olx*Nx-1,Nx)
47 jGl(J,bj) = 1+mod(myYGlobalLo-1+(bj-1)*sNy+J+Oly*Ny-1,Ny)
48
49 #ifdef ALLOW_TIMEAVE
50 C Initialize averages to zero
51 DO bj = myByLo(myThid), myByHi(myThid)
52 DO bi = myBxLo(myThid), myBxHi(myThid)
53 CALL TIMEAVE_RESET(FUtave ,1,bi,bj,myThid)
54 CALL TIMEAVE_RESET(FVtave ,1,bi,bj,myThid)
55 CALL TIMEAVE_RESET(EmPmRtave,1,bi,bj,myThid)
56 CALL TIMEAVE_RESET(QNETtave ,1,bi,bj,myThid)
57 CALL TIMEAVE_RESET(QSWtave ,1,bi,bj,myThid)
58 CALL TIMEAVE_RESET(UICEtave ,1,bi,bj,myThid)
59 CALL TIMEAVE_RESET(VICEtave ,1,bi,bj,myThid)
60 CALL TIMEAVE_RESET(HEFFtave ,1,bi,bj,myThid)
61 CALL TIMEAVE_RESET(AREAtave ,1,bi,bj,myThid)
62 DO k=1,Nr
63 SEAICE_TimeAve(k,bi,bj)=ZERO
64 ENDDO
65 ENDDO
66 ENDDO
67 #endif /* ALLOW_TIMEAVE */
68
69 cph(
70 cph make sure TAF sees proper initialisation
71 cph to avoid partial recomputation issues
72 DO bj=myByLo(myThid),myByHi(myThid)
73 DO bi=myBxLo(myThid),myBxHi(myThid)
74 c
75 DO K=1,3
76 DO J=1-OLy,sNy+OLy
77 DO I=1-OLx,sNx+OLx
78 HEFF(I,J,k,bi,bj)=ZERO
79 AREA(I,J,k,bi,bj)=ZERO
80 UICE(I,J,k,bi,bj)=ZERO
81 VICE(I,J,k,bi,bj)=ZERO
82 ENDDO
83 ENDDO
84 ENDDO
85 c
86 DO J=1-OLy,sNy+OLy
87 DO I=1-OLx,sNx+OLx
88 HSNOW(I,J,bi,bj)=ZERO
89 ZETA(I,J,bi,bj)=ZERO
90 ENDDO
91 ENDDO
92 c
93 ENDDO
94 ENDDO
95 cph)
96
97 C_JZ- The following is for a rotated model grid (for the UW regional arctic model)
98 C For each tile ...
99 DO bj = myByLo(myThid), myByHi(myThid)
100 DO bi = myBxLo(myThid), myBxHi(myThid)
101
102 C-- "Global" index (place holder)
103 jG = myYGlobalLo + (bj-1)*sNy
104 iG = myXGlobalLo + (bi-1)*sNx
105
106 C-- First find coordinate of tile corner (meaning outer corner of halo)
107 xG0 = thetaMin
108 C Find the X-coordinate of the outer grid-line of the "real" tile
109 DO i=1, iG-1
110 xG0 = xG0 + delX(i)
111 ENDDO
112 C Back-step to the outer grid-line of the "halo" region
113 DO i=1, Olx
114 xG0 = xG0 - delX( 1+mod(Olx*Nx-1+iG-i,Nx) )
115 ENDDO
116 C Find the Y-coordinate of the outer grid-line of the "real" tile
117 yG0 = phiMin
118 DO j=1, jG-1
119 yG0 = yG0 + delY(j)
120 ENDDO
121 C Back-step to the outer grid-line of the "halo" region
122 DO j=1, Oly
123 yG0 = yG0 - delY( 1+mod(Oly*Ny-1+jG-j,Ny) )
124 ENDDO
125
126 C-- Calculate coordinates of cell corners for N+1 grid-lines
127 DO J=1-Oly,sNy+Oly +1
128 xGloc(1-Olx,J) = xG0
129 DO I=1-Olx,sNx+Olx
130 c xGloc(I+1,J) = xGloc(I,J) + delX(1+mod(Nx-1+iG-1+i,Nx))
131 xGloc(I+1,J) = xGloc(I,J) + delX( iGl(I,bi) )
132 ENDDO
133 ENDDO
134 DO I=1-Olx,sNx+Olx +1
135 yGloc(I,1-Oly) = yG0
136 DO J=1-Oly,sNy+Oly
137 c yGloc(I,J+1) = yGloc(I,J) + delY(1+mod(Ny-1+jG-1+j,Ny))
138 yGloc(I,J+1) = yGloc(I,J) + delY( jGl(J,bj) )
139 ENDDO
140 ENDDO
141 DO J=1-Oly,sNy+Oly
142 DO I=1-Olx,sNx+Olx
143 CSUICE(i,j,bi,bj) =cos(yGloc(I,J)*deg2rad)
144 SINEICE(i,j,bi,bj) =sin(yGloc(I,J)*deg2rad)
145 TNGICE(i,j,bi,bj)=SINEICE(i,j,bi,bj)/CSUICE(i,j,bi,bj)
146 ENDDO
147 ENDDO
148
149 C-- Calculate [xC,yC], coordinates of cell centers by averaging
150 DO J=1-Oly,sNy+Oly
151 DO I=1-Olx,sNx+Olx
152 CSTICE(i,j,bi,bj) =cos(0.25*(
153 & yGloc(I,J)+yGloc(I+1,J)+yGloc(I,J+1)+yGloc(I+1,J+1))*deg2rad)
154 SINEICE(i,j,bi,bj) =sin(0.25*(
155 & yGloc(I,J)+yGloc(I+1,J)+yGloc(I,J+1)+yGloc(I+1,J+1))*deg2rad)
156 TNGTICE(i,j,bi,bj)=SINEICE(i,j,bi,bj)/CSTICE(i,j,bi,bj)
157 ENDDO
158 ENDDO
159
160 ENDDO ! bi
161 ENDDO ! bj
162
163 C-- Initialize grid info
164 DO bj=myByLo(myThid),myByHi(myThid)
165 DO bi=myBxLo(myThid),myBxHi(myThid)
166 DO J=1-OLy,sNy+OLy
167 DO I=1-OLx,sNx+OLx
168 SINEICE(i,j,bi,bj)=sin(yG(I,J,bi,bj)*deg2rad)
169 DXTICE(i,j,bi,bj)=dxF(i,j,bi,bj)/CSTICE(i,j,bi,bj)
170 DXUICE(i,j,bi,bj)=dxV(i,j,bi,bj)/CSUICE(i,j,bi,bj)
171 DYTICE(i,j,bi,bj)=dyF(i,j,bi,bj)
172 DYUICE(i,j,bi,bj)=dyU(i,j,bi,bj)
173 ENDDO
174 ENDDO
175 DO j=1-OLy,sNy+OLy
176 DO i=1-OLx,sNx+OLx
177 HEFFM(i,j,bi,bj)=ONE
178 IF (_hFacC(i,j,1,bi,bj).eq.0.) HEFFM(i,j,bi,bj)=ZERO
179 ENDDO
180 ENDDO
181 DO J=2-OLy,sNy+OLy
182 DO I=2-OLx,sNx+OLx
183 UVM(i,j,bi,bj)=ZERO
184 mask_uice=HEFFM(I,J, bi,bj)+HEFFM(I-1,J-1,bi,bj)
185 & +HEFFM(I,J-1,bi,bj)+HEFFM(I-1,J, bi,bj)
186 IF(mask_uice.GT.3.5) UVM(I,J,bi,bj)=ONE
187 ENDDO
188 ENDDO
189
190 #ifdef ALLOW_EXCH2
191 C-- Special stuff for cubed sphere: assume grid is rectangular and
192 C set UV mask to zero except for Arctic and Antarctic cube faces.
193 IF (useCubedSphereExchange) THEN
194 DO J=1-OLy,sNy+OLy
195 DO I=1-OLx,sNx+OLx
196 CSTICE(i,j,bi,bj) = ONE
197 CSUICE(i,j,bi,bj) = ONE
198 TNGTICE(i,j,bi,bj)= ZERO
199 TNGICE(i,j,bi,bj) = ZERO
200 DXTICE(i,j,bi,bj) = dxF(i,j,bi,bj)
201 DXUICE(i,j,bi,bj) = dxV(i,j,bi,bj)
202 ENDDO
203 ENDDO
204 myTile = W2_myTileList(bi)
205 IF ( exch2_myFace(myTile) .EQ. 1 .OR.
206 & exch2_myFace(myTile) .EQ. 2 .OR.
207 & exch2_myFace(myTile) .EQ. 4 .OR.
208 & exch2_myFace(myTile) .EQ. 5 ) THEN
209 DO J=1-OLy,sNy+OLy
210 DO I=1-OLx,sNx+OLx
211 UVM(i,j,bi,bj)=ZERO
212 ENDDO
213 ENDDO
214 ELSEIF ( exch2_isWedge(myTile) .EQ. 1 ) THEN
215 I=1
216 DO J=1-OLy,sNy+OLy
217 UVM(i,j,bi,bj)=ZERO
218 ENDDO
219 ELSEIF ( exch2_isSedge(myTile) .EQ. 1 ) THEN
220 J=1
221 DO I=1-OLx,sNx+OLx
222 UVM(i,j,bi,bj)=ZERO
223 ENDDO
224 ENDIF
225 ENDIF
226 #endif /* ALLOW_EXCH2 */
227
228 DO j=1-OLy,sNy+OLy
229 DO i=1-OLx,sNx+OLx
230 TICE(I,J,bi,bj)=273.0 _d 0
231 #ifdef SEAICE_MULTILEVEL
232 DO k=1,MULTDIM
233 TICES(I,J,k,bi,bj)=273.0 _d 0
234 ENDDO
235 #endif
236 UICEC(I,J,bi,bj)=ZERO
237 VICEC(I,J,bi,bj)=ZERO
238 AMASS(I,J,bi,bj)=1000.0 _d 0
239 ENDDO
240 ENDDO
241 ENDDO
242 ENDDO
243
244 C-- Update overlap regions
245 _EXCH_XY_R8(UVM, myThid)
246
247 myIter=0
248 CALL PLOT_FIELD_XYRL( CSTICE , 'Current CSTICE ' ,
249 & myIter, myThid )
250 CALL PLOT_FIELD_XYRL( CSUICE , 'Current CSUICE ' ,
251 & myIter, myThid )
252 CALL PLOT_FIELD_XYRL( TNGTICE , 'Current TNGTICE ' ,
253 & myIter, myThid )
254 CALL PLOT_FIELD_XYRL( TNGICE , 'Current TNGICE ' ,
255 & myIter, myThid )
256 CALL PLOT_FIELD_XYRL( SINEICE , 'Current SINEICE ' ,
257 & myIter, myThid )
258 CALL PLOT_FIELD_XYRL( DXTICE , 'Current DXTICE ' ,
259 & myIter, myThid )
260 CALL PLOT_FIELD_XYRL( DXUICE , 'Current DXUICE ' ,
261 & myIter, myThid )
262 CALL PLOT_FIELD_XYRL( DYTICE , 'Current DYTICE ' ,
263 & myIter, myThid )
264 CALL PLOT_FIELD_XYRL( DYUICE , 'Current DYUICE ' ,
265 & myIter, myThid )
266 CALL PLOT_FIELD_XYRL( HEFFM , 'Current HEFFM ' ,
267 & myIter, myThid )
268 CALL PLOT_FIELD_XYRL( UVM , 'Current UVM ' ,
269 & myIter, myThid )
270
271 C-- Set model variables to initial/restart conditions
272 IF ( nIter0 .NE. 0 ) 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 HSNOW(I,J,bi,bj)=0.2*HEFFM(i,j,bi,bj)
283 YNEG(I,J,bi,bj)=ZERO
284 TMIX(I,J,bi,bj)=TICE(I,J,bi,bj)
285 DO k=1,3
286 HEFF(I,J,k,bi,bj)=SEAICE_initialHEFF*HEFFM(i,j,bi,bj)
287 UICE(I,J,k,bi,bj)=ZERO
288 VICE(I,J,k,bi,bj)=ZERO
289 ENDDO
290 ENDDO
291 ENDDO
292 ENDDO
293 ENDDO
294
295 C-- Read initial sea-ice thickness from file if available.
296 IF ( HeffFile .NE. ' ' ) THEN
297 _BEGIN_MASTER( myThid )
298 CALL READ_FLD_XY_RL( HeffFile, ' ', ZETA, 0, myThid )
299 _END_MASTER(myThid)
300 _EXCH_XY_R8(ZETA,myThid)
301 DO bj=myByLo(myThid),myByHi(myThid)
302 DO bi=myBxLo(myThid),myBxHi(myThid)
303 DO j=1-OLy,sNy+OLy
304 DO i=1-OLx,sNx+OLx
305 DO k=1,3
306 HEFF(I,J,k,bi,bj) = MAX(ZETA(i,j,bi,bj),ZERO)
307 ENDDO
308 ENDDO
309 ENDDO
310 ENDDO
311 ENDDO
312 ENDIF
313
314 DO bj=myByLo(myThid),myByHi(myThid)
315 DO bi=myBxLo(myThid),myBxHi(myThid)
316 DO j=1-OLy,sNy+OLy
317 DO i=1-OLx,sNx+OLx
318 DO k=1,3
319 IF(HEFF(I,J,k,bi,bj).GT.ZERO) AREA(I,J,k,bi,bj)=ONE
320 ENDDO
321 ENDDO
322 ENDDO
323 ENDDO
324 ENDDO
325
326 ENDIF
327
328 C--- Complete initialization
329 DO bj=myByLo(myThid),myByHi(myThid)
330 DO bi=myBxLo(myThid),myBxHi(myThid)
331 DO j=1-OLy,sNy+OLy
332 DO i=1-OLx,sNx+OLx
333 ZETA(I,J,bi,bj)=HEFF(I,J,1,bi,bj)*(1.0 _d 11)
334 ETA(I,J,bi,bj)=ZETA(I,J,bi,bj)/4.0 _d 0
335 ENDDO
336 ENDDO
337 ENDDO
338 ENDDO
339
340 #endif /* ALLOW_SEAICE */
341
342 RETURN
343 END

  ViewVC Help
Powered by ViewVC 1.1.22