/[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.5 - (show annotations) (download)
Wed Oct 19 01:28:45 2011 UTC (12 years, 7 months ago) by dfer
Branch: MAIN
CVS Tags: checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g
Changes since 1.4: +23 -3 lines
Include potential density as new coordinate (Thanks to David Munday)

1 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_output.F,v 1.4 2009/12/28 02:37:52 jmc 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
43 C !LOCAL VARIABLES:
44 C == Local variables ==
45 LOGICAL DIFFERENT_MULTIPLE
46 EXTERNAL DIFFERENT_MULTIPLE
47 CHARACTER*(MAX_LEN_MBUF) suff
48 #ifdef ALLOW_TIMEAVE
49 INTEGER bi, bj
50 #endif
51 #ifdef ALLOW_MNC
52 CHARACTER*(1) pf
53 #endif
54
55 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
56
57 #ifdef ALLOW_MNC
58 IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
59 pf(1:1) = 'D'
60 ELSE
61 pf(1:1) = 'R'
62 ENDIF
63 #endif /* ALLOW_MNC */
64
65 IF ( DIFFERENT_MULTIPLE(layers_diagFreq,myTime,deltaTClock)
66 & .AND. myIter .GT. 0 ) THEN
67
68 IF ( layers_MDSIO ) THEN
69 WRITE(suff,'(I10.10)') myIter
70 #ifdef LAYERS_UFLUX
71 CALL WRITE_FLD_3D_RL( 'layers_UFlux.', suff, Nlayers,
72 & layers_UFlux, myIter, myThid )
73 #ifdef LAYERS_THICKNESS
74 CALL WRITE_FLD_3D_RL( 'layers_HU.', suff, Nlayers,
75 & layers_HU, myIter, myThid )
76 #endif /* LAYERS_THICKNESS */
77 #endif /* LAYERS_UFLUX */
78 #ifdef LAYERS_VFLUX
79 CALL WRITE_FLD_3D_RL( 'layers_VFlux.', suff, Nlayers,
80 & layers_VFlux, myIter, myThid )
81 #ifdef LAYERS_THICKNESS
82 CALL WRITE_FLD_3D_RL( 'layers_HV.', suff, Nlayers,
83 & layers_HV, myIter, myThid )
84 #endif /* LAYERS_THICKNESS */
85 #endif /* LAYERS_VFLUX */
86 #ifdef LAYERS_PRHO_REF
87 CALL WRITE_FLD_3D_RL( 'layers_prho.', suff, Nr,
88 & prho, myIter, myThid )
89 #endif /* LAYERS_PRHO_REF */
90 ENDIF
91 ENDIF
92
93 #ifdef ALLOW_MNC
94 #ifdef LAYERS_MNC
95 IF ( layers_MNC) THEN
96 C Do MNC output...
97 C But how?
98 ENDIF
99 #endif /* LAYERS_MNC */
100 #endif /* ALLOW_MNC */
101
102 #ifdef ALLOW_TIMEAVE
103 IF ( layers_taveFreq.GT.0. ) THEN
104
105 c set arrays to zero if first timestep
106 IF ( myIter.EQ.nIter0 ) THEN
107 DO bj = myByLo(myThid), myByHi(myThid)
108 DO bi = myBxLo(myThid), myBxHi(myThid)
109 #ifdef LAYERS_UFLUX
110 CALL TIMEAVE_RESET(layers_UFlux_T,Nlayers,bi,bj,myThid)
111 #ifdef LAYERS_THICKNESS
112 CALL TIMEAVE_RESET(layers_HU_T,Nlayers,bi,bj,myThid)
113 #endif /* LAYERS_THICKNESS */
114 #endif /* LAYERS_UFLUX */
115
116 #ifdef LAYERS_VFLUX
117 CALL TIMEAVE_RESET(layers_VFlux_T,Nlayers,bi,bj,myThid)
118 #ifdef LAYERS_THICKNESS
119 CALL TIMEAVE_RESET(layers_HV_T,Nlayers,bi,bj,myThid)
120 #endif /* LAYERS_THICKNESS */
121 #endif /* LAYERS_VFLUX */
122 layers_TimeAve(bi,bj) = 0.
123 ENDDO
124 ENDDO
125
126 C Dump files and restart average computation if needed
127 ELSEIF (
128 & DIFFERENT_MULTIPLE(layers_taveFreq,myTime,deltaTClock)
129 & ) THEN
130
131 C Normalize by integrated time
132 DO bj = myByLo(myThid), myByHi(myThid)
133 DO bi = myBxLo(myThid), myBxHi(myThid)
134
135 #ifdef LAYERS_UFLUX
136 CALL TIMEAVE_NORMALIZE( layers_UFlux_T, layers_timeave,
137 & Nlayers, bi, bj, myThid )
138 #ifdef LAYERS_THICKNESS
139 CALL TIMEAVE_NORMALIZE( layers_HU_T, layers_timeave,
140 & Nlayers, bi, bj, myThid )
141 #endif /* LAYERS_THICKNESS */
142 #endif /* LAYERS_UFLUX */
143
144 #ifdef LAYERS_VFLUX
145 CALL TIMEAVE_NORMALIZE( layers_VFlux_T, layers_timeave,
146 & Nlayers, bi, bj, myThid )
147 #ifdef LAYERS_THICKNESS
148 CALL TIMEAVE_NORMALIZE( layers_HV_T, layers_timeave,
149 & Nlayers, bi, bj, myThid )
150 #endif /* LAYERS_THICKNESS */
151 #endif /* LAYERS_VFLUX */
152
153 #ifdef LAYERS_PRHO_REF
154 CALL TIMEAVE_NORMALIZE( prho_tave, layers_timeave,
155 & Nr, bi, bj, myThid )
156 #endif /* LAYERS_PRHO_REF */
157
158 ENDDO
159 ENDDO
160
161 IF ( layers_MDSIO ) THEN
162 WRITE(suff,'(I10.10)') myIter
163 #ifdef LAYERS_UFLUX
164 CALL WRITE_FLD_3D_RL( 'layers_UFlux-tave.', suff, Nlayers,
165 & layers_UFlux_T, myIter, myThid )
166 #ifdef LAYERS_THICKNESS
167 CALL WRITE_FLD_3D_RL( 'layers_HU-tave.', suff, Nlayers,
168 & layers_HU_T, myIter, myThid )
169 #endif /* LAYERS_THICKNESS */
170 #endif /* LAYERS_UFLUX */
171 #ifdef LAYERS_VFLUX
172 CALL WRITE_FLD_3D_RL( 'layers_VFlux-tave.', suff, Nlayers,
173 & layers_VFlux_T, myIter, myThid )
174 #ifdef LAYERS_THICKNESS
175 CALL WRITE_FLD_3D_RL( 'layers_HV-tave.', suff, Nlayers,
176 & layers_HV_T, myIter, myThid )
177 #endif /* LAYERS_THICKNESS */
178 #endif /* LAYERS_VFLUX */
179
180 #ifdef LAYERS_PRHO_REF
181 CALL WRITE_FLD_3D_RL( 'layers_prho-tave.', suff, Nr,
182 & prho_tave, myIter, myThid )
183 #endif /* LAYERS_PRHO_REF */
184
185 ENDIF
186
187 #ifdef ALLOW_MNC
188 C Do MNC output.
189 #endif
190
191 C Reset averages to zero
192 DO bj = myByLo(myThid), myByHi(myThid)
193 DO bi = myBxLo(myThid), myBxHi(myThid)
194 #ifdef LAYERS_UFLUX
195 CALL TIMEAVE_RESET(layers_UFlux_T,Nlayers,bi,bj,myThid)
196 #ifdef LAYERS_THICKNESS
197 CALL TIMEAVE_RESET(layers_HU_T,Nlayers,bi,bj,myThid)
198 #endif /* LAYERS_THICKNESS */
199 #endif /* LAYERS_UFLUX */
200
201 #ifdef LAYERS_VFLUX
202 CALL TIMEAVE_RESET(layers_VFlux_T,Nlayers,bi,bj,myThid)
203 #ifdef LAYERS_THICKNESS
204 CALL TIMEAVE_RESET(layers_HV_T,Nlayers,bi,bj,myThid)
205 #endif /* LAYERS_THICKNESS */
206 #endif /* LAYERS_VFLUX */
207
208 #ifdef LAYERS_PRHO_REF
209 CALL TIMEAVE_RESET(prho_tave,Nr,bi,bj,myThid)
210 #endif /* LAYERS_PRHO_REF */
211
212 layers_TimeAve(bi,bj) = 0.
213 ENDDO
214 ENDDO
215
216 C-- end of bloc: if time is a multiple of layers_taveFreq
217 ENDIF
218
219 ENDIF
220 #endif /* ALLOW_TIMEAVE */
221
222 #endif /* ALLOW_LAYERS */
223
224 RETURN
225 END

  ViewVC Help
Powered by ViewVC 1.1.22