/[MITgcm]/MITgcm/pkg/shap_filt/shap_filt_readparms.F
ViewVC logotype

Diff of /MITgcm/pkg/shap_filt/shap_filt_readparms.F

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

revision 1.4 by adcroft, Fri Aug 17 19:09:23 2001 UTC revision 1.4.2.2 by jmc, Wed Jul 31 21:25:02 2002 UTC
# Line 22  C     === Routine arguments === Line 22  C     === Routine arguments ===
22    
23  #ifdef ALLOW_SHAP_FILT  #ifdef ALLOW_SHAP_FILT
24    
25        NAMELIST /SHAP_PARM01/ Shap_funct,        NAMELIST /SHAP_PARM01/
26       &          nShapT,  nShapTrPhys,       & Shap_funct, shap_filt_uvStar, shap_filt_TrStagg,
27       &          nShapUV, nShapUVPhys,       &    nShapT,  nShapTrPhys, Shap_Trtau, Shap_TrLength,
28       &          Shap_Trtau, Shap_TrLength,       &    nShapUV, nShapUVPhys, Shap_uvtau, Shap_uvLength
      &          Shap_uvtau, Shap_uvLength  
29    
30  C     === Local variables ===  C     === Local variables ===
31  C     msgBuf      - Informational/error meesage buffer  C     msgBuf      - Informational/error meesage buffer
# Line 51  c     SHAPIsOn=.TRUE. Line 50  c     SHAPIsOn=.TRUE.
50    
51  C--   Default flags and values for Shapiro Filter  C--   Default flags and values for Shapiro Filter
52        Shap_funct = 2        Shap_funct = 2
53          shap_filt_uvStar =.TRUE.
54          shap_filt_TrStagg=.TRUE.
55        nShapT = 0        nShapT = 0
56        nShapUV = 0        nShapUV = 0
57        nShapTrPhys = 0        nShapTrPhys = 0
# Line 68  C--   Read parameters from open data fil Line 69  C--   Read parameters from open data fil
69        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
70       &                    SQUEEZE_RIGHT , 1)       &                    SQUEEZE_RIGHT , 1)
71    
       write(*,*) 'Shap_funct, nShap_Tr,UV _Phys=',  
      &    Shap_funct, nShapT, nShapUV, nShapTrPhys, nShapUVPhys  
       write(*,*) 'Shap_Trtau,Shap_uvtau=',Shap_Trtau,Shap_uvtau  
72  C--   Close the open data file  C--   Close the open data file
73        CLOSE(iUnit)        CLOSE(iUnit)
74    
75    C- print out some kee parameters :
76          CALL WRITE_0D_I( Shap_funct, INDEX_NONE,
77         & 'Shap_funct =',
78         & '   /* select Shapiro filter function */')  
79          CALL WRITE_0D_I( nShapT , INDEX_NONE,
80         & 'nShapTr =',
81         & '   /* power of Shapiro filter for Tracers */')  
82          CALL WRITE_0D_I( nShapUV, INDEX_NONE,
83         & 'nShapUV =',
84         & '   /* power of Shapiro filter for momentum */')  
85    
86          IF (Shap_funct.EQ.2) THEN
87           CALL WRITE_0D_I( nShapTrPhys, INDEX_NONE,
88         & 'nShapTrPhys =',
89         & '   /* power of physical-space filter (Tracer) */')  
90           CALL WRITE_0D_I( nShapUVPhys, INDEX_NONE,
91         & 'nShapUVPhys =',
92         & '   /* power of physical-space filter (Momentum) */')  
93          ENDIF
94    
95          CALL WRITE_0D_R8( Shap_Trtau, INDEX_NONE,
96         & 'Shap_Trtau =',
97         & '   /* time scale of Shapiro filter (Tracer) */')  
98          CALL WRITE_0D_R8( Shap_TrLength, INDEX_NONE,
99         & 'Shap_TrLength =',
100         & '   /* Length scale of Shapiro filter (Tracer) */')  
101          CALL WRITE_0D_R8( Shap_uvtau, INDEX_NONE,
102         & 'Shap_uvtau =',
103         & '   /* time scale of Shapiro filter (Momentum) */')  
104          CALL WRITE_0D_R8( Shap_uvLength, INDEX_NONE,
105         & 'Shap_uvLength =',
106         & '   /* Length scale of Shapiro filter (Momentum) */')  
107    
108        _END_MASTER(myThid)        _END_MASTER(myThid)
109    
110  C--   Everyone else must wait for the parameters to be loaded  C--   Everyone else must wait for the parameters to be loaded
111        _BARRIER        _BARRIER
112    
113    C--   Check the parameters :
114    
115           IF ( .NOT.shap_filt_uvStar ) THEN
116    
117    C- Notes: applying the filter at the end of the time step (after SOLVE_FOR_P)
118    C    affects the barotropic flow divergence ; this might not be consistent
119    C    with some option of the code.
120    
121            IF ( rigidLid ) THEN
122             WRITE(msgBuf,'(2A)') 'SHAP_FILT with rigidLid ',
123         &                         'needs shap_filt_uvStar=.true.'
124             CALL PRINT_ERROR( msgBuf , 1)
125             STOP 'ABNORMAL END: S/R SHAP_FILT_READPARMS'
126            ELSEIF ( .NOT.exactConserv ) THEN
127             WRITE(msgBuf,'(2A)') 'S/R SHAP_FILT_READPARMS: WARNING <<< ',
128         &    'applying Filter after SOLVE_FOR_P (shap_filt_uvStar=FALSE)'
129             CALL PRINT_MESSAGE(msgBuf, errorMessageUnit, SQUEEZE_RIGHT,1)
130             WRITE(msgBuf,'(2A)') 'S/R SHAP_FILT_READPARMS: WARNING <<< ',
131         &    'requires to recompute Eta after ==> turn on exactConserv '
132             CALL PRINT_MESSAGE(msgBuf, errorMessageUnit, SQUEEZE_RIGHT,1)
133            ENDIF
134    
135           ENDIF                    
136    
137  #endif /* ALLOW_SHAP_FILT */  #endif /* ALLOW_SHAP_FILT */
138        RETURN        RETURN
139        END        END

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

  ViewVC Help
Powered by ViewVC 1.1.22