1 |
torge |
1.23 |
C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_fixed.F,v 1.22 2012/12/17 15:06:02 mlosch 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 |
jmc |
1.12 |
_BEGIN_MASTER(myThid) |
71 |
heimbach |
1.3 |
#ifdef SHORTWAVE_HEATING |
72 |
jmc |
1.12 |
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 |
heimbach |
1.3 |
U swfracba, |
79 |
jmc |
1.12 |
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 |
jmc |
1.18 |
|
86 |
jmc |
1.19 |
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 |
gforget |
1.16 |
C-- SItracer specifications for basic tracers |
102 |
|
|
#ifdef ALLOW_SITRACER |
103 |
jmc |
1.19 |
_BEGIN_MASTER(myThid) |
104 |
gforget |
1.16 |
DO iTracer = 1, SItrNumInUse |
105 |
jmc |
1.19 |
C "ice concentration" tracer that should remain .EQ.1. |
106 |
|
|
IF (SItrName(iTracer).EQ.'one') THEN |
107 |
gforget |
1.16 |
SItrFromOcean0(iTracer) =ONE |
108 |
|
|
SItrFromFlood0(iTracer) =ONE |
109 |
|
|
SItrExpand0(iTracer) =ONE |
110 |
|
|
SItrFromOceanFrac(iTracer) =ZERO |
111 |
|
|
SItrFromFloodFrac(iTracer) =ZERO |
112 |
jmc |
1.19 |
ENDIF |
113 |
|
|
C age tracer: no age in ocean, or effect from ice cover changes |
114 |
|
|
IF (SItrName(iTracer).EQ.'age') THEN |
115 |
gforget |
1.16 |
SItrFromOcean0(iTracer) =ZERO |
116 |
|
|
SItrFromFlood0(iTracer) =ZERO |
117 |
|
|
SItrExpand0(iTracer) =ZERO |
118 |
|
|
SItrFromOceanFrac(iTracer) =ZERO |
119 |
|
|
SItrFromFloodFrac(iTracer) =ZERO |
120 |
jmc |
1.19 |
ENDIf |
121 |
|
|
C salinity tracer: |
122 |
|
|
IF (SItrName(iTracer).EQ.'salinity') THEN |
123 |
gforget |
1.16 |
SItrMate(iTracer) ='HEFF' |
124 |
|
|
SItrExpand0(iTracer) =ZERO |
125 |
jmc |
1.19 |
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 |
gforget |
1.16 |
SItrMate(iTracer) ='AREA' |
133 |
|
|
SItrFromOcean0(iTracer) =ZERO |
134 |
|
|
SItrFromFlood0(iTracer) =ZERO |
135 |
|
|
SItrExpand0(iTracer) =ZERO |
136 |
|
|
SItrFromOceanFrac(iTracer) =ZERO |
137 |
|
|
SItrFromFloodFrac(iTracer) =ZERO |
138 |
jmc |
1.19 |
ENDIF |
139 |
gforget |
1.16 |
ENDDO |
140 |
jmc |
1.19 |
_END_MASTER(myThid) |
141 |
gforget |
1.16 |
#endif |
142 |
|
|
|
143 |
jmc |
1.19 |
C-- all threads wait for master to finish initialisation of shared params |
144 |
|
|
_BARRIER |
145 |
gforget |
1.14 |
|
146 |
mlosch |
1.4 |
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 |
jmc |
1.11 |
#ifndef SEAICE_CGRID |
162 |
mlosch |
1.4 |
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 |
jmc |
1.11 |
#endif /* SEAICE_CGRID */ |
171 |
mlosch |
1.4 |
ENDDO |
172 |
|
|
ENDDO |
173 |
|
|
|
174 |
heimbach |
1.20 |
#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 |
torge |
1.23 |
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 |
heimbach |
1.20 |
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 = ', |
206 |
|
|
& Hlimit(0:nITD-1),Hlimit(nITD) |
207 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
208 |
|
|
& SQUEEZE_RIGHT , myThid) |
209 |
|
|
|
210 |
|
|
#endif |
211 |
jmc |
1.6 |
C coefficients for metric terms |
212 |
mlosch |
1.4 |
DO bj=myByLo(myThid),myByHi(myThid) |
213 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
214 |
mlosch |
1.10 |
#ifdef SEAICE_CGRID |
215 |
mlosch |
1.4 |
DO j=1-OLy,sNy+OLy |
216 |
|
|
DO i=1-OLx,sNx+OLx |
217 |
|
|
k1AtC(I,J,bi,bj) = 0.0 _d 0 |
218 |
|
|
k1AtZ(I,J,bi,bj) = 0.0 _d 0 |
219 |
|
|
k2AtC(I,J,bi,bj) = 0.0 _d 0 |
220 |
|
|
k2AtZ(I,J,bi,bj) = 0.0 _d 0 |
221 |
|
|
ENDDO |
222 |
|
|
ENDDO |
223 |
|
|
IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN |
224 |
jmc |
1.6 |
C This is the only case where tan(phi) is not zero. In this case |
225 |
|
|
C C and U points, and Z and V points have the same phi, so that we |
226 |
mlosch |
1.4 |
C only need a copy here. Do not use tan(YC) and tan(YG), because these |
227 |
jmc |
1.6 |
C can be the geographical coordinates and not the correct grid |
228 |
mlosch |
1.4 |
C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0) |
229 |
|
|
DO j=1-OLy,sNy+OLy |
230 |
|
|
DO i=1-OLx,sNx+OLx |
231 |
|
|
k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
232 |
|
|
k2AtZ(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere |
233 |
|
|
ENDDO |
234 |
|
|
ENDDO |
235 |
|
|
ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN |
236 |
|
|
C compute metric term coefficients from finite difference approximation |
237 |
|
|
DO j=1-OLy,sNy+OLy |
238 |
|
|
DO i=1-OLx,sNx+OLx-1 |
239 |
|
|
k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj) |
240 |
|
|
& * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) ) |
241 |
|
|
& * _recip_dxF(I,J,bi,bj) |
242 |
|
|
ENDDO |
243 |
|
|
ENDDO |
244 |
|
|
DO j=1-OLy,sNy+OLy |
245 |
|
|
DO i=1-OLx+1,sNx+OLx |
246 |
jmc |
1.6 |
k1AtZ(I,J,bi,bj) = _recip_dyU(I,J,bi,bj) |
247 |
mlosch |
1.4 |
& * ( _dyC(I,J,bi,bj) - _dyC(I-1,J,bi,bj) ) |
248 |
|
|
& * _recip_dxV(I,J,bi,bj) |
249 |
|
|
ENDDO |
250 |
|
|
ENDDO |
251 |
|
|
DO j=1-OLy,sNy+OLy-1 |
252 |
|
|
DO i=1-OLx,sNx+OLx |
253 |
jmc |
1.6 |
k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj) |
254 |
mlosch |
1.4 |
& * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) ) |
255 |
|
|
& * _recip_dyF(I,J,bi,bj) |
256 |
|
|
ENDDO |
257 |
|
|
ENDDO |
258 |
|
|
DO j=1-OLy+1,sNy+OLy |
259 |
|
|
DO i=1-OLx,sNx+OLx |
260 |
mlosch |
1.9 |
k2AtZ(I,J,bi,bj) = _recip_dxV(I,J,bi,bj) |
261 |
mlosch |
1.4 |
& * ( _dxC(I,J,bi,bj) - _dxC(I,J-1,bi,bj) ) |
262 |
|
|
& * _recip_dyU(I,J,bi,bj) |
263 |
|
|
ENDDO |
264 |
|
|
ENDDO |
265 |
|
|
ENDIF |
266 |
mlosch |
1.10 |
#else /* not SEAICE_CGRID */ |
267 |
|
|
DO j=1-OLy,sNy+OLy |
268 |
|
|
DO i=1-OLx,sNx+OLx |
269 |
|
|
k1AtC(I,J,bi,bj) = 0.0 _d 0 |
270 |
|
|
k1AtU(I,J,bi,bj) = 0.0 _d 0 |
271 |
|
|
k1AtV(I,J,bi,bj) = 0.0 _d 0 |
272 |
|
|
k2AtC(I,J,bi,bj) = 0.0 _d 0 |
273 |
|
|
k2AtU(I,J,bi,bj) = 0.0 _d 0 |
274 |
|
|
k2AtV(I,J,bi,bj) = 0.0 _d 0 |
275 |
|
|
ENDDO |
276 |
|
|
ENDDO |
277 |
|
|
IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN |
278 |
|
|
C This is the only case where tan(phi) is not zero. In this case |
279 |
|
|
C C and U points, and Z and V points have the same phi, so that we |
280 |
|
|
C only need a copy here. Do not use tan(YC) and tan(YG), because these |
281 |
|
|
C can be the geographical coordinates and not the correct grid |
282 |
|
|
C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0) |
283 |
|
|
DO j=1-OLy,sNy+OLy |
284 |
|
|
DO i=1-OLx,sNx+OLx |
285 |
|
|
k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
286 |
|
|
k2AtU(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
287 |
|
|
k2AtV(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere |
288 |
|
|
ENDDO |
289 |
|
|
ENDDO |
290 |
|
|
ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN |
291 |
|
|
C compute metric term coefficients from finite difference approximation |
292 |
|
|
DO j=1-OLy,sNy+OLy |
293 |
|
|
DO i=1-OLx,sNx+OLx-1 |
294 |
|
|
k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj) |
295 |
|
|
& * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) ) |
296 |
|
|
& * _recip_dxF(I,J,bi,bj) |
297 |
|
|
ENDDO |
298 |
|
|
ENDDO |
299 |
|
|
DO j=1-OLy,sNy+OLy |
300 |
|
|
DO i=1-OLx+1,sNx+OLx |
301 |
|
|
k1AtU(I,J,bi,bj) = _recip_dyG(I,J,bi,bj) |
302 |
|
|
& * ( _dyF(I,J,bi,bj) - _dyF(I-1,J,bi,bj) ) |
303 |
|
|
& * _recip_dxC(I,J,bi,bj) |
304 |
|
|
ENDDO |
305 |
|
|
ENDDO |
306 |
|
|
DO j=1-OLy,sNy+OLy |
307 |
|
|
DO i=1-OLx,sNx+OLx-1 |
308 |
|
|
k1AtV(I,J,bi,bj) = _recip_dyC(I,J,bi,bj) |
309 |
|
|
& * ( _dyU(I+1,J,bi,bj) - _dyU(I,J,bi,bj) ) |
310 |
|
|
& * _recip_dxG(I,J,bi,bj) |
311 |
|
|
ENDDO |
312 |
|
|
ENDDO |
313 |
|
|
DO j=1-OLy,sNy+OLy-1 |
314 |
|
|
DO i=1-OLx,sNx+OLx |
315 |
|
|
k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj) |
316 |
|
|
& * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) ) |
317 |
|
|
& * _recip_dyF(I,J,bi,bj) |
318 |
|
|
ENDDO |
319 |
|
|
ENDDO |
320 |
|
|
DO j=1-OLy,sNy+OLy-1 |
321 |
|
|
DO i=1-OLx,sNx+OLx |
322 |
|
|
k2AtU(I,J,bi,bj) = _recip_dxC(I,J,bi,bj) |
323 |
|
|
& * ( _dxV(I,J+1,bi,bj) - _dxV(I,J,bi,bj) ) |
324 |
|
|
& * _recip_dyG(I,J,bi,bj) |
325 |
|
|
ENDDO |
326 |
|
|
ENDDO |
327 |
|
|
DO j=1-OLy+1,sNy+OLy |
328 |
|
|
DO i=1-OLx,sNx+OLx |
329 |
|
|
k2AtV(I,J,bi,bj) = _recip_dxG(I,J,bi,bj) |
330 |
|
|
& * ( _dxF(I,J,bi,bj) - _dxF(I,J-1,bi,bj) ) |
331 |
|
|
& * _recip_dyC(I,J,bi,bj) |
332 |
|
|
ENDDO |
333 |
|
|
ENDDO |
334 |
|
|
ENDIF |
335 |
|
|
#endif /* not SEAICE_CGRID */ |
336 |
mlosch |
1.4 |
ENDDO |
337 |
|
|
ENDDO |
338 |
|
|
|
339 |
|
|
#ifndef SEAICE_CGRID |
340 |
|
|
C-- Choose a proxy level for geostrophic velocity, |
341 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
342 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
343 |
|
|
DO j=1-OLy,sNy+OLy |
344 |
|
|
DO i=1-OLx,sNx+OLx |
345 |
|
|
KGEO(i,j,bi,bj) = 0 |
346 |
|
|
ENDDO |
347 |
|
|
ENDDO |
348 |
|
|
DO j=1-OLy,sNy+OLy |
349 |
|
|
DO i=1-OLx,sNx+OLx |
350 |
|
|
#ifdef SEAICE_BICE_STRESS |
351 |
|
|
KGEO(i,j,bi,bj) = 1 |
352 |
|
|
#else /* SEAICE_BICE_STRESS */ |
353 |
|
|
IF (klowc(i,j,bi,bj) .LT. 2) THEN |
354 |
|
|
KGEO(i,j,bi,bj) = 1 |
355 |
|
|
ELSE |
356 |
|
|
KGEO(i,j,bi,bj) = 2 |
357 |
|
|
DO WHILE ( abs(rC(KGEO(i,j,bi,bj))) .LT. 50.0 _d 0 .AND. |
358 |
|
|
& KGEO(i,j,bi,bj) .LT. (klowc(i,j,bi,bj)-1) ) |
359 |
|
|
KGEO(i,j,bi,bj) = KGEO(i,j,bi,bj) + 1 |
360 |
|
|
ENDDO |
361 |
|
|
ENDIF |
362 |
|
|
#endif /* SEAICE_BICE_STRESS */ |
363 |
|
|
ENDDO |
364 |
|
|
ENDDO |
365 |
|
|
ENDDO |
366 |
|
|
ENDDO |
367 |
|
|
#endif /* SEAICE_CGRID */ |
368 |
|
|
|
369 |
mlosch |
1.21 |
#ifdef SEAICE_ALLOW_JFNK |
370 |
|
|
C initialise some diagnostic counters for the JFNK solver |
371 |
|
|
totalNewtonIters = 0 |
372 |
|
|
totalNewtonFails = 0 |
373 |
|
|
totalKrylovIters = 0 |
374 |
|
|
totalKrylovFails = 0 |
375 |
|
|
totalJFNKtimeSteps = 0 |
376 |
mlosch |
1.22 |
C this cannot be done here, because globalArea is only defined |
377 |
|
|
C after S/R PACKAGES_INIT_FIXED, so we move it to S/R SEAICE_INIT_VARIA |
378 |
|
|
CML CALL SEAICE_MAP2VEC( nVec, rAw, rAs, |
379 |
|
|
CML & scalarProductMetric, .TRUE., myThid ) |
380 |
|
|
CML DO bj=myByLo(myThid),myByHi(myThid) |
381 |
|
|
CML DO bi=myBxLo(myThid),myBxHi(myThid) |
382 |
|
|
CML DO i=1,nVec |
383 |
|
|
CML scalarProductMetric(i,1,bi,bj) = |
384 |
|
|
CML & scalarProductMetric(i,1,bi,bj)/globalArea |
385 |
|
|
CML ENDDO |
386 |
|
|
CML ENDDO |
387 |
|
|
CML ENDDO |
388 |
mlosch |
1.21 |
#endif /* SEAICE_ALLOW_JFNK */ |
389 |
|
|
|
390 |
mlosch |
1.4 |
#ifdef ALLOW_DIAGNOSTICS |
391 |
|
|
IF ( useDiagnostics ) THEN |
392 |
|
|
CALL SEAICE_DIAGNOSTICS_INIT( myThid ) |
393 |
|
|
ENDIF |
394 |
|
|
#endif |
395 |
|
|
|
396 |
gforget |
1.13 |
C-- Summarise pkg/seaice configuration |
397 |
|
|
CALL SEAICE_SUMMARY( myThid ) |
398 |
jmc |
1.18 |
|
399 |
heimbach |
1.1 |
RETURN |
400 |
|
|
END |