1 |
C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_output.F,v 1.5 2011/10/19 01:28:45 dfer 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 |
IF (LAYER_nb .EQ. 3) THEN |
88 |
CALL WRITE_FLD_3D_RL( 'layers_prho.', suff, Nr, |
89 |
& prho, myIter, myThid ) |
90 |
ENDIF |
91 |
#endif /* LAYERS_PRHO_REF */ |
92 |
ENDIF |
93 |
ENDIF |
94 |
|
95 |
#ifdef ALLOW_MNC |
96 |
#ifdef LAYERS_MNC |
97 |
IF ( layers_MNC) THEN |
98 |
C Do MNC output... |
99 |
C But how? |
100 |
ENDIF |
101 |
#endif /* LAYERS_MNC */ |
102 |
#endif /* ALLOW_MNC */ |
103 |
|
104 |
#ifdef ALLOW_TIMEAVE |
105 |
IF ( layers_taveFreq.GT.0. ) THEN |
106 |
|
107 |
c set arrays to zero if first timestep |
108 |
IF ( myIter.EQ.nIter0 ) THEN |
109 |
DO bj = myByLo(myThid), myByHi(myThid) |
110 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
111 |
#ifdef LAYERS_UFLUX |
112 |
CALL TIMEAVE_RESET(layers_UFlux_T,Nlayers,bi,bj,myThid) |
113 |
#ifdef LAYERS_THICKNESS |
114 |
CALL TIMEAVE_RESET(layers_HU_T,Nlayers,bi,bj,myThid) |
115 |
#endif /* LAYERS_THICKNESS */ |
116 |
#endif /* LAYERS_UFLUX */ |
117 |
|
118 |
#ifdef LAYERS_VFLUX |
119 |
CALL TIMEAVE_RESET(layers_VFlux_T,Nlayers,bi,bj,myThid) |
120 |
#ifdef LAYERS_THICKNESS |
121 |
CALL TIMEAVE_RESET(layers_HV_T,Nlayers,bi,bj,myThid) |
122 |
#endif /* LAYERS_THICKNESS */ |
123 |
#endif /* LAYERS_VFLUX */ |
124 |
layers_TimeAve(bi,bj) = 0. |
125 |
ENDDO |
126 |
ENDDO |
127 |
|
128 |
C Dump files and restart average computation if needed |
129 |
ELSEIF ( |
130 |
& DIFFERENT_MULTIPLE(layers_taveFreq,myTime,deltaTClock) |
131 |
& ) THEN |
132 |
|
133 |
C Normalize by integrated time |
134 |
DO bj = myByLo(myThid), myByHi(myThid) |
135 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
136 |
|
137 |
#ifdef LAYERS_UFLUX |
138 |
CALL TIMEAVE_NORMALIZE( layers_UFlux_T, layers_timeave, |
139 |
& Nlayers, bi, bj, myThid ) |
140 |
#ifdef LAYERS_THICKNESS |
141 |
CALL TIMEAVE_NORMALIZE( layers_HU_T, layers_timeave, |
142 |
& Nlayers, bi, bj, myThid ) |
143 |
#endif /* LAYERS_THICKNESS */ |
144 |
#endif /* LAYERS_UFLUX */ |
145 |
|
146 |
#ifdef LAYERS_VFLUX |
147 |
CALL TIMEAVE_NORMALIZE( layers_VFlux_T, layers_timeave, |
148 |
& Nlayers, bi, bj, myThid ) |
149 |
#ifdef LAYERS_THICKNESS |
150 |
CALL TIMEAVE_NORMALIZE( layers_HV_T, layers_timeave, |
151 |
& Nlayers, bi, bj, myThid ) |
152 |
#endif /* LAYERS_THICKNESS */ |
153 |
#endif /* LAYERS_VFLUX */ |
154 |
|
155 |
#ifdef LAYERS_PRHO_REF |
156 |
IF (LAYER_nb .EQ. 3) THEN |
157 |
CALL TIMEAVE_NORMALIZE( prho_tave, layers_timeave, |
158 |
& Nr, bi, bj, myThid ) |
159 |
ENDIF |
160 |
#endif /* LAYERS_PRHO_REF */ |
161 |
|
162 |
ENDDO |
163 |
ENDDO |
164 |
|
165 |
IF ( layers_MDSIO ) THEN |
166 |
WRITE(suff,'(I10.10)') myIter |
167 |
#ifdef LAYERS_UFLUX |
168 |
CALL WRITE_FLD_3D_RL( 'layers_UFlux-tave.', suff, Nlayers, |
169 |
& layers_UFlux_T, myIter, myThid ) |
170 |
#ifdef LAYERS_THICKNESS |
171 |
CALL WRITE_FLD_3D_RL( 'layers_HU-tave.', suff, Nlayers, |
172 |
& layers_HU_T, myIter, myThid ) |
173 |
#endif /* LAYERS_THICKNESS */ |
174 |
#endif /* LAYERS_UFLUX */ |
175 |
#ifdef LAYERS_VFLUX |
176 |
CALL WRITE_FLD_3D_RL( 'layers_VFlux-tave.', suff, Nlayers, |
177 |
& layers_VFlux_T, myIter, myThid ) |
178 |
#ifdef LAYERS_THICKNESS |
179 |
CALL WRITE_FLD_3D_RL( 'layers_HV-tave.', suff, Nlayers, |
180 |
& layers_HV_T, myIter, myThid ) |
181 |
#endif /* LAYERS_THICKNESS */ |
182 |
#endif /* LAYERS_VFLUX */ |
183 |
|
184 |
#ifdef LAYERS_PRHO_REF |
185 |
IF (LAYER_nb .EQ. 3) THEN |
186 |
CALL WRITE_FLD_3D_RL( 'layers_prho-tave.', suff, Nr, |
187 |
& prho_tave, myIter, myThid ) |
188 |
ENDIF |
189 |
#endif /* LAYERS_PRHO_REF */ |
190 |
|
191 |
ENDIF |
192 |
|
193 |
#ifdef ALLOW_MNC |
194 |
C Do MNC output. |
195 |
#endif |
196 |
|
197 |
C Reset averages to zero |
198 |
DO bj = myByLo(myThid), myByHi(myThid) |
199 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
200 |
#ifdef LAYERS_UFLUX |
201 |
CALL TIMEAVE_RESET(layers_UFlux_T,Nlayers,bi,bj,myThid) |
202 |
#ifdef LAYERS_THICKNESS |
203 |
CALL TIMEAVE_RESET(layers_HU_T,Nlayers,bi,bj,myThid) |
204 |
#endif /* LAYERS_THICKNESS */ |
205 |
#endif /* LAYERS_UFLUX */ |
206 |
|
207 |
#ifdef LAYERS_VFLUX |
208 |
CALL TIMEAVE_RESET(layers_VFlux_T,Nlayers,bi,bj,myThid) |
209 |
#ifdef LAYERS_THICKNESS |
210 |
CALL TIMEAVE_RESET(layers_HV_T,Nlayers,bi,bj,myThid) |
211 |
#endif /* LAYERS_THICKNESS */ |
212 |
#endif /* LAYERS_VFLUX */ |
213 |
|
214 |
#ifdef LAYERS_PRHO_REF |
215 |
IF (LAYER_nb .EQ. 3) THEN |
216 |
CALL TIMEAVE_RESET(prho_tave,Nr,bi,bj,myThid) |
217 |
ENDIF |
218 |
#endif /* LAYERS_PRHO_REF */ |
219 |
|
220 |
layers_TimeAve(bi,bj) = 0. |
221 |
ENDDO |
222 |
ENDDO |
223 |
|
224 |
C-- end of bloc: if time is a multiple of layers_taveFreq |
225 |
ENDIF |
226 |
|
227 |
ENDIF |
228 |
#endif /* ALLOW_TIMEAVE */ |
229 |
|
230 |
#endif /* ALLOW_LAYERS */ |
231 |
|
232 |
RETURN |
233 |
END |