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

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

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


Revision 1.24 - (show annotations) (download)
Sat May 4 17:32:04 2013 UTC (11 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64h
Changes since 1.23: +5 -6 lines
change to f77 syntax (to enable to compile with g77 & pgf77)

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_fixed.F,v 1.23 2013/05/03 19:43:09 torge Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5
6 CStartOfInterface
7 SUBROUTINE SEAICE_INIT_FIXED( myThid )
8 C *==========================================================*
9 C | SUBROUTINE SEAICE_INIT_FIXED
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 "FFIELDS.h"
21 #include "SEAICE_SIZE.h"
22 #include "SEAICE_PARAMS.h"
23 #include "SEAICE.h"
24 #include "SEAICE_TRACER.h"
25
26 C === Routine arguments ===
27 C myThid - Thread no. that called this routine.
28 INTEGER myThid
29 CEndOfInterface
30
31 C === Local variables ===
32 #ifdef SEAICE_ITD
33 C msgBuf :: Informational/error message buffer
34 CHARACTER*(MAX_LEN_MBUF) msgBuf
35 CHARACTER*15 HlimitMsgFormat
36 C k - loop counter for ITD categories
37 INTEGER k
38 #endif
39 C i,j,k,bi,bj - Loop counters
40
41 INTEGER i, j, bi, bj
42 INTEGER kSurface
43 #ifdef ALLOW_SITRACER
44 INTEGER iTracer
45 #endif
46 #ifndef SEAICE_CGRID
47 _RS mask_uice
48 #endif
49 #ifdef SHORTWAVE_HEATING
50 cif Helper variable for determining the fraction of sw radiation
51 cif penetrating the model shallowest layer
52 _RL dummyTime
53 _RL swfracba(2)
54 _RL tmpFac
55 #endif /* SHORTWAVE_HEATING */
56
57 IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
58 kSurface = Nr
59 ELSE
60 kSurface = 1
61 ENDIF
62
63 C Initialize MNC variable information for SEAICE
64 IF ( useMNC .AND.
65 & (seaice_tave_mnc.OR.seaice_dump_mnc.OR.SEAICE_mon_mnc)
66 & ) THEN
67 CALL SEAICE_MNC_INIT( myThid )
68 ENDIF
69
70 _BEGIN_MASTER(myThid)
71 #ifdef SHORTWAVE_HEATING
72 tmpFac = -1.0
73 dummyTime = 1.0
74 swfracba(1) = ABS(rF(1))
75 swfracba(2) = ABS(rF(2))
76 CALL SWFRAC(
77 I 2, tmpFac,
78 U swfracba,
79 I dummyTime, 0, myThid )
80 SWFracB = swfracba(2)
81 #else /* SHORTWAVE_HEATING */
82 SWFracB = 0. _d 0
83 #endif /* SHORTWAVE_HEATING */
84 _END_MASTER(myThid)
85
86 C-- Set mcPheePiston coeff (if still unset)
87 _BEGIN_MASTER(myThid)
88 IF ( SEAICE_mcPheePiston.EQ.UNSET_RL ) THEN
89 IF ( SEAICE_availHeatFrac.NE.UNSET_RL ) THEN
90 SEAICE_mcPheePiston = SEAICE_availHeatFrac
91 & * drF(kSurface)/SEAICE_deltaTtherm
92 ELSE
93 SEAICE_mcPheePiston = MCPHEE_TAPER_FAC
94 & * STANTON_NUMBER * USTAR_BASE
95 SEAICE_mcPheePiston = MIN( SEAICE_mcPheePiston,
96 & drF(kSurface)/SEAICE_deltaTtherm )
97 ENDIF
98 ENDIF
99 _END_MASTER(myThid)
100
101 C-- SItracer specifications for basic tracers
102 #ifdef ALLOW_SITRACER
103 _BEGIN_MASTER(myThid)
104 DO iTracer = 1, SItrNumInUse
105 C "ice concentration" tracer that should remain .EQ.1.
106 IF (SItrName(iTracer).EQ.'one') THEN
107 SItrFromOcean0(iTracer) =ONE
108 SItrFromFlood0(iTracer) =ONE
109 SItrExpand0(iTracer) =ONE
110 SItrFromOceanFrac(iTracer) =ZERO
111 SItrFromFloodFrac(iTracer) =ZERO
112 ENDIF
113 C age tracer: no age in ocean, or effect from ice cover changes
114 IF (SItrName(iTracer).EQ.'age') THEN
115 SItrFromOcean0(iTracer) =ZERO
116 SItrFromFlood0(iTracer) =ZERO
117 SItrExpand0(iTracer) =ZERO
118 SItrFromOceanFrac(iTracer) =ZERO
119 SItrFromFloodFrac(iTracer) =ZERO
120 ENDIf
121 C salinity tracer:
122 IF (SItrName(iTracer).EQ.'salinity') THEN
123 SItrMate(iTracer) ='HEFF'
124 SItrExpand0(iTracer) =ZERO
125 IF ( SEAICE_salinityTracer ) THEN
126 SEAICE_salt0 = ZERO
127 SEAICE_saltFrac = ZERO
128 ENDIF
129 ENDIF
130 C simple, made up, ice surface roughness index prototype
131 IF (SItrName(iTracer).EQ.'ridge') THEN
132 SItrMate(iTracer) ='AREA'
133 SItrFromOcean0(iTracer) =ZERO
134 SItrFromFlood0(iTracer) =ZERO
135 SItrExpand0(iTracer) =ZERO
136 SItrFromOceanFrac(iTracer) =ZERO
137 SItrFromFloodFrac(iTracer) =ZERO
138 ENDIF
139 ENDDO
140 _END_MASTER(myThid)
141 #endif
142
143 C-- all threads wait for master to finish initialisation of shared params
144 _BARRIER
145
146 C-- Initialize grid info
147 DO bj=myByLo(myThid),myByHi(myThid)
148 DO bi=myBxLo(myThid),myBxHi(myThid)
149 DO j=1-OLy,sNy+OLy
150 DO i=1-OLx,sNx+OLx
151 HEFFM(i,j,bi,bj) = 0. _d 0
152 ENDDO
153 ENDDO
154 DO j=1-OLy,sNy+OLy
155 DO i=1-OLx,sNx+OLx
156 HEFFM(i,j,bi,bj)= 1. _d 0
157 IF (_hFacC(i,j,kSurface,bi,bj).eq.0.)
158 & HEFFM(i,j,bi,bj)= 0. _d 0
159 ENDDO
160 ENDDO
161 #ifndef SEAICE_CGRID
162 DO j=1-OLy+1,sNy+OLy
163 DO i=1-OLx+1,sNx+OLx
164 UVM(i,j,bi,bj)=0. _d 0
165 mask_uice=HEFFM(i,j, bi,bj)+HEFFM(i-1,j-1,bi,bj)
166 & +HEFFM(i,j-1,bi,bj)+HEFFM(i-1,j, bi,bj)
167 IF(mask_uice.GT.3.5 _d 0) UVM(i,j,bi,bj)=1. _d 0
168 ENDDO
169 ENDDO
170 #endif /* SEAICE_CGRID */
171 ENDDO
172 ENDDO
173
174 #ifdef SEAICE_ITD
175 C use Equ. 22 of Lipscomb et al. (2001, JGR) to generate ice
176 C thickness category limits:
177 C - dependends on given number of categories nITD
178 C - choose between original parameters of Lipscomb et al. (2001):
179 C c1=3.0/N, c2=15*c1, c3=3.0
180 C or emphasize thin end of ITD (in order to enhance ice growth):
181 C c1=1.5/N, c2=42*c1, c3=3.3
182 C -> HINT: set parameters c1, c2 and c3 in seaice_readparms.F
183 C
184 Hlimit(0) = 0.0
185 Hlimit_c1=Hlimit_c1/float(nITD)
186 Hlimit_c2=Hlimit_c2*Hlimit_c1
187 IF (nITD .gt. 1) THEN
188 DO k=1,nITD-1
189 Hlimit(k) = Hlimit(k-1)
190 & + Hlimit_c1
191 & + Hlimit_c2
192 & *( 1.+tanh( Hlimit_c3 *(float(k-1)/float(nITD) -1. ) ) )
193 ENDDO
194 ENDIF
195 C thickest category is unlimited
196 Hlimit(nITD) = 999.9
197
198 WRITE(msgBuf,'(A,I2,A)')
199 & ' SEAICE_INIT_FIXED: ', nITD,
200 & ' sea ice thickness categories'
201 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
202 & SQUEEZE_RIGHT , myThid)
203 WRITE(HlimitMsgFormat,'(A,I2,A)') '(A,',nITD,'F6.2,F6.1)'
204 WRITE(msgBuf,HlimitMsgFormat)
205 & ' SEAICE_INIT_FIXED: Hlimit = ', (Hlimit(k),k=0,nITD)
206 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
207 & SQUEEZE_RIGHT , myThid)
208
209 #endif
210 C coefficients for metric terms
211 DO bj=myByLo(myThid),myByHi(myThid)
212 DO bi=myBxLo(myThid),myBxHi(myThid)
213 #ifdef SEAICE_CGRID
214 DO j=1-OLy,sNy+OLy
215 DO i=1-OLx,sNx+OLx
216 k1AtC(I,J,bi,bj) = 0.0 _d 0
217 k1AtZ(I,J,bi,bj) = 0.0 _d 0
218 k2AtC(I,J,bi,bj) = 0.0 _d 0
219 k2AtZ(I,J,bi,bj) = 0.0 _d 0
220 ENDDO
221 ENDDO
222 IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN
223 C This is the only case where tan(phi) is not zero. In this case
224 C C and U points, and Z and V points have the same phi, so that we
225 C only need a copy here. Do not use tan(YC) and tan(YG), because these
226 C can be the geographical coordinates and not the correct grid
227 C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0)
228 DO j=1-OLy,sNy+OLy
229 DO i=1-OLx,sNx+OLx
230 k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere
231 k2AtZ(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere
232 ENDDO
233 ENDDO
234 ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN
235 C compute metric term coefficients from finite difference approximation
236 DO j=1-OLy,sNy+OLy
237 DO i=1-OLx,sNx+OLx-1
238 k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj)
239 & * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) )
240 & * _recip_dxF(I,J,bi,bj)
241 ENDDO
242 ENDDO
243 DO j=1-OLy,sNy+OLy
244 DO i=1-OLx+1,sNx+OLx
245 k1AtZ(I,J,bi,bj) = _recip_dyU(I,J,bi,bj)
246 & * ( _dyC(I,J,bi,bj) - _dyC(I-1,J,bi,bj) )
247 & * _recip_dxV(I,J,bi,bj)
248 ENDDO
249 ENDDO
250 DO j=1-OLy,sNy+OLy-1
251 DO i=1-OLx,sNx+OLx
252 k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj)
253 & * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) )
254 & * _recip_dyF(I,J,bi,bj)
255 ENDDO
256 ENDDO
257 DO j=1-OLy+1,sNy+OLy
258 DO i=1-OLx,sNx+OLx
259 k2AtZ(I,J,bi,bj) = _recip_dxV(I,J,bi,bj)
260 & * ( _dxC(I,J,bi,bj) - _dxC(I,J-1,bi,bj) )
261 & * _recip_dyU(I,J,bi,bj)
262 ENDDO
263 ENDDO
264 ENDIF
265 #else /* not SEAICE_CGRID */
266 DO j=1-OLy,sNy+OLy
267 DO i=1-OLx,sNx+OLx
268 k1AtC(I,J,bi,bj) = 0.0 _d 0
269 k1AtU(I,J,bi,bj) = 0.0 _d 0
270 k1AtV(I,J,bi,bj) = 0.0 _d 0
271 k2AtC(I,J,bi,bj) = 0.0 _d 0
272 k2AtU(I,J,bi,bj) = 0.0 _d 0
273 k2AtV(I,J,bi,bj) = 0.0 _d 0
274 ENDDO
275 ENDDO
276 IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN
277 C This is the only case where tan(phi) is not zero. In this case
278 C C and U points, and Z and V points have the same phi, so that we
279 C only need a copy here. Do not use tan(YC) and tan(YG), because these
280 C can be the geographical coordinates and not the correct grid
281 C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0)
282 DO j=1-OLy,sNy+OLy
283 DO i=1-OLx,sNx+OLx
284 k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere
285 k2AtU(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere
286 k2AtV(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere
287 ENDDO
288 ENDDO
289 ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN
290 C compute metric term coefficients from finite difference approximation
291 DO j=1-OLy,sNy+OLy
292 DO i=1-OLx,sNx+OLx-1
293 k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj)
294 & * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) )
295 & * _recip_dxF(I,J,bi,bj)
296 ENDDO
297 ENDDO
298 DO j=1-OLy,sNy+OLy
299 DO i=1-OLx+1,sNx+OLx
300 k1AtU(I,J,bi,bj) = _recip_dyG(I,J,bi,bj)
301 & * ( _dyF(I,J,bi,bj) - _dyF(I-1,J,bi,bj) )
302 & * _recip_dxC(I,J,bi,bj)
303 ENDDO
304 ENDDO
305 DO j=1-OLy,sNy+OLy
306 DO i=1-OLx,sNx+OLx-1
307 k1AtV(I,J,bi,bj) = _recip_dyC(I,J,bi,bj)
308 & * ( _dyU(I+1,J,bi,bj) - _dyU(I,J,bi,bj) )
309 & * _recip_dxG(I,J,bi,bj)
310 ENDDO
311 ENDDO
312 DO j=1-OLy,sNy+OLy-1
313 DO i=1-OLx,sNx+OLx
314 k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj)
315 & * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) )
316 & * _recip_dyF(I,J,bi,bj)
317 ENDDO
318 ENDDO
319 DO j=1-OLy,sNy+OLy-1
320 DO i=1-OLx,sNx+OLx
321 k2AtU(I,J,bi,bj) = _recip_dxC(I,J,bi,bj)
322 & * ( _dxV(I,J+1,bi,bj) - _dxV(I,J,bi,bj) )
323 & * _recip_dyG(I,J,bi,bj)
324 ENDDO
325 ENDDO
326 DO j=1-OLy+1,sNy+OLy
327 DO i=1-OLx,sNx+OLx
328 k2AtV(I,J,bi,bj) = _recip_dxG(I,J,bi,bj)
329 & * ( _dxF(I,J,bi,bj) - _dxF(I,J-1,bi,bj) )
330 & * _recip_dyC(I,J,bi,bj)
331 ENDDO
332 ENDDO
333 ENDIF
334 #endif /* not SEAICE_CGRID */
335 ENDDO
336 ENDDO
337
338 #ifndef SEAICE_CGRID
339 C-- Choose a proxy level for geostrophic velocity,
340 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 KGEO(i,j,bi,bj) = 0
345 ENDDO
346 ENDDO
347 DO j=1-OLy,sNy+OLy
348 DO i=1-OLx,sNx+OLx
349 #ifdef SEAICE_BICE_STRESS
350 KGEO(i,j,bi,bj) = 1
351 #else /* SEAICE_BICE_STRESS */
352 IF (klowc(i,j,bi,bj) .LT. 2) THEN
353 KGEO(i,j,bi,bj) = 1
354 ELSE
355 KGEO(i,j,bi,bj) = 2
356 DO WHILE ( abs(rC(KGEO(i,j,bi,bj))) .LT. 50.0 _d 0 .AND.
357 & KGEO(i,j,bi,bj) .LT. (klowc(i,j,bi,bj)-1) )
358 KGEO(i,j,bi,bj) = KGEO(i,j,bi,bj) + 1
359 ENDDO
360 ENDIF
361 #endif /* SEAICE_BICE_STRESS */
362 ENDDO
363 ENDDO
364 ENDDO
365 ENDDO
366 #endif /* SEAICE_CGRID */
367
368 #ifdef SEAICE_ALLOW_JFNK
369 C initialise some diagnostic counters for the JFNK solver
370 totalNewtonIters = 0
371 totalNewtonFails = 0
372 totalKrylovIters = 0
373 totalKrylovFails = 0
374 totalJFNKtimeSteps = 0
375 C this cannot be done here, because globalArea is only defined
376 C after S/R PACKAGES_INIT_FIXED, so we move it to S/R SEAICE_INIT_VARIA
377 CML CALL SEAICE_MAP2VEC( nVec, rAw, rAs,
378 CML & scalarProductMetric, .TRUE., myThid )
379 CML DO bj=myByLo(myThid),myByHi(myThid)
380 CML DO bi=myBxLo(myThid),myBxHi(myThid)
381 CML DO i=1,nVec
382 CML scalarProductMetric(i,1,bi,bj) =
383 CML & scalarProductMetric(i,1,bi,bj)/globalArea
384 CML ENDDO
385 CML ENDDO
386 CML ENDDO
387 #endif /* SEAICE_ALLOW_JFNK */
388
389 #ifdef ALLOW_DIAGNOSTICS
390 IF ( useDiagnostics ) THEN
391 CALL SEAICE_DIAGNOSTICS_INIT( myThid )
392 ENDIF
393 #endif
394
395 C-- Summarise pkg/seaice configuration
396 CALL SEAICE_SUMMARY( myThid )
397
398 RETURN
399 END

  ViewVC Help
Powered by ViewVC 1.1.22