1 |
C $Header: /u/gcmpack/MITgcm/pkg/diagnostics/diagnostics_readparms.F,v 1.28 2011/06/27 22:23:09 jmc Exp $ |
2 |
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 |
|
16 |
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 |
#include "DIAGSTATS_REGIONS.h" |
24 |
|
25 |
C !INPUT PARAMETERS: |
26 |
INTEGER myThid |
27 |
CEOP |
28 |
|
29 |
C !LOCAL VARIABLES: |
30 |
C ldimLoc :: Max Number of Lists (in data.diagnostics) |
31 |
C kdimLoc :: Max Number of Levels (in data.diagnostics) |
32 |
C fdimLoc :: Max Number of Fields (in data.diagnostics) |
33 |
C frequency :: Frequency (in s) of Output (ouput every "frequency" second) |
34 |
C timePhase :: phase (in s) within the "frequency" period to write output |
35 |
C averagingFreq :: frequency (in s) for periodic averaging interval |
36 |
C averagingPhase :: phase (in s) for periodic averaging interval |
37 |
C repeatCycle :: number of averaging intervals in 1 cycle |
38 |
C missing_value :: missing value for real-type fields in output file |
39 |
C missing_value_int :: missing value for integers in output |
40 |
C levels :: List Output Levels |
41 |
C fields :: List Output Fields |
42 |
C fileName :: List Output Filename |
43 |
C-- for regional-statistics |
44 |
C set_regMask(n) :: region-mask set-index that define the region "n" |
45 |
C val_regMask(n) :: corresponding mask value of region "n" in the region-mask |
46 |
C-- per level statistics output: |
47 |
C stat_freq :: Frequency (in s) of statistics output |
48 |
C stat_phase :: phase (in s) to write statistics output |
49 |
C stat_region :: List of statistics output Regions |
50 |
C stat_fields :: List of statistics output Fields |
51 |
C stat_fName :: List of statistics output Filename |
52 |
INTEGER ldimLoc, kdimLoc, fdimLoc, rdimLoc |
53 |
PARAMETER ( ldimLoc = 2*numLists ) |
54 |
PARAMETER ( kdimLoc = 2*numLevels ) |
55 |
PARAMETER ( fdimLoc = 2*numperList ) |
56 |
PARAMETER ( rdimLoc = nRegions+21 ) |
57 |
_RL frequency(ldimLoc), timePhase(ldimLoc) |
58 |
_RL averagingFreq(ldimLoc), averagingPhase(ldimLoc) |
59 |
INTEGER repeatCycle(ldimLoc) |
60 |
_RL missing_value(ldimLoc) |
61 |
INTEGER missing_value_int(ldimLoc) |
62 |
_RL levels(kdimLoc,ldimLoc) |
63 |
_RL stat_freq(ldimLoc), stat_phase(ldimLoc) |
64 |
CHARACTER*8 fields(fdimLoc,ldimLoc) |
65 |
CHARACTER*8 stat_fields(fdimLoc,ldimLoc) |
66 |
CHARACTER*80 fileName(ldimLoc), blkFilName |
67 |
CHARACTER*80 stat_fname(ldimLoc) |
68 |
CHARACTER*8 fileFlags(ldimLoc) |
69 |
CHARACTER*8 blk8c |
70 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
71 |
CHARACTER*12 suffix |
72 |
INTEGER stat_region(rdimLoc,ldimLoc) |
73 |
INTEGER set_regMask(rdimLoc) |
74 |
_RS val_regMask(rdimLoc) |
75 |
INTEGER ku, stdUnit |
76 |
INTEGER j,k,l,n,m,nf |
77 |
INTEGER iLen, regionCount |
78 |
INTEGER ILNBLNK |
79 |
EXTERNAL ILNBLNK |
80 |
|
81 |
C-- full level output: |
82 |
NAMELIST / DIAGNOSTICS_LIST / |
83 |
& frequency, timePhase, |
84 |
& averagingFreq, averagingPhase, repeatCycle, |
85 |
& missing_value, missing_value_int, |
86 |
& levels, fields, fileName, fileFlags, |
87 |
& dumpAtLast, diag_mnc, useMissingValue, |
88 |
& diag_pickup_read, diag_pickup_write, |
89 |
& diag_pickup_read_mnc, diag_pickup_write_mnc |
90 |
|
91 |
C-- per level statistics output: |
92 |
NAMELIST / DIAG_STATIS_PARMS / |
93 |
& stat_freq, stat_phase, stat_region, stat_fields, |
94 |
& stat_fname, diagSt_mnc, |
95 |
& set_regMask, val_regMask, |
96 |
& diagSt_regMaskFile, nSetRegMskFile |
97 |
|
98 |
C Initialize and Read Diagnostics Namelist |
99 |
_BEGIN_MASTER(myThid) |
100 |
|
101 |
blk8c = ' ' |
102 |
DO k=1,LEN(blkFilName) |
103 |
blkFilName(k:k) = ' ' |
104 |
ENDDO |
105 |
|
106 |
DO l = 1,ldimLoc |
107 |
frequency(l) = 0. |
108 |
timePhase(l) = UNSET_RL |
109 |
averagingFreq(l) = 0. |
110 |
averagingPhase(l)= 0. |
111 |
repeatCycle(l) = 0 |
112 |
fileName(l) = blkFilName |
113 |
missing_value(l) = UNSET_RL |
114 |
missing_value_int(l) = UNSET_I |
115 |
fileFlags(l) = blk8c |
116 |
DO k = 1,kdimLoc |
117 |
levels(k,l) = UNSET_RL |
118 |
ENDDO |
119 |
DO m = 1,fdimLoc |
120 |
fields(m,l) = blk8c |
121 |
ENDDO |
122 |
ENDDO |
123 |
diagLoc_ioUnit = 0 |
124 |
settingDiags = .FALSE. |
125 |
dumpAtLast = .FALSE. |
126 |
diag_mnc = useMNC |
127 |
useMissingValue = .FALSE. |
128 |
diag_pickup_read = .FALSE. |
129 |
diag_pickup_write = .FALSE. |
130 |
diag_pickup_read_mnc = .FALSE. |
131 |
diag_pickup_write_mnc = .FALSE. |
132 |
|
133 |
diagSt_regMaskFile = ' ' |
134 |
nSetRegMskFile = 0 |
135 |
DO k = 1,rdimLoc |
136 |
set_regMask(k) = 0 |
137 |
val_regMask(k) = 0. |
138 |
ENDDO |
139 |
DO l = 1,ldimLoc |
140 |
stat_freq(l) = 0. |
141 |
stat_phase(l) = UNSET_RL |
142 |
stat_fname(l) = blkFilName |
143 |
DO k = 1,rdimLoc |
144 |
stat_region(k,l) = UNSET_I |
145 |
ENDDO |
146 |
DO m = 1,fdimLoc |
147 |
stat_fields(m,l) = blk8c |
148 |
ENDDO |
149 |
ENDDO |
150 |
|
151 |
WRITE(msgBuf,'(2A)') |
152 |
& ' DIAGNOSTICS_READPARMS: opening data.diagnostics' |
153 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,SQUEEZE_RIGHT,1) |
154 |
|
155 |
CALL OPEN_COPY_DATA_FILE('data.diagnostics', |
156 |
& 'DIAGNOSTICS_READPARMS', ku, myThid ) |
157 |
|
158 |
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
159 |
& ' read namelist "diagnostics_list": start' |
160 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
161 |
& SQUEEZE_RIGHT , 1) |
162 |
READ (ku,NML=diagnostics_list) |
163 |
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
164 |
& ' read namelist "diagnostics_list": OK' |
165 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
166 |
& SQUEEZE_RIGHT , 1) |
167 |
|
168 |
C- set default for statistics output according to the main flag |
169 |
diag_mnc = diag_mnc .AND. useMNC |
170 |
diagSt_mnc = diag_mnc |
171 |
|
172 |
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
173 |
& ' read namelist "DIAG_STATIS_PARMS": start' |
174 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
175 |
& SQUEEZE_RIGHT , 1) |
176 |
READ (ku,NML=DIAG_STATIS_PARMS) |
177 |
WRITE(msgBuf,'(2A)') 'S/R DIAGNOSTICS_READPARMS,', |
178 |
& ' read namelist "DIAG_STATIS_PARMS": OK' |
179 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
180 |
& SQUEEZE_RIGHT , 1) |
181 |
|
182 |
CLOSE (ku) |
183 |
|
184 |
C Initialise DIAG_SELECT common block (except pointers) |
185 |
nlists = 0 |
186 |
DO n = 1,numLists |
187 |
freq(n) = 0. |
188 |
phase(n) = 0. |
189 |
averageFreq(n) = 0. |
190 |
averagePhase(n) = 0. |
191 |
averageCycle(n) = 1 |
192 |
nlevels(n) = 0 |
193 |
nfields(n) = 0 |
194 |
fnames(n) = blkFilName |
195 |
misvalFlt(n) = UNSET_RL |
196 |
misvalInt(n) = UNSET_I |
197 |
DO k = 1,numLevels |
198 |
levs(k,n) = 0 |
199 |
ENDDO |
200 |
DO m = 1,numperList |
201 |
flds(m,n) = blk8c |
202 |
ENDDO |
203 |
fflags(n) = blk8c |
204 |
ENDDO |
205 |
|
206 |
C useMNC is confusing (can be T at this point & turned off later, whereas |
207 |
C for all other pkgs, model stops if use${PKG}= T with #undef ALLOW_${PKG}) |
208 |
#ifndef ALLOW_MNC |
209 |
C Fix to avoid running without getting any output: |
210 |
diag_mnc = .FALSE. |
211 |
diagSt_mnc = .FALSE. |
212 |
#endif |
213 |
|
214 |
C Fill Diagnostics Common Block with Namelist Info |
215 |
diagSt_mnc = diagSt_mnc .AND. useMNC |
216 |
diag_mdsio = (.NOT. diag_mnc) .OR. outputTypesInclusive |
217 |
diag_pickup_read_mnc = diag_pickup_read_mnc .AND. diag_mnc |
218 |
diag_pickup_write_mnc = diag_pickup_write_mnc .AND. diag_mnc |
219 |
diag_pickup_read_mdsio = |
220 |
& diag_pickup_read .AND. (.NOT. diag_pickup_read_mnc) |
221 |
diag_pickup_write_mdsio = diag_pickup_write .AND. |
222 |
& ((.NOT. diag_pickup_write_mnc) .OR. outputTypesInclusive) |
223 |
diagSt_ascii = (.NOT. diagSt_mnc) .OR. outputTypesInclusive |
224 |
|
225 |
DO l = 1,ldimLoc |
226 |
iLen = ILNBLNK(fileName(l)) |
227 |
C- Only lists with non-empty file name (iLen>0) are considered |
228 |
IF ( iLen.GE.1 .AND. nlists.LT.numLists ) THEN |
229 |
n = nlists + 1 |
230 |
freq(n) = frequency(l) |
231 |
IF ( timePhase(l).NE. UNSET_RL ) THEN |
232 |
phase(n) = timePhase(l) |
233 |
ELSEIF ( frequency(l) .LT. 0. ) THEN |
234 |
phase(n) = -0.5 _d 0 * frequency(l) |
235 |
ENDIF |
236 |
IF ( averagingFreq(l).GT.0. .AND. repeatCycle(l).GT.1 ) THEN |
237 |
averageFreq(n) = averagingFreq(l) |
238 |
averagePhase(n) = averagingPhase(l) |
239 |
averageCycle(n) = repeatCycle(l) |
240 |
ELSEIF (averagingFreq(l).NE.0. .OR. repeatCycle(l).NE.0) THEN |
241 |
WRITE(msgBuf,'(2A,F18.6,I4)') 'DIAGNOSTICS_READPARMS: ', |
242 |
& 'unvalid Average-Freq & Cycle:', |
243 |
& averagingFreq(l), repeatCycle(l) |
244 |
CALL PRINT_ERROR( msgBuf , myThid ) |
245 |
WRITE(msgBuf,'(2A,I4,2A)') 'DIAGNOSTICS_READPARMS: ', |
246 |
& ' for list l=', l, ', fileName: ', fileName(l) |
247 |
CALL PRINT_ERROR( msgBuf , myThid ) |
248 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
249 |
ELSEIF ( frequency(l) .EQ. 0. ) THEN |
250 |
averageFreq(n) = nTimeSteps*deltaTClock |
251 |
averagePhase(n) = phase(n) |
252 |
ELSEIF ( frequency(l) .GT. 0. ) THEN |
253 |
averageFreq(n) = frequency(l) |
254 |
averagePhase(n) = phase(n) |
255 |
ENDIF |
256 |
IF ( missing_value(l) .NE. UNSET_RL ) |
257 |
& misvalFlt(n) = missing_value(l) |
258 |
IF ( missing_value_int(l) .NE. UNSET_I ) |
259 |
& misvalInt(n) = missing_value_int(l) |
260 |
fnames(n) = fileName (l) |
261 |
fflags(n) = fileFlags(l) |
262 |
nlevels(n) = 0 |
263 |
IF ( levels(1,l).NE.UNSET_RL ) THEN |
264 |
DO k=1,kdimLoc |
265 |
IF ( levels(k,l).NE.UNSET_RL .AND. |
266 |
& nlevels(n).LT.numLevels ) THEN |
267 |
nlevels(n) = nlevels(n) + 1 |
268 |
levs(nlevels(n),n) = levels(k,l) |
269 |
ELSEIF ( levels(k,l).NE.UNSET_RL ) THEN |
270 |
WRITE(msgBuf,'(2A,I4)') 'DIAGNOSTICS_READPARMS: ', |
271 |
& 'Exceed Max.Num. of Levels numLevels=', numLevels |
272 |
CALL PRINT_ERROR( msgBuf , myThid ) |
273 |
WRITE(msgBuf,'(2A,I4,A,F8.0)') 'DIAGNOSTICS_READPARMS: ', |
274 |
& 'when trying to add level(k=', k, ' )=', levels(k,l) |
275 |
CALL PRINT_ERROR( msgBuf , myThid ) |
276 |
WRITE(msgBuf,'(2A,I4,2A)') 'DIAGNOSTICS_READPARMS: ', |
277 |
& ' for list l=', l, ', fileName: ', fileName(l) |
278 |
CALL PRINT_ERROR( msgBuf , myThid ) |
279 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
280 |
ENDIF |
281 |
ENDDO |
282 |
ELSE |
283 |
C- will set levels later, once the Nb of levels of each diag is known |
284 |
nlevels(n) = -1 |
285 |
ENDIF |
286 |
nfields(n) = 0 |
287 |
DO m=1,fdimLoc |
288 |
IF ( fields(m,l).NE.blk8c .AND. |
289 |
& nfields(n).LT.numperList ) THEN |
290 |
nfields(n) = nfields(n) + 1 |
291 |
flds(nfields(n),n) = fields(m,l) |
292 |
ELSEIF ( fields(m,l).NE.blk8c ) THEN |
293 |
WRITE(msgBuf,'(2A,I4)') 'DIAGNOSTICS_READPARMS: ', |
294 |
& 'Exceed Max.Num. of Fields/list numperList=', numperList |
295 |
CALL PRINT_ERROR( msgBuf , myThid ) |
296 |
WRITE(msgBuf,'(2A,I4,3A,I4,2A)') 'DIAGNOSTICS_READPARMS: ', |
297 |
& 'when trying to add field (m=', m, ' ): ',fields(m,l) |
298 |
CALL PRINT_ERROR( msgBuf , myThid ) |
299 |
WRITE(msgBuf,'(2A,I4,2A)') 'DIAGNOSTICS_READPARMS: ', |
300 |
& ' in list l=', l, ', fileName: ', fileName(l) |
301 |
CALL PRINT_ERROR( msgBuf , myThid ) |
302 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
303 |
ENDIF |
304 |
ENDDO |
305 |
nlists = nlists + 1 |
306 |
c write(6,*) 'list summary:',n,nfields(n),nlevels(n) |
307 |
ELSEIF ( iLen.GE.1 ) THEN |
308 |
WRITE(msgBuf,'(2A,I4)') 'DIAGNOSTICS_READPARMS: ', |
309 |
& 'Exceed Max.Num. of list numLists=', numLists |
310 |
CALL PRINT_ERROR( msgBuf , myThid ) |
311 |
WRITE(msgBuf,'(2A,I4)') 'DIAGNOSTICS_READPARMS: ', |
312 |
& 'when trying to add list l=', l |
313 |
CALL PRINT_ERROR( msgBuf , myThid ) |
314 |
WRITE(msgBuf,'(2A,F18.6,2A)') 'DIAGNOSTICS_READPARMS: ', |
315 |
& ' Frq=', frequency(l), ', fileName: ', fileName(l) |
316 |
CALL PRINT_ERROR( msgBuf , myThid ) |
317 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
318 |
ENDIF |
319 |
ENDDO |
320 |
|
321 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
322 |
|
323 |
C- Initialise DIAG_STATS_REGMASK common block (except the mask) |
324 |
nSetRegMask = 0 |
325 |
DO j = 0,nRegions |
326 |
diagSt_kRegMsk(j) = 0 |
327 |
diagSt_vRegMsk(j) = 0. |
328 |
ENDDO |
329 |
C Global statistics (region # 0) |
330 |
diagSt_kRegMsk(0) = 1 |
331 |
|
332 |
C- Initialise DIAG_STATIS common block (except pointers) |
333 |
diagSt_nbLists = 0 |
334 |
DO n = 1,numLists |
335 |
diagSt_freq(n) = 0. |
336 |
diagSt_phase(n) = 0. |
337 |
diagSt_nbFlds(n) = 0 |
338 |
diagSt_ioUnit(n) = 0 |
339 |
diagSt_Fname(n) = blkFilName |
340 |
DO j = 0,nRegions |
341 |
diagSt_region(j,n) = 0 |
342 |
ENDDO |
343 |
DO m = 1,numperList |
344 |
diagSt_Flds(m,n) = blk8c |
345 |
ENDDO |
346 |
ENDDO |
347 |
|
348 |
C Fill Diagnostics Common Block with Namelist Info |
349 |
diagSt_ascii = (.NOT. diagSt_mnc) .OR. outputTypesInclusive |
350 |
|
351 |
C- Region mask correspondence table: |
352 |
C note: this table should be build when regions are defined ; |
353 |
C for now, simpler just to read it from namelist in data.diagnostics |
354 |
j = 0 |
355 |
DO k = 1,rdimLoc |
356 |
IF ( set_regMask(k).NE.0 .OR. val_regMask(k).NE.0. ) THEN |
357 |
j = j+1 |
358 |
IF ( j.LE.nRegions ) THEN |
359 |
diagSt_kRegMsk(j) = set_regMask(k) |
360 |
diagSt_vRegMsk(j) = val_regMask(k) |
361 |
ENDIF |
362 |
ENDIF |
363 |
ENDDO |
364 |
IF ( j.GT.nRegions ) THEN |
365 |
WRITE(msgBuf,'(2A,I4,A)') 'DIAGNOSTICS_READPARMS: ', |
366 |
& 'set_regMask & val_regMask lists assume at least',j,' regions' |
367 |
CALL PRINT_ERROR( msgBuf , myThid ) |
368 |
WRITE(msgBuf,'(2A)') 'DIAGNOSTICS_READPARMS: ', |
369 |
& 'Need to increase "nRegions" in DIAGNOSTICS_SIZE.h' |
370 |
CALL PRINT_ERROR( msgBuf , myThid ) |
371 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
372 |
ENDIF |
373 |
|
374 |
DO l = 1,ldimLoc |
375 |
iLen = ILNBLNK(stat_fname(l)) |
376 |
C- Only lists with non-empty file name (iLen>0) are considered |
377 |
IF ( iLen.GE.1 .AND. diagSt_nbLists.LT.numLists)THEN |
378 |
n = diagSt_nbLists + 1 |
379 |
diagSt_freq(n) = stat_freq(l) |
380 |
IF ( stat_phase(l).NE. UNSET_RL ) THEN |
381 |
diagSt_phase(n) = stat_phase(l) |
382 |
ELSEIF ( stat_freq(l) .LT. 0. ) THEN |
383 |
diagSt_phase(n) = -0.5 _d 0 * stat_freq(l) |
384 |
ENDIF |
385 |
diagSt_Fname(n) = stat_fname(l) |
386 |
regionCount = 0 |
387 |
DO k=1,rdimLoc |
388 |
j = stat_region(k,l) |
389 |
IF ( j.NE.UNSET_I .AND. j.GE.0 .AND. j.LE.nRegions ) THEN |
390 |
IF ( diagSt_region(j,n).EQ.0 ) THEN |
391 |
diagSt_region(j,n) = 1 |
392 |
regionCount = regionCount + 1 |
393 |
ELSE |
394 |
WRITE(msgBuf,'(2A,I4,2A)') |
395 |
& 'DIAGNOSTICS_READPARMS:', |
396 |
& ' in list l=', l, ', stat_fname: ', stat_fname(l) |
397 |
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
398 |
& SQUEEZE_RIGHT , myThid ) |
399 |
WRITE(msgBuf,'(A,I4,A)') |
400 |
& 'DIAGNOSTICS_READPARMS: region=',j, |
401 |
& ' can only be selected once => ignore 2nd selection' |
402 |
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
403 |
& SQUEEZE_RIGHT , myThid ) |
404 |
ENDIF |
405 |
ELSEIF ( j.NE.UNSET_I ) THEN |
406 |
WRITE(msgBuf,'(A,I4,A,I4,2A)') |
407 |
& 'DIAGNOSTICS_READPARMS: region=',j, |
408 |
& ' in list l=', l, ', stat_fname: ', stat_fname(l) |
409 |
CALL PRINT_ERROR( msgBuf , myThid ) |
410 |
WRITE(msgBuf,'(2A,I4,A,I4,2A)') |
411 |
& 'DIAGNOSTICS_READPARMS: ==> exceed Max.Nb of regions', |
412 |
& '(=',nRegions,' )' |
413 |
CALL PRINT_ERROR( msgBuf , myThid ) |
414 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
415 |
ENDIF |
416 |
ENDDO |
417 |
IF ( regionCount.EQ.0 ) THEN |
418 |
C- no region selected => default is Global statistics (region Id: 0) |
419 |
diagSt_region(0,n) = 1 |
420 |
ENDIF |
421 |
diagSt_nbFlds(n) = 0 |
422 |
DO m=1,fdimLoc |
423 |
IF ( stat_fields(m,l).NE.blk8c .AND. |
424 |
& diagSt_nbFlds(n).LT.numperList ) THEN |
425 |
diagSt_nbFlds(n) = diagSt_nbFlds(n) + 1 |
426 |
diagSt_Flds(diagSt_nbFlds(n),n) = stat_fields(m,l) |
427 |
ELSEIF ( stat_fields(m,l).NE.blk8c ) THEN |
428 |
WRITE(msgBuf,'(2A,I4)') 'DIAGNOSTICS_READPARMS: ', |
429 |
& 'Exceed Max.Num. of Fields/list numperList=', numperList |
430 |
CALL PRINT_ERROR( msgBuf , myThid ) |
431 |
WRITE(msgBuf,'(2A,I4,3A,I4,2A)') 'DIAGNOSTICS_READPARMS: ', |
432 |
& 'when trying to add stat_field (m=', m, |
433 |
& ' ): ',stat_fields(m,l) |
434 |
CALL PRINT_ERROR( msgBuf , myThid ) |
435 |
WRITE(msgBuf,'(2A,I4,2A)') 'DIAGNOSTICS_READPARMS: ', |
436 |
& ' in list l=', l, ', stat_fname: ', stat_fname(l) |
437 |
CALL PRINT_ERROR( msgBuf , myThid ) |
438 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
439 |
ENDIF |
440 |
ENDDO |
441 |
diagSt_nbLists = diagSt_nbLists + 1 |
442 |
c write(6,*) 'stat-list summary:',n,diagSt_nbFlds(n),regionCount |
443 |
ELSEIF ( iLen.GE.1 ) THEN |
444 |
WRITE(msgBuf,'(2A,I4)') 'DIAGNOSTICS_READPARMS: ', |
445 |
& 'Exceed Max.Num. of list numLists=', numLists |
446 |
CALL PRINT_ERROR( msgBuf , myThid ) |
447 |
WRITE(msgBuf,'(2A,I4)') 'DIAGNOSTICS_READPARMS: ', |
448 |
& 'when trying to add stat_list l=', l |
449 |
CALL PRINT_ERROR( msgBuf , myThid ) |
450 |
WRITE(msgBuf,'(2A,F18.6,2A)') 'DIAGNOSTICS_READPARMS: ', |
451 |
& ' Frq=', stat_freq(l), ', stat_fname: ', stat_fname(l) |
452 |
CALL PRINT_ERROR( msgBuf , myThid ) |
453 |
STOP 'ABNORMAL END: S/R DIAGNOSTICS_READPARMS' |
454 |
ENDIF |
455 |
ENDDO |
456 |
|
457 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
458 |
C Echo History List Data Structure |
459 |
stdUnit = standardMessageUnit |
460 |
WRITE(msgBuf,'(A)') |
461 |
& ' DIAGNOSTICS_READPARMS: global parameter summary:' |
462 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
463 |
CALL WRITE_0D_L( dumpAtLast, INDEX_NONE, |
464 |
& ' dumpAtLast =',' /* always write time-ave diags at the end */') |
465 |
CALL WRITE_0D_L( diag_mnc, INDEX_NONE, |
466 |
& ' diag_mnc =', ' /* write NetCDF output files */') |
467 |
CALL WRITE_0D_L( useMissingValue, INDEX_NONE, |
468 |
& ' useMissingValue =', ' /* put MissingValue where mask = 0 */') |
469 |
WRITE(msgBuf,'(A)') |
470 |
& '-----------------------------------------------------' |
471 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
472 |
WRITE(msgBuf,'(A)') |
473 |
& ' DIAGNOSTICS_READPARMS: active diagnostics summary:' |
474 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
475 |
WRITE(msgBuf,'(A)') |
476 |
& '-----------------------------------------------------' |
477 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
478 |
DO n = 1,nlists |
479 |
WRITE(msgBuf,'(2A)') 'Creating Output Stream: ', fnames(n) |
480 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
481 |
WRITE(msgBuf,'(2(A,F18.6))') 'Output Frequency:', freq(n), |
482 |
& ' ; Phase: ', phase(n) |
483 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
484 |
WRITE(msgBuf,'(2(A,F18.6),A,I4)') |
485 |
& ' Averaging Freq.:', averageFreq(n), |
486 |
& ' , Phase: ', averagePhase(n), ' , Cycle:', averageCycle(n) |
487 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
488 |
IF ( fflags(n).EQ.blk8c ) THEN |
489 |
WRITE(msgBuf,'(A,1PE20.12,A,I12,3A)') |
490 |
& ' missing value:', misvalFlt(n), |
491 |
& ' ; for integers:', misvalInt(n) |
492 |
ELSE |
493 |
WRITE(msgBuf,'(A,1PE20.12,A,I12,3A)') |
494 |
& ' missing value:', misvalFlt(n), |
495 |
& ' ; for integers:', misvalInt(n), |
496 |
& ' ; F-Flags="', fflags(n),'"' |
497 |
ENDIF |
498 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
499 |
IF ( nlevels(n).EQ.-1 .AND. fflags(n)(2:2).EQ.'I' ) THEN |
500 |
WRITE(msgBuf,'(A)') ' Cumulate all Levels (to be set later)' |
501 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
502 |
ELSEIF ( nlevels(n).EQ.-1 ) THEN |
503 |
WRITE(msgBuf,'(A,A)') ' Levels: ','will be set later' |
504 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
505 |
ELSEIF ( fflags(n)(2:2).EQ.'P' ) THEN |
506 |
DO l=1,nlevels(n),10 |
507 |
m = MIN(nlevels(n),l+9) |
508 |
WRITE(msgBuf,'(A,1P10E10.3)')' interp: ', (levs(k,n),k=l,m) |
509 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
510 |
ENDDO |
511 |
ELSE |
512 |
suffix = ' Levels: ' |
513 |
IF ( fflags(n)(2:2).EQ.'I' ) suffix = ' Sum Levels:' |
514 |
DO l=1,nlevels(n),20 |
515 |
m = MIN(nlevels(n),l+19) |
516 |
WRITE(msgBuf,'(A,20F5.0)') suffix, (levs(k,n),k=l,m) |
517 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
518 |
ENDDO |
519 |
ENDIF |
520 |
DO nf = 1,nfields(n),10 |
521 |
m = MIN(nfields(n),nf+9) |
522 |
WRITE(msgBuf,'(21A)') ' Fields: ',(' ',flds(l,n),l=nf,m) |
523 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
524 |
ENDDO |
525 |
ENDDO |
526 |
WRITE(msgBuf,'(A)') |
527 |
& '-----------------------------------------------------' |
528 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
529 |
WRITE(msgBuf,'(A)') |
530 |
& ' DIAGNOSTICS_READPARMS: statistics diags. summary:' |
531 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
532 |
DO n = 1,diagSt_nbLists |
533 |
WRITE(msgBuf,'(2A)') 'Creating Stats. Output Stream: ', |
534 |
& diagSt_Fname(n) |
535 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
536 |
WRITE(msgBuf,'(2(A,F18.6))') 'Output Frequency:', |
537 |
& diagSt_freq(n), ' ; Phase: ', diagSt_phase(n) |
538 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
539 |
WRITE(msgBuf,'(A)') ' Regions: ' |
540 |
l = 10 |
541 |
DO j=0,nRegions |
542 |
IF ( diagSt_region(j,n).GE.1 ) THEN |
543 |
l = l+3 |
544 |
IF (l.LE.MAX_LEN_MBUF) WRITE(msgBuf(l-2:l),'(I3)') j |
545 |
ENDIF |
546 |
ENDDO |
547 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
548 |
DO nf = 1,diagSt_nbFlds(n),10 |
549 |
m = MIN(diagSt_nbFlds(n),nf+9) |
550 |
WRITE(msgBuf,'(21A)') ' Fields: ', |
551 |
& (' ',diagSt_Flds(l,n),l=nf,m) |
552 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
553 |
ENDDO |
554 |
ENDDO |
555 |
WRITE(msgBuf,'(A)') |
556 |
& '-----------------------------------------------------' |
557 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
558 |
WRITE(msgBuf,'(A)') |
559 |
CALL PRINT_MESSAGE( msgBuf, stdUnit,SQUEEZE_RIGHT, myThid) |
560 |
|
561 |
_END_MASTER(myThid) |
562 |
|
563 |
C-- Everyone else must wait for the parameters to be loaded |
564 |
_BARRIER |
565 |
|
566 |
RETURN |
567 |
END |