/[MITgcm]/MITgcm/pkg/layers/layers_calc.F
ViewVC logotype

Contents of /MITgcm/pkg/layers/layers_calc.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.32 - (show annotations) (download)
Fri Mar 24 23:38:56 2017 UTC (7 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.31: +15 -11 lines
use new S/R RW_GET_SUFFIX to get file suffix (according to "rwSuffixType")

1 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_calc.F,v 1.31 2015/06/16 21:43:10 rpa Exp $
2 C $Name: $
3
4 #include "LAYERS_OPTIONS.h"
5 #ifdef ALLOW_GMREDI
6 #include "GMREDI_OPTIONS.h"
7 #endif
8
9 CBOP 0
10 C !ROUTINE: LAYERS_CALC
11
12 C !INTERFACE:
13 SUBROUTINE LAYERS_CALC(
14 I myTime, myIter, myThid )
15
16 C !DESCRIPTION:
17 C ===================================================================
18 C Calculate the transport in isopycnal layers.
19 C This was the meat of the LAYERS package, which
20 C has been moved to S/R LAYERS_FLUXCALC.F
21 C ===================================================================
22
23 C !USES:
24 IMPLICIT NONE
25 #include "SIZE.h"
26 #include "EEPARAMS.h"
27 #include "PARAMS.h"
28 #include "GRID.h"
29 #include "DYNVARS.h"
30 #include "LAYERS_SIZE.h"
31 #include "LAYERS.h"
32 #ifdef ALLOW_GMREDI
33 # include "GMREDI.h"
34 #endif
35
36 C !INPUT PARAMETERS:
37 C myTime :: Current time in simulation
38 C myIter :: Current iteration number
39 C myThid :: my Thread Id number
40 _RL myTime
41 INTEGER myIter
42 INTEGER myThid
43 CEOP
44
45 #ifdef ALLOW_LAYERS
46 C !FUNCTIONS:
47 LOGICAL DIFFERENT_MULTIPLE
48 EXTERNAL DIFFERENT_MULTIPLE
49
50 C !LOCAL VARIABLES:
51 C -- 3D Layers fields. The vertical dimension in these fields is Nlayers,
52 C i.e. the isopycnal coordinate.
53 C layers_UH :: U integrated over layer (m^2/s)
54 C layers_VH :: V integrated over layer (m^2/s)
55 C layers_Hw :: Layer thickness at the U point (m)
56 C layers_Hs :: Layer thickness at the V point (m)
57 C layers_PIw :: 1 if layer exists, 0 otherwise
58 C layers_PIs :: 1 if layer exists, 0 otherwise
59 C layers_U :: mean zonal velocity in layer (only if layer exists) (m/s)
60 C layers_V :: mean meridional velocity in layer (only if layer exists) (m/s)
61 #ifdef LAYERS_UFLUX
62 _RL layers_UH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
63 # ifdef LAYERS_THICKNESS
64 _RL layers_Hw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
65 _RL layers_PIw(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
66 _RL layers_U (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
67 # endif /* LAYERS_THICKNESS */
68 #endif /* LAYERS_UFLUX */
69 #ifdef LAYERS_VFLUX
70 _RL layers_VH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
71 # ifdef LAYERS_THICKNESS
72 _RL layers_Hs (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
73 _RL layers_PIs(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
74 _RL layers_V (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
75 # endif /* LAYERS_THICKNESS */
76 #endif /* LAYERS_VFLUX */
77 #ifdef LAYERS_PRHO_REF
78 _RL prho(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
79 _RL rhoShift
80 INTEGER i, j, k
81 #endif
82 C -- other local variables:
83 C bi, bj :: tile indices
84 C i,j :: horizontal indices
85 C iLa :: layer-type index
86 C k :: vertical index for model grid
87 INTEGER bi, bj, iLa
88 CHARACTER*(10) sufx
89 CHARACTER*(13) suff
90
91 #ifdef LAYERS_THERMODYNAMICS
92 INTEGER iTracer
93 #endif
94 #ifdef ALLOW_DIAGNOSTICS
95 CHARACTER*8 diagName
96 #endif
97 c#ifdef ALLOW_MNC
98 c CHARACTER*(1) pf
99 c#endif
100
101 #ifndef LAYERS_UFLUX
102 _RL layers_UH(1)
103 #endif
104 #ifndef LAYERS_VFLUX
105 _RL layers_VH(1)
106 #endif
107 #if !(defined LAYERS_THICKNESS) || !(defined LAYERS_UFLUX)
108 _RL layers_Hw(1), layers_PIw(1), layers_U(1)
109 #endif
110 #if !(defined LAYERS_THICKNESS) || !(defined LAYERS_VFLUX)
111 _RL layers_Hs(1), layers_PIs(1), layers_V(1)
112 #endif
113
114 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
115
116 IF ( myIter.EQ.nIter0 ) RETURN
117
118 #ifdef LAYERS_THERMODYNAMICS
119 CALL LAYERS_CALC_RHS(myThid)
120 #endif
121
122 DO iLa=1,layers_maxNum
123
124 IF ( layers_num(iLa) .EQ. 1 ) THEN
125 CALL LAYERS_FLUXCALC( uVel,vVel,theta,iLa,
126 & layers_UH, layers_VH,
127 & layers_Hw, layers_Hs,
128 & layers_PIw,layers_PIs,
129 & layers_U, layers_V,
130 & myThid )
131 #ifdef LAYERS_THERMODYNAMICS
132 CALL LAYERS_DIAPYCNAL( theta,iLa,
133 & layers_TtendSurf,
134 & layers_TtendDiffh, layers_TtendDiffr,
135 & layers_TtendAdvh, layers_TtendAdvr,
136 & layers_Ttendtot,
137 & layers_StendSurf,
138 & layers_StendDiffh, layers_StendDiffr,
139 & layers_StendAdvh, layers_StendAdvr,
140 & layers_Stendtot,
141 & layers_Hc, layers_PIc,
142 & myThid)
143 #endif
144 ELSEIF ( layers_num(iLa) .EQ. 2 ) THEN
145 CALL LAYERS_FLUXCALC( uVel,vVel,salt,iLa,
146 & layers_UH, layers_VH,
147 & layers_Hw, layers_Hs,
148 & layers_PIw,layers_PIs,
149 & layers_U, layers_V,
150 & myThid )
151 #ifdef LAYERS_THERMODYNAMICS
152 CALL LAYERS_DIAPYCNAL( salt,iLa,
153 & layers_TtendSurf,
154 & layers_TtendDiffh, layers_TtendDiffr,
155 & layers_TtendAdvh, layers_TtendAdvr,
156 & layers_Ttendtot,
157 & layers_StendSurf,
158 & layers_StendDiffh, layers_StendDiffr,
159 & layers_StendAdvh, layers_StendAdvr,
160 & layers_Stendtot,
161 & layers_Hc, layers_PIc,
162 & myThid)
163 #endif
164 ELSEIF ( layers_num(iLa) .EQ. 3 ) THEN
165 #ifdef LAYERS_PRHO_REF
166 C For layers_num(iLa) = 3, calculate the potential density (minus 1000)
167 C referenced to the model level given by layers_krho.
168 rhoShift = rhoConst - 1000. _d 0
169 DO bj=myByLo(myThid),myByHi(myThid)
170 DO bi=myBxLo(myThid),myBxHi(myThid)
171 C Initialise layers variable prho:
172 DO k=1,Nr
173 DO j=1-OLy,sNy+OLy
174 DO i=1-OLx,sNx+OLx
175 prho(i,j,k,bi,bj) = 0. _d 0
176 ENDDO
177 ENDDO
178 ENDDO
179 DO k = 1,Nr
180 CALL FIND_RHO_2D( 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
181 & layers_krho(iLa),
182 & theta(1-OLx,1-OLy,k,bi,bj),
183 & salt(1-OLx,1-OLy,k,bi,bj),
184 & prho(1-OLx,1-OLy,k,bi,bj),
185 & k, bi, bj, myThid )
186 #ifdef LAYERS_THERMODYNAMICS
187 C -- it might be more memory efficient not to store alpha and beta
188 C but to multiply the fluxes in place here
189 CALL FIND_ALPHA( bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
190 & k, layers_krho(iLa),
191 & layers_alpha(1-OLx,1-OLy,k,bi,bj), myThid )
192 CALL FIND_BETA( bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
193 & k, layers_krho(iLa),
194 & layers_beta(1-OLx,1-OLy,k,bi,bj), myThid )
195 #endif /* LAYERS_THERMODYNAMICS */
196 DO j = 1-OLy,sNy+OLy
197 DO i = 1-OLx,sNx+OLx
198 prho(i,j,k,bi,bj) = prho(i,j,k,bi,bj) + rhoShift
199 ENDDO
200 ENDDO
201 ENDDO
202 ENDDO
203 ENDDO
204 CALL LAYERS_FLUXCALC( uVel,vVel, prho, iLa,
205 & layers_UH, layers_VH,
206 & layers_Hw, layers_Hs,
207 & layers_PIw,layers_PIs,
208 & layers_U, layers_V,
209 & myThid )
210 #ifdef LAYERS_THERMODYNAMICS
211 CALL LAYERS_DIAPYCNAL( prho,iLa,
212 & layers_TtendSurf,
213 & layers_TtendDiffh, layers_TtendDiffr,
214 & layers_TtendAdvh, layers_TtendAdvr,
215 & layers_Ttendtot,
216 & layers_StendSurf,
217 & layers_StendDiffh, layers_StendDiffr,
218 & layers_StendAdvh, layers_StendAdvr,
219 & layers_Stendtot,
220 & layers_Hc, layers_PIc,
221 & myThid)
222 #endif
223 #endif /* LAYERS_PRHO_REF */
224 ENDIF
225
226 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
227 C-- Direct Snap-shot output
228 IF ( DIFFERENT_MULTIPLE(layers_diagFreq,myTime,deltaTClock)
229 & .AND. layers_num(iLa).NE.0 ) THEN
230
231 IF ( layers_MDSIO ) THEN
232 IF ( rwSuffixType.EQ.0 ) THEN
233 WRITE(suff,'(I2.2,A1,I10.10)') iLa, '.', myIter
234 ELSE
235 CALL RW_GET_SUFFIX( sufx, myTime, myIter, myThid )
236 WRITE(suff,'(I2.2,A1,A)') iLa, '.', sufx
237 ENDIF
238 #ifdef LAYERS_UFLUX
239 CALL WRITE_FLD_3D_RL( 'layers_UH.', suff, Nlayers,
240 & layers_UH, myIter, myThid )
241 #ifdef LAYERS_THICKNESS
242 CALL WRITE_FLD_3D_RL( 'layers_Hw.', suff, Nlayers,
243 & layers_Hw, myIter, myThid )
244 #endif /* LAYERS_THICKNESS */
245 #endif /* LAYERS_UFLUX */
246 #ifdef LAYERS_VFLUX
247 CALL WRITE_FLD_3D_RL( 'layers_VH.', suff, Nlayers,
248 & layers_VH, myIter, myThid )
249 #ifdef LAYERS_THICKNESS
250 CALL WRITE_FLD_3D_RL( 'layers_Hs.', suff, Nlayers,
251 & layers_Hs, myIter, myThid )
252 #endif /* LAYERS_THICKNESS */
253 #endif /* LAYERS_VFLUX */
254 #ifdef LAYERS_PRHO_REF
255 IF ( layers_num(iLa).EQ.3 ) THEN
256 CALL WRITE_FLD_3D_RL( 'layers_prho.', suff, Nr,
257 & prho, myIter, myThid )
258 ENDIF
259 #endif /* LAYERS_PRHO_REF */
260
261 #ifdef LAYERS_THERMODYNAMICS
262 CALL WRITE_FLD_3D_RL( 'layers_Ttottend.', suff, 2*Nr,
263 & layers_tottend, myIter, myThid )
264 #ifdef SHORTWAVE_HEATING
265 CALL WRITE_FLD_3D_RL( 'layers_sw.', suff, Nr,
266 & layers_sw, myIter, myThid )
267 #endif /* LAYERS_SHORTWAVE */
268 CALL WRITE_FLD_3D_RL( 'layers_surfflux.', suff, 2,
269 & layers_surfflux, myIter, myThid )
270 CALL WRITE_FLD_3D_RL( 'layers_dfx.', suff, 2*Nr,
271 & layers_dfx, myIter, myThid )
272 CALL WRITE_FLD_3D_RL( 'layers_dfy.', suff, 2*Nr,
273 & layers_dfy, myIter, myThid )
274 CALL WRITE_FLD_3D_RL( 'layers_dfr.', suff, 2*Nr,
275 & layers_dfr, myIter, myThid )
276 CALL WRITE_FLD_3D_RL( 'layers_afx.', suff, 2*Nr,
277 & layers_afx, myIter, myThid )
278 CALL WRITE_FLD_3D_RL( 'layers_afy.', suff, 2*Nr,
279 & layers_afy, myIter, myThid )
280 CALL WRITE_FLD_3D_RL( 'layers_afr.', suff, 2*Nr,
281 & layers_afr, myIter, myThid )
282 #endif /* LAYERS_THERMODYNAMICS */
283 ENDIF
284
285 c#ifdef ALLOW_MNC
286 c#ifdef LAYERS_MNC
287 c IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
288 c pf(1:1) = 'D'
289 c ELSE
290 c pf(1:1) = 'R'
291 c ENDIF
292 c IF ( layers_MNC) THEN
293 C Do MNC output... But how?
294 c ENDIF
295 c#endif /* LAYERS_MNC */
296 c#endif /* ALLOW_MNC */
297
298 ENDIF
299
300 #ifdef ALLOW_DIAGNOSTICS
301 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
302 C-- Fill-in diagnostics
303 IF ( useDiagnostics .AND. layers_num(iLa).NE.0 ) THEN
304
305 #ifdef LAYERS_PRHO_REF
306 IF ( layers_num(iLa).EQ.3 ) THEN
307 WRITE(diagName,'(A4,I1,A3)') 'LaTr',iLa,layers_name(iLa)
308 CALL DIAGNOSTICS_FILL( prho,
309 & diagName, 0, Nr, 0, 1, 1, myThid )
310 ENDIF
311 #endif /* LAYERS_PRHO_REF */
312
313 #ifdef LAYERS_UFLUX
314 WRITE(diagName,'(A4,I1,A3)') 'LaUH',iLa,layers_name(iLa)
315 CALL DIAGNOSTICS_FILL( layers_UH,
316 & diagName,0,Nlayers, 0, 1, 1, myThid )
317 # ifdef LAYERS_THICKNESS
318 WRITE(diagName,'(A4,I1,A3)') 'LaHw',iLa,layers_name(iLa)
319 CALL DIAGNOSTICS_FILL( layers_Hw,
320 & diagName,0,Nlayers, 0, 1, 1, myThid )
321 WRITE(diagName,'(A4,I1,A3)') 'LaPw',iLa,layers_name(iLa)
322 CALL DIAGNOSTICS_FILL( layers_PIw,
323 & diagName,0,Nlayers, 0, 1, 1, myThid )
324 WRITE(diagName,'(A4,I1,A3)') 'LaUa',iLa,layers_name(iLa)
325 CALL DIAGNOSTICS_FILL( layers_U,
326 & diagName,0,Nlayers, 0, 1, 1, myThid )
327 # endif
328 #endif /* LAYERS_UFLUX */
329
330 #ifdef LAYERS_VFLUX
331 WRITE(diagName,'(A4,I1,A3)') 'LaVH',iLa,layers_name(iLa)
332 CALL DIAGNOSTICS_FILL( layers_VH,
333 & diagName,0,Nlayers, 0, 1, 1, myThid )
334 # ifdef LAYERS_THICKNESS
335 WRITE(diagName,'(A4,I1,A3)') 'LaHs',iLa,layers_name(iLa)
336 CALL DIAGNOSTICS_FILL( layers_Hs,
337 & diagName,0,Nlayers, 0, 1, 1, myThid )
338 WRITE(diagName,'(A4,I1,A3)') 'LaPs',iLa,layers_name(iLa)
339 CALL DIAGNOSTICS_FILL( layers_PIs,
340 & diagName,0,Nlayers, 0, 1, 1, myThid )
341 WRITE(diagName,'(A4,I1,A3)') 'LaVa',iLa,layers_name(iLa)
342 CALL DIAGNOSTICS_FILL( layers_V,
343 & diagName,0,Nlayers, 0, 1, 1, myThid )
344 # endif
345 #endif /* LAYERS_VFLUX */
346
347 #ifdef LAYERS_THERMODYNAMICS
348 WRITE(diagName,'(A4,I1,A3)') 'LaHc',iLa,layers_name(iLa)
349 CALL DIAGNOSTICS_FILL( layers_Hc,
350 & diagName,0,Nlayers, 0, 1, 1, myThid )
351 WRITE(diagName,'(A4,I1,A3)') 'LaPc',iLa,layers_name(iLa)
352 CALL DIAGNOSTICS_FILL( layers_PIc,
353 & diagName,0,Nlayers, 0, 1, 1, myThid )
354 WRITE(diagName,'(A4,I1,A3)') 'LaTs',iLa,layers_name(iLa)
355 CALL DIAGNOSTICS_FILL( layers_TtendSurf,
356 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
357 WRITE(diagName,'(A4,I1,A3)') 'LaTh',iLa,layers_name(iLa)
358 CALL DIAGNOSTICS_FILL( layers_TtendDiffh,
359 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
360 WRITE(diagName,'(A4,I1,A3)') 'LaTz',iLa,layers_name(iLa)
361 CALL DIAGNOSTICS_FILL( layers_TtendDiffr,
362 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
363 WRITE(diagName,'(A4,I1,A3)') 'LTha',iLa,layers_name(iLa)
364 CALL DIAGNOSTICS_FILL( layers_TtendAdvh,
365 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
366 WRITE(diagName,'(A4,I1,A3)') 'LTza',iLa,layers_name(iLa)
367 CALL DIAGNOSTICS_FILL( layers_TtendAdvr,
368 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
369 WRITE(diagName,'(A4,I1,A3)') 'LTto',iLa,layers_name(iLa)
370 CALL DIAGNOSTICS_FILL( layers_Ttendtot,
371 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
372 WRITE(diagName,'(A4,I1,A3)') 'LaSs',iLa,layers_name(iLa)
373 CALL DIAGNOSTICS_FILL( layers_StendSurf,
374 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
375 WRITE(diagName,'(A4,I1,A3)') 'LaSh',iLa,layers_name(iLa)
376 CALL DIAGNOSTICS_FILL( layers_StendDiffh,
377 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
378 WRITE(diagName,'(A4,I1,A3)') 'LaSz',iLa,layers_name(iLa)
379 CALL DIAGNOSTICS_FILL( layers_StendDiffr,
380 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
381 WRITE(diagName,'(A4,I1,A3)') 'LSha',iLa,layers_name(iLa)
382 CALL DIAGNOSTICS_FILL( layers_StendAdvh,
383 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
384 WRITE(diagName,'(A4,I1,A3)') 'LSza',iLa,layers_name(iLa)
385 CALL DIAGNOSTICS_FILL( layers_StendAdvr,
386 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
387 WRITE(diagName,'(A4,I1,A3)') 'LSto',iLa,layers_name(iLa)
388 CALL DIAGNOSTICS_FILL( layers_Stendtot,
389 & diagName,0,Nlayers-1, 0, 1, 1, myThid )
390 #endif /* LAYERS_THERMODYNAMICS */
391
392 ENDIF
393 #endif /* ALLOW_DIAGNOSTICS */
394
395 #ifdef ALLOW_TIMEAVE
396 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
397 C-- Time-average
398 cgf layers_maxNum loop and dimension would be needed for
399 cgf the following and tave output to work beyond iLa.EQ.1
400 IF ( layers_taveFreq.GT.0. .AND. iLa.EQ.1 ) THEN
401 C --- The tile loops
402 DO bj=myByLo(myThid),myByHi(myThid)
403 DO bi=myBxLo(myThid),myBxHi(myThid)
404
405 #ifdef LAYERS_UFLUX
406 CALL TIMEAVE_CUMULATE( layers_UH_T, layers_UH, Nlayers,
407 & deltaTClock, bi, bj, myThid )
408 #ifdef LAYERS_THICKNESS
409 CALL TIMEAVE_CUMULATE( layers_Hw_T, layers_Hw, Nlayers,
410 & deltaTClock, bi, bj, myThid )
411 #endif /* LAYERS_THICKNESS */
412 #endif /* LAYERS_UFLUX */
413 #ifdef LAYERS_VFLUX
414 CALL TIMEAVE_CUMULATE( layers_VH_T, layers_VH, Nlayers,
415 & deltaTClock, bi, bj, myThid )
416 #ifdef LAYERS_THICKNESS
417 CALL TIMEAVE_CUMULATE( layers_Hs_T, layers_Hs, Nlayers,
418 & deltaTClock, bi, bj, myThid )
419 #endif /* LAYERS_THICKNESS */
420 #endif /* LAYERS_VFLUX */
421
422 #ifdef LAYERS_PRHO_REF
423 IF ( layers_num(iLa) .EQ. 3 )
424 & CALL TIMEAVE_CUMULATE( prho_tave, prho, Nr,
425 & deltaTClock, bi, bj, myThid )
426 #endif /* LAYERS_PRHO_REF */
427
428 layers_TimeAve(bi,bj)=layers_TimeAve(bi,bj)+deltaTClock
429
430 C --- End bi,bj loop
431 ENDDO
432 ENDDO
433 ENDIF
434 #endif /* ALLOW_TIMEAVE */
435
436 ENDDO !DO iLa=1,layers_maxNum
437
438 #ifdef LAYERS_THERMODYNAMICS
439 C-- Reset temporary flux arrays to zero
440 DO bj = myByLo(myThid), myByHi(myThid)
441 DO bi = myBxLo(myThid), myBxHi(myThid)
442 DO iTracer = 1,2
443 DO j=1-OLy,sNy+OLy
444 DO i=1-OLx,sNx+OLx
445 layers_surfflux(i,j,1,iTracer,bi,bj) = 0. _d 0
446 ENDDO
447 ENDDO
448 DO k=1,Nr
449 DO j=1-OLy,sNy+OLy
450 DO i=1-OLx,sNx+OLx
451 layers_dfx (i,j,k,iTracer,bi,bj) = 0. _d 0
452 layers_dfy (i,j,k,iTracer,bi,bj) = 0. _d 0
453 layers_dfr (i,j,k,iTracer,bi,bj) = 0. _d 0
454 layers_afx (i,j,k,iTracer,bi,bj) = 0. _d 0
455 layers_afy (i,j,k,iTracer,bi,bj) = 0. _d 0
456 layers_afr (i,j,k,iTracer,bi,bj) = 0. _d 0
457 layers_tottend (i,j,k,iTracer,bi,bj) = 0. _d 0
458 #ifdef SHORTWAVE_HEATING
459 layers_sw (i,j,k,1 ,bi,bj) = 0. _d 0
460 #endif /* SHORTWAVE_HEATING */
461 ENDDO
462 ENDDO
463 ENDDO
464 ENDDO
465 ENDDO
466 ENDDO
467 #endif /* LAYERS_THERMODYNAMICS */
468
469 #endif /* ALLOW_LAYERS */
470
471 RETURN
472 END

  ViewVC Help
Powered by ViewVC 1.1.22