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

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

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


Revision 1.4 - (hide annotations) (download)
Mon Dec 28 02:37:52 2009 UTC (14 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63a, checkpoint63b, checkpoint63c, checkpoint63, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x
Changes since 1.3: +56 -46 lines
- use standard RW pkg S/R to write the output
  and remove local write_fld version (which had a problem).

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

  ViewVC Help
Powered by ViewVC 1.1.22