1 |
dimitri |
1.1 |
C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_fixed.F,v 1.19 2012/03/11 13:41:38 jmc 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 |
|
|
C i,j,k,bi,bj - Loop counters |
33 |
|
|
|
34 |
|
|
INTEGER i, j, bi, bj |
35 |
|
|
INTEGER kSurface |
36 |
|
|
#ifdef ALLOW_SITRACER |
37 |
|
|
INTEGER iTracer |
38 |
|
|
#endif |
39 |
|
|
#ifndef SEAICE_CGRID |
40 |
|
|
_RS mask_uice |
41 |
|
|
#endif |
42 |
|
|
#ifdef SHORTWAVE_HEATING |
43 |
|
|
cif Helper variable for determining the fraction of sw radiation |
44 |
|
|
cif penetrating the model shallowest layer |
45 |
|
|
_RL dummyTime |
46 |
|
|
_RL swfracba(2) |
47 |
|
|
_RL tmpFac |
48 |
|
|
#endif /* SHORTWAVE_HEATING */ |
49 |
|
|
|
50 |
|
|
IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN |
51 |
|
|
kSurface = Nr |
52 |
|
|
ELSE |
53 |
|
|
kSurface = 1 |
54 |
|
|
ENDIF |
55 |
|
|
|
56 |
|
|
C Initialize MNC variable information for SEAICE |
57 |
|
|
IF ( useMNC .AND. |
58 |
|
|
& (seaice_tave_mnc.OR.seaice_dump_mnc.OR.SEAICE_mon_mnc) |
59 |
|
|
& ) THEN |
60 |
|
|
CALL SEAICE_MNC_INIT( myThid ) |
61 |
|
|
ENDIF |
62 |
|
|
|
63 |
|
|
_BEGIN_MASTER(myThid) |
64 |
|
|
#ifdef SHORTWAVE_HEATING |
65 |
|
|
tmpFac = -1.0 |
66 |
|
|
dummyTime = 1.0 |
67 |
|
|
swfracba(1) = ABS(rF(1)) |
68 |
|
|
swfracba(2) = ABS(rF(2)) |
69 |
|
|
CALL SWFRAC( |
70 |
|
|
I 2, tmpFac, |
71 |
|
|
U swfracba, |
72 |
|
|
I dummyTime, 0, myThid ) |
73 |
|
|
SWFracB = swfracba(2) |
74 |
|
|
#else /* SHORTWAVE_HEATING */ |
75 |
|
|
SWFracB = 0. _d 0 |
76 |
|
|
#endif /* SHORTWAVE_HEATING */ |
77 |
|
|
_END_MASTER(myThid) |
78 |
|
|
|
79 |
|
|
C-- Set mcPheePiston coeff (if still unset) |
80 |
|
|
_BEGIN_MASTER(myThid) |
81 |
|
|
IF ( SEAICE_mcPheePiston.EQ.UNSET_RL ) THEN |
82 |
|
|
IF ( SEAICE_availHeatFrac.NE.UNSET_RL ) THEN |
83 |
|
|
SEAICE_mcPheePiston = SEAICE_availHeatFrac |
84 |
|
|
& * drF(kSurface)/SEAICE_deltaTtherm |
85 |
|
|
ELSE |
86 |
|
|
SEAICE_mcPheePiston = MCPHEE_TAPER_FAC |
87 |
|
|
& * STANTON_NUMBER * USTAR_BASE |
88 |
|
|
SEAICE_mcPheePiston = MIN( SEAICE_mcPheePiston, |
89 |
|
|
& drF(kSurface)/SEAICE_deltaTtherm ) |
90 |
|
|
ENDIF |
91 |
|
|
ENDIF |
92 |
|
|
_END_MASTER(myThid) |
93 |
|
|
|
94 |
|
|
C-- SItracer specifications for basic tracers |
95 |
|
|
#ifdef ALLOW_SITRACER |
96 |
|
|
_BEGIN_MASTER(myThid) |
97 |
|
|
DO iTracer = 1, SItrNumInUse |
98 |
|
|
C "ice concentration" tracer that should remain .EQ.1. |
99 |
|
|
IF (SItrName(iTracer).EQ.'one') THEN |
100 |
|
|
SItrFromOcean0(iTracer) =ONE |
101 |
|
|
SItrFromFlood0(iTracer) =ONE |
102 |
|
|
SItrExpand0(iTracer) =ONE |
103 |
|
|
SItrFromOceanFrac(iTracer) =ZERO |
104 |
|
|
SItrFromFloodFrac(iTracer) =ZERO |
105 |
|
|
ENDIF |
106 |
|
|
C age tracer: no age in ocean, or effect from ice cover changes |
107 |
|
|
IF (SItrName(iTracer).EQ.'age') THEN |
108 |
|
|
SItrFromOcean0(iTracer) =ZERO |
109 |
|
|
SItrFromFlood0(iTracer) =ZERO |
110 |
|
|
SItrExpand0(iTracer) =ZERO |
111 |
|
|
SItrFromOceanFrac(iTracer) =ZERO |
112 |
|
|
SItrFromFloodFrac(iTracer) =ZERO |
113 |
|
|
ENDIf |
114 |
|
|
C salinity tracer: |
115 |
|
|
IF (SItrName(iTracer).EQ.'salinity') THEN |
116 |
|
|
SItrMate(iTracer) ='HEFF' |
117 |
|
|
SItrExpand0(iTracer) =ZERO |
118 |
|
|
IF ( SEAICE_salinityTracer ) THEN |
119 |
|
|
SEAICE_salt0 = ZERO |
120 |
|
|
SEAICE_saltFrac = ZERO |
121 |
|
|
ENDIF |
122 |
|
|
ENDIF |
123 |
|
|
C simple, made up, ice surface roughness index prototype |
124 |
|
|
IF (SItrName(iTracer).EQ.'ridge') THEN |
125 |
|
|
SItrMate(iTracer) ='AREA' |
126 |
|
|
SItrFromOcean0(iTracer) =ZERO |
127 |
|
|
SItrFromFlood0(iTracer) =ZERO |
128 |
|
|
SItrExpand0(iTracer) =ZERO |
129 |
|
|
SItrFromOceanFrac(iTracer) =ZERO |
130 |
|
|
SItrFromFloodFrac(iTracer) =ZERO |
131 |
|
|
ENDIF |
132 |
|
|
ENDDO |
133 |
|
|
_END_MASTER(myThid) |
134 |
|
|
#endif |
135 |
|
|
|
136 |
|
|
C-- all threads wait for master to finish initialisation of shared params |
137 |
|
|
_BARRIER |
138 |
|
|
|
139 |
|
|
C-- Initialize grid info |
140 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
141 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
142 |
|
|
DO j=1-OLy,sNy+OLy |
143 |
|
|
DO i=1-OLx,sNx+OLx |
144 |
|
|
HEFFM(i,j,bi,bj) = 0. _d 0 |
145 |
|
|
ENDDO |
146 |
|
|
ENDDO |
147 |
|
|
DO j=1-OLy,sNy+OLy |
148 |
|
|
DO i=1-OLx,sNx+OLx |
149 |
|
|
HEFFM(i,j,bi,bj)= 1. _d 0 |
150 |
|
|
IF (_hFacC(i,j,kSurface,bi,bj).eq.0.) |
151 |
|
|
& HEFFM(i,j,bi,bj)= 0. _d 0 |
152 |
|
|
ENDDO |
153 |
|
|
ENDDO |
154 |
|
|
#ifndef SEAICE_CGRID |
155 |
|
|
DO j=1-OLy+1,sNy+OLy |
156 |
|
|
DO i=1-OLx+1,sNx+OLx |
157 |
|
|
UVM(i,j,bi,bj)=0. _d 0 |
158 |
|
|
mask_uice=HEFFM(i,j, bi,bj)+HEFFM(i-1,j-1,bi,bj) |
159 |
|
|
& +HEFFM(i,j-1,bi,bj)+HEFFM(i-1,j, bi,bj) |
160 |
|
|
IF(mask_uice.GT.3.5 _d 0) UVM(i,j,bi,bj)=1. _d 0 |
161 |
|
|
ENDDO |
162 |
|
|
ENDDO |
163 |
|
|
#endif /* SEAICE_CGRID */ |
164 |
|
|
ENDDO |
165 |
|
|
ENDDO |
166 |
|
|
|
167 |
|
|
C coefficients for metric terms |
168 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
169 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
170 |
|
|
#ifdef SEAICE_CGRID |
171 |
|
|
DO j=1-OLy,sNy+OLy |
172 |
|
|
DO i=1-OLx,sNx+OLx |
173 |
|
|
k1AtC(I,J,bi,bj) = 0.0 _d 0 |
174 |
|
|
k1AtZ(I,J,bi,bj) = 0.0 _d 0 |
175 |
|
|
k2AtC(I,J,bi,bj) = 0.0 _d 0 |
176 |
|
|
k2AtZ(I,J,bi,bj) = 0.0 _d 0 |
177 |
|
|
ENDDO |
178 |
|
|
ENDDO |
179 |
|
|
IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN |
180 |
|
|
C This is the only case where tan(phi) is not zero. In this case |
181 |
|
|
C C and U points, and Z and V points have the same phi, so that we |
182 |
|
|
C only need a copy here. Do not use tan(YC) and tan(YG), because these |
183 |
|
|
C can be the geographical coordinates and not the correct grid |
184 |
|
|
C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0) |
185 |
|
|
DO j=1-OLy,sNy+OLy |
186 |
|
|
DO i=1-OLx,sNx+OLx |
187 |
|
|
k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
188 |
|
|
k2AtZ(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere |
189 |
|
|
ENDDO |
190 |
|
|
ENDDO |
191 |
|
|
ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN |
192 |
|
|
C compute metric term coefficients from finite difference approximation |
193 |
|
|
DO j=1-OLy,sNy+OLy |
194 |
|
|
DO i=1-OLx,sNx+OLx-1 |
195 |
|
|
k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj) |
196 |
|
|
& * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) ) |
197 |
|
|
& * _recip_dxF(I,J,bi,bj) |
198 |
|
|
ENDDO |
199 |
|
|
ENDDO |
200 |
|
|
DO j=1-OLy,sNy+OLy |
201 |
|
|
DO i=1-OLx+1,sNx+OLx |
202 |
|
|
k1AtZ(I,J,bi,bj) = _recip_dyU(I,J,bi,bj) |
203 |
|
|
& * ( _dyC(I,J,bi,bj) - _dyC(I-1,J,bi,bj) ) |
204 |
|
|
& * _recip_dxV(I,J,bi,bj) |
205 |
|
|
ENDDO |
206 |
|
|
ENDDO |
207 |
|
|
DO j=1-OLy,sNy+OLy-1 |
208 |
|
|
DO i=1-OLx,sNx+OLx |
209 |
|
|
k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj) |
210 |
|
|
& * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) ) |
211 |
|
|
& * _recip_dyF(I,J,bi,bj) |
212 |
|
|
ENDDO |
213 |
|
|
ENDDO |
214 |
|
|
DO j=1-OLy+1,sNy+OLy |
215 |
|
|
DO i=1-OLx,sNx+OLx |
216 |
|
|
k2AtZ(I,J,bi,bj) = _recip_dxV(I,J,bi,bj) |
217 |
|
|
& * ( _dxC(I,J,bi,bj) - _dxC(I,J-1,bi,bj) ) |
218 |
|
|
& * _recip_dyU(I,J,bi,bj) |
219 |
|
|
ENDDO |
220 |
|
|
ENDDO |
221 |
|
|
ENDIF |
222 |
|
|
#else /* not SEAICE_CGRID */ |
223 |
|
|
DO j=1-OLy,sNy+OLy |
224 |
|
|
DO i=1-OLx,sNx+OLx |
225 |
|
|
k1AtC(I,J,bi,bj) = 0.0 _d 0 |
226 |
|
|
k1AtU(I,J,bi,bj) = 0.0 _d 0 |
227 |
|
|
k1AtV(I,J,bi,bj) = 0.0 _d 0 |
228 |
|
|
k2AtC(I,J,bi,bj) = 0.0 _d 0 |
229 |
|
|
k2AtU(I,J,bi,bj) = 0.0 _d 0 |
230 |
|
|
k2AtV(I,J,bi,bj) = 0.0 _d 0 |
231 |
|
|
ENDDO |
232 |
|
|
ENDDO |
233 |
|
|
IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN |
234 |
|
|
C This is the only case where tan(phi) is not zero. In this case |
235 |
|
|
C C and U points, and Z and V points have the same phi, so that we |
236 |
|
|
C only need a copy here. Do not use tan(YC) and tan(YG), because these |
237 |
|
|
C can be the geographical coordinates and not the correct grid |
238 |
|
|
C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0) |
239 |
|
|
DO j=1-OLy,sNy+OLy |
240 |
|
|
DO i=1-OLx,sNx+OLx |
241 |
|
|
k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
242 |
|
|
k2AtU(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
243 |
|
|
k2AtV(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere |
244 |
|
|
ENDDO |
245 |
|
|
ENDDO |
246 |
|
|
ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN |
247 |
|
|
C compute metric term coefficients from finite difference approximation |
248 |
|
|
DO j=1-OLy,sNy+OLy |
249 |
|
|
DO i=1-OLx,sNx+OLx-1 |
250 |
|
|
k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj) |
251 |
|
|
& * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) ) |
252 |
|
|
& * _recip_dxF(I,J,bi,bj) |
253 |
|
|
ENDDO |
254 |
|
|
ENDDO |
255 |
|
|
DO j=1-OLy,sNy+OLy |
256 |
|
|
DO i=1-OLx+1,sNx+OLx |
257 |
|
|
k1AtU(I,J,bi,bj) = _recip_dyG(I,J,bi,bj) |
258 |
|
|
& * ( _dyF(I,J,bi,bj) - _dyF(I-1,J,bi,bj) ) |
259 |
|
|
& * _recip_dxC(I,J,bi,bj) |
260 |
|
|
ENDDO |
261 |
|
|
ENDDO |
262 |
|
|
DO j=1-OLy,sNy+OLy |
263 |
|
|
DO i=1-OLx,sNx+OLx-1 |
264 |
|
|
k1AtV(I,J,bi,bj) = _recip_dyC(I,J,bi,bj) |
265 |
|
|
& * ( _dyU(I+1,J,bi,bj) - _dyU(I,J,bi,bj) ) |
266 |
|
|
& * _recip_dxG(I,J,bi,bj) |
267 |
|
|
ENDDO |
268 |
|
|
ENDDO |
269 |
|
|
DO j=1-OLy,sNy+OLy-1 |
270 |
|
|
DO i=1-OLx,sNx+OLx |
271 |
|
|
k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj) |
272 |
|
|
& * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) ) |
273 |
|
|
& * _recip_dyF(I,J,bi,bj) |
274 |
|
|
ENDDO |
275 |
|
|
ENDDO |
276 |
|
|
DO j=1-OLy,sNy+OLy-1 |
277 |
|
|
DO i=1-OLx,sNx+OLx |
278 |
|
|
k2AtU(I,J,bi,bj) = _recip_dxC(I,J,bi,bj) |
279 |
|
|
& * ( _dxV(I,J+1,bi,bj) - _dxV(I,J,bi,bj) ) |
280 |
|
|
& * _recip_dyG(I,J,bi,bj) |
281 |
|
|
ENDDO |
282 |
|
|
ENDDO |
283 |
|
|
DO j=1-OLy+1,sNy+OLy |
284 |
|
|
DO i=1-OLx,sNx+OLx |
285 |
|
|
k2AtV(I,J,bi,bj) = _recip_dxG(I,J,bi,bj) |
286 |
|
|
& * ( _dxF(I,J,bi,bj) - _dxF(I,J-1,bi,bj) ) |
287 |
|
|
& * _recip_dyC(I,J,bi,bj) |
288 |
|
|
ENDDO |
289 |
|
|
ENDDO |
290 |
|
|
ENDIF |
291 |
|
|
#endif /* not SEAICE_CGRID */ |
292 |
|
|
ENDDO |
293 |
|
|
ENDDO |
294 |
|
|
|
295 |
|
|
#ifndef SEAICE_CGRID |
296 |
|
|
C-- Choose a proxy level for geostrophic velocity, |
297 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
298 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
299 |
|
|
DO j=1-OLy,sNy+OLy |
300 |
|
|
DO i=1-OLx,sNx+OLx |
301 |
|
|
KGEO(i,j,bi,bj) = 0 |
302 |
|
|
ENDDO |
303 |
|
|
ENDDO |
304 |
|
|
DO j=1-OLy,sNy+OLy |
305 |
|
|
DO i=1-OLx,sNx+OLx |
306 |
|
|
#ifdef SEAICE_BICE_STRESS |
307 |
|
|
KGEO(i,j,bi,bj) = 1 |
308 |
|
|
#else /* SEAICE_BICE_STRESS */ |
309 |
|
|
IF (klowc(i,j,bi,bj) .LT. 2) THEN |
310 |
|
|
KGEO(i,j,bi,bj) = 1 |
311 |
|
|
ELSE |
312 |
|
|
KGEO(i,j,bi,bj) = 2 |
313 |
|
|
DO WHILE ( abs(rC(KGEO(i,j,bi,bj))) .LT. 50.0 _d 0 .AND. |
314 |
|
|
& KGEO(i,j,bi,bj) .LT. (klowc(i,j,bi,bj)-1) ) |
315 |
|
|
KGEO(i,j,bi,bj) = KGEO(i,j,bi,bj) + 1 |
316 |
|
|
ENDDO |
317 |
|
|
ENDIF |
318 |
|
|
#endif /* SEAICE_BICE_STRESS */ |
319 |
|
|
ENDDO |
320 |
|
|
ENDDO |
321 |
|
|
ENDDO |
322 |
|
|
ENDDO |
323 |
|
|
#endif /* SEAICE_CGRID */ |
324 |
|
|
|
325 |
|
|
#ifdef ALLOW_DIAGNOSTICS |
326 |
|
|
IF ( useDiagnostics ) THEN |
327 |
|
|
CALL SEAICE_DIAGNOSTICS_INIT( myThid ) |
328 |
|
|
ENDIF |
329 |
|
|
#endif |
330 |
|
|
|
331 |
|
|
C-- Summarise pkg/seaice configuration |
332 |
|
|
CALL SEAICE_SUMMARY( myThid ) |
333 |
|
|
|
334 |
|
|
RETURN |
335 |
|
|
END |