1 |
gforget |
1.15 |
C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_init_fixed.F,v 1.14 2012/02/09 03:42:32 gforget 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 |
|
|
#include "SEAICE_PARAMS.h" |
22 |
heimbach |
1.1 |
#include "SEAICE.h" |
23 |
|
|
|
24 |
|
|
C === Routine arguments === |
25 |
|
|
C myThid - Thread no. that called this routine. |
26 |
|
|
INTEGER myThid |
27 |
|
|
CEndOfInterface |
28 |
jmc |
1.6 |
|
29 |
heimbach |
1.1 |
C === Local variables === |
30 |
|
|
C i,j,k,bi,bj - Loop counters |
31 |
|
|
|
32 |
jmc |
1.7 |
INTEGER i, j, bi, bj |
33 |
mlosch |
1.4 |
INTEGER kSurface |
34 |
jmc |
1.6 |
#ifndef SEAICE_CGRID |
35 |
mlosch |
1.4 |
_RS mask_uice |
36 |
jmc |
1.6 |
#endif |
37 |
jmc |
1.12 |
#ifdef SHORTWAVE_HEATING |
38 |
heimbach |
1.3 |
cif Helper variable for determining the fraction of sw radiation |
39 |
jmc |
1.8 |
cif penetrating the model shallowest layer |
40 |
heimbach |
1.3 |
_RL dummyTime |
41 |
|
|
_RL swfracba(2) |
42 |
jmc |
1.12 |
_RL tmpFac |
43 |
|
|
#endif /* SHORTWAVE_HEATING */ |
44 |
heimbach |
1.1 |
|
45 |
mlosch |
1.4 |
IF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN |
46 |
|
|
kSurface = Nr |
47 |
|
|
ELSE |
48 |
|
|
kSurface = 1 |
49 |
heimbach |
1.1 |
ENDIF |
50 |
|
|
|
51 |
jmc |
1.6 |
C Initialize MNC variable information for SEAICE |
52 |
|
|
IF ( useMNC .AND. |
53 |
|
|
& (seaice_tave_mnc.OR.seaice_dump_mnc.OR.SEAICE_mon_mnc) |
54 |
|
|
& ) THEN |
55 |
|
|
CALL SEAICE_MNC_INIT( myThid ) |
56 |
|
|
ENDIF |
57 |
|
|
|
58 |
jmc |
1.12 |
_BEGIN_MASTER(myThid) |
59 |
heimbach |
1.3 |
#ifdef SHORTWAVE_HEATING |
60 |
jmc |
1.12 |
tmpFac = -1.0 |
61 |
|
|
dummyTime = 1.0 |
62 |
|
|
swfracba(1) = ABS(rF(1)) |
63 |
|
|
swfracba(2) = ABS(rF(2)) |
64 |
|
|
CALL SWFRAC( |
65 |
|
|
I 2, tmpFac, |
66 |
heimbach |
1.3 |
U swfracba, |
67 |
jmc |
1.12 |
I dummyTime, 0, myThid ) |
68 |
|
|
SWFracB = swfracba(2) |
69 |
|
|
#else /* SHORTWAVE_HEATING */ |
70 |
|
|
SWFracB = 0. _d 0 |
71 |
|
|
#endif /* SHORTWAVE_HEATING */ |
72 |
|
|
_END_MASTER(myThid) |
73 |
gforget |
1.13 |
|
74 |
|
|
|
75 |
|
|
C-- efficiency of ocean-ice turbulent flux ... |
76 |
|
|
if (SEAICEturbFluxFormula.EQ.1) then |
77 |
|
|
c ... can be specified as fractions between 0 and 1 |
78 |
|
|
IF ( SEAICE_availHeatFrac .EQ. UNSET_RL ) |
79 |
|
|
& SEAICE_availHeatFrac = ONE |
80 |
|
|
IF ( SEAICE_availHeatFracFrz .EQ. UNSET_RL ) |
81 |
|
|
& SEAICE_availHeatFracFrz = SEAICE_availHeatFrac |
82 |
|
|
elseif (SEAICEturbFluxFormula.EQ.2) then |
83 |
|
|
c ... can be specified as time scales (>SEAICE_deltaTtherm) |
84 |
|
|
IF ( SEAICE_gamma_t .EQ. UNSET_RL ) |
85 |
|
|
& SEAICE_gamma_t=SEAICE_deltaTtherm |
86 |
|
|
IF ( SEAICE_gamma_t_frz .EQ. UNSET_RL ) |
87 |
|
|
& SEAICE_gamma_t_frz=SEAICE_gamma_t |
88 |
|
|
IF ( SEAICE_gamma_t .LE. SEAICE_deltaTtherm ) THEN |
89 |
|
|
SEAICE_availHeatFracFrz = 1. _d 0 |
90 |
|
|
ELSE |
91 |
|
|
SEAICE_availHeatFrac = SEAICE_deltaTtherm/SEAICE_gamma_t |
92 |
|
|
ENDIF |
93 |
|
|
IF ( SEAICE_gamma_t_frz .LE. SEAICE_deltaTtherm ) THEN |
94 |
|
|
SEAICE_availHeatFracFrz = 1. _d 0 |
95 |
|
|
ELSE |
96 |
|
|
SEAICE_availHeatFracFrz = |
97 |
|
|
& SEAICE_deltaTtherm/SEAICE_gamma_t_frz |
98 |
|
|
ENDIF |
99 |
|
|
elseif ( (SEAICEturbFluxFormula.EQ.3).OR. |
100 |
|
|
& (SEAICEturbFluxFormula.EQ.4) ) then |
101 |
|
|
c ... is hard-coded (after McPhee) |
102 |
|
|
SEAICE_availHeatFrac= MCPHEE_TAPER_FAC * STANTON_NUMBER * |
103 |
|
|
& USTAR_BASE / dRf(kSurface) * SEAICE_deltaTtherm |
104 |
|
|
SEAICE_availHeatFracFrz=ZERO |
105 |
|
|
endif |
106 |
gforget |
1.15 |
C-- tapering of ocean-ice turbulent flux as AREA increases |
107 |
|
|
if ( (SEAICEturbFluxFormula.EQ.3).OR. |
108 |
|
|
& (SEAICEturbFluxFormula.EQ.4) ) then |
109 |
|
|
c hard-coded (after McPhee) |
110 |
|
|
SEAICE_availHeatTaper = |
111 |
|
|
& (MCPHEE_TAPER_FAC-ONE)/MCPHEE_TAPER_FAC |
112 |
|
|
elseif (SEAICE_availHeatTaper.EQ.UNSET_RL) then |
113 |
|
|
SEAICE_availHeatTaper = ZERO |
114 |
|
|
endif |
115 |
heimbach |
1.3 |
|
116 |
gforget |
1.14 |
C-- convert SEAICE_doOpenWaterGrowth/Melt logical switch to numerical facOpenGrow/Melt |
117 |
|
|
facOpenGrow=0. _d 0 |
118 |
|
|
facOpenMelt=0. _d 0 |
119 |
|
|
if (SEAICE_doOpenWaterGrowth) facOpenGrow=1. _d 0 |
120 |
|
|
if (SEAICE_doOpenWaterMelt) facOpenMelt=1. _d 0 |
121 |
|
|
|
122 |
mlosch |
1.4 |
C-- Initialize grid info |
123 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
124 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
125 |
|
|
DO j=1-OLy,sNy+OLy |
126 |
|
|
DO i=1-OLx,sNx+OLx |
127 |
|
|
HEFFM(i,j,bi,bj) = 0. _d 0 |
128 |
|
|
ENDDO |
129 |
|
|
ENDDO |
130 |
|
|
DO j=1-OLy,sNy+OLy |
131 |
|
|
DO i=1-OLx,sNx+OLx |
132 |
|
|
HEFFM(i,j,bi,bj)= 1. _d 0 |
133 |
|
|
IF (_hFacC(i,j,kSurface,bi,bj).eq.0.) |
134 |
|
|
& HEFFM(i,j,bi,bj)= 0. _d 0 |
135 |
|
|
ENDDO |
136 |
|
|
ENDDO |
137 |
jmc |
1.11 |
#ifndef SEAICE_CGRID |
138 |
mlosch |
1.4 |
DO j=1-OLy+1,sNy+OLy |
139 |
|
|
DO i=1-OLx+1,sNx+OLx |
140 |
|
|
UVM(i,j,bi,bj)=0. _d 0 |
141 |
|
|
mask_uice=HEFFM(i,j, bi,bj)+HEFFM(i-1,j-1,bi,bj) |
142 |
|
|
& +HEFFM(i,j-1,bi,bj)+HEFFM(i-1,j, bi,bj) |
143 |
|
|
IF(mask_uice.GT.3.5 _d 0) UVM(i,j,bi,bj)=1. _d 0 |
144 |
|
|
ENDDO |
145 |
|
|
ENDDO |
146 |
jmc |
1.11 |
#endif /* SEAICE_CGRID */ |
147 |
mlosch |
1.4 |
ENDDO |
148 |
|
|
ENDDO |
149 |
|
|
|
150 |
jmc |
1.6 |
C coefficients for metric terms |
151 |
mlosch |
1.4 |
DO bj=myByLo(myThid),myByHi(myThid) |
152 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
153 |
mlosch |
1.10 |
#ifdef SEAICE_CGRID |
154 |
mlosch |
1.4 |
DO j=1-OLy,sNy+OLy |
155 |
|
|
DO i=1-OLx,sNx+OLx |
156 |
|
|
k1AtC(I,J,bi,bj) = 0.0 _d 0 |
157 |
|
|
k1AtZ(I,J,bi,bj) = 0.0 _d 0 |
158 |
|
|
k2AtC(I,J,bi,bj) = 0.0 _d 0 |
159 |
|
|
k2AtZ(I,J,bi,bj) = 0.0 _d 0 |
160 |
|
|
ENDDO |
161 |
|
|
ENDDO |
162 |
|
|
IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN |
163 |
jmc |
1.6 |
C This is the only case where tan(phi) is not zero. In this case |
164 |
|
|
C C and U points, and Z and V points have the same phi, so that we |
165 |
mlosch |
1.4 |
C only need a copy here. Do not use tan(YC) and tan(YG), because these |
166 |
jmc |
1.6 |
C can be the geographical coordinates and not the correct grid |
167 |
mlosch |
1.4 |
C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0) |
168 |
|
|
DO j=1-OLy,sNy+OLy |
169 |
|
|
DO i=1-OLx,sNx+OLx |
170 |
|
|
k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
171 |
|
|
k2AtZ(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere |
172 |
|
|
ENDDO |
173 |
|
|
ENDDO |
174 |
|
|
ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN |
175 |
|
|
C compute metric term coefficients from finite difference approximation |
176 |
|
|
DO j=1-OLy,sNy+OLy |
177 |
|
|
DO i=1-OLx,sNx+OLx-1 |
178 |
|
|
k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj) |
179 |
|
|
& * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) ) |
180 |
|
|
& * _recip_dxF(I,J,bi,bj) |
181 |
|
|
ENDDO |
182 |
|
|
ENDDO |
183 |
|
|
DO j=1-OLy,sNy+OLy |
184 |
|
|
DO i=1-OLx+1,sNx+OLx |
185 |
jmc |
1.6 |
k1AtZ(I,J,bi,bj) = _recip_dyU(I,J,bi,bj) |
186 |
mlosch |
1.4 |
& * ( _dyC(I,J,bi,bj) - _dyC(I-1,J,bi,bj) ) |
187 |
|
|
& * _recip_dxV(I,J,bi,bj) |
188 |
|
|
ENDDO |
189 |
|
|
ENDDO |
190 |
|
|
DO j=1-OLy,sNy+OLy-1 |
191 |
|
|
DO i=1-OLx,sNx+OLx |
192 |
jmc |
1.6 |
k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj) |
193 |
mlosch |
1.4 |
& * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) ) |
194 |
|
|
& * _recip_dyF(I,J,bi,bj) |
195 |
|
|
ENDDO |
196 |
|
|
ENDDO |
197 |
|
|
DO j=1-OLy+1,sNy+OLy |
198 |
|
|
DO i=1-OLx,sNx+OLx |
199 |
mlosch |
1.9 |
k2AtZ(I,J,bi,bj) = _recip_dxV(I,J,bi,bj) |
200 |
mlosch |
1.4 |
& * ( _dxC(I,J,bi,bj) - _dxC(I,J-1,bi,bj) ) |
201 |
|
|
& * _recip_dyU(I,J,bi,bj) |
202 |
|
|
ENDDO |
203 |
|
|
ENDDO |
204 |
|
|
ENDIF |
205 |
mlosch |
1.10 |
#else /* not SEAICE_CGRID */ |
206 |
|
|
DO j=1-OLy,sNy+OLy |
207 |
|
|
DO i=1-OLx,sNx+OLx |
208 |
|
|
k1AtC(I,J,bi,bj) = 0.0 _d 0 |
209 |
|
|
k1AtU(I,J,bi,bj) = 0.0 _d 0 |
210 |
|
|
k1AtV(I,J,bi,bj) = 0.0 _d 0 |
211 |
|
|
k2AtC(I,J,bi,bj) = 0.0 _d 0 |
212 |
|
|
k2AtU(I,J,bi,bj) = 0.0 _d 0 |
213 |
|
|
k2AtV(I,J,bi,bj) = 0.0 _d 0 |
214 |
|
|
ENDDO |
215 |
|
|
ENDDO |
216 |
|
|
IF ( usingSphericalPolarGrid .AND. SEAICEuseMetricTerms ) THEN |
217 |
|
|
C This is the only case where tan(phi) is not zero. In this case |
218 |
|
|
C C and U points, and Z and V points have the same phi, so that we |
219 |
|
|
C only need a copy here. Do not use tan(YC) and tan(YG), because these |
220 |
|
|
C can be the geographical coordinates and not the correct grid |
221 |
|
|
C coordinates when the grid is rotated (phi/theta/psiEuler .NE. 0) |
222 |
|
|
DO j=1-OLy,sNy+OLy |
223 |
|
|
DO i=1-OLx,sNx+OLx |
224 |
|
|
k2AtC(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
225 |
|
|
k2AtU(I,J,bi,bj) = - _tanPhiAtU(I,J,bi,bj)*recip_rSphere |
226 |
|
|
k2AtV(I,J,bi,bj) = - _tanPhiAtV(I,J,bi,bj)*recip_rSphere |
227 |
|
|
ENDDO |
228 |
|
|
ENDDO |
229 |
|
|
ELSEIF ( usingCurvilinearGrid .AND. SEAICEuseMetricTerms ) THEN |
230 |
|
|
C compute metric term coefficients from finite difference approximation |
231 |
|
|
DO j=1-OLy,sNy+OLy |
232 |
|
|
DO i=1-OLx,sNx+OLx-1 |
233 |
|
|
k1AtC(I,J,bi,bj) = _recip_dyF(I,J,bi,bj) |
234 |
|
|
& * ( _dyG(I+1,J,bi,bj) - _dyG(I,J,bi,bj) ) |
235 |
|
|
& * _recip_dxF(I,J,bi,bj) |
236 |
|
|
ENDDO |
237 |
|
|
ENDDO |
238 |
|
|
DO j=1-OLy,sNy+OLy |
239 |
|
|
DO i=1-OLx+1,sNx+OLx |
240 |
|
|
k1AtU(I,J,bi,bj) = _recip_dyG(I,J,bi,bj) |
241 |
|
|
& * ( _dyF(I,J,bi,bj) - _dyF(I-1,J,bi,bj) ) |
242 |
|
|
& * _recip_dxC(I,J,bi,bj) |
243 |
|
|
ENDDO |
244 |
|
|
ENDDO |
245 |
|
|
DO j=1-OLy,sNy+OLy |
246 |
|
|
DO i=1-OLx,sNx+OLx-1 |
247 |
|
|
k1AtV(I,J,bi,bj) = _recip_dyC(I,J,bi,bj) |
248 |
|
|
& * ( _dyU(I+1,J,bi,bj) - _dyU(I,J,bi,bj) ) |
249 |
|
|
& * _recip_dxG(I,J,bi,bj) |
250 |
|
|
ENDDO |
251 |
|
|
ENDDO |
252 |
|
|
DO j=1-OLy,sNy+OLy-1 |
253 |
|
|
DO i=1-OLx,sNx+OLx |
254 |
|
|
k2AtC(I,J,bi,bj) = _recip_dxF(I,J,bi,bj) |
255 |
|
|
& * ( _dxG(I,J+1,bi,bj) - _dxG(I,J,bi,bj) ) |
256 |
|
|
& * _recip_dyF(I,J,bi,bj) |
257 |
|
|
ENDDO |
258 |
|
|
ENDDO |
259 |
|
|
DO j=1-OLy,sNy+OLy-1 |
260 |
|
|
DO i=1-OLx,sNx+OLx |
261 |
|
|
k2AtU(I,J,bi,bj) = _recip_dxC(I,J,bi,bj) |
262 |
|
|
& * ( _dxV(I,J+1,bi,bj) - _dxV(I,J,bi,bj) ) |
263 |
|
|
& * _recip_dyG(I,J,bi,bj) |
264 |
|
|
ENDDO |
265 |
|
|
ENDDO |
266 |
|
|
DO j=1-OLy+1,sNy+OLy |
267 |
|
|
DO i=1-OLx,sNx+OLx |
268 |
|
|
k2AtV(I,J,bi,bj) = _recip_dxG(I,J,bi,bj) |
269 |
|
|
& * ( _dxF(I,J,bi,bj) - _dxF(I,J-1,bi,bj) ) |
270 |
|
|
& * _recip_dyC(I,J,bi,bj) |
271 |
|
|
ENDDO |
272 |
|
|
ENDDO |
273 |
|
|
ENDIF |
274 |
|
|
#endif /* not SEAICE_CGRID */ |
275 |
mlosch |
1.4 |
ENDDO |
276 |
|
|
ENDDO |
277 |
|
|
|
278 |
|
|
#ifndef SEAICE_CGRID |
279 |
|
|
C-- Choose a proxy level for geostrophic velocity, |
280 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
281 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
282 |
|
|
DO j=1-OLy,sNy+OLy |
283 |
|
|
DO i=1-OLx,sNx+OLx |
284 |
|
|
KGEO(i,j,bi,bj) = 0 |
285 |
|
|
ENDDO |
286 |
|
|
ENDDO |
287 |
|
|
DO j=1-OLy,sNy+OLy |
288 |
|
|
DO i=1-OLx,sNx+OLx |
289 |
|
|
#ifdef SEAICE_BICE_STRESS |
290 |
|
|
KGEO(i,j,bi,bj) = 1 |
291 |
|
|
#else /* SEAICE_BICE_STRESS */ |
292 |
|
|
IF (klowc(i,j,bi,bj) .LT. 2) THEN |
293 |
|
|
KGEO(i,j,bi,bj) = 1 |
294 |
|
|
ELSE |
295 |
|
|
KGEO(i,j,bi,bj) = 2 |
296 |
|
|
DO WHILE ( abs(rC(KGEO(i,j,bi,bj))) .LT. 50.0 _d 0 .AND. |
297 |
|
|
& KGEO(i,j,bi,bj) .LT. (klowc(i,j,bi,bj)-1) ) |
298 |
|
|
KGEO(i,j,bi,bj) = KGEO(i,j,bi,bj) + 1 |
299 |
|
|
ENDDO |
300 |
|
|
ENDIF |
301 |
|
|
#endif /* SEAICE_BICE_STRESS */ |
302 |
|
|
ENDDO |
303 |
|
|
ENDDO |
304 |
|
|
ENDDO |
305 |
|
|
ENDDO |
306 |
|
|
#endif /* SEAICE_CGRID */ |
307 |
|
|
|
308 |
|
|
#ifdef ALLOW_DIAGNOSTICS |
309 |
|
|
IF ( useDiagnostics ) THEN |
310 |
|
|
CALL SEAICE_DIAGNOSTICS_INIT( myThid ) |
311 |
|
|
ENDIF |
312 |
|
|
#endif |
313 |
|
|
|
314 |
gforget |
1.13 |
C-- Summarise pkg/seaice configuration |
315 |
|
|
CALL SEAICE_SUMMARY( myThid ) |
316 |
|
|
|
317 |
heimbach |
1.1 |
RETURN |
318 |
|
|
END |