/[MITgcm]/MITgcm/pkg/diagnostics/diagnostics_check.F
ViewVC logotype

Diff of /MITgcm/pkg/diagnostics/diagnostics_check.F

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

revision 1.4 by jmc, Fri May 20 07:28:49 2005 UTC revision 1.8 by jmc, Tue Feb 5 15:31:19 2008 UTC
# Line 12  C     !INTERFACE: Line 12  C     !INTERFACE:
12    
13  C     !DESCRIPTION:  C     !DESCRIPTION:
14  C     Check option and parameter consistency  C     Check option and parameter consistency
15          
16  C     !USES:  C     !USES:
17        IMPLICIT NONE        IMPLICIT NONE
18  #include "SIZE.h"  #include "SIZE.h"
19  #include "EEPARAMS.h"  #include "EEPARAMS.h"
20  #include "PARAMS.h"  #include "PARAMS.h"
21    #include "GRID.h"
22  #include "DIAGNOSTICS_SIZE.h"  #include "DIAGNOSTICS_SIZE.h"
23  #include "DIAGNOSTICS.h"  #include "DIAGNOSTICS.h"
24    
# Line 27  CEOP Line 28  CEOP
28    
29  C     !LOCAL VARIABLES:  C     !LOCAL VARIABLES:
30        CHARACTER*(MAX_LEN_MBUF) msgBuf        CHARACTER*(MAX_LEN_MBUF) msgBuf
31        INTEGER k,l,n,m        INTEGER ld,md,nd
32          INTEGER k,m
33          INTEGER jpoint1, ipoint1, jpoint2, ipoint2
34          _RL     margin
35    
36        _BEGIN_MASTER(myThid)        _BEGIN_MASTER(myThid)
37    
# Line 67  C-    stop if trying to use part of the Line 71  C-    stop if trying to use part of the
71          STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'          STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'
72        ENDIF        ENDIF
73  #endif /* DIAGNOSTICS_HAS_PICKUP */  #endif /* DIAGNOSTICS_HAS_PICKUP */
74          
75  C-    File names:  C-    File names:
76        DO n = 2,nlists        DO ld = 2,nlists
77         DO m = 1,n-1         DO m = 1,ld-1
78          IF ( fnames(n).EQ.fnames(m) ) THEN          IF ( fnames(ld).EQ.fnames(m) ) THEN
79           WRITE(msgBuf,'(2A)') 'DIAGNOSTICS_CHECK: ',           WRITE(msgBuf,'(2A)') 'DIAGNOSTICS_CHECK: ',
80       &            'found 2 identical filenames:'       &            'found 2 identical filenames:'
81           CALL PRINT_ERROR( msgBuf , myThid )           CALL PRINT_ERROR( msgBuf , myThid )
82           WRITE(msgBuf,'(2A,I3,2A)') 'DIAGNOSTICS_CHECK: ',           WRITE(msgBuf,'(2A,I5,2A)') 'DIAGNOSTICS_CHECK: ',
83       &    '1rst (m=', m, ' ): ', fnames(m)       &    '1rst (m=', m, ' ): ', fnames(m)
84           CALL PRINT_ERROR( msgBuf , myThid )           CALL PRINT_ERROR( msgBuf , myThid )
85           WRITE(msgBuf,'(2A,I3,2A)') 'DIAGNOSTICS_CHECK: ',           WRITE(msgBuf,'(2A,I5,2A)') 'DIAGNOSTICS_CHECK: ',
86       &    ' 2nd (n=', n, ' ): ', fnames(n)       &    ' 2nd (n=', ld, ' ): ', fnames(ld)
87           CALL PRINT_ERROR( msgBuf , myThid )           CALL PRINT_ERROR( msgBuf , myThid )
88           STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'           STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'
89          ENDIF          ENDIF
90         ENDDO         ENDDO
91        ENDDO        ENDDO
92    
93        DO n = 2,diagSt_nbLists        DO ld = 2,diagSt_nbLists
94         DO m = 1,n-1         DO m = 1,ld-1
95          IF ( diagSt_Fname(n).EQ.diagSt_Fname(m) ) THEN          IF ( diagSt_Fname(ld).EQ.diagSt_Fname(m) ) THEN
96           WRITE(msgBuf,'(2A)') 'DIAGNOSTICS_CHECK: ',           WRITE(msgBuf,'(2A)') 'DIAGNOSTICS_CHECK: ',
97       &            'found 2 identical stat_fname:'       &            'found 2 identical stat_fname:'
98           CALL PRINT_ERROR( msgBuf , myThid )           CALL PRINT_ERROR( msgBuf , myThid )
99           WRITE(msgBuf,'(2A,I3,2A)') 'DIAGNOSTICS_CHECK: ',           WRITE(msgBuf,'(2A,I5,2A)') 'DIAGNOSTICS_CHECK: ',
100       &    '1rst (m=', m, ' ): ', diagSt_Fname(m)       &    '1rst (m=', m, ' ): ', diagSt_Fname(m)
101           CALL PRINT_ERROR( msgBuf , myThid )           CALL PRINT_ERROR( msgBuf , myThid )
102           WRITE(msgBuf,'(2A,I3,2A)') 'DIAGNOSTICS_CHECK: ',           WRITE(msgBuf,'(2A,I5,2A)') 'DIAGNOSTICS_CHECK: ',
103       &    ' 2nd (n=', n, ' ): ', diagSt_Fname(n)       &    ' 2nd (n=', ld, ' ): ', diagSt_Fname(ld)
104           CALL PRINT_ERROR( msgBuf , myThid )           CALL PRINT_ERROR( msgBuf , myThid )
105           STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'           STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'
106          ENDIF          ENDIF
# Line 105  C-    File names: Line 109  C-    File names:
109    
110  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
111  C-    Check for field that appears 2 times (or more) with differents frequency:  C-    Check for field that appears 2 times (or more) with differents frequency:
112    C     disable this checking since now diagnostics pkg can handle this case.
113    
114        DO n = 2,nlists  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
115         DO m = 1,n-1  
116          IF ( freq(m).NE.freq(n) .OR. phase(m).NE.phase(n) ) THEN  C--   Vertical Interpolation: check for compatibility:
117  C--     once the SWITCH_ONOFF is changed to only turns ON diag with <0 freq  C     better to stop here, rather much later, when trying to write output
118  C                 and CLRDIAG is changed to turns OFF diag with <0 freq,        DO ld = 1,nlists
119  C       then we can allow 1 diag to be used with 2 differents <0 freq.         IF ( fflags(ld)(2:2).EQ.'P' ) THEN
120  C       and this would become:          IF ( fluidIsAir ) THEN
121  c       IF (  ( freq(m).GT.0. .OR. freq(n).GT.0. )  C-    check that interpolated levels are >0 & fall within the domain +/- X %
122  c    &   .AND.( freq(m).NE.freq(n) .OR. phase(m).NE.phase(n) )  C      (needs p>0 for p^kappa ; here take a 10 % margin)
123  c    &     ) THEN            margin = rkSign*(rF(Nr+1)-rF(1))*0.1 _d 0
124           DO k = 1,nActive(n)            DO k=1,nlevels(ld)
125            DO l = 1,nActive(m)             IF ( levs(k,ld)-MAX(rF(1),rF(Nr+1)).GT.margin
126             IF ( flds(k,n).EQ.flds(l,m) ) THEN       &     .OR. levs(k,ld)-MIN(rF(1),rF(Nr+1)).LT.-margin
127              WRITE(msgBuf,'(4A)') 'DIAGNOSTICS_CHECK: ',       &     .OR. levs(k,ld).LE.0. ) THEN
128       &       'field : ',flds(k,n),' use 2 different freq. :'  
129              CALL PRINT_ERROR( msgBuf , myThid )              WRITE(msgBuf,'(2A,I5,2A)') 'DIAGNOSTICS_CHECK: ',
130              WRITE(msgBuf,'(2A,I3,A,2F17.6,2A)') 'DIAGNOSTICS_CHECK: ',       &       'Vertical Interp. for list l=', ld,
131       &       '1rst (m=', m, ' ) freq,phase=', freq(m),phase(m),       &       ', filename: ', fnames(ld)
132       &       ' file:',fnames(m)              CALL PRINT_ERROR( msgBuf , myThid )
133              CALL PRINT_ERROR( msgBuf , myThid )              WRITE(msgBuf,'(2A,I4,3(A,F16.8))') 'DIAGNOSTICS_CHECK: ',
134              WRITE(msgBuf,'(2A,I3,A,2F17.6,2A)') 'DIAGNOSTICS_CHECK: ',       &       ' lev(k=', k, ') p=', levs(k,ld),
135       &       ' 2nd (n=', n, ' ) freq,phase=', freq(n),phase(n),       &       ' not in the domain:',rF(1),' :',rF(Nr+1)
      &       ' file:',fnames(n)  
136              CALL PRINT_ERROR( msgBuf , myThid )              CALL PRINT_ERROR( msgBuf , myThid )
137              STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'              STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'
138             ENDIF             ENDIF
139            ENDDO            ENDDO
140           ENDDO          ELSE
141    C-    p^kappa interpolation: meaningfull only if Atmosphere & P-coordiante
142              WRITE(msgBuf,'(2A)') 'DIAGNOSTICS_CHECK: ',
143         &       'INTERP_VERT not allowed in this config'
144              CALL PRINT_ERROR( msgBuf , myThid )
145               WRITE(msgBuf,'(2A,I5,2A)') 'DIAGNOSTICS_CHECK: ',
146         &       ' for list l=', ld, ', filename: ', fnames(ld)
147              CALL PRINT_ERROR( msgBuf , myThid )
148              STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'
149          ENDIF          ENDIF
150         ENDDO          IF (select_rStar.GT.0) THEN
151        ENDDO  C-    If nonlinear free surf is active, need averaged pressures
152             DO md = 1,nfields(ld)
153        DO n = 2,diagSt_nbLists            nd = jdiag(md,ld)
154         DO m = 1,n-1            CALL DIAGNOSTICS_GET_POINTERS( 'RSURF   ', ld,
155          IF    ( diagSt_freq(m) .NE. diagSt_freq(n) .OR.       &                                   jpoint1, ipoint1, myThid )
156       &          diagSt_phase(m).NE.diagSt_phase(n) ) THEN            IF ( useFIZHI .AND.
157  c       IF (  ( diagSt_freq(m).GT.0. .OR. diagSt_freq(n).GT.0. )       &          gdiag(nd)(10:10) .EQ. 'L') THEN
158  c    &   .AND.( diagSt_freq(m) .NE. diagSt_freq(n) .OR.             CALL DIAGNOSTICS_GET_POINTERS('FIZPRES ', ld,
159  c    &          diagSt_phase(m).NE.diagSt_phase(n) )       &                                   jpoint2, ipoint2, myThid )
160  c    &     ) THEN            ELSE
161           DO k = 1,diagSt_nbActv(n)             CALL DIAGNOSTICS_GET_POINTERS('RCENTER ', ld,
162            DO l = 1,diagSt_nbActv(m)       &                                   jpoint2, ipoint2, myThid )
163             IF ( diagSt_Flds(k,n).EQ.diagSt_Flds(l,m) ) THEN            ENDIF
164              WRITE(msgBuf,'(4A)') 'DIAGNOSTICS_CHECK: ',            IF ( ipoint1.EQ.0 .OR. ipoint2.EQ.0 ) THEN
165       &       'field : ',diagSt_Flds(k,n),' use 2 different stat_freq.:'              WRITE(msgBuf,'(2A,I5)') 'DIAGNOSTICS_CHECK: ',
166              CALL PRINT_ERROR( msgBuf , myThid )       &      'to interpolate diags from output list:', ld
167              WRITE(msgBuf,'(2A,I3,A,2F17.6,2A)') 'DIAGNOSTICS_CHECK: ',              CALL PRINT_ERROR( msgBuf , myThid )
168       &       '1rst (m=', m, ' ) freq,phase=', diagSt_freq(m),              IF ( ipoint1.EQ.0 .AND. jpoint1.EQ.0 ) THEN
169       &       diagSt_phase(m), ' file:', diagSt_Fname(m)                WRITE(msgBuf,'(3A)') 'DIAGNOSTICS_CHECK: ',
170              CALL PRINT_ERROR( msgBuf , myThid )       &        'needs to turn ON surface pressure diagnostic "RSURF   "'
171              WRITE(msgBuf,'(2A,I3,A,2F17.6,2A)') 'DIAGNOSTICS_CHECK: ',                CALL PRINT_ERROR( msgBuf , myThid )
172       &       ' 2nd (n=', n, ' ) freq,phase=', diagSt_freq(n),              ELSEIF ( ipoint1.EQ.0 ) THEN
173       &       diagSt_phase(n), ' file:', diagSt_Fname(n)                WRITE(msgBuf,'(3A)') 'DIAGNOSTICS_CHECK: ',
174              CALL PRINT_ERROR( msgBuf , myThid )       &        'needs surface pressure diagnostic "RSURF   " ',
175         &        'with same output time'
176                  CALL PRINT_ERROR( msgBuf , myThid )
177                ENDIF
178                IF ( ipoint2.EQ.0 .AND. jpoint2.EQ.0 ) THEN
179                  WRITE(msgBuf,'(3A)') 'DIAGNOSTICS_CHECK: ',
180         &        'needs to turn ON  3-D pressure diagnostic "RCENTER "'
181                  CALL PRINT_ERROR( msgBuf , myThid )
182                ELSEIF ( ipoint2.EQ.0 ) THEN
183                  WRITE(msgBuf,'(3A)') 'DIAGNOSTICS_CHECK: ',
184         &        'needs  3-D pressure diagnostic "RCENTER " ',
185         &        'with same output time'
186                  CALL PRINT_ERROR( msgBuf , myThid )
187                ENDIF
188              STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'              STOP 'ABNORMAL END: S/R DIAGNOSTICS_CHECK'
189             ENDIF            ENDIF
           ENDDO  
190           ENDDO           ENDDO
191          ENDIF          ENDIF
192         ENDDO         ENDIF
193        ENDDO        ENDDO
194    
 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|  
   
195        _END_MASTER(myThid)        _END_MASTER(myThid)
196    
197        RETURN        RETURN

Legend:
Removed from v.1.4  
changed lines
  Added in v.1.8

  ViewVC Help
Powered by ViewVC 1.1.22