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

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

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


Revision 1.10 - (show annotations) (download)
Thu Oct 18 13:04:59 2012 UTC (11 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.9: +54 -55 lines
- fix time-ave output (broken since layers_maxNum addition, last Sept).
- snapshot: output all type of layer (add layer type to file name)

1 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_output.F,v 1.9 2012/10/17 18:49:15 rpa Exp $
2 C $Name: $
3
4 #include "LAYERS_OPTIONS.h"
5
6 CBOP 0
7 C !ROUTINE: LAYERS_OUTPUT
8
9 C !INTERFACE:
10 SUBROUTINE LAYERS_OUTPUT( myTime, myIter, myThid )
11
12 C !DESCRIPTION: \bv
13 C *==========================================================*
14 C | SUBROUTINE LAYERS_OUTPUT
15 C | o general routine for LAYERS output
16 C *==========================================================*
17 C | write time-average & snap-shot output
18 C *==========================================================*
19 C \ev
20
21 C !USES:
22 IMPLICIT NONE
23
24 C === Global variables ===
25 #include "SIZE.h"
26 #include "EEPARAMS.h"
27 #include "PARAMS.h"
28 #include "LAYERS_SIZE.h"
29 #include "LAYERS.h"
30
31 C !INPUT PARAMETERS:
32 C == Routine arguments ==
33 C myTime :: Current time of simulation ( s )
34 C myIter :: Iteration number
35 C myThid :: my Thread Id number
36 _RL myTime
37 INTEGER myIter
38 INTEGER myThid
39 CEOP
40
41 #ifdef ALLOW_LAYERS
42 C !LOCAL VARIABLES:
43 C == Local variables ==
44 LOGICAL DIFFERENT_MULTIPLE
45 EXTERNAL DIFFERENT_MULTIPLE
46 CHARACTER*(MAX_LEN_MBUF) suff
47 INTEGER iLa
48 #ifdef ALLOW_TIMEAVE
49 INTEGER bi, bj
50 #endif
51 c#ifdef ALLOW_MNC
52 c CHARACTER*(1) pf
53 c#endif
54
55 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
56
57 #ifdef ALLOW_DIAGNOSTICS
58 IF ( useDiagnostics )
59 & CALL LAYERS_DIAGNOSTICS_FILL( myThid )
60 #endif
61
62 c#ifdef ALLOW_MNC
63 c IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
64 c pf(1:1) = 'D'
65 c ELSE
66 c pf(1:1) = 'R'
67 c ENDIF
68 c#endif /* ALLOW_MNC */
69
70 IF ( DIFFERENT_MULTIPLE(layers_diagFreq,myTime,deltaTClock)
71 & .AND. myIter .GT. 0 ) THEN
72
73 DO iLa=1,layers_maxNum
74
75 IF ( layers_MDSIO ) THEN
76 WRITE(suff,'(I2.2,A1,I10.10)') iLa, '.', myIter
77 #ifdef LAYERS_UFLUX
78 CALL WRITE_FLD_3D_RL( 'layers_UH.', suff, Nlayers,
79 & layers_UH(1-OLx,1-OLy,1,1,1,iLa),
80 & myIter, myThid )
81 #ifdef LAYERS_THICKNESS
82 CALL WRITE_FLD_3D_RL( 'layers_Hw.', suff, Nlayers,
83 & layers_Hw(1-OLx,1-OLy,1,1,1,iLa),
84 & myIter, myThid )
85 #endif /* LAYERS_THICKNESS */
86 #endif /* LAYERS_UFLUX */
87 #ifdef LAYERS_VFLUX
88 CALL WRITE_FLD_3D_RL( 'layers_VH.', suff, Nlayers,
89 & layers_VH(1-OLx,1-OLy,1,1,1,iLa),
90 & myIter, myThid )
91 #ifdef LAYERS_THICKNESS
92 CALL WRITE_FLD_3D_RL( 'layers_Hs.', suff, Nlayers,
93 & layers_Hs(1-OLx,1-OLy,1,1,1,iLa),
94 & myIter, myThid )
95 #endif /* LAYERS_THICKNESS */
96 #endif /* LAYERS_VFLUX */
97 #ifdef LAYERS_PRHO_REF
98 IF ( layers_num(1).EQ.3 ) THEN
99 CALL WRITE_FLD_3D_RL( 'layers_prho.', suff, Nr,
100 & prho(1-OLx,1-OLy,1,1,1,iLa),
101 & myIter, myThid )
102 ENDIF
103 #endif /* LAYERS_PRHO_REF */
104 ENDIF
105
106 c#ifdef ALLOW_MNC
107 c#ifdef LAYERS_MNC
108 c IF ( layers_MNC) THEN
109 C Do MNC output...
110 C But how?
111 c ENDIF
112 c#endif /* LAYERS_MNC */
113 c#endif /* ALLOW_MNC */
114
115 ENDDO
116 ENDIF
117
118 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
119
120 #ifdef ALLOW_TIMEAVE
121 IF ( layers_taveFreq.GT.0. ) THEN
122 cgf layers_maxNum loop and dimension would be needed for
123 cgf the following and tave output to work beyond iLa=1
124 c DO iLa=1,layers_maxNum
125 iLa=1
126
127 c set arrays to zero if first timestep
128 IF ( myIter.EQ.nIter0 ) THEN
129 DO bj = myByLo(myThid), myByHi(myThid)
130 DO bi = myBxLo(myThid), myBxHi(myThid)
131 #ifdef LAYERS_UFLUX
132 CALL TIMEAVE_RESET(layers_UH_T,Nlayers,bi,bj,myThid)
133 #ifdef LAYERS_THICKNESS
134 CALL TIMEAVE_RESET(layers_Hw_T,Nlayers,bi,bj,myThid)
135 #endif /* LAYERS_THICKNESS */
136 #endif /* LAYERS_UFLUX */
137
138 #ifdef LAYERS_VFLUX
139 CALL TIMEAVE_RESET(layers_VH_T,Nlayers,bi,bj,myThid)
140 #ifdef LAYERS_THICKNESS
141 CALL TIMEAVE_RESET(layers_Hs_T,Nlayers,bi,bj,myThid)
142 #endif /* LAYERS_THICKNESS */
143 #endif /* LAYERS_VFLUX */
144 layers_TimeAve(bi,bj) = 0.
145 ENDDO
146 ENDDO
147
148 C Dump files and restart average computation if needed
149 ELSEIF (
150 & DIFFERENT_MULTIPLE(layers_taveFreq,myTime,deltaTClock)
151 & ) THEN
152
153 C Normalize by integrated time
154 DO bj = myByLo(myThid), myByHi(myThid)
155 DO bi = myBxLo(myThid), myBxHi(myThid)
156
157 #ifdef LAYERS_UFLUX
158 CALL TIMEAVE_NORMALIZE( layers_UH_T, layers_timeave,
159 & Nlayers, bi, bj, myThid )
160 #ifdef LAYERS_THICKNESS
161 CALL TIMEAVE_NORMALIZE( layers_Hw_T, layers_timeave,
162 & Nlayers, bi, bj, myThid )
163 #endif /* LAYERS_THICKNESS */
164 #endif /* LAYERS_UFLUX */
165
166 #ifdef LAYERS_VFLUX
167 CALL TIMEAVE_NORMALIZE( layers_VH_T, layers_timeave,
168 & Nlayers, bi, bj, myThid )
169 #ifdef LAYERS_THICKNESS
170 CALL TIMEAVE_NORMALIZE( layers_Hs_T, layers_timeave,
171 & Nlayers, bi, bj, myThid )
172 #endif /* LAYERS_THICKNESS */
173 #endif /* LAYERS_VFLUX */
174
175 #ifdef LAYERS_PRHO_REF
176 IF ( layers_num(1).EQ.3 )
177 & CALL TIMEAVE_NORMALIZE( prho_tave, layers_timeave,
178 & Nr, bi, bj, myThid )
179 #endif /* LAYERS_PRHO_REF */
180
181 ENDDO
182 ENDDO
183
184 IF ( layers_MDSIO ) THEN
185 WRITE(suff,'(I10.10)') myIter
186 #ifdef LAYERS_UFLUX
187 CALL WRITE_FLD_3D_RL( 'layers_UH-tave.', suff, Nlayers,
188 & layers_UH_T, myIter, myThid )
189 #ifdef LAYERS_THICKNESS
190 CALL WRITE_FLD_3D_RL( 'layers_Hw-tave.', suff, Nlayers,
191 & layers_Hw_T, myIter, myThid )
192 #endif /* LAYERS_THICKNESS */
193 #endif /* LAYERS_UFLUX */
194 #ifdef LAYERS_VFLUX
195 CALL WRITE_FLD_3D_RL( 'layers_VH-tave.', suff, Nlayers,
196 & layers_VH_T, myIter, myThid )
197 #ifdef LAYERS_THICKNESS
198 CALL WRITE_FLD_3D_RL( 'layers_Hs-tave.', suff, Nlayers,
199 & layers_Hs_T, myIter, myThid )
200 #endif /* LAYERS_THICKNESS */
201 #endif /* LAYERS_VFLUX */
202
203 #ifdef LAYERS_PRHO_REF
204 IF ( layers_num(1).EQ.3 )
205 & CALL WRITE_FLD_3D_RL( 'layers_prho-tave.', suff, Nr,
206 & prho_tave, myIter, myThid )
207 #endif /* LAYERS_PRHO_REF */
208
209 ENDIF
210
211 c#ifdef ALLOW_MNC
212 C Do MNC output.
213 c#endif
214
215 C Reset averages to zero
216 DO bj = myByLo(myThid), myByHi(myThid)
217 DO bi = myBxLo(myThid), myBxHi(myThid)
218 #ifdef LAYERS_UFLUX
219 CALL TIMEAVE_RESET(layers_UH_T,Nlayers,bi,bj,myThid)
220 #ifdef LAYERS_THICKNESS
221 CALL TIMEAVE_RESET(layers_Hw_T,Nlayers,bi,bj,myThid)
222 #endif /* LAYERS_THICKNESS */
223 #endif /* LAYERS_UFLUX */
224
225 #ifdef LAYERS_VFLUX
226 CALL TIMEAVE_RESET(layers_VH_T,Nlayers,bi,bj,myThid)
227 #ifdef LAYERS_THICKNESS
228 CALL TIMEAVE_RESET(layers_Hs_T,Nlayers,bi,bj,myThid)
229 #endif /* LAYERS_THICKNESS */
230 #endif /* LAYERS_VFLUX */
231
232 #ifdef LAYERS_PRHO_REF
233 IF (layers_num(1) .EQ. 3) THEN
234 CALL TIMEAVE_RESET(prho_tave,Nr,bi,bj,myThid)
235 ENDIF
236 #endif /* LAYERS_PRHO_REF */
237
238 layers_TimeAve(bi,bj) = 0.
239 ENDDO
240 ENDDO
241
242 C-- end of bloc: if time is a multiple of layers_taveFreq
243 ENDIF
244
245 ENDIF
246 #endif /* ALLOW_TIMEAVE */
247
248 #endif /* ALLOW_LAYERS */
249
250 RETURN
251 END

  ViewVC Help
Powered by ViewVC 1.1.22