/[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.2 - (hide annotations) (download)
Sun Sep 20 20:54:48 2009 UTC (14 years, 9 months ago) by rpa
Branch: MAIN
CVS Tags: checkpoint61v
Changes since 1.1: +3 -3 lines
Fixed bugs, volume is now conserved

1 rpa 1.2 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_output.F,v 1.1 2009/09/16 21:25:47 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     IF ( DIFFERENT_MULTIPLE(dumpFreq,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     C Dump files and restart average computation if needed
89     IF ( DIFFERENT_MULTIPLE(taveFreq,myTime,deltaTClock)
90 rpa 1.2 & .AND. myIter .GT. 0 ) THEN
91 rpa 1.1
92     C Normalize by integrated time
93     DO bj = myByLo(myThid), myByHi(myThid)
94     DO bi = myBxLo(myThid), myBxHi(myThid)
95    
96     #ifdef LAYERS_UFLUX
97     CALL TIMEAVE_NORMALIZ(layers_UFlux_T,layers_timeave,Nlayers,
98     & bi,bj,myThid)
99     #ifdef LAYERS_THICKNESS
100     CALL TIMEAVE_NORMALIZ(layers_HU_T,layers_timeave,Nlayers,
101     & bi,bj,myThid)
102     #endif /* LAYERS_THICKNESS */
103     #endif /* LAYERS_UFLUX */
104    
105     #ifdef LAYERS_VFLUX
106     CALL TIMEAVE_NORMALIZ(layers_VFlux_T,layers_timeave,Nlayers,
107     & bi,bj,myThid)
108     #ifdef LAYERS_THICKNESS
109     CALL TIMEAVE_NORMALIZ(layers_HV_T,layers_timeave,Nlayers,
110     & bi,bj,myThid)
111     #endif /* LAYERS_THICKNESS */
112     #endif /* LAYERS_VFLUX */
113    
114     ENDDO
115     ENDDO
116    
117     IF ( layers_MDSIO ) THEN
118     WRITE(suff,'(I10.10)') myIter
119     #ifdef LAYERS_UFLUX
120     CALL WRITE_FLD_XYG_RL(
121     & 'layers_UFlux-tave.',suff,layers_UFlux_T,myIter,myThid)
122     #ifdef LAYERS_THICKNESS
123     CALL WRITE_FLD_XYG_RL(
124     & 'layers_HU-tave.',suff,layers_HU_T,myIter,myThid)
125     #endif /* LAYERS_THICKNESS */
126     #endif /* LAYERS_UFLUX */
127     #ifdef LAYERS_VFLUX
128     CALL WRITE_FLD_XYG_RL(
129     & 'layers_VFlux-tave.',suff,layers_VFlux_T,myIter,myThid)
130     #ifdef LAYERS_THICKNESS
131     CALL WRITE_FLD_XYG_RL(
132     & 'layers_HV-tave.',suff,layers_HV_T,myIter,myThid)
133     #endif /* LAYERS_THICKNESS */
134     #endif /* LAYERS_VFLUX */
135     ENDIF
136    
137     #ifdef ALLOW_MNC
138     C Do MNC output.
139     #endif
140    
141     C Reset averages to zero
142     DO bj = myByLo(myThid), myByHi(myThid)
143     DO bi = myBxLo(myThid), myBxHi(myThid)
144     #ifdef LAYERS_UFLUX
145     CALL TIMEAVE_RESET(layers_UFlux_T,Nlayers,bi,bj,myThid)
146     #ifdef LAYERS_THICKNESS
147     CALL TIMEAVE_RESET(layers_HU_T,Nlayers,bi,bj,myThid)
148     #endif /* LAYERS_THICKNESS */
149     #endif /* LAYERS_UFLUX */
150    
151     #ifdef LAYERS_VFLUX
152     CALL TIMEAVE_RESET(layers_VFlux_T,Nlayers,bi,bj,myThid)
153     #ifdef LAYERS_THICKNESS
154     CALL TIMEAVE_RESET(layers_HV_T,Nlayers,bi,bj,myThid)
155     #endif /* LAYERS_THICKNESS */
156     #endif /* LAYERS_VFLUX */
157     DO k=1,Nlayers
158     layers_TimeAve(k,bi,bj)=0.
159     ENDDO
160     ENDDO
161     ENDDO
162    
163     ENDIF
164    
165    
166     #endif /* ALLOW_TIMEAVE */
167    
168     #endif /* ALLOW_LAYERS */
169    
170     RETURN
171     END

  ViewVC Help
Powered by ViewVC 1.1.22