1 |
jmc |
1.10 |
C $Header: /u/gcmpack/MITgcm/pkg/diagnostics/diagnostics_readparms.F,v 1.9 2005/05/14 20:45:28 jmc Exp $ |
2 |
jmc |
1.1 |
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "DIAG_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
7 |
|
|
CBOP 0 |
8 |
|
|
C !ROUTINE: DIAGNOSTICS_READPARMS |
9 |
|
|
|
10 |
|
|
C !INTERFACE: |
11 |
|
|
SUBROUTINE DIAGNOSTICS_READPARMS(myThid) |
12 |
|
|
|
13 |
|
|
C !DESCRIPTION: |
14 |
|
|
C Read Diagnostics Namelists to specify output sequence. |
15 |
jmc |
1.3 |
|
16 |
jmc |
1.1 |
C !USES: |
17 |
|
|
IMPLICIT NONE |
18 |
|
|
#include "SIZE.h" |
19 |
|
|
#include "EEPARAMS.h" |
20 |
|
|
#include "PARAMS.h" |
21 |
|
|
#include "DIAGNOSTICS_SIZE.h" |
22 |
|
|
#include "DIAGNOSTICS.h" |
23 |
|
|
|
24 |
|
|
C !INPUT PARAMETERS: |
25 |
|
|
INTEGER myThid |
26 |
|
|
CEOP |
27 |
|
|
|
28 |
|
|
C !LOCAL VARIABLES: |
29 |
|
|
C ldimLoc :: Max Number of Lists |
30 |
|
|
C kdimLoc :: Max Number of Levels |
31 |
|
|
C fdimLoc :: Max Number of Fields |
32 |
jmc |
1.10 |
C frequency :: Frequency (in s) of Output (ouput every "frequency" second) |
33 |
|
|
C timePhase :: phase (in s) within the "frequency" period to write output |
34 |
jmc |
1.1 |
C levels :: List Output Levels |
35 |
|
|
C fields :: List Output Fields |
36 |
|
|
C filename :: List Output Filename |
37 |
jmc |
1.10 |
C-- per level statistics output: |
38 |
|
|
C stat_freq :: Frequency (in s) of statistics output |
39 |
|
|
C stat_phase :: phase (in s) to write statistics output |
40 |
|
|
C stat_region :: List of statistics output Regions |
41 |
|
|
C stat_fields :: List of statistics output Fields |
42 |
|
|
C stat_fname :: List of statistics output Filename |
43 |
|
|
INTEGER ldimLoc, kdimLoc, fdimLoc, rdimLoc |
44 |
jmc |
1.1 |
PARAMETER ( ldimLoc = 2*numlists ) |
45 |
|
|
PARAMETER ( kdimLoc = 2*numLevels ) |
46 |
|
|
PARAMETER ( fdimLoc = 2*numperlist ) |
47 |
jmc |
1.10 |
PARAMETER ( rdimLoc = 2*nRegions+1 ) |
48 |
jmc |
1.9 |
_RL frequency(ldimLoc), timePhase(ldimLoc) |
49 |
jmc |
1.1 |
_RL levels(kdimLoc,ldimLoc) |
50 |
jmc |
1.10 |
_RL stat_freq(ldimLoc), stat_phase(ldimLoc) |
51 |
jmc |
1.1 |
CHARACTER*8 fields(fdimLoc,ldimLoc) |
52 |
jmc |
1.10 |
CHARACTER*8 stat_fields(fdimLoc,ldimLoc) |
53 |
jmc |
1.6 |
CHARACTER*80 filename(ldimLoc), blkFilName |
54 |
jmc |
1.10 |
CHARACTER*80 stat_fname(ldimLoc) |
55 |
edhill |
1.7 |
CHARACTER*8 fileflags(ldimLoc) |
56 |
jmc |
1.1 |
CHARACTER*8 blk8c |
57 |
|
|
CHARACTER*(MAX_LEN_MBUF) msgBuf |
58 |
jmc |
1.10 |
INTEGER stat_region(rdimLoc,ldimLoc) |
59 |
jmc |
1.1 |
INTEGER ku, stdUnit |
60 |
jmc |
1.10 |
INTEGER j,k,l,n,m |
61 |
|
|
INTEGER iL, regionCount |
62 |
jmc |
1.1 |
_RL undef, getcon |
63 |
jmc |
1.6 |
INTEGER ILNBLNK |
64 |
|
|
EXTERNAL ILNBLNK |
65 |
jmc |
1.1 |
|
66 |
jmc |
1.10 |
C-- full level output: |
67 |
jmc |
1.1 |
NAMELIST / diagnostics_list / |
68 |
jmc |
1.9 |
& frequency, timePhase, levels, fields, filename, fileflags, |
69 |
edhill |
1.5 |
& diag_mnc, |
70 |
|
|
& diag_pickup_read, diag_pickup_write, |
71 |
|
|
& diag_pickup_read_mnc, diag_pickup_write_mnc |
72 |
jmc |
1.1 |
|
73 |
jmc |
1.10 |
C-- per level statistics output: |
74 |
|
|
NAMELIST / DIAG_STATIS_PARMS / |
75 |
|
|
& stat_freq, stat_phase, stat_region, stat_fields, |
76 |
|
|
& stat_fname, |
77 |
|
|
& diagSt_mnc |
78 |
|
|
|
79 |
jmc |
1.1 |
C Initialize and Read Diagnostics Namelist |
80 |
|
|
_BEGIN_MASTER(myThid) |
81 |
|
|
|
82 |
|
|
undef = getcon('UNDEF') |
83 |
|
|
blk8c = ' ' |
84 |
jmc |
1.6 |
DO k=1,LEN(blkFilName) |
85 |
|
|
blkFilName(k:k) = ' ' |
86 |
|
|
ENDDO |
87 |
jmc |
1.1 |
|
88 |
|
|
DO l = 1,ldimLoc |
89 |
molod |
1.8 |
frequency(l) = 0. |
90 |
jmc |
1.9 |
timePhase(l) = UNSET_RL |
91 |
jmc |
1.6 |
filename (l) = blkFilName |
92 |
edhill |
1.7 |
C eight spaces: 12345678 |
93 |
|
|
fileflags(l)(1:8) = ' ' |
94 |
jmc |
1.1 |
DO k = 1,kdimLoc |
95 |
|
|
levels (k,l) = undef |
96 |
|
|
ENDDO |
97 |
|
|
DO m = 1,fdimLoc |
98 |
|
|
fields (m,l) = blk8c |
99 |
|
|
ENDDO |
100 |
|
|
ENDDO |
101 |
jmc |
1.4 |
diag_mnc = useMNC |
102 |
edhill |
1.5 |
diag_pickup_read = .FALSE. |
103 |
|
|
diag_pickup_write = .FALSE. |
104 |
|
|
diag_pickup_read_mnc = .FALSE. |
105 |
|
|
diag_pickup_write_mnc = .FALSE. |
106 |
jmc |
1.1 |
|
107 |
jmc |
1.10 |
DO l = 1,ldimLoc |
108 |
|
|
stat_freq(l) = 0. |
109 |
|
|
stat_phase(l) = UNSET_RL |
110 |
|
|
stat_fname(l) = blkFilName |
111 |
|
|
DO k = 1,rdimLoc |
112 |
|
|
stat_region(k,l) = UNSET_I |
113 |
|
|
ENDDO |
114 |
|
|
DO m = 1,fdimLoc |
115 |
|
|
stat_fields(m,l) = blk8c |
116 |
|
|
ENDDO |
117 |
|
|
ENDDO |
118 |
|
|
|
119 |
|
|
WRITE(msgBuf,'(2A)') |
120 |
jmc |
1.1 |
& ' DIAGNOSTICS_READPARMS: opening data.diagnostics' |
121 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,SQUEEZE_RIGHT,1) |
122 |
jmc |
1.3 |
|
123 |
|
|
CALL OPEN_COPY_DATA_FILE('data.diagnostics', |
124 |
jmc |
1.1 |
& 'DIAGNOSTICS_READPARMS', ku, myThid ) |
125 |
jmc |
1.10 |
|
126 |
|
|
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
127 |
|
|
& ' read namelist "diagnostics_list": start' |
128 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
129 |
|
|
& SQUEEZE_RIGHT , 1) |
130 |
jmc |
1.1 |
READ (ku,NML=diagnostics_list) |
131 |
jmc |
1.10 |
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
132 |
|
|
& ' read namelist "diagnostics_list": OK' |
133 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
134 |
|
|
& SQUEEZE_RIGHT , 1) |
135 |
|
|
|
136 |
|
|
C- set default for statistics output according to the main flag |
137 |
|
|
diag_mnc = diag_mnc .AND. useMNC |
138 |
|
|
diagSt_mnc = diag_mnc |
139 |
|
|
|
140 |
|
|
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
141 |
|
|
& ' read namelist "DIAG_STATIS_PARMS": start' |
142 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
143 |
|
|
& SQUEEZE_RIGHT , 1) |
144 |
|
|
c STOP 'before reading namelist: DIAG_STATIS_PARMS' |
145 |
|
|
READ (ku,NML=DIAG_STATIS_PARMS) |
146 |
|
|
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
147 |
|
|
& ' read namelist "DIAG_STATIS_PARMS": OK' |
148 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
149 |
|
|
& SQUEEZE_RIGHT , 1) |
150 |
|
|
|
151 |
jmc |
1.1 |
CLOSE (ku) |
152 |
|
|
|
153 |
|
|
C Initialise diag_choices common block |
154 |
|
|
nlists = 0 |
155 |
|
|
DO n = 1,numlists |
156 |
molod |
1.8 |
freq(n) = 0. |
157 |
jmc |
1.9 |
phase(n) = 0. |
158 |
jmc |
1.1 |
nlevels(n) = 0 |
159 |
|
|
nfields(n) = 0 |
160 |
jmc |
1.6 |
fnames(n) = blkFilName |
161 |
jmc |
1.1 |
DO k = 1,numLevels |
162 |
|
|
levs(k,n) = 0 |
163 |
|
|
ENDDO |
164 |
|
|
DO m = 1,numperlist |
165 |
jmc |
1.6 |
flds(m,n) = blk8c |
166 |
jmc |
1.1 |
jdiag(m,n) = 0 |
167 |
|
|
ENDDO |
168 |
|
|
ENDDO |
169 |
|
|
|
170 |
jmc |
1.10 |
C useMNC is confusing (can be T at this point & turned off later, whereas |
171 |
|
|
C for all other pkgs, model stops if use${PKG}= T with #undef ALLOW_${PKG}) |
172 |
|
|
#ifndef ALLOW_MNC |
173 |
|
|
C Fix to avoid running without getting any output: |
174 |
|
|
diag_mnc = .FALSE. |
175 |
|
|
diagSt_mnc = .FALSE. |
176 |
|
|
#endif |
177 |
|
|
|
178 |
jmc |
1.1 |
C Fill Diagnostics Common Block with Namelist Info |
179 |
jmc |
1.10 |
diagSt_mnc = diagSt_mnc .AND. useMNC |
180 |
jmc |
1.4 |
diag_mdsio = (.NOT. diag_mnc) .OR. outputTypesInclusive |
181 |
edhill |
1.5 |
diag_pickup_read_mnc = diag_pickup_read_mnc .AND. diag_mnc |
182 |
|
|
diag_pickup_write_mnc = diag_pickup_write_mnc .AND. diag_mnc |
183 |
jmc |
1.10 |
diag_pickup_read_mdsio = |
184 |
edhill |
1.5 |
& diag_pickup_read .AND. (.NOT. diag_pickup_read_mnc) |
185 |
|
|
diag_pickup_write_mdsio = diag_pickup_write .AND. |
186 |
|
|
& ((.NOT. diag_pickup_write_mnc) .OR. outputTypesInclusive) |
187 |
jmc |
1.10 |
diagSt_ascii = (.NOT. diagSt_mnc) .OR. outputTypesInclusive |
188 |
jmc |
1.3 |
|
189 |
jmc |
1.1 |
DO l = 1,ldimLoc |
190 |
jmc |
1.6 |
iL = ILNBLNK(filename(l)) |
191 |
molod |
1.8 |
IF ( frequency(l).NE.0. .AND. iL.EQ.0 ) THEN |
192 |
jmc |
1.10 |
WRITE(msgBuf,'(2A,I3,A,F17.6)') 'DIAGNOSTICS_READPARMS: ', |
193 |
jmc |
1.6 |
& 'Empty File-name ! (list l=', l, ' ), freq:',frequency(l) |
194 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
195 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
196 |
|
|
ENDIF |
197 |
molod |
1.8 |
IF ( frequency(l).NE.0. .AND. nlists.LT.numlists ) THEN |
198 |
jmc |
1.1 |
n = nlists + 1 |
199 |
|
|
freq(n) = frequency(l) |
200 |
jmc |
1.9 |
IF ( timePhase(l).NE. UNSET_RL ) THEN |
201 |
|
|
phase(n) = timePhase(l) |
202 |
|
|
ELSEIF ( frequency(l) .LT. 0. ) THEN |
203 |
|
|
phase(n) = -0.5 _d 0 * frequency(l) |
204 |
|
|
ENDIF |
205 |
jmc |
1.1 |
fnames(n) = filename (l) |
206 |
edhill |
1.7 |
fflags(n) = fileflags(l) |
207 |
jmc |
1.1 |
nlevels(n) = 0 |
208 |
|
|
IF ( levels(1,l).NE.undef ) THEN |
209 |
|
|
DO k=1,kdimLoc |
210 |
jmc |
1.3 |
IF ( levels(k,l).NE.undef .AND. |
211 |
jmc |
1.1 |
& nlevels(n).LT.numLevels ) THEN |
212 |
|
|
nlevels(n) = nlevels(n) + 1 |
213 |
|
|
levs(nlevels(n),n) = levels(k,l) |
214 |
|
|
ELSEIF ( levels(k,l).NE.undef ) THEN |
215 |
|
|
WRITE(msgBuf,'(2A,I3)') 'DIAGNOSTICS_READPARMS: ', |
216 |
|
|
& 'Exceed Max.Num. of Levels numLevels=', numLevels |
217 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
218 |
|
|
WRITE(msgBuf,'(2A,I3,A,F3.0)') 'DIAGNOSTICS_READPARMS: ', |
219 |
jmc |
1.3 |
& 'when trying to add level(k=', k, ' )=', levels(k,l) |
220 |
jmc |
1.1 |
CALL PRINT_ERROR( msgBuf , myThid ) |
221 |
|
|
WRITE(msgBuf,'(2A,I3,2A)') 'DIAGNOSTICS_READPARMS: ', |
222 |
|
|
& ' for list l=', l, ', filename: ', filename(l) |
223 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
224 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
225 |
|
|
ENDIF |
226 |
|
|
ENDDO |
227 |
|
|
ELSE |
228 |
jmc |
1.3 |
C- will set levels later, once the Nb of levels of each diag is known |
229 |
|
|
nlevels(n) = -1 |
230 |
jmc |
1.1 |
ENDIF |
231 |
|
|
nfields(n) = 0 |
232 |
|
|
DO m=1,fdimLoc |
233 |
jmc |
1.3 |
IF ( fields(m,l).NE.blk8c .AND. |
234 |
jmc |
1.1 |
& nfields(n).LT.numperlist ) THEN |
235 |
|
|
nfields(n) = nfields(n) + 1 |
236 |
|
|
flds(nfields(n),n) = fields(m,l) |
237 |
jmc |
1.2 |
ELSEIF ( fields(m,l).NE.blk8c ) THEN |
238 |
jmc |
1.1 |
WRITE(msgBuf,'(2A,I3)') 'DIAGNOSTICS_READPARMS: ', |
239 |
|
|
& 'Exceed Max.Num. of Fields/list numperlist=', numperlist |
240 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
241 |
|
|
WRITE(msgBuf,'(2A,I3,3A,I3,2A)') 'DIAGNOSTICS_READPARMS: ', |
242 |
|
|
& 'when trying to add field (m=', m, ' ): ',fields(m,l) |
243 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
244 |
|
|
WRITE(msgBuf,'(2A,I3,2A)') 'DIAGNOSTICS_READPARMS: ', |
245 |
|
|
& ' in list l=', l, ', filename: ', filename(l) |
246 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
247 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
248 |
|
|
ENDIF |
249 |
|
|
ENDDO |
250 |
|
|
nlists = nlists + 1 |
251 |
jmc |
1.2 |
c write(6,*) 'list summary:',n,nfields(n),nlevels(n) |
252 |
molod |
1.8 |
ELSEIF ( frequency(l).NE.0. ) THEN |
253 |
jmc |
1.1 |
WRITE(msgBuf,'(2A,I3)') 'DIAGNOSTICS_READPARMS: ', |
254 |
|
|
& 'Exceed Max.Num. of list numlists=', numlists |
255 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
256 |
|
|
WRITE(msgBuf,'(2A,I3)') 'DIAGNOSTICS_READPARMS: ', |
257 |
|
|
& 'when trying to add list l=', l |
258 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
259 |
jmc |
1.10 |
WRITE(msgBuf,'(2A,F17.6,2A)') 'DIAGNOSTICS_READPARMS: ', |
260 |
jmc |
1.1 |
& ' Frq=', frequency(l), ', filename: ', filename(l) |
261 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
262 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
263 |
|
|
ENDIF |
264 |
|
|
ENDDO |
265 |
|
|
|
266 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
267 |
jmc |
1.10 |
|
268 |
|
|
C Initialise DIAG_STATIS common block |
269 |
|
|
diagSt_nbLists = 0 |
270 |
|
|
DO n = 1,numlists |
271 |
|
|
diagSt_freq(n) = 0. |
272 |
|
|
diagSt_phase(n) = 0. |
273 |
|
|
diagSt_nbFlds(n) = 0 |
274 |
|
|
diagSt_ioUnit(n) = 0 |
275 |
|
|
diagSt_Fname(n) = blkFilName |
276 |
|
|
DO j = 0,nRegions |
277 |
|
|
diagSt_region(j,n) = 0 |
278 |
|
|
ENDDO |
279 |
|
|
DO m = 1,numperlist |
280 |
|
|
diagSt_Flds(m,n) = blk8c |
281 |
|
|
jSdiag(m,n) = 0 |
282 |
|
|
ENDDO |
283 |
|
|
ENDDO |
284 |
|
|
|
285 |
|
|
C Fill Diagnostics Common Block with Namelist Info |
286 |
|
|
diagSt_ascii = (.NOT. diagSt_mnc) .OR. outputTypesInclusive |
287 |
|
|
|
288 |
|
|
DO l = 1,ldimLoc |
289 |
|
|
iL = ILNBLNK(stat_fname(l)) |
290 |
|
|
IF ( stat_freq(l).NE.0. .AND. iL.EQ.0 ) THEN |
291 |
|
|
WRITE(msgBuf,'(2A,I3,A,F17.6)') 'DIAGNOSTICS_READPARMS: ', |
292 |
|
|
& 'Empty File-name ! (list l=', l, ' ), stat_freq:',stat_freq(l) |
293 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
294 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
295 |
|
|
ENDIF |
296 |
|
|
IF ( stat_freq(l).NE.0. .AND. diagSt_nbLists.LT.numlists ) THEN |
297 |
|
|
n = diagSt_nbLists + 1 |
298 |
|
|
diagSt_freq(n) = stat_freq(l) |
299 |
|
|
IF ( stat_phase(l).NE. UNSET_RL ) THEN |
300 |
|
|
diagSt_phase(n) = stat_phase(l) |
301 |
|
|
ELSEIF ( stat_freq(l) .LT. 0. ) THEN |
302 |
|
|
diagSt_phase(n) = -0.5 _d 0 * stat_freq(l) |
303 |
|
|
ENDIF |
304 |
|
|
diagSt_Fname(n) = stat_fname(l) |
305 |
|
|
regionCount = 0 |
306 |
|
|
DO k=1,rdimLoc |
307 |
|
|
j = stat_region(k,l) |
308 |
|
|
IF ( j.NE.UNSET_I .AND. j.GE.0 .AND. j.LE.nRegions ) THEN |
309 |
|
|
diagSt_region(j,n) = 1 |
310 |
|
|
regionCount = regionCount + 1 |
311 |
|
|
ELSEIF ( j.NE.UNSET_I ) THEN |
312 |
|
|
WRITE(msgBuf,'(A,I3,A,I3,2A)') |
313 |
|
|
& 'DIAGNOSTICS_READPARMS: region=',j, |
314 |
|
|
& ' in list l=', l, ', stat_fname: ', stat_fname(l) |
315 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
316 |
|
|
WRITE(msgBuf,'(2A,I3,A,I3,2A)') |
317 |
|
|
& 'DIAGNOSTICS_READPARMS: ==> exceed Max.Nb of regions', |
318 |
|
|
& '(=',nRegions,' )' |
319 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
320 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
321 |
|
|
ENDIF |
322 |
|
|
ENDDO |
323 |
|
|
IF ( regionCount.EQ.0 ) THEN |
324 |
|
|
C- no region selected => default is Global statistics (region Id: 0) |
325 |
|
|
diagSt_region(0,n) = 1 |
326 |
|
|
ENDIF |
327 |
|
|
diagSt_nbFlds(n) = 0 |
328 |
|
|
DO m=1,fdimLoc |
329 |
|
|
IF ( stat_fields(m,l).NE.blk8c .AND. |
330 |
|
|
& diagSt_nbFlds(n).LT.numperlist ) THEN |
331 |
|
|
diagSt_nbFlds(n) = diagSt_nbFlds(n) + 1 |
332 |
|
|
diagSt_Flds(diagSt_nbFlds(n),n) = stat_fields(m,l) |
333 |
|
|
ELSEIF ( stat_fields(m,l).NE.blk8c ) THEN |
334 |
|
|
WRITE(msgBuf,'(2A,I3)') 'DIAGNOSTICS_READPARMS: ', |
335 |
|
|
& 'Exceed Max.Num. of Fields/list numperlist=', numperlist |
336 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
337 |
|
|
WRITE(msgBuf,'(2A,I3,3A,I3,2A)') 'DIAGNOSTICS_READPARMS: ', |
338 |
|
|
& 'when trying to add stat_field (m=', m, |
339 |
|
|
& ' ): ',stat_fields(m,l) |
340 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
341 |
|
|
WRITE(msgBuf,'(2A,I3,2A)') 'DIAGNOSTICS_READPARMS: ', |
342 |
|
|
& ' in list l=', l, ', stat_fname: ', stat_fname(l) |
343 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
344 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
345 |
|
|
ENDIF |
346 |
|
|
ENDDO |
347 |
|
|
diagSt_nbLists = diagSt_nbLists + 1 |
348 |
|
|
c write(6,*) 'stat-list summary:',n,diagSt_nbFlds(n),regionCount |
349 |
|
|
ELSEIF ( stat_freq(l).NE.0. ) THEN |
350 |
|
|
WRITE(msgBuf,'(2A,I3)') 'DIAGNOSTICS_READPARMS: ', |
351 |
|
|
& 'Exceed Max.Num. of list numlists=', numlists |
352 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
353 |
|
|
WRITE(msgBuf,'(2A,I3)') 'DIAGNOSTICS_READPARMS: ', |
354 |
|
|
& 'when trying to add stat_list l=', l |
355 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
356 |
|
|
WRITE(msgBuf,'(2A,F17.6,2A)') 'DIAGNOSTICS_READPARMS: ', |
357 |
|
|
& ' Frq=', stat_freq(l), ', stat_fname: ', stat_fname(l) |
358 |
|
|
CALL PRINT_ERROR( msgBuf , myThid ) |
359 |
|
|
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
360 |
|
|
ENDIF |
361 |
|
|
ENDDO |
362 |
|
|
|
363 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
364 |
jmc |
1.1 |
C Echo History List Data Structure |
365 |
|
|
stdUnit = standardMessageUnit |
366 |
|
|
WRITE(msgBuf,'(A)') |
367 |
|
|
& '-----------------------------------------------------' |
368 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
369 |
|
|
WRITE(msgBuf,'(A)') |
370 |
|
|
& ' DIAGNOSTICS_READPARMS: active diagnostics summary:' |
371 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
372 |
|
|
WRITE(msgBuf,'(A)') |
373 |
|
|
& '-----------------------------------------------------' |
374 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
375 |
|
|
DO n = 1,nlists |
376 |
|
|
WRITE(msgBuf,'(2a)') 'Creating Output Stream: ',fnames(n) |
377 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
378 |
jmc |
1.9 |
WRITE(msgBuf,'(2(A,F17.6))') 'Frequency : ',freq(n), |
379 |
|
|
& ' ; Phase: ', phase(n) |
380 |
jmc |
1.1 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
381 |
jmc |
1.3 |
IF ( nlevels(n).EQ.-1 ) THEN |
382 |
|
|
WRITE(msgBuf,'(A,A)') ' Levels: ','will be set later' |
383 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
384 |
|
|
ELSE |
385 |
|
|
DO l=1,nlevels(n),20 |
386 |
jmc |
1.1 |
m = MIN(nlevels(n),l+19) |
387 |
|
|
WRITE(msgBuf,'(A,20F5.0)') ' Levels: ', (levs(k,n),k=l,m) |
388 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
389 |
jmc |
1.3 |
ENDDO |
390 |
|
|
ENDIF |
391 |
jmc |
1.1 |
WRITE(msgBuf,*) 'Fields: ',(' ',flds(l,n),l=1,nfields(n)) |
392 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
393 |
|
|
ENDDO |
394 |
|
|
WRITE(msgBuf,'(A)') |
395 |
|
|
& '-----------------------------------------------------' |
396 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
397 |
|
|
WRITE(msgBuf,'(A)') |
398 |
jmc |
1.10 |
& ' DIAGNOSTICS_READPARMS: statistics diags. summary:' |
399 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
400 |
|
|
DO n = 1,diagSt_nbLists |
401 |
|
|
WRITE(msgBuf,'(2a)') 'Creating Stats. Output Stream: ', |
402 |
|
|
& diagSt_Fname(n) |
403 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
404 |
|
|
WRITE(msgBuf,'(2(A,F17.6))') 'Frequency : ',diagSt_freq(n), |
405 |
|
|
& ' ; Phase: ', diagSt_phase(n) |
406 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
407 |
|
|
WRITE(msgBuf,'(A)') ' Regions : ' |
408 |
|
|
l = 12 |
409 |
|
|
DO j=0,nRegions |
410 |
|
|
IF ( diagSt_region(j,n).GE.1 ) THEN |
411 |
|
|
IF (l+3.LE.MAX_LEN_MBUF) WRITE(msgBuf,'(A,I3)') msgBuf(1:l),j |
412 |
|
|
l = l+3 |
413 |
|
|
ENDIF |
414 |
|
|
ENDDO |
415 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
416 |
|
|
WRITE(msgBuf,*) 'Fields: ', |
417 |
|
|
& (' ',diagSt_Flds(l,n),l=1,diagSt_nbFlds(n)) |
418 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
419 |
|
|
ENDDO |
420 |
|
|
WRITE(msgBuf,'(A)') |
421 |
|
|
& '-----------------------------------------------------' |
422 |
|
|
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
423 |
|
|
WRITE(msgBuf,'(A)') |
424 |
jmc |
1.1 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
425 |
|
|
|
426 |
|
|
_END_MASTER(myThid) |
427 |
|
|
|
428 |
|
|
RETURN |
429 |
|
|
END |