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

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

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


Revision 1.18 - (hide annotations) (download)
Thu Oct 18 19:51:14 2012 UTC (11 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.17: +121 -38 lines
do snap-shot output and fill diagnostics from inside iLa loop (layers_calc.F);
this allows to remove "layers_maxNum" dimension from all layers arrays

1 jmc 1.18 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_calc.F,v 1.17 2012/10/18 13:55:30 jmc Exp $
2 rpa 1.1 C $Name: $
3    
4     #include "LAYERS_OPTIONS.h"
5 jmc 1.6 #ifdef ALLOW_GMREDI
6 dfer 1.5 #include "GMREDI_OPTIONS.h"
7 jmc 1.6 #endif
8 rpa 1.1
9 jmc 1.4 CBOP 0
10     C !ROUTINE: LAYERS_CALC
11 rpa 1.1
12 jmc 1.4 C !INTERFACE:
13 rpa 1.1 SUBROUTINE LAYERS_CALC(
14 jmc 1.4 I myTime, myIter, myThid )
15 rpa 1.1
16 jmc 1.4 C !DESCRIPTION:
17 rpa 1.1 C ===================================================================
18     C Calculate the transport in isopycnal layers.
19 gforget 1.13 C This was the meat of the LAYERS package, which
20     C has been moved to S/R LAYERS_FLUXCALC.F
21 rpa 1.1 C ===================================================================
22    
23 jmc 1.4 C !USES:
24 rpa 1.1 IMPLICIT NONE
25     #include "SIZE.h"
26 jmc 1.4 #include "EEPARAMS.h"
27     #include "PARAMS.h"
28 rpa 1.1 #include "GRID.h"
29     #include "DYNVARS.h"
30     #include "LAYERS_SIZE.h"
31     #include "LAYERS.h"
32 jmc 1.6 #ifdef ALLOW_GMREDI
33     # include "GMREDI.h"
34 dfer 1.5 #endif
35 rpa 1.1
36 jmc 1.4 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 rpa 1.1
45     #ifdef ALLOW_LAYERS
46 jmc 1.18 C !FUNCTIONS:
47     LOGICAL DIFFERENT_MULTIPLE
48     EXTERNAL DIFFERENT_MULTIPLE
49 rpa 1.1
50 jmc 1.4 C !LOCAL VARIABLES:
51     C bi, bj :: tile indices
52 rpa 1.1 C i,j :: horizontal indices
53 gforget 1.14 C iLa :: layer coordinate index
54 rpa 1.1 C k :: vertical index for model grid
55 gforget 1.14 INTEGER bi, bj, iLa
56 jmc 1.18 CHARACTER*(MAX_LEN_MBUF) suff
57 gforget 1.13 #ifdef LAYERS_PRHO_REF
58     INTEGER i, j, k
59 jmc 1.6 #endif
60 jmc 1.18 #ifdef ALLOW_DIAGNOSTICS
61     CHARACTER*8 diagName
62     #endif
63     c#ifdef ALLOW_MNC
64     c CHARACTER*(1) pf
65     c#endif
66 rpa 1.1
67 jmc 1.17 #ifndef LAYERS_UFLUX
68 jmc 1.18 _RL layers_UH(1)
69 jmc 1.17 #endif
70     #ifndef LAYERS_VFLUX
71 jmc 1.18 _RL layers_VH(1)
72 jmc 1.17 #endif
73     #if !(defined LAYERS_THICKNESS) || !(defined LAYERS_UFLUX)
74 jmc 1.18 _RL layers_Hw(1), layers_PIw(1), layers_U(1)
75 jmc 1.17 #endif
76     #if !(defined LAYERS_THICKNESS) || !(defined LAYERS_VFLUX)
77 jmc 1.18 _RL layers_Hs(1), layers_PIs(1), layers_V(1)
78 jmc 1.17 #endif
79    
80 jmc 1.4 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
81    
82 jmc 1.18 IF ( myIter.EQ.nIter0 ) RETURN
83    
84 gforget 1.14 DO iLa=1,layers_maxNum
85 jmc 1.16
86     IF ( layers_num(iLa) .EQ. 1 ) THEN
87     CALL LAYERS_FLUXCALC( uVel,vVel,theta,iLa,
88 jmc 1.18 & layers_UH, layers_VH,
89     & layers_Hw, layers_Hs,
90     & layers_PIw,layers_PIs,
91     & layers_U, layers_V,
92 jmc 1.16 & myThid )
93     ELSEIF ( layers_num(iLa) .EQ. 2 ) THEN
94     CALL LAYERS_FLUXCALC( uVel,vVel,salt,iLa,
95 jmc 1.18 & layers_UH, layers_VH,
96     & layers_Hw, layers_Hs,
97     & layers_PIw,layers_PIs,
98     & layers_U, layers_V,
99 jmc 1.16 & myThid )
100     ELSEIF ( layers_num(iLa) .EQ. 3 ) THEN
101 dfer 1.10 #ifdef LAYERS_PRHO_REF
102 gforget 1.14 C For layers_num(iLa) = 3, calculate the potential density referenced to
103     C the model level given by layers_krho.
104 jmc 1.16 DO bj=myByLo(myThid),myByHi(myThid)
105     DO bi=myBxLo(myThid),myBxHi(myThid)
106     DO k = 1,Nr
107     CALL FIND_RHO_2D( 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
108     & layers_krho(iLa),
109     & theta(1-OLx,1-OLy,k,bi,bj),
110     & salt(1-OLx,1-OLy,k,bi,bj),
111 jmc 1.18 & prho(1-OLx,1-OLy,k,bi,bj),
112 jmc 1.16 & k, bi, bj, myThid )
113     DO j = 1-OLy,sNy+OLy
114     DO i = 1-OLx,sNx+OLx
115 jmc 1.18 prho(i,j,k,bi,bj) = rhoConst + prho(i,j,k,bi,bj)
116 jmc 1.16 ENDDO
117     ENDDO
118     ENDDO
119 jmc 1.11 ENDDO
120     ENDDO
121 jmc 1.18 CALL LAYERS_FLUXCALC( uVel,vVel, prho, iLa,
122     & layers_UH, layers_VH,
123     & layers_Hw, layers_Hs,
124     & layers_PIw,layers_PIs,
125     & layers_U, layers_V,
126 jmc 1.16 & myThid )
127 gforget 1.14 #endif
128 jmc 1.16 ENDIF
129 rpa 1.1
130 jmc 1.18 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
131     C-- Direct Snap-shot output
132     IF ( DIFFERENT_MULTIPLE(layers_diagFreq,myTime,deltaTClock)
133     & .AND. layers_num(iLa).NE.0 ) THEN
134    
135     IF ( layers_MDSIO ) THEN
136     WRITE(suff,'(I2.2,A1,I10.10)') iLa, '.', myIter
137     #ifdef LAYERS_UFLUX
138     CALL WRITE_FLD_3D_RL( 'layers_UH.', suff, Nlayers,
139     & layers_UH, myIter, myThid )
140     #ifdef LAYERS_THICKNESS
141     CALL WRITE_FLD_3D_RL( 'layers_Hw.', suff, Nlayers,
142     & layers_Hw, myIter, myThid )
143     #endif /* LAYERS_THICKNESS */
144     #endif /* LAYERS_UFLUX */
145     #ifdef LAYERS_VFLUX
146     CALL WRITE_FLD_3D_RL( 'layers_VH.', suff, Nlayers,
147     & layers_VH, myIter, myThid )
148     #ifdef LAYERS_THICKNESS
149     CALL WRITE_FLD_3D_RL( 'layers_Hs.', suff, Nlayers,
150     & layers_Hs, myIter, myThid )
151     #endif /* LAYERS_THICKNESS */
152     #endif /* LAYERS_VFLUX */
153     #ifdef LAYERS_PRHO_REF
154     IF ( layers_num(1).EQ.3 ) THEN
155     CALL WRITE_FLD_3D_RL( 'layers_prho.', suff, Nr,
156     & prho, myIter, myThid )
157     ENDIF
158     #endif /* LAYERS_PRHO_REF */
159     ENDIF
160    
161     c#ifdef ALLOW_MNC
162     c#ifdef LAYERS_MNC
163     c IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
164     c pf(1:1) = 'D'
165     c ELSE
166     c pf(1:1) = 'R'
167     c ENDIF
168     c IF ( layers_MNC) THEN
169     C Do MNC output... But how?
170     c ENDIF
171     c#endif /* LAYERS_MNC */
172     c#endif /* ALLOW_MNC */
173    
174     ENDIF
175    
176     #ifdef ALLOW_DIAGNOSTICS
177     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
178     C-- Fill-in diagnostics
179     IF ( useDiagnostics .AND. layers_num(iLa).NE.0 ) THEN
180    
181     #ifdef LAYERS_UFLUX
182     WRITE(diagName,'(A4,I1,A3)') 'LaUH',iLa,layers_name(iLa)
183     CALL DIAGNOSTICS_FILL( layers_UH,
184     & diagName,0,Nlayers, 0, 1, 1, myThid )
185     # ifdef LAYERS_THICKNESS
186     WRITE(diagName,'(A4,I1,A3)') 'LaHw',iLa,layers_name(iLa)
187     CALL DIAGNOSTICS_FILL( layers_Hw,
188     & diagName,0,Nlayers, 0, 1, 1, myThid )
189     WRITE(diagName,'(A4,I1,A3)') 'LaPw',iLa,layers_name(iLa)
190     CALL DIAGNOSTICS_FILL( layers_PIw,
191     & diagName,0,Nlayers, 0, 1, 1, myThid )
192     WRITE(diagName,'(A4,I1,A3)') 'LaUa',iLa,layers_name(iLa)
193     CALL DIAGNOSTICS_FILL( layers_U,
194     & diagName,0,Nlayers, 0, 1, 1, myThid )
195     # endif
196     #endif /* LAYERS_UFLUX */
197    
198     #ifdef LAYERS_VFLUX
199     WRITE(diagName,'(A4,I1,A3)') 'LaVH',iLa,layers_name(iLa)
200     CALL DIAGNOSTICS_FILL( layers_VH,
201     & diagName,0,Nlayers, 0, 1, 1, myThid )
202     # ifdef LAYERS_THICKNESS
203     WRITE(diagName,'(A4,I1,A3)') 'LaHs',iLa,layers_name(iLa)
204     CALL DIAGNOSTICS_FILL( layers_Hs,
205     & diagName,0,Nlayers, 0, 1, 1, myThid )
206     WRITE(diagName,'(A4,I1,A3)') 'LaPs',iLa,layers_name(iLa)
207     CALL DIAGNOSTICS_FILL( layers_PIs,
208     & diagName,0,Nlayers, 0, 1, 1, myThid )
209     WRITE(diagName,'(A4,I1,A3)') 'LaVa',iLa,layers_name(iLa)
210     CALL DIAGNOSTICS_FILL( layers_V,
211     & diagName,0,Nlayers, 0, 1, 1, myThid )
212     # endif
213     #endif /* LAYERS_VFLUX */
214    
215     ENDIF
216     #endif /* ALLOW_DIAGNOSTICS */
217    
218 rpa 1.1 #ifdef ALLOW_TIMEAVE
219 jmc 1.18 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
220 rpa 1.1 C-- Time-average
221 gforget 1.14 cgf layers_maxNum loop and dimension would be needed for
222     cgf the following and tave output to work beyond iLa.EQ.1
223 jmc 1.16 IF ( layers_taveFreq.GT.0. .AND. iLa.EQ.1 ) THEN
224 gforget 1.13 C --- The tile loops
225 jmc 1.16 DO bj=myByLo(myThid),myByHi(myThid)
226     DO bi=myBxLo(myThid),myBxHi(myThid)
227 rpa 1.1
228     #ifdef LAYERS_UFLUX
229 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_UH_T, layers_UH, Nlayers,
230     & deltaTclock, bi, bj, myThid )
231 rpa 1.1 #ifdef LAYERS_THICKNESS
232 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_Hw_T, layers_Hw, Nlayers,
233     & deltaTclock, bi, bj, myThid )
234     CALL TIMEAVE_CUMULATE( layers_PIw_T, layers_PIw, Nlayers,
235     & deltaTclock, bi, bj, myThid )
236     CALL TIMEAVE_CUMULATE( layers_U_T, layers_U, Nlayers,
237     & deltaTclock, bi, bj, myThid )
238 rpa 1.1 #endif /* LAYERS_THICKNESS */
239     #endif /* LAYERS_UFLUX */
240     #ifdef LAYERS_VFLUX
241 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_VH_T, layers_VH, Nlayers,
242     & deltaTclock, bi, bj, myThid )
243 rpa 1.1 #ifdef LAYERS_THICKNESS
244 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_Hs_T, layers_Hs, Nlayers,
245     & deltaTclock, bi, bj, myThid )
246     CALL TIMEAVE_CUMULATE( layers_PIs_T, layers_PIs, Nlayers,
247     & deltaTclock, bi, bj, myThid )
248     CALL TIMEAVE_CUMULATE( layers_V_T, layers_V, Nlayers,
249     & deltaTclock, bi, bj, myThid )
250 rpa 1.1 #endif /* LAYERS_THICKNESS */
251     #endif /* LAYERS_VFLUX */
252    
253 dfer 1.10 #ifdef LAYERS_PRHO_REF
254 jmc 1.16 IF ( layers_num(iLa) .EQ. 3 )
255     & CALL TIMEAVE_CUMULATE( prho_tave, prho, Nr,
256     & deltaTclock, bi, bj, myThid )
257 dfer 1.10 #endif /* LAYERS_PRHO_REF */
258    
259 jmc 1.16 layers_TimeAve(bi,bj)=layers_TimeAve(bi,bj)+deltaTclock
260 rpa 1.1
261     C --- End bi,bj loop
262 jmc 1.16 ENDDO
263     ENDDO
264     ENDIF
265 gforget 1.13 #endif /* ALLOW_TIMEAVE */
266 rpa 1.1
267 jmc 1.16 ENDDO !DO iLa=1,layers_maxNum
268 gforget 1.14
269 rpa 1.1 #endif /* ALLOW_LAYERS */
270    
271     RETURN
272     END

  ViewVC Help
Powered by ViewVC 1.1.22