/[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.16 - (show annotations) (download)
Thu Feb 16 01:28:58 2012 UTC (12 years, 3 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint63j
Changes since 1.15: +50 -1 lines
- SItracer : added run time params, and allow coupling to ocn salinity.

-> seaice_tracer_phys.F :
	use new run time params
	remove some logic that is now better placed in init & check.
	when SEAICE_salinityTracer is set, which implies SIsal0&SIsalFrac=0,
	activate salinity tracer coupling with ocn (i.e. saltFlux & saltPlumeFlux)
	This feature used to work but is due for further re-testing.
-> seaice_readparms.F :
	init and read new run time params from SEAICE_TRACER.h
-> seaice_init_fixed.F :
	set SITR params accordingly for 'age', salinity', 'one', ridge'
	if SEAICE_salinityTracer is set, then reset to SIsal0&SIsalFrac=0.
-> seaice_check.F :
	stop if ALLOW_SITRACER is undef and SEAICE_salinityTracer or SEAICE_ageTracer run time is set
	prevent erroneous spec. of SItrFromOceanFrac/SItrFromFloodFrac (must be in [0. 1.])
	restrict use of SItrFromOceanFrac/SItrFromFloodFrac to salinity case (for now at least)
	+ a bit of re-organization

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

  ViewVC Help
Powered by ViewVC 1.1.22