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

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

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


Revision 1.8 - (show annotations) (download)
Sun Aug 3 03:12:36 2003 UTC (20 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint52l_pre, hrcube4, hrcube5, checkpoint52d_pre, checkpoint52j_pre, checkpoint51o_pre, checkpoint51l_post, checkpoint52l_post, checkpoint52k_post, checkpoint53, checkpoint52, checkpoint52f_post, checkpoint51f_post, checkpoint51t_post, checkpoint51n_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint51s_post, checkpoint51j_post, checkpoint52e_pre, checkpoint52e_post, checkpoint51n_pre, checkpoint53d_post, checkpoint52b_pre, checkpoint51l_pre, checkpoint52m_post, checkpoint51q_post, checkpoint52b_post, checkpoint52c_post, checkpoint51h_pre, checkpoint52f_pre, checkpoint53c_post, branchpoint-genmake2, checkpoint51r_post, checkpoint51i_post, checkpoint53a_post, checkpoint52d_post, checkpoint52a_pre, checkpoint52i_post, checkpoint51i_pre, checkpoint52h_pre, checkpoint53f_post, checkpoint52j_post, branch-netcdf, checkpoint52n_post, checkpoint53b_pre, checkpoint51e_post, checkpoint51o_post, checkpoint51f_pre, checkpoint53b_post, checkpoint52a_post, checkpoint51g_post, ecco_c52_e35, checkpoint51m_post, checkpoint53d_pre, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-genmake2, branch-nonh, tg2-branch, netcdf-sm0, checkpoint51n_branch
Changes since 1.7: +7 -3 lines
add diagnostic (instantaneous) of Shapiro Filter effect on T,S & UV

1 C $Header: /u/gcmpack/MITgcm/pkg/shap_filt/shap_filt_readparms.F,v 1.7 2002/07/31 21:01:25 jmc Exp $
2 C $Name: $
3
4 #include "SHAP_FILT_OPTIONS.h"
5
6 SUBROUTINE SHAP_FILT_READPARMS( myThid )
7 C /==========================================================\
8 C | SUBROUTINE SHAP_FILT_READPARMS |
9 C | o Routine to initialize Shapiro Filter parameters |
10 C |==========================================================|
11 C \==========================================================/
12 IMPLICIT NONE
13
14 C === Global variables ===
15 #include "SIZE.h"
16 #include "EEPARAMS.h"
17 #include "PARAMS.h"
18 #include "SHAP_FILT.h"
19
20 C === Routine arguments ===
21 INTEGER myThid
22
23 #ifdef ALLOW_SHAP_FILT
24
25 NAMELIST /SHAP_PARM01/
26 & Shap_funct, shap_filt_uvStar, shap_filt_TrStagg,
27 & nShapT, nShapTrPhys, Shap_Trtau, Shap_TrLength,
28 & nShapUV, nShapUVPhys, Shap_uvtau, Shap_uvLength,
29 & Shap_noSlip, Shap_diagFreq
30
31 C === Local variables ===
32 C msgBuf - Informational/error meesage buffer
33 C iUnit - Work variable for IO unit number
34 CHARACTER*(MAX_LEN_MBUF) msgBuf
35 INTEGER iUnit
36
37 C-- SHAP_FILT_READPARMS has been called so we know that
38 C the package is active.
39 c SHAPIsOn=.TRUE.
40
41 _BEGIN_MASTER(myThid)
42
43 WRITE(msgBuf,'(A)') ' SHAP_FILT_READPARMS: opening data.shap'
44 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
45 & SQUEEZE_RIGHT , 1)
46
47 CALL OPEN_COPY_DATA_FILE(
48 I 'data.shap', 'SHAP_FILT_READPARMS',
49 O iUnit,
50 I myThid )
51
52 C-- Default flags and values for Shapiro Filter
53 Shap_funct = 2
54 shap_filt_uvStar =.TRUE.
55 shap_filt_TrStagg=.TRUE.
56 nShapT = 0
57 nShapUV = 0
58 nShapTrPhys = 0
59 nShapUVPhys = 0
60 Shap_Trtau = deltaTtracer
61 Shap_TrLength = 0.
62 Shap_uvtau = deltaTMom
63 Shap_TrLength = 0.
64 Shap_noSlip = 0.
65 Shap_diagFreq = diagFreq
66
67 C-- Read parameters from open data file
68 READ(UNIT=iUnit,NML=SHAP_PARM01)
69
70 WRITE(msgBuf,'(A)')
71 & ' SHAP_FILT_READPARMS: finished reading data.shap'
72 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
73 & SQUEEZE_RIGHT , 1)
74
75 C-- Close the open data file
76 CLOSE(iUnit)
77
78 C- print out some kee parameters :
79 CALL WRITE_0D_I( Shap_funct, INDEX_NONE,
80 & 'Shap_funct =',
81 & ' /* select Shapiro filter function */')
82 CALL WRITE_0D_I( nShapT , INDEX_NONE,
83 & 'nShapTr =',
84 & ' /* power of Shapiro filter for Tracers */')
85 CALL WRITE_0D_I( nShapUV, INDEX_NONE,
86 & 'nShapUV =',
87 & ' /* power of Shapiro filter for momentum */')
88
89 IF (Shap_funct.EQ.2) THEN
90 CALL WRITE_0D_I( nShapTrPhys, INDEX_NONE,
91 & 'nShapTrPhys =',
92 & ' /* power of physical-space filter (Tracer) */')
93 CALL WRITE_0D_I( nShapUVPhys, INDEX_NONE,
94 & 'nShapUVPhys =',
95 & ' /* power of physical-space filter (Momentum) */')
96 ENDIF
97
98 CALL WRITE_0D_R8( Shap_Trtau, INDEX_NONE,
99 & 'Shap_Trtau =',
100 & ' /* time scale of Shapiro filter (Tracer) */')
101 CALL WRITE_0D_R8( Shap_TrLength, INDEX_NONE,
102 & 'Shap_TrLength =',
103 & ' /* Length scale of Shapiro filter (Tracer) */')
104 CALL WRITE_0D_R8( Shap_uvtau, INDEX_NONE,
105 & 'Shap_uvtau =',
106 & ' /* time scale of Shapiro filter (Momentum) */')
107 CALL WRITE_0D_R8( Shap_uvLength, INDEX_NONE,
108 & 'Shap_uvLength =',
109 & ' /* Length scale of Shapiro filter (Momentum) */')
110 CALL WRITE_0D_R8( Shap_noSlip, INDEX_NONE,
111 & 'Shap_noSlip =',
112 & ' /* No-slip parameter (0=Free-slip ; 1=No-slip)*/')
113 CALL WRITE_0D_R8( Shap_diagFreq, INDEX_NONE,
114 & 'Shap_diagFreq =',
115 & ' /* Frequency^-1 for diagnostic output (s)*/')
116
117 _END_MASTER(myThid)
118
119 C-- Everyone else must wait for the parameters to be loaded
120 _BARRIER
121
122 C-- Check the Options :
123 #ifndef USE_OLD_SHAPIRO_FILTERS
124 #ifdef NO_SLIP_SHAP
125 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
126 WRITE(msgBuf,'(2A)') 'SHAP_FILT: CPP-option NO_SLIP_SHAP',
127 & ' only in OLD_SHAPIRO S/R ;'
128 CALL PRINT_ERROR( msgBuf , 1)
129 WRITE(msgBuf,'(2A)') ' ==> use parameter Shap_noSlip=1. ',
130 & '(in "data.shap") instead'
131 CALL PRINT_ERROR( msgBuf , 1)
132 STOP 'ABNORMAL END: S/R SHAP_FILT_READPARMS'
133 #endif
134 #endif
135
136 C-- Check the parameters :
137
138 IF ( .NOT.shap_filt_uvStar ) THEN
139
140 C- Notes: applying the filter at the end of the time step (after SOLVE_FOR_P)
141 C affects the barotropic flow divergence ; this might not be consistent
142 C with some option of the code.
143
144 IF ( rigidLid ) THEN
145 WRITE(msgBuf,'(2A)') 'SHAP_FILT with rigidLid ',
146 & 'needs shap_filt_uvStar=.true.'
147 CALL PRINT_ERROR( msgBuf , 1)
148 STOP 'ABNORMAL END: S/R SHAP_FILT_READPARMS'
149 ELSEIF ( .NOT.exactConserv ) THEN
150 WRITE(msgBuf,'(2A)') 'S/R SHAP_FILT_READPARMS: WARNING <<< ',
151 & 'applying Filter after SOLVE_FOR_P (shap_filt_uvStar=FALSE)'
152 CALL PRINT_MESSAGE(msgBuf, errorMessageUnit, SQUEEZE_RIGHT,1)
153 WRITE(msgBuf,'(2A)') 'S/R SHAP_FILT_READPARMS: WARNING <<< ',
154 & 'requires to recompute Eta after ==> turn on exactConserv '
155 CALL PRINT_MESSAGE(msgBuf, errorMessageUnit, SQUEEZE_RIGHT,1)
156 ENDIF
157
158 ENDIF
159
160 #endif /* ALLOW_SHAP_FILT */
161 RETURN
162 END

  ViewVC Help
Powered by ViewVC 1.1.22