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

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

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


Revision 1.25 - (hide annotations) (download)
Thu May 30 14:07:19 2013 UTC (11 years ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64t, checkpoint64i, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n
Changes since 1.24: +5 -1 lines
add Adams-Bashforth2 time discretization for ice dynamics;
so far only for JFNK-solver
 - requires an additional time level that is stored in the pickup
 - you can start with AB2 from a pickup without this time level with
   pickupStrictlyMatch = .false.

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

  ViewVC Help
Powered by ViewVC 1.1.22