1 |
C $Header: /u/gcmpack/MITgcm/pkg/thsice/thsice_advdiff.F,v 1.5 2007/04/16 22:38:24 heimbach Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "THSICE_OPTIONS.h" |
5 |
#ifdef ALLOW_GENERIC_ADVDIFF |
6 |
# include "GAD_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
CBOP |
10 |
C !ROUTINE: THSICE_ADVDIFF |
11 |
|
12 |
C !INTERFACE: ========================================================== |
13 |
SUBROUTINE THSICE_ADVDIFF( |
14 |
U uIce, vIce, |
15 |
I bi, bj, myTime, myIter, myThid ) |
16 |
|
17 |
C !DESCRIPTION: \bv |
18 |
C *===========================================================* |
19 |
C | SUBROUTINE THSICE_ADVDIFF |
20 |
C | o driver for different advection routines |
21 |
C | calls an adaption of gad_advection to call different |
22 |
C | advection routines of pkg/generic_advdiff |
23 |
C *===========================================================* |
24 |
C \ev |
25 |
|
26 |
C !USES: =============================================================== |
27 |
IMPLICIT NONE |
28 |
|
29 |
C === Global variables === |
30 |
C oceFWfx :: fresh water flux to the ocean [kg/m^2/s] |
31 |
C oceSflx :: salt flux to the ocean [psu.kg/m^2/s] (~g/m^2/s) |
32 |
C oceQnet :: heat flux to the ocean [W/m^2] |
33 |
|
34 |
#include "SIZE.h" |
35 |
#include "EEPARAMS.h" |
36 |
#include "PARAMS.h" |
37 |
#include "GRID.h" |
38 |
#include "THSICE_SIZE.h" |
39 |
#include "THSICE_PARAMS.h" |
40 |
#include "THSICE_VARS.h" |
41 |
#include "THSICE_2DYN.h" |
42 |
#ifdef ALLOW_GENERIC_ADVDIFF |
43 |
# include "GAD.h" |
44 |
#endif |
45 |
#ifdef ALLOW_AUTODIFF_TAMC |
46 |
# include "tamc.h" |
47 |
# include "tamc_keys.h" |
48 |
#endif |
49 |
|
50 |
C !INPUT PARAMETERS: =================================================== |
51 |
C === Routine arguments === |
52 |
C uIce/vIce :: ice velocity on C-grid [m/s] |
53 |
C bi,bj :: Tile indices |
54 |
C myTime :: Current time in simulation (s) |
55 |
C myIter :: Current iteration number |
56 |
C myThid :: My Thread Id. number |
57 |
_RL uIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
58 |
_RL vIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
59 |
INTEGER bi,bj |
60 |
_RL myTime |
61 |
INTEGER myIter |
62 |
INTEGER myThid |
63 |
CEndOfInterface |
64 |
|
65 |
#ifdef ALLOW_THSICE |
66 |
C !LOCAL VARIABLES: ==================================================== |
67 |
C === Local variables === |
68 |
C i,j, :: Loop counters |
69 |
C uTrans :: sea-ice area transport, x direction |
70 |
C vTrans :: sea-ice area transport, y direction |
71 |
C uTrIce :: sea-ice volume transport, x direction |
72 |
C vTrIce :: sea-ice volume transport, y direction |
73 |
C afx :: horizontal advective flux, x direction |
74 |
C afy :: horizontal advective flux, y direction |
75 |
C iceFrc :: (new) sea-ice fraction |
76 |
C iceVol :: temporary array used in advection S/R |
77 |
C oldVol :: (old) sea-ice volume |
78 |
C msgBuf :: Informational/error meesage buffer |
79 |
INTEGER i, j |
80 |
LOGICAL thSIce_multiDimAdv |
81 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
82 |
|
83 |
_RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
84 |
_RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
85 |
_RL uTrIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
86 |
_RL vTrIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
87 |
_RL afx (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
88 |
_RL afy (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
89 |
_RS maskOce (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
90 |
_RL iceFrc (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
91 |
_RL iceVol (1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
92 |
_RL oldVol (1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
93 |
_RL r_minArea |
94 |
_RL meanCellArea, areaEpsil, vol_Epsil |
95 |
#ifdef ALLOW_DIAGNOSTICS |
96 |
CHARACTER*8 diagName |
97 |
CHARACTER*4 THSICE_DIAG_SUFX, diagSufx |
98 |
EXTERNAL THSICE_DIAG_SUFX |
99 |
#endif |
100 |
#ifdef ALLOW_DBUG_THSICE |
101 |
_RL tmpVar, sumVar1, sumVar2 |
102 |
#endif |
103 |
LOGICAL dBugFlag |
104 |
#include "THSICE_DEBUG.h" |
105 |
CEOP |
106 |
|
107 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
108 |
|
109 |
#ifdef ALLOW_AUTODIFF_TAMC |
110 |
act1 = bi - myBxLo(myThid) |
111 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
112 |
act2 = bj - myByLo(myThid) |
113 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
114 |
act3 = myThid - 1 |
115 |
max3 = nTx*nTy |
116 |
act4 = ikey_dynamics - 1 |
117 |
iicekey = (act1 + 1) + act2*max1 |
118 |
& + act3*max1*max2 |
119 |
& + act4*max1*max2*max3 |
120 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
121 |
|
122 |
C areaEpsil, vol_Epsil are 2 small numbers for ice area & ice volume: |
123 |
C if ice area (=ice fraction * grid-cell area) or ice volume (= effective |
124 |
C thickness * grid-cell area) are too small (i.e.: < areaEpsil,vol_Epsil) |
125 |
C will assume that ice is gone, and will loose mass or energy. |
126 |
C However, if areaEpsil,vol_Epsil are much smaller than minimun ice area |
127 |
C (iceMaskMin*rAc) and minimum ice volume (iceMaskMin*hIceMin*rAc), |
128 |
C good chance that this will never happen within 1 time step. |
129 |
|
130 |
dBugFlag = debugLevel.GE.debLevB |
131 |
|
132 |
C- definitively not an accurate computation of mean grid-cell area; |
133 |
C but what matter here is just to have the right order of magnitude. |
134 |
meanCellArea = Nx*Ny |
135 |
meanCellArea = globalArea / meanCellArea |
136 |
areaEpsil = 1. _d -10 * meanCellArea |
137 |
vol_Epsil = 1. _d -15 * meanCellArea |
138 |
|
139 |
r_minArea = 0. _d 0 |
140 |
IF ( iceMaskMin.GT.0. _d 0 ) r_minArea = 1. _d 0 / iceMaskMin |
141 |
|
142 |
thSIce_multiDimAdv = .TRUE. |
143 |
#ifdef ALLOW_GENERIC_ADVDIFF |
144 |
IF ( thSIceAdvScheme.EQ.ENUM_CENTERED_2ND |
145 |
& .OR.thSIceAdvScheme.EQ.ENUM_UPWIND_3RD |
146 |
& .OR.thSIceAdvScheme.EQ.ENUM_CENTERED_4TH ) THEN |
147 |
thSIce_multiDimAdv = .FALSE. |
148 |
ENDIF |
149 |
#endif /* ALLOW_GENERIC_ADVDIFF */ |
150 |
|
151 |
C-- Initialisation (+ build oceanic mask) |
152 |
DO j=1-OLy,sNy+OLy |
153 |
DO i=1-OLx,sNx+OLx |
154 |
maskOce(i,j) = 0. _d 0 |
155 |
IF ( hOceMxL(i,j,bi,bj).GT.0. ) maskOce(i,j) = 1. |
156 |
iceVol(i,j) = 0. _d 0 |
157 |
uTrans(i,j) = 0. _d 0 |
158 |
vTrans(i,j) = 0. _d 0 |
159 |
uTrIce(i,j) = 0. _d 0 |
160 |
vTrIce(i,j) = 0. _d 0 |
161 |
oceFWfx(i,j,bi,bj) = 0. _d 0 |
162 |
oceSflx(i,j,bi,bj) = 0. _d 0 |
163 |
oceQnet(i,j,bi,bj) = 0. _d 0 |
164 |
ENDDO |
165 |
ENDDO |
166 |
|
167 |
#ifdef ALLOW_AUTODIFF_TAMC |
168 |
CADJ STORE iceHeight(i,j,bi,bj) |
169 |
CADJ & = comlev1_bibj, key=iicekey, byte=isbyte |
170 |
CADJ STORE snowHeight(i,j,bi,bj) |
171 |
CADJ & = comlev1_bibj, key=iicekey, byte=isbyte |
172 |
#endif |
173 |
IF ( thSIce_diffK .GT. 0. ) THEN |
174 |
CALL THSICE_DIFFUSION( |
175 |
I maskOce, |
176 |
U uIce, vIce, |
177 |
I bi, bj, myTime, myIter, myThid ) |
178 |
ENDIF |
179 |
|
180 |
IF ( thSIce_multiDimAdv ) THEN |
181 |
|
182 |
C- Calculate ice transports through tracer cell faces. |
183 |
DO j=1-Oly,sNy+Oly |
184 |
DO i=1-Olx+1,sNx+Olx |
185 |
uTrIce(i,j) = uIce(i,j)*_dyG(i,j,bi,bj) |
186 |
& *maskOce(i-1,j)*maskOce(i,j) |
187 |
ENDDO |
188 |
ENDDO |
189 |
DO j=1-Oly+1,sNy+Oly |
190 |
DO i=1-Olx,sNx+Olx |
191 |
vTrIce(i,j) = vIce(i,j)*_dxG(i,j,bi,bj) |
192 |
& *maskOce(i,j-1)*maskOce(i,j) |
193 |
ENDDO |
194 |
ENDDO |
195 |
|
196 |
C-- Fractional area |
197 |
DO j=1-Oly,sNy+Oly |
198 |
DO i=1-Olx,sNx+Olx |
199 |
iceFrc(i,j) = iceMask(i,j,bi,bj) |
200 |
ENDDO |
201 |
ENDDO |
202 |
#ifdef ALLOW_AUTODIFF_TAMC |
203 |
CADJ STORE icevol(i,j) = comlev1_bibj, key=iicekey, byte=isbyte |
204 |
CADJ STORE utrice(i,j) = comlev1_bibj, key=iicekey, byte=isbyte |
205 |
CADJ STORE vtrice(i,j) = comlev1_bibj, key=iicekey, byte=isbyte |
206 |
#endif |
207 |
CALL THSICE_ADVECTION( |
208 |
I GAD_SI_FRAC, thSIceAdvScheme, .TRUE., |
209 |
I uTrIce, vTrIce, maskOce, thSIce_deltaT, areaEpsil, |
210 |
U iceVol, iceFrc, |
211 |
O uTrans, vTrans, |
212 |
I bi, bj, myTime, myIter, myThid ) |
213 |
|
214 |
C-- Snow thickness |
215 |
DO j=1-Oly,sNy+Oly |
216 |
DO i=1-Olx,sNx+Olx |
217 |
iceVol(i,j) = iceMask(i,j,bi,bj)*rA(i,j,bi,bj) |
218 |
ENDDO |
219 |
ENDDO |
220 |
CALL THSICE_ADVECTION( |
221 |
I GAD_SI_HSNOW, thSIceAdvScheme, .FALSE., |
222 |
I uTrans, vTrans, maskOce, thSIce_deltaT, areaEpsil, |
223 |
U iceVol, snowHeight(1-Olx,1-Oly,bi,bj), |
224 |
O afx, afy, |
225 |
I bi, bj, myTime, myIter, myThid ) |
226 |
|
227 |
#ifdef ALLOW_AUTODIFF_TAMC |
228 |
CADJ STORE iceHeight(i,j,bi,bj) = comlev1_bibj, key=iicekey, byte=isbyte |
229 |
CADJ STORE iceMask(i,j,bi,bj) = comlev1_bibj, key=iicekey, byte=isbyte |
230 |
#endif |
231 |
C-- sea-ice Thickness |
232 |
DO j=1-Oly,sNy+Oly |
233 |
DO i=1-Olx,sNx+Olx |
234 |
iceVol(i,j) = iceMask(i,j,bi,bj)*rA(i,j,bi,bj) |
235 |
oldVol(i,j) = iceVol(i,j)*iceHeight(i,j,bi,bj) |
236 |
ENDDO |
237 |
ENDDO |
238 |
CALL THSICE_ADVECTION( |
239 |
I GAD_SI_HICE, thSIceAdvScheme, .FALSE., |
240 |
I uTrans, vTrans, maskOce, thSIce_deltaT, areaEpsil, |
241 |
U iceVol, iceHeight(1-Olx,1-Oly,bi,bj), |
242 |
O uTrIce, vTrIce, |
243 |
I bi, bj, myTime, myIter, myThid ) |
244 |
|
245 |
#ifdef ALLOW_AUTODIFF_TAMC |
246 |
CADJ STORE qice2(i,j,bi,bj) = comlev1_bibj, key=iicekey, byte=isbyte |
247 |
CADJ STORE utrice(i,j) = comlev1_bibj, key=iicekey, byte=isbyte |
248 |
CADJ STORE vtrice(i,j) = comlev1_bibj, key=iicekey, byte=isbyte |
249 |
#endif |
250 |
|
251 |
#ifdef ALLOW_DBUG_THSICE |
252 |
IF ( dBugFlag ) THEN |
253 |
sumVar1 = 0. |
254 |
sumVar2 = 0. |
255 |
DO j=1,sNy |
256 |
DO i=1,sNx |
257 |
C- Check that updated iceVol = iceFrc*rA |
258 |
tmpVar = ABS(iceVol(i,j)-iceFrc(i,j)*rA(i,j,bi,bj)) |
259 |
IF ( tmpVar.GT.0. ) THEN |
260 |
sumVar1 = sumVar1 + 1. |
261 |
sumVar2 = sumVar2 + tmpVar |
262 |
ENDIF |
263 |
IF ( tmpVar.GT.vol_Epsil ) THEN |
264 |
WRITE(6,'(A,2I4,2I2,I12)') 'ARE_ADV: ij,bij,it=', |
265 |
& i,j,bi,bj,myIter |
266 |
WRITE(6,'(2(A,1P2E14.6))') 'ARE_ADV: iceVol,iceFrc*rA=', |
267 |
& iceVol(i,j),iceFrc(i,j)*rA(i,j,bi,bj), |
268 |
& ' , diff=', tmpVar |
269 |
ENDIF |
270 |
IF ( dBug(i,j,bi,bj) ) THEN |
271 |
WRITE(6,'(A,2I4,2I2,I12)') 'ICE_ADV: ij,bij,it=', |
272 |
& i,j,bi,bj,myIter |
273 |
WRITE(6,'(2(A,1P2E14.6))') |
274 |
& 'ICE_ADV: uIce=', uIce(i,j), uIce(i+1,j), |
275 |
& ' , vIce=', vIce(i,j), vIce(i,j+1) |
276 |
WRITE(6,'(2(A,1P2E14.6))') |
277 |
c & 'ICE_ADV: heff_b,a=', HEFF(i,j,2,bi,bj),HEFF(i,j,1,bi,bj) |
278 |
c WRITE(6,'(A,1P4E14.6)') 'ICE_ADV: mFx=', gFld(i,j) |
279 |
ENDIF |
280 |
ENDDO |
281 |
ENDDO |
282 |
IF ( sumVar2.GT.vol_Epsil ) |
283 |
& WRITE(6,'(A,2I2,I10,A,I4,1P2E14.6)') 'ARE_ADV: bij,it:', |
284 |
& bi,bj,myIter, ' ; Npts,aveDiff,Epsil=', |
285 |
& INT(sumVar1),sumVar2/sumVar1,vol_Epsil |
286 |
ENDIF |
287 |
#endif |
288 |
#ifdef ALLOW_DIAGNOSTICS |
289 |
C-- Diagnosse advective fluxes (ice-fraction, snow & ice thickness): |
290 |
IF ( useDiagnostics ) THEN |
291 |
diagSufx = THSICE_DIAG_SUFX( GAD_SI_FRAC, myThid ) |
292 |
diagName = 'ADVx'//diagSufx |
293 |
CALL DIAGNOSTICS_FILL( uTrans, diagName, 1,1,2,bi,bj, myThid ) |
294 |
diagName = 'ADVy'//diagSufx |
295 |
CALL DIAGNOSTICS_FILL( vTrans, diagName, 1,1,2,bi,bj, myThid ) |
296 |
|
297 |
diagSufx = THSICE_DIAG_SUFX( GAD_SI_HSNOW, myThid ) |
298 |
diagName = 'ADVx'//diagSufx |
299 |
CALL DIAGNOSTICS_FILL( afx, diagName, 1,1,2,bi,bj, myThid ) |
300 |
diagName = 'ADVy'//diagSufx |
301 |
CALL DIAGNOSTICS_FILL( afy, diagName, 1,1,2,bi,bj, myThid ) |
302 |
|
303 |
diagSufx = THSICE_DIAG_SUFX( GAD_SI_HICE, myThid ) |
304 |
diagName = 'ADVx'//diagSufx |
305 |
CALL DIAGNOSTICS_FILL( uTrIce, diagName, 1,1,2,bi,bj, myThid ) |
306 |
diagName = 'ADVy'//diagSufx |
307 |
CALL DIAGNOSTICS_FILL( vTrIce, diagName, 1,1,2,bi,bj, myThid ) |
308 |
ENDIF |
309 |
#endif |
310 |
|
311 |
C-- Enthalpy in layer 1 |
312 |
DO j=1-Oly,sNy+Oly |
313 |
DO i=1-Olx,sNx+Olx |
314 |
iceVol(i,j) = oldVol(i,j) |
315 |
ENDDO |
316 |
ENDDO |
317 |
CALL THSICE_ADVECTION( |
318 |
I GAD_SI_QICE1, thSIceAdvScheme, .FALSE., |
319 |
I uTrIce, vTrIce, maskOce, thSIce_deltaT, vol_Epsil, |
320 |
U iceVol, Qice1(1-Olx,1-Oly,bi,bj), |
321 |
O afx, afy, |
322 |
I bi, bj, myTime, myIter, myThid ) |
323 |
#ifdef ALLOW_DBUG_THSICE |
324 |
IF ( dBugFlag ) THEN |
325 |
DO j=1,sNy |
326 |
DO i=1,sNx |
327 |
IF ( dBug(i,j,bi,bj) ) THEN |
328 |
c WRITE(6,'(A,1P4E14.6)') 'ICE_ADV: Qice1_b,a=', |
329 |
c & Qice1(i,j,bi,bj), |
330 |
c & ( iceFld(i,j) + thSIce_deltaT * gFld(i,j) |
331 |
c & ) * recip_heff(i,j) |
332 |
c WRITE(6,'(A,1P4E14.6)') 'ICE_ADV: q1Fx=', gFld(i,j) |
333 |
ENDIF |
334 |
ENDDO |
335 |
ENDDO |
336 |
ENDIF |
337 |
#endif |
338 |
#ifdef ALLOW_DIAGNOSTICS |
339 |
IF ( useDiagnostics ) THEN |
340 |
diagSufx = THSICE_DIAG_SUFX( GAD_SI_QICE1, myThid ) |
341 |
diagName = 'ADVx'//diagSufx |
342 |
CALL DIAGNOSTICS_FILL( afx, diagName, 1,1,2,bi,bj, myThid ) |
343 |
diagName = 'ADVy'//diagSufx |
344 |
CALL DIAGNOSTICS_FILL( afy, diagName, 1,1,2,bi,bj, myThid ) |
345 |
ENDIF |
346 |
#endif |
347 |
|
348 |
C-- Enthalpy in layer 2 |
349 |
DO j=1-Oly,sNy+Oly |
350 |
DO i=1-Olx,sNx+Olx |
351 |
iceVol(i,j) = oldVol(i,j) |
352 |
ENDDO |
353 |
ENDDO |
354 |
CALL THSICE_ADVECTION( |
355 |
I GAD_SI_QICE2, thSIceAdvScheme, .FALSE., |
356 |
I uTrIce, vTrIce, maskOce, thSIce_deltaT, vol_Epsil, |
357 |
U iceVol, Qice2(1-Olx,1-Oly,bi,bj), |
358 |
O afx, afy, |
359 |
I bi, bj, myTime, myIter, myThid ) |
360 |
#ifdef ALLOW_DBUG_THSICE |
361 |
IF ( dBugFlag ) THEN |
362 |
sumVar1 = 0. |
363 |
sumVar2 = 0. |
364 |
DO j=1,sNy |
365 |
DO i=1,sNx |
366 |
C- Check that updated iceVol = Hic*Frc*rA |
367 |
tmpVar = ABS(iceVol(i,j) |
368 |
& -iceHeight(i,j,bi,bj)*iceFrc(i,j)*rA(i,j,bi,bj)) |
369 |
IF ( tmpVar.GT.0. ) THEN |
370 |
sumVar1 = sumVar1 + 1. |
371 |
sumVar2 = sumVar2 + tmpVar |
372 |
ENDIF |
373 |
IF ( tmpVar.GT.vol_Epsil ) THEN |
374 |
WRITE(6,'(A,2I4,2I2,I12)') 'VOL_ADV: ij,bij,it=', |
375 |
& i,j,bi,bj,myIter |
376 |
WRITE(6,'(2(A,1P2E14.6))') 'VOL_ADV: iceVol,Hic*Frc*rA=', |
377 |
& iceVol(i,j),iceHeight(i,j,bi,bj)*iceFrc(i,j)*rA(i,j,bi,bj), |
378 |
& ' , diff=', tmpVar |
379 |
ENDIF |
380 |
IF ( dBug(i,j,bi,bj) ) THEN |
381 |
c WRITE(6,'(A,1P4E14.6)') 'ICE_ADV: Qice2_b,a=', |
382 |
c & Qice2(i,j,bi,bj), |
383 |
c & ( iceFld(i,j) + thSIce_deltaT * gFld(i,j) |
384 |
c & ) * recip_heff(i,j) |
385 |
c WRITE(6,'(A,1P4E14.6)') 'ICE_ADV: q2Fx=', gFld(i,j) |
386 |
ENDIF |
387 |
ENDDO |
388 |
ENDDO |
389 |
IF ( sumVar2.GT.vol_Epsil ) |
390 |
& WRITE(6,'(A,2I2,I10,A,I4,1P2E14.6)') 'VOL_ADV: bij,it:', |
391 |
& bi,bj,myIter, ' ; Npts,aveDiff,Epsil=', |
392 |
& INT(sumVar1),sumVar2/sumVar1,vol_Epsil |
393 |
ENDIF |
394 |
#endif |
395 |
#ifdef ALLOW_DIAGNOSTICS |
396 |
IF ( useDiagnostics ) THEN |
397 |
diagSufx = THSICE_DIAG_SUFX( GAD_SI_QICE2, myThid ) |
398 |
diagName = 'ADVx'//diagSufx |
399 |
CALL DIAGNOSTICS_FILL( afx, diagName, 1,1,2,bi,bj, myThid ) |
400 |
diagName = 'ADVy'//diagSufx |
401 |
CALL DIAGNOSTICS_FILL( afy, diagName, 1,1,2,bi,bj, myThid ) |
402 |
ENDIF |
403 |
#endif |
404 |
|
405 |
#ifdef ALLOW_AUTODIFF_TAMC |
406 |
CADJ STORE iceHeight(:,:,bi,bj) = |
407 |
CADJ & comlev1_bibj, key=iicekey, byte=isbyte |
408 |
CADJ STORE snowHeight(:,:,bi,bj) = |
409 |
CADJ & comlev1_bibj, key=iicekey, byte=isbyte |
410 |
CADJ STORE iceFrc(:,:) = |
411 |
CADJ & comlev1_bibj, key=iicekey, byte=isbyte |
412 |
#endif |
413 |
|
414 |
C-- Update Ice Fraction: ensure that fraction is > iceMaskMin & < 1 |
415 |
C and adjust Ice thickness and snow thickness accordingly |
416 |
DO j=1,sNy |
417 |
DO i=1,sNx |
418 |
IF ( iceFrc(i,j) .GT. 1. _d 0 ) THEN |
419 |
iceMask(i,j,bi,bj) = 1. _d 0 |
420 |
iceHeight(i,j,bi,bj) = iceHeight(i,j,bi,bj) *iceFrc(i,j) |
421 |
snowHeight(i,j,bi,bj) = snowHeight(i,j,bi,bj)*iceFrc(i,j) |
422 |
ELSEIF ( iceFrc(i,j) .LT. iceMaskMin ) THEN |
423 |
iceMask(i,j,bi,bj) = iceMaskMin |
424 |
iceHeight(i,j,bi,bj) = iceHeight(i,j,bi,bj) |
425 |
& *iceFrc(i,j)*r_minArea |
426 |
snowHeight(i,j,bi,bj) = snowHeight(i,j,bi,bj) |
427 |
& *iceFrc(i,j)*r_minArea |
428 |
ELSE |
429 |
iceMask(i,j,bi,bj) = iceFrc(i,j) |
430 |
ENDIF |
431 |
ENDDO |
432 |
ENDDO |
433 |
C- adjust sea-ice state if not enough ice. |
434 |
DO j=1,sNy |
435 |
DO i=1,sNx |
436 |
IF ( iceHeight(i,j,bi,bj).LT.hIceMin ) THEN |
437 |
C- Not enough ice, melt the tiny amount of snow & ice: |
438 |
C and return fresh-water, salt & energy to the ocean (flx > 0 = into ocean) |
439 |
C- - Note: using 1rst.Order Upwind, I can get the same results as when |
440 |
C using seaice_advdiff (with SEAICEadvScheme=1) providing I comment |
441 |
C out the following lines (and then loose conservation). |
442 |
C- - |
443 |
oceFWfx(i,j,bi,bj) = ( rhos*snowHeight(i,j,bi,bj) |
444 |
& +rhoi*iceHeight(i,j,bi,bj) |
445 |
& )*iceMask(i,j,bi,bj)/thSIce_deltaT |
446 |
oceSflx(i,j,bi,bj) = rhoi*iceHeight(i,j,bi,bj)*saltIce |
447 |
& *iceMask(i,j,bi,bj)/thSIce_deltaT |
448 |
oceQnet(i,j,bi,bj) = -( rhos*snowHeight(i,j,bi,bj)*qsnow |
449 |
& +rhoi*iceHeight(i,j,bi,bj) |
450 |
& *( Qice1(i,j,bi,bj) |
451 |
& +Qice2(i,j,bi,bj) )*0.5 _d 0 |
452 |
& )*iceMask(i,j,bi,bj)/thSIce_deltaT |
453 |
C- - |
454 |
c flx2oc (i,j) = flx2oc (i,j) + |
455 |
c frw2oc (i,j) = frw2oc (i,j) + |
456 |
c fsalt (i,j) = fsalt (i,j) + |
457 |
iceMask (i,j,bi,bj) = 0. _d 0 |
458 |
iceHeight (i,j,bi,bj) = 0. _d 0 |
459 |
snowHeight(i,j,bi,bj) = 0. _d 0 |
460 |
Qice1 (i,j,bi,bj) = 0. _d 0 |
461 |
Qice2 (i,j,bi,bj) = 0. _d 0 |
462 |
snowAge (i,j,bi,bj) = 0. _d 0 |
463 |
ENDIF |
464 |
ENDDO |
465 |
ENDDO |
466 |
|
467 |
#ifdef ALLOW_DBUG_THSICE |
468 |
IF ( dBugFlag ) THEN |
469 |
DO j=1,sNy |
470 |
DO i=1,sNx |
471 |
IF ( dBug(i,j,bi,bj) ) THEN |
472 |
WRITE(6,'(2(A,1P2E14.6))') |
473 |
c & 'ICE_ADV: area_b,a=', AREA(i,j,2,bi,bj),AREA(i,j,1,bi,bj) |
474 |
c WRITE(6,'(A,1P4E14.6)') 'ICE_ADV: mFx=', gFld(i,j) |
475 |
ENDIF |
476 |
ENDDO |
477 |
ENDDO |
478 |
ENDIF |
479 |
#endif |
480 |
|
481 |
ELSE |
482 |
C--- if not multiDimAdvection |
483 |
|
484 |
WRITE(msgBuf,'(2A)') 'S/R THSICE_ADVDIFF: ', |
485 |
& 'traditional advection/diffusion not yet implemented' |
486 |
CALL PRINT_ERROR( msgBuf , myThid) |
487 |
WRITE(msgBuf,'(2A)') ' ', |
488 |
& 'for ThSice variable Qice1, Qice2, SnowHeight. Sorry!' |
489 |
CALL PRINT_ERROR( msgBuf , myThid) |
490 |
STOP 'ABNORMAL: END: S/R THSICE_ADVDIFF' |
491 |
|
492 |
C--- end if multiDimAdvection |
493 |
ENDIF |
494 |
|
495 |
#endif /* ALLOW_THSICE */ |
496 |
|
497 |
RETURN |
498 |
END |