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