1 |
dfer |
1.3 |
C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_output.F,v 1.2 2009/09/20 20:54:48 rpa Exp $ |
2 |
rpa |
1.1 |
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "LAYERS_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
SUBROUTINE LAYERS_OUTPUT( myTime, myIter, myThid ) |
7 |
|
|
|
8 |
|
|
C !DESCRIPTION: \bv |
9 |
|
|
C *==========================================================* |
10 |
|
|
C | SUBROUTINE LAYERS_OUTPUT |
11 |
|
|
C | o general routine for LAYERS output |
12 |
|
|
C *==========================================================* |
13 |
|
|
C | write time-average & snap-shot output |
14 |
|
|
C *==========================================================* |
15 |
|
|
C \ev |
16 |
|
|
|
17 |
|
|
C !USES: |
18 |
|
|
IMPLICIT NONE |
19 |
|
|
|
20 |
|
|
C === Global variables === |
21 |
|
|
#include "SIZE.h" |
22 |
|
|
#include "EEPARAMS.h" |
23 |
|
|
#include "PARAMS.h" |
24 |
|
|
#include "LAYERS_SIZE.h" |
25 |
|
|
#include "LAYERS.h" |
26 |
|
|
|
27 |
|
|
C !INPUT PARAMETERS: |
28 |
|
|
C == Routine arguments == |
29 |
|
|
C myTime :: Current time of simulation ( s ) |
30 |
|
|
C myIter :: Iteration number |
31 |
|
|
C myThid :: my Thread Id number |
32 |
|
|
_RL myTime |
33 |
|
|
INTEGER myIter |
34 |
|
|
INTEGER myThid |
35 |
|
|
CEOP |
36 |
|
|
|
37 |
|
|
#ifdef ALLOW_LAYERS |
38 |
|
|
|
39 |
|
|
C !LOCAL VARIABLES: |
40 |
|
|
C == Local variables == |
41 |
|
|
LOGICAL DIFFERENT_MULTIPLE |
42 |
|
|
EXTERNAL DIFFERENT_MULTIPLE |
43 |
|
|
INTEGER bi, bj, K |
44 |
|
|
CHARACTER*(MAX_LEN_MBUF) suff |
45 |
|
|
CHARACTER*(1) pf |
46 |
|
|
|
47 |
|
|
IF ( writeBinaryPrec .EQ. precFloat64 ) THEN |
48 |
|
|
pf(1:1) = 'D' |
49 |
|
|
ELSE |
50 |
|
|
pf(1:1) = 'R' |
51 |
|
|
ENDIF |
52 |
|
|
|
53 |
dfer |
1.3 |
IF ( DIFFERENT_MULTIPLE(layers_diagFreq,myTime,deltaTClock) |
54 |
rpa |
1.2 |
& .AND. myIter .GT. 0 ) THEN |
55 |
rpa |
1.1 |
|
56 |
|
|
IF ( layers_MDSIO ) THEN |
57 |
|
|
WRITE(suff,'(I10.10)') myIter |
58 |
|
|
#ifdef LAYERS_UFLUX |
59 |
|
|
CALL WRITE_FLD_XYG_RL( |
60 |
|
|
& 'layers_UFlux.',suff,layers_UFlux,myIter,myThid) |
61 |
|
|
#ifdef LAYERS_THICKNESS |
62 |
|
|
CALL WRITE_FLD_XYG_RL( |
63 |
|
|
& 'layers_HU.',suff,layers_HU,myIter,myThid) |
64 |
|
|
#endif /* LAYERS_THICKNESS */ |
65 |
|
|
#endif /* LAYESR_UFLUX */ |
66 |
|
|
#ifdef LAYERS_VFLUX |
67 |
|
|
CALL WRITE_FLD_XYG_RL( |
68 |
|
|
& 'layers_VFlux.',suff,layers_VFlux,myIter,myThid) |
69 |
|
|
#ifdef LAYERS_THICKNESS |
70 |
|
|
CALL WRITE_FLD_XYG_RL( |
71 |
|
|
& 'layers_HV.',suff,layers_HV,myIter,myThid) |
72 |
|
|
#endif /* LAYERS_THICKNESS */ |
73 |
|
|
#endif /* LAYESR_VFLUX */ |
74 |
|
|
ENDIF |
75 |
|
|
ENDIF |
76 |
|
|
|
77 |
|
|
#ifdef ALLOW_MNC |
78 |
|
|
#ifdef LAYERS_MNC |
79 |
|
|
IF ( layers_MNC) THEN |
80 |
|
|
C Do MNC output... |
81 |
|
|
C But how? |
82 |
|
|
ENDIF |
83 |
|
|
#endif /* LAYERS_MNC */ |
84 |
|
|
#endif /* ALLOW_MNC */ |
85 |
|
|
|
86 |
|
|
#ifdef ALLOW_TIMEAVE |
87 |
|
|
|
88 |
dfer |
1.3 |
c set arrays to zero if first timestep |
89 |
|
|
IF ( myIter.EQ.nIter0 ) THEN |
90 |
|
|
DO bj = myByLo(myThid), myByHi(myThid) |
91 |
|
|
DO bi = myBxLo(myThid), myBxHi(myThid) |
92 |
|
|
#ifdef LAYERS_UFLUX |
93 |
|
|
CALL TIMEAVE_RESET(layers_UFlux_T,Nlayers,bi,bj,myThid) |
94 |
|
|
#ifdef LAYERS_THICKNESS |
95 |
|
|
CALL TIMEAVE_RESET(layers_HU_T,Nlayers,bi,bj,myThid) |
96 |
|
|
#endif /* LAYERS_THICKNESS */ |
97 |
|
|
#endif /* LAYERS_UFLUX */ |
98 |
|
|
|
99 |
|
|
#ifdef LAYERS_VFLUX |
100 |
|
|
CALL TIMEAVE_RESET(layers_VFlux_T,Nlayers,bi,bj,myThid) |
101 |
|
|
#ifdef LAYERS_THICKNESS |
102 |
|
|
CALL TIMEAVE_RESET(layers_HV_T,Nlayers,bi,bj,myThid) |
103 |
|
|
#endif /* LAYERS_THICKNESS */ |
104 |
|
|
#endif /* LAYERS_VFLUX */ |
105 |
|
|
DO k=1,Nlayers |
106 |
|
|
layers_TimeAve(k,bi,bj)=0. |
107 |
|
|
ENDDO |
108 |
|
|
ENDDO |
109 |
|
|
ENDDO |
110 |
|
|
|
111 |
rpa |
1.1 |
C Dump files and restart average computation if needed |
112 |
dfer |
1.3 |
ELSEIF ( |
113 |
|
|
& DIFFERENT_MULTIPLE(layers_taveFreq,myTime,deltaTClock) |
114 |
|
|
& ) THEN |
115 |
rpa |
1.1 |
|
116 |
|
|
C Normalize by integrated time |
117 |
|
|
DO bj = myByLo(myThid), myByHi(myThid) |
118 |
|
|
DO bi = myBxLo(myThid), myBxHi(myThid) |
119 |
|
|
|
120 |
|
|
#ifdef LAYERS_UFLUX |
121 |
|
|
CALL TIMEAVE_NORMALIZ(layers_UFlux_T,layers_timeave,Nlayers, |
122 |
|
|
& bi,bj,myThid) |
123 |
|
|
#ifdef LAYERS_THICKNESS |
124 |
|
|
CALL TIMEAVE_NORMALIZ(layers_HU_T,layers_timeave,Nlayers, |
125 |
|
|
& bi,bj,myThid) |
126 |
|
|
#endif /* LAYERS_THICKNESS */ |
127 |
|
|
#endif /* LAYERS_UFLUX */ |
128 |
|
|
|
129 |
|
|
#ifdef LAYERS_VFLUX |
130 |
|
|
CALL TIMEAVE_NORMALIZ(layers_VFlux_T,layers_timeave,Nlayers, |
131 |
|
|
& bi,bj,myThid) |
132 |
|
|
#ifdef LAYERS_THICKNESS |
133 |
|
|
CALL TIMEAVE_NORMALIZ(layers_HV_T,layers_timeave,Nlayers, |
134 |
|
|
& bi,bj,myThid) |
135 |
|
|
#endif /* LAYERS_THICKNESS */ |
136 |
|
|
#endif /* LAYERS_VFLUX */ |
137 |
|
|
|
138 |
|
|
ENDDO |
139 |
|
|
ENDDO |
140 |
|
|
|
141 |
|
|
IF ( layers_MDSIO ) THEN |
142 |
|
|
WRITE(suff,'(I10.10)') myIter |
143 |
|
|
#ifdef LAYERS_UFLUX |
144 |
|
|
CALL WRITE_FLD_XYG_RL( |
145 |
|
|
& 'layers_UFlux-tave.',suff,layers_UFlux_T,myIter,myThid) |
146 |
|
|
#ifdef LAYERS_THICKNESS |
147 |
|
|
CALL WRITE_FLD_XYG_RL( |
148 |
|
|
& 'layers_HU-tave.',suff,layers_HU_T,myIter,myThid) |
149 |
|
|
#endif /* LAYERS_THICKNESS */ |
150 |
|
|
#endif /* LAYERS_UFLUX */ |
151 |
|
|
#ifdef LAYERS_VFLUX |
152 |
|
|
CALL WRITE_FLD_XYG_RL( |
153 |
|
|
& 'layers_VFlux-tave.',suff,layers_VFlux_T,myIter,myThid) |
154 |
|
|
#ifdef LAYERS_THICKNESS |
155 |
|
|
CALL WRITE_FLD_XYG_RL( |
156 |
|
|
& 'layers_HV-tave.',suff,layers_HV_T,myIter,myThid) |
157 |
|
|
#endif /* LAYERS_THICKNESS */ |
158 |
|
|
#endif /* LAYERS_VFLUX */ |
159 |
|
|
ENDIF |
160 |
|
|
|
161 |
|
|
#ifdef ALLOW_MNC |
162 |
|
|
C Do MNC output. |
163 |
|
|
#endif |
164 |
|
|
|
165 |
|
|
C Reset averages to zero |
166 |
|
|
DO bj = myByLo(myThid), myByHi(myThid) |
167 |
|
|
DO bi = myBxLo(myThid), myBxHi(myThid) |
168 |
|
|
#ifdef LAYERS_UFLUX |
169 |
|
|
CALL TIMEAVE_RESET(layers_UFlux_T,Nlayers,bi,bj,myThid) |
170 |
|
|
#ifdef LAYERS_THICKNESS |
171 |
|
|
CALL TIMEAVE_RESET(layers_HU_T,Nlayers,bi,bj,myThid) |
172 |
|
|
#endif /* LAYERS_THICKNESS */ |
173 |
|
|
#endif /* LAYERS_UFLUX */ |
174 |
|
|
|
175 |
|
|
#ifdef LAYERS_VFLUX |
176 |
|
|
CALL TIMEAVE_RESET(layers_VFlux_T,Nlayers,bi,bj,myThid) |
177 |
|
|
#ifdef LAYERS_THICKNESS |
178 |
|
|
CALL TIMEAVE_RESET(layers_HV_T,Nlayers,bi,bj,myThid) |
179 |
|
|
#endif /* LAYERS_THICKNESS */ |
180 |
|
|
#endif /* LAYERS_VFLUX */ |
181 |
|
|
DO k=1,Nlayers |
182 |
|
|
layers_TimeAve(k,bi,bj)=0. |
183 |
|
|
ENDDO |
184 |
|
|
ENDDO |
185 |
|
|
ENDDO |
186 |
|
|
|
187 |
|
|
ENDIF |
188 |
|
|
|
189 |
|
|
|
190 |
|
|
#endif /* ALLOW_TIMEAVE */ |
191 |
|
|
|
192 |
|
|
#endif /* ALLOW_LAYERS */ |
193 |
|
|
|
194 |
|
|
RETURN |
195 |
|
|
END |