1 |
C $Header: /u/gcmpack/MITgcm/model/src/config_summary.F,v 1.50 2004/06/09 14:03:35 adcroft Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: CONFIG_SUMMARY |
8 |
C !INTERFACE: |
9 |
SUBROUTINE CONFIG_SUMMARY( myThid ) |
10 |
C !DESCRIPTION: \bv |
11 |
C *=========================================================* |
12 |
C | SUBROUTINE CONFIG_SUMMARY |
13 |
C | o Summarize model parameter settings. |
14 |
C *=========================================================* |
15 |
C | This routine writes a tabulated summary of the kernel |
16 |
C | model configuration. Information describes all the |
17 |
C | parameter setting in force and the meaning and units of |
18 |
C | those parameters. Individal packages report a similar |
19 |
C | table for each package using the same format as employed |
20 |
C | here. If parameters are missing or incorrectly described |
21 |
C | or dimensioned please contact support@mitgcm.org |
22 |
C *=========================================================* |
23 |
C \ev |
24 |
|
25 |
C !USES: |
26 |
IMPLICIT NONE |
27 |
C === Global variables === |
28 |
#include "SIZE.h" |
29 |
#include "EEPARAMS.h" |
30 |
#include "PARAMS.h" |
31 |
#include "EOS.h" |
32 |
#include "GRID.h" |
33 |
#include "DYNVARS.h" |
34 |
|
35 |
C !INPUT/OUTPUT PARAMETERS: |
36 |
C == Routine arguments == |
37 |
C myThid - Number of this instance of CONFIG_SUMMARY |
38 |
INTEGER myThid |
39 |
CEndOfInterface |
40 |
|
41 |
C !LOCAL VARIABLES: |
42 |
C == Local variables == |
43 |
C msgBuf :: Temp. for building output string. |
44 |
C I,J,K :: Loop counters. |
45 |
C bi,bj :: Tile loop counters. |
46 |
C xcoord :: Temps. for building lists of values for uni-dimensionally |
47 |
C ycoord :: varying parameters. |
48 |
C zcoord :: |
49 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
50 |
INTEGER I,J,K |
51 |
INTEGER bi, bj |
52 |
_RL xcoord(Nx) |
53 |
_RL ycoord(Ny) |
54 |
_RL rcoord(Nr+1) |
55 |
INTEGER coordLine |
56 |
INTEGER tileLine |
57 |
CEOP |
58 |
|
59 |
|
60 |
_BARRIER |
61 |
_BEGIN_MASTER(myThid) |
62 |
|
63 |
WRITE(msgBuf,'(A)') |
64 |
&'// =======================================================' |
65 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
66 |
& SQUEEZE_RIGHT , 1) |
67 |
WRITE(msgBuf,'(A)') '// Model configuration' |
68 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
69 |
& SQUEEZE_RIGHT , 1) |
70 |
WRITE(msgBuf,'(A)') |
71 |
&'// =======================================================' |
72 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
73 |
& SQUEEZE_RIGHT , 1) |
74 |
|
75 |
WRITE(msgBuf,'(A)') '// ' |
76 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
77 |
& SQUEEZE_RIGHT , 1) |
78 |
WRITE(msgBuf,'(A)') |
79 |
& '// "Physical" paramters ( PARM01 in namelist ) ' |
80 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
81 |
& SQUEEZE_RIGHT , 1) |
82 |
WRITE(msgBuf,'(A)') '// ' |
83 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
84 |
& SQUEEZE_RIGHT , 1) |
85 |
WRITE(msgBuf,'(A,A40)') 'buoyancyRelation = ', buoyancyRelation |
86 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
87 |
& SQUEEZE_RIGHT , 1) |
88 |
CALL WRITE_1D_R8( tRef, Nr, INDEX_K,'tRef =', |
89 |
&' /* Reference temperature profile ( oC or oK ) */') |
90 |
CALL WRITE_1D_R8( sRef, Nr, INDEX_K,'sRef =', |
91 |
&' /* Reference salinity profile ( ppt ) */') |
92 |
CALL WRITE_0D_R8( viscAh, INDEX_NONE,'viscAh =', |
93 |
&' /* Lateral eddy viscosity ( m^2/s ) */') |
94 |
CALL WRITE_0D_R8( viscAhMax, INDEX_NONE,'viscAhMax =', |
95 |
&' /* Maximum lateral eddy viscosity ( m^2/s ) */') |
96 |
CALL WRITE_0D_R8( viscAhGrid, INDEX_NONE,'viscAhGrid =', |
97 |
&' /* Grid dependent lateral eddy viscosity ( non-dim. ) */') |
98 |
CALL WRITE_0D_R8( viscC2leith, INDEX_NONE,'viscC2leith =', |
99 |
&' /* Leith harmonic viscosity factor ( non-dom. ) */') |
100 |
CALL WRITE_0D_R8( viscA4, INDEX_NONE,'viscA4 =', |
101 |
&' /* Lateral biharmonic viscosity ( m^4/s ) */') |
102 |
CALL WRITE_0D_R8( viscA4Max, INDEX_NONE,'viscA4Max =', |
103 |
&' /* Maximum biharmonic viscosity ( m^2/s ) */') |
104 |
CALL WRITE_0D_R8( viscA4Grid, INDEX_NONE,'viscA4Grid =', |
105 |
&' /* Grid dependent biharmonic viscosity ( non-dim. ) */') |
106 |
CALL WRITE_0D_R8( viscC4leith, INDEX_NONE,'viscC4leith =', |
107 |
&' /* Leith biharmonic viscosity factor ( non-dom. ) */') |
108 |
CALL WRITE_0D_L( no_slip_sides, INDEX_NONE, |
109 |
& 'no_slip_sides =', ' /* Viscous BCs: No-slip sides */') |
110 |
c IF ( viscAz .NE. UNSET_RL ) THEN |
111 |
c CALL WRITE_0D_R8( viscAz, INDEX_NONE,'viscAz =', |
112 |
c & ' /* Vertical eddy viscosity ( m^2/s ) */') |
113 |
c ENDIF |
114 |
c IF ( viscAp .NE. UNSET_RL ) THEN |
115 |
c CALL WRITE_0D_R8( viscAp, INDEX_NONE,'viscAp =', |
116 |
c & ' /* Vertical eddy viscosity ( Pa^2/s ) */') |
117 |
c ENDIF |
118 |
CALL WRITE_0D_R8( viscAr, INDEX_NONE,'viscAr =', |
119 |
&' /* Vertical eddy viscosity ( units of r^2/s ) */') |
120 |
CALL WRITE_0D_R8( diffKhT, INDEX_NONE,'diffKhT =', |
121 |
&' /* Laplacian diffusion of heat laterally ( m^2/s ) */') |
122 |
CALL WRITE_0D_R8( diffK4T, INDEX_NONE,'diffK4T =', |
123 |
&' /* Bihaarmonic diffusion of heat laterally ( m^4/s ) */') |
124 |
c CALL WRITE_0D_R8( diffKzT, INDEX_NONE,'diffKzT =', |
125 |
c &' /* Laplacian diffusion of heat vertically ( m^2/s ) */') |
126 |
CALL WRITE_0D_R8( diffKrT, INDEX_NONE,'diffKrT =', |
127 |
&' /* Laplacian diffusion of heat vertically ( m^2/s ) */') |
128 |
CALL WRITE_0D_R8( diffKhS, INDEX_NONE,'diffKhS =', |
129 |
&' /* Laplacian diffusion of salt laterally ( m^2/s ) */') |
130 |
CALL WRITE_0D_R8( diffK4S, INDEX_NONE,'diffK4S =', |
131 |
&' /* Bihaarmonic diffusion of salt laterally ( m^4/s ) */') |
132 |
c CALL WRITE_0D_R8( diffKzS, INDEX_NONE,'diffKzS =', |
133 |
c &' /* Laplacian diffusion of salt vertically ( m^2/s ) */') |
134 |
CALL WRITE_0D_R8( diffKrS, INDEX_NONE,'diffKrS =', |
135 |
&' /* Laplacian diffusion of salt vertically ( m^2/s ) */') |
136 |
CALL WRITE_0D_R8( diffKrBL79surf, INDEX_NONE,'diffKrBL79surf =', |
137 |
&' /* Surface diffusion for Bryan and Lewis 1979 ( m^2/s ) */') |
138 |
CALL WRITE_0D_R8( diffKrBL79deep, INDEX_NONE,'diffKrBL79deep =', |
139 |
&' /* Deep diffusion for Bryan and Lewis 1979 ( m^2/s ) */') |
140 |
CALL WRITE_0D_R8( diffKrBL79scl, INDEX_NONE,'diffKrBL79scl =', |
141 |
&' /* Depth scale for Bryan and Lewis 1979 ( m ) */') |
142 |
CALL WRITE_0D_R8( diffKrBL79Ho, INDEX_NONE,'diffKrBL79Ho =', |
143 |
&' /* Turning depth for Bryan and Lewis 1979 ( m ) */') |
144 |
WRITE(msgBuf,'(2A)') ' Equation of State : eosType = ', eosType |
145 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
146 |
& SQUEEZE_RIGHT , 1) |
147 |
CALL WRITE_0D_R8( tAlpha, INDEX_NONE,'tAlpha =', |
148 |
&' /* Linear EOS thermal expansion coefficient ( 1/degree ) */') |
149 |
CALL WRITE_0D_R8( sBeta, INDEX_NONE,'sBeta =', |
150 |
&' /* Linear EOS haline contraction coefficient ( 1/ppt ) */') |
151 |
IF ( eosType .EQ. 'POLY3' ) THEN |
152 |
WRITE(msgBuf,'(A)') |
153 |
& '// Polynomial EQS parameters ( from POLY3.COEFFS ) ' |
154 |
DO K = 1, Nr |
155 |
WRITE(msgBuf,'(I3,13F8.3)') |
156 |
& K,eosRefT(K),eosRefS(K),eosSig0(K), (eosC(I,K),I=1,9) |
157 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
158 |
& SQUEEZE_RIGHT , 1) |
159 |
ENDDO |
160 |
ENDIF |
161 |
IF ( buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN |
162 |
CALL WRITE_0D_R8( atm_Rd, INDEX_NONE, 'atm_Rd =', |
163 |
& ' /* gas constant for dry air ( J/kg/K ) */') |
164 |
CALL WRITE_0D_R8( atm_Cp, INDEX_NONE, 'atm_Cp =', |
165 |
& ' /* specific heat (Cp) of dry air ( J/kg/K ) */') |
166 |
CALL WRITE_0D_R8( atm_kappa, INDEX_NONE, 'atm_kappa =', |
167 |
& ' /* kappa (=Rd/Cp ) of dry air */') |
168 |
CALL WRITE_0D_R8( atm_Rq, INDEX_NONE, 'atm_Rq =', |
169 |
& ' /* water vap. specific vol. anomaly relative to dry air */') |
170 |
CALL WRITE_0D_R8( atm_Po, INDEX_NONE, 'atm_Po =', |
171 |
& ' /* standard reference pressure ( Pa ) */') |
172 |
CALL WRITE_0D_I( integr_GeoPot, INDEX_NONE, 'integr_GeoPot =', |
173 |
& ' /* select how the geopotential is integrated */') |
174 |
CALL WRITE_0D_I( selectFindRoSurf, INDEX_NONE, |
175 |
& 'selectFindRoSurf=', |
176 |
& ' /* select how Surf.Ref. pressure is defined */') |
177 |
ENDIF |
178 |
CALL WRITE_0D_R8( rhonil, INDEX_NONE,'rhonil =', |
179 |
&' /* Reference density ( kg/m^3 ) */') |
180 |
CALL WRITE_0D_R8( rhoConst, INDEX_NONE,'rhoConst =', |
181 |
&' /* Reference density ( kg/m^3 ) */') |
182 |
CALL WRITE_0D_R8( rhoConstFresh, INDEX_NONE,'rhoConstFresh =', |
183 |
&' /* Reference density ( kg/m^3 ) */') |
184 |
CALL WRITE_0D_R8( gravity, INDEX_NONE,'gravity =', |
185 |
&' /* Gravitational acceleration ( m/s^2 ) */') |
186 |
CALL WRITE_0D_R8( gBaro, INDEX_NONE,'gBaro =', |
187 |
&' /* Barotropic gravity ( m/s^2 ) */') |
188 |
CALL WRITE_0D_R8(rotationPeriod,INDEX_NONE,'rotationPeriod =', |
189 |
&' /* Rotation Period ( s ) */') |
190 |
CALL WRITE_0D_R8( omega, INDEX_NONE,'omega =', |
191 |
&' /* Angular velocity ( rad/s ) */') |
192 |
CALL WRITE_0D_R8( f0, INDEX_NONE,'f0 =', |
193 |
&' /* Reference coriolis parameter ( 1/s ) */') |
194 |
CALL WRITE_0D_R8( beta, INDEX_NONE,'beta =', |
195 |
&' /* Beta ( 1/(m.s) ) */') |
196 |
|
197 |
CALL WRITE_0D_R8( freeSurfFac, INDEX_NONE,'freeSurfFac =', |
198 |
&' /* Implicit free surface factor */') |
199 |
CALL WRITE_0D_L( implicitFreeSurface, INDEX_NONE, |
200 |
& 'implicitFreeSurface =', |
201 |
&' /* Implicit free surface on/off flag */') |
202 |
CALL WRITE_0D_L( rigidLid, INDEX_NONE, |
203 |
& 'rigidLid =', |
204 |
&' /* Rigid lid on/off flag */') |
205 |
CALL WRITE_0D_R8( implicSurfPress, INDEX_NONE, |
206 |
&'implicSurfPress =', |
207 |
&' /* Surface Pressure implicit factor (0-1)*/') |
208 |
CALL WRITE_0D_R8( implicDiv2Dflow, INDEX_NONE, |
209 |
&'implicDiv2Dflow =', |
210 |
&' /* Barot. Flow Div. implicit factor (0-1)*/') |
211 |
CALL WRITE_0D_L( exactConserv, INDEX_NONE, |
212 |
&'exactConserv =', |
213 |
&' /* Exact Volume Conservation on/off flag*/') |
214 |
CALL WRITE_0D_L( uniformLin_PhiSurf, INDEX_NONE, |
215 |
&'uniformLin_PhiSurf =', |
216 |
&' /* use uniform Bo_surf on/off flag*/') |
217 |
CALL WRITE_0D_I( nonlinFreeSurf, INDEX_NONE, |
218 |
&'nonlinFreeSurf =', |
219 |
&' /* Non-linear Free Surf. options (-1,0,1,2,3)*/') |
220 |
WRITE(msgBuf,'(2A)') ' -1,0= Off ; 1,2,3= On,', |
221 |
& ' 2=+rescale gU,gV, 3=+update cg2d solv.' |
222 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
223 |
& SQUEEZE_RIGHT , 1) |
224 |
CALL WRITE_0D_R8( hFacInf, INDEX_NONE, |
225 |
&'hFacInf =', |
226 |
&' /* lower threshold for hFac (nonlinFreeSurf only)*/') |
227 |
CALL WRITE_0D_R8( hFacSup, INDEX_NONE, |
228 |
&'hFacSup =', |
229 |
&' /* upper threshold for hFac (nonlinFreeSurf only)*/') |
230 |
CALL WRITE_0D_I( select_rStar, INDEX_NONE, |
231 |
&'select_rStar =', |
232 |
&' /* r* Coordinate options (not yet implemented)*/') |
233 |
CALL WRITE_0D_L( useRealFreshWaterFlux, INDEX_NONE, |
234 |
&'useRealFreshWaterFlux =', |
235 |
&' /* Real Fresh Water Flux on/off flag*/') |
236 |
IF (useRealFreshWaterFlux .AND. nonlinFreeSurf.GT.0) THEN |
237 |
CALL WRITE_0D_R8( temp_EvPrRn, INDEX_NONE, |
238 |
&'temp_EvPrRn =', |
239 |
&' /* Temp. of Evap/Prec/R (UNSET=use local T)(oC)*/') |
240 |
CALL WRITE_0D_R8( salt_EvPrRn, INDEX_NONE, |
241 |
&'salt_EvPrRn =', |
242 |
&' /* Salin. of Evap/Prec/R (UNSET=use local S)(ppt)*/') |
243 |
CALL WRITE_0D_R8( trac_EvPrRn, INDEX_NONE, |
244 |
&'trac_EvPrRn =', |
245 |
&' /* Tracer in Evap/Prec/R (UNSET=use local Tr)*/') |
246 |
ELSE |
247 |
CALL WRITE_0D_R8( convertFW2Salt, INDEX_NONE, |
248 |
&'convertFW2Salt =', |
249 |
&' /* convert F.W. Flux to Salt Flux (-1=use local S)(ppt)*/') |
250 |
ENDIF |
251 |
|
252 |
CALL WRITE_0D_L( nonHydrostatic, INDEX_NONE, |
253 |
& 'nonHydrostatic =', ' /* Non-Hydrostatic on/off flag */') |
254 |
CALL WRITE_0D_L( momStepping, INDEX_NONE, |
255 |
& 'momStepping =', ' /* Momentum equation on/off flag */') |
256 |
CALL WRITE_0D_L( momAdvection, INDEX_NONE, |
257 |
& 'momAdvection =', ' /* Momentum advection on/off flag */') |
258 |
CALL WRITE_0D_L( momViscosity, INDEX_NONE, |
259 |
& 'momViscosity =', ' /* Momentum viscosity on/off flag */') |
260 |
CALL WRITE_0D_L( momImplVertAdv, INDEX_NONE, 'momImplVertAdv =', |
261 |
& '/* Momentum implicit vert. advection on/off*/') |
262 |
CALL WRITE_0D_L( implicitViscosity, INDEX_NONE, |
263 |
& 'implicitViscosity =', ' /* Implicit viscosity on/off flag */') |
264 |
CALL WRITE_0D_L( useCoriolis, INDEX_NONE, |
265 |
& 'useCoriolis =', ' /* Coriolis on/off flag */') |
266 |
CALL WRITE_0D_L( useCDscheme, INDEX_NONE, |
267 |
& 'useCDscheme =', ' /* CD scheme on/off flag */') |
268 |
CALL WRITE_0D_L( useJamartWetPoints, INDEX_NONE, |
269 |
& 'useJamartWetPoints=',' /* Coriolis WetPoints method flag */') |
270 |
CALL WRITE_0D_L( useJamartMomAdv, INDEX_NONE, |
271 |
& 'useJamartMomAdv=',' /* V.I. Non-linear terms Jamart flag */') |
272 |
CALL WRITE_0D_L( SadournyCoriolis, INDEX_NONE, |
273 |
& 'SadournyCoriolis=',' /* Sadourny Coriolis discr. flag */') |
274 |
CALL WRITE_0D_L( upwindVorticity, INDEX_NONE, |
275 |
& 'upwindVorticity=',' /* Upwind bias vorticity flag */') |
276 |
CALL WRITE_0D_L( useAbsVorticity, INDEX_NONE, |
277 |
& 'useAbsVorticity=',' /* Work with f+zeta in Coriolis */') |
278 |
CALL WRITE_0D_L( highOrderVorticity, INDEX_NONE, |
279 |
& 'highOrderVorticity=',' /* High order interp. of vort. flag */') |
280 |
CALL WRITE_0D_L( momForcing, INDEX_NONE, |
281 |
& 'momForcing =', ' /* Momentum forcing on/off flag */') |
282 |
CALL WRITE_0D_L( momPressureForcing, INDEX_NONE, |
283 |
& 'momPressureForcing =', |
284 |
& ' /* Momentum pressure term on/off flag */') |
285 |
CALL WRITE_0D_L( staggerTimeStep, INDEX_NONE, |
286 |
& 'staggerTimeStep =', |
287 |
&' /* Stagger time stepping on/off flag */') |
288 |
CALL WRITE_0D_L( multiDimAdvection, INDEX_NONE, |
289 |
& 'multiDimAdvection =', |
290 |
&' /* enable/disable Multi-Dim Advection */') |
291 |
CALL WRITE_0D_L( implicitDiffusion, INDEX_NONE, |
292 |
& 'implicitDiffusion =','/* Implicit Diffusion on/off flag */') |
293 |
CALL WRITE_0D_L( tempStepping, INDEX_NONE, |
294 |
& 'tempStepping =', ' /* Temperature equation on/off flag */') |
295 |
CALL WRITE_0D_L( tempAdvection, INDEX_NONE, |
296 |
& 'tempAdvection=', ' /* Temperature advection on/off flag */') |
297 |
CALL WRITE_0D_L( tempImplVertAdv,INDEX_NONE,'tempImplVertAdv =', |
298 |
& '/* Temp. implicit vert. advection on/off */') |
299 |
CALL WRITE_0D_L( tempForcing, INDEX_NONE, |
300 |
& 'tempForcing =', ' /* Temperature forcing on/off flag */') |
301 |
CALL WRITE_0D_L( saltStepping, INDEX_NONE, |
302 |
& 'saltStepping =', ' /* Salinity equation on/off flag */') |
303 |
CALL WRITE_0D_L( saltAdvection, INDEX_NONE, |
304 |
& 'saltAdvection=', ' /* Salinity advection on/off flag */') |
305 |
CALL WRITE_0D_L( saltImplVertAdv,INDEX_NONE,'saltImplVertAdv =', |
306 |
& '/* Sali. implicit vert. advection on/off */') |
307 |
CALL WRITE_0D_L( saltForcing, INDEX_NONE, |
308 |
& 'saltForcing =', ' /* Salinity forcing on/off flag */') |
309 |
WRITE(msgBuf,'(A)') '// ' |
310 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
311 |
& SQUEEZE_RIGHT , 1) |
312 |
|
313 |
WRITE(msgBuf,'(A)') |
314 |
& '// Elliptic solver(s) paramters ( PARM02 in namelist ) ' |
315 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
316 |
& SQUEEZE_RIGHT , 1) |
317 |
WRITE(msgBuf,'(A)') '// ' |
318 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
319 |
& SQUEEZE_RIGHT , 1) |
320 |
CALL WRITE_0D_I( cg2dMaxIters, INDEX_NONE,'cg2dMaxIters =', |
321 |
&' /* Upper limit on 2d con. grad iterations */') |
322 |
CALL WRITE_0D_I( cg2dChkResFreq, INDEX_NONE,'cg2dChkResFreq =', |
323 |
&' /* 2d con. grad convergence test frequency */') |
324 |
CALL WRITE_0D_R8( cg2dTargetResidual, INDEX_NONE, |
325 |
& 'cg2dTargetResidual =', |
326 |
&' /* 2d con. grad target residual */') |
327 |
|
328 |
WRITE(msgBuf,'(A)') '// ' |
329 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
330 |
& SQUEEZE_RIGHT , 1) |
331 |
WRITE(msgBuf,'(A)') |
332 |
& '// Time stepping paramters ( PARM03 in namelist ) ' |
333 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
334 |
& SQUEEZE_RIGHT , 1) |
335 |
WRITE(msgBuf,'(A)') '// ' |
336 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
337 |
& SQUEEZE_RIGHT , 1) |
338 |
CALL WRITE_0D_I( nIter0, INDEX_NONE,'nIter0 =', |
339 |
&' /* Base timestep number */') |
340 |
CALL WRITE_0D_I( nTimeSteps, INDEX_NONE,'nTimeSteps =', |
341 |
&' /* Number of timesteps */') |
342 |
CALL WRITE_0D_R8( deltaTmom, INDEX_NONE,'deltatTmom =', |
343 |
&' /* Momentum equation timestep ( s ) */') |
344 |
CALL WRITE_0D_R8( deltaTfreesurf,INDEX_NONE,'deltaTfreesurf =', |
345 |
&' /* FreeSurface equation timestep ( s ) */') |
346 |
CALL WRITE_0D_R8( deltaTtracer, INDEX_NONE,'deltatTtracer =', |
347 |
&' /* Tracer equation timestep ( s ) */') |
348 |
CALL WRITE_0D_R8( deltaTClock, INDEX_NONE,'deltatTClock =', |
349 |
&' /* Model clock timestep ( s ) */') |
350 |
CALL WRITE_0D_R8( cAdjFreq, INDEX_NONE,'cAdjFreq =', |
351 |
&' /* Convective adjustment interval ( s ) */') |
352 |
CALL WRITE_0D_L( forcing_In_AB,INDEX_NONE,'forcing_In_AB =', |
353 |
&' /* put T,S Forcing in Adams-Bash. stepping */') |
354 |
CALL WRITE_0D_R8( abeps, INDEX_NONE,'abeps =', |
355 |
&' /* Adams-Bashforth stabilizing weight */') |
356 |
IF (useCDscheme) THEN |
357 |
CALL WRITE_0D_R8( tauCD, INDEX_NONE,'tauCD =', |
358 |
&' /* CD coupling time-scale ( s ) */') |
359 |
CALL WRITE_0D_R8( rCD, INDEX_NONE,'rCD =', |
360 |
&' /* Normalised CD coupling parameter */') |
361 |
ENDIF |
362 |
CALL WRITE_0D_R8( startTime, INDEX_NONE,'startTime =', |
363 |
&' /* Run start time ( s ). */') |
364 |
CALL WRITE_0D_R8( endTime, INDEX_NONE,'endTime =', |
365 |
&' /* Integration ending time ( s ). */') |
366 |
CALL WRITE_0D_R8( pChkPtFreq, INDEX_NONE,'pChkPtFreq =', |
367 |
&' /* Permanent restart/checkpoint file interval ( s ). */') |
368 |
CALL WRITE_0D_R8( chkPtFreq, INDEX_NONE,'chkPtFreq =', |
369 |
&' /* Rolling restart/checkpoint file interval ( s ). */') |
370 |
CALL WRITE_0D_R8( dumpFreq, INDEX_NONE,'dumpFreq =', |
371 |
&' /* Model state write out interval ( s ). */') |
372 |
CALL WRITE_0D_R8( externForcingPeriod, INDEX_NONE, |
373 |
& 'externForcingPeriod =', ' /* forcing period (s) */') |
374 |
CALL WRITE_0D_R8( externForcingCycle, INDEX_NONE, |
375 |
& 'externForcingCycle =', ' /* period of the cyle (s). */') |
376 |
CALL WRITE_0D_R8( tauThetaClimRelax, INDEX_NONE, |
377 |
& 'tauThetaClimRelax =', ' /* relaxation time scale (s) */') |
378 |
CALL WRITE_0D_R8( tauSaltClimRelax, INDEX_NONE, |
379 |
& 'tauSaltClimRelax =', ' /* relaxation time scale (s) */') |
380 |
CALL WRITE_0D_R8( latBandClimRelax, INDEX_NONE, |
381 |
& 'latBandClimRelax =', ' /* max. Lat. where relaxation */') |
382 |
WRITE(msgBuf,'(A)') '// ' |
383 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
384 |
& SQUEEZE_RIGHT , 1) |
385 |
WRITE(msgBuf,'(A)') |
386 |
& '// Gridding paramters ( PARM04 in namelist ) ' |
387 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
388 |
& SQUEEZE_RIGHT , 1) |
389 |
WRITE(msgBuf,'(A)') '// ' |
390 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
391 |
& SQUEEZE_RIGHT , 1) |
392 |
CALL WRITE_0D_L( usingCartesianGrid, INDEX_NONE, |
393 |
& 'usingCartesianGrid =', |
394 |
&' /* Cartesian coordinates flag ( True / False ) */') |
395 |
CALL WRITE_0D_L( usingSphericalPolarGrid, INDEX_NONE, |
396 |
& 'usingSphericalPolarGrid =', |
397 |
&' /* Spherical coordinates flag ( True / False ) */') |
398 |
CALL WRITE_0D_L( groundAtK1, INDEX_NONE, 'groundAtK1 =', |
399 |
&' /* Lower Boundary (ground) at the surface(k=1) ( T / F ) */') |
400 |
CALL WRITE_0D_R8( Ro_SeaLevel, INDEX_NONE,'Ro_SeaLevel =', |
401 |
&' /* r(1) ( units of r ) */') |
402 |
CALL WRITE_0D_R8( rkFac, INDEX_NONE,'rkFac =', |
403 |
&' /* minus Vertical index orientation */') |
404 |
CALL WRITE_0D_R8( horiVertRatio, INDEX_NONE,'horiVertRatio =', |
405 |
&' /* Ratio on units : Horiz - Vertical */') |
406 |
c CALL WRITE_1D_R8( delZ,Nr, INDEX_K,'delZ = ', |
407 |
c &' /* W spacing ( m ) */') |
408 |
c CALL WRITE_1D_R8( delP,Nr, INDEX_K,'delP = ', |
409 |
c &' /* W spacing ( Pa ) */') |
410 |
c CALL WRITE_1D_R8( delR,Nr, INDEX_K,'delR = ', |
411 |
c &' /* W spacing ( units of r ) */') |
412 |
CALL WRITE_1D_R8( drC,Nr, INDEX_K,'drC = ', |
413 |
&' /* C spacing ( units of r ) */') |
414 |
CALL WRITE_1D_R8( drF,Nr, INDEX_K,'drF = ', |
415 |
&' /* W spacing ( units of r ) */') |
416 |
CALL WRITE_1D_R8( delX, Nx, INDEX_I,'delX = ', |
417 |
&' /* U spacing ( m - cartesian, degrees - spherical ) */') |
418 |
CALL WRITE_1D_R8( delY, Ny, INDEX_J,'delY = ', |
419 |
&' /* V spacing ( m - cartesian, degrees - spherical ) */') |
420 |
CALL WRITE_0D_R8( phiMin, INDEX_NONE,'phiMin = ', |
421 |
&' /* South edge (ignored - cartesian, degrees - spherical ) */') |
422 |
CALL WRITE_0D_R8( thetaMin, INDEX_NONE,'thetaMin = ', |
423 |
&' /* West edge ( ignored - cartesian, degrees - spherical ) */') |
424 |
CALL WRITE_0D_R8( rSphere, INDEX_NONE,'rSphere = ', |
425 |
&' /* Radius ( ignored - cartesian, m - spherical ) */') |
426 |
DO bi=1,nSx |
427 |
DO I=1,sNx |
428 |
xcoord((bi-1)*sNx+I) = xC(I,1,bi,1) |
429 |
ENDDO |
430 |
ENDDO |
431 |
CALL WRITE_1D_R8( xcoord, sNx*nSx, INDEX_I,'xcoord = ', |
432 |
&' /* P-point X coord ( m - cartesian, degrees - spherical ) */') |
433 |
DO bj=1,nSy |
434 |
DO J=1,sNy |
435 |
ycoord((bj-1)*sNy+J) = yC(1,J,1,bj) |
436 |
ENDDO |
437 |
ENDDO |
438 |
CALL WRITE_1D_R8( ycoord, sNy*nSy, INDEX_J,'ycoord = ', |
439 |
&' /* P-point Y coord ( m - cartesian, degrees - spherical ) */') |
440 |
DO K=1,Nr |
441 |
rcoord(K) = rC(K) |
442 |
ENDDO |
443 |
CALL WRITE_1D_R8( rcoord, Nr, INDEX_K,'rcoord = ', |
444 |
&' /* P-point R coordinate ( units of r ) */') |
445 |
DO K=1,Nr+1 |
446 |
rcoord(K) = rF(K) |
447 |
ENDDO |
448 |
CALL WRITE_1D_R8( rcoord, Nr+1, INDEX_K,'rF = ', |
449 |
&' /* W-Interf. R coordinate ( units of r ) */') |
450 |
|
451 |
C Grid along selected grid lines |
452 |
coordLine = 1 |
453 |
tileLine = 1 |
454 |
CALL WRITE_XY_XLINE_RS( dxF, coordLine, tileLine, |
455 |
I 'dxF','( m - cartesian, degrees - spherical )') |
456 |
CALL WRITE_XY_YLINE_RS( dxF, coordLine, tileLine, |
457 |
I 'dxF','( m - cartesian, degrees - spherical )') |
458 |
CALL WRITE_XY_XLINE_RS( dyF, coordLine, tileLine, |
459 |
I 'dyF','( m - cartesian, degrees - spherical )') |
460 |
CALL WRITE_XY_YLINE_RS( dyF, coordLine, tileLine, |
461 |
I 'dyF','( m - cartesian, degrees - spherical )') |
462 |
CALL WRITE_XY_XLINE_RS( dxG, coordLine, tileLine, |
463 |
I 'dxG','( m - cartesian, degrees - spherical )') |
464 |
CALL WRITE_XY_YLINE_RS( dxG, coordLine, tileLine, |
465 |
I 'dxG','( m - cartesian, degrees - spherical )') |
466 |
CALL WRITE_XY_XLINE_RS( dyG, coordLine, tileLine, |
467 |
I 'dyG','( m - cartesian, degrees - spherical )') |
468 |
CALL WRITE_XY_YLINE_RS( dyG, coordLine, tileLine, |
469 |
I 'dyG','( m - cartesian, degrees - spherical )') |
470 |
CALL WRITE_XY_XLINE_RS( dxC, coordLine, tileLine, |
471 |
I 'dxC','( m - cartesian, degrees - spherical )') |
472 |
CALL WRITE_XY_YLINE_RS( dxC, coordLine, tileLine, |
473 |
I 'dxC','( m - cartesian, degrees - spherical )') |
474 |
CALL WRITE_XY_XLINE_RS( dyC, coordLine, tileLine, |
475 |
I 'dyC','( m - cartesian, degrees - spherical )') |
476 |
CALL WRITE_XY_YLINE_RS( dyC, coordLine, tileLine, |
477 |
I 'dyC','( m - cartesian, degrees - spherical )') |
478 |
CALL WRITE_XY_XLINE_RS( dxV, coordLine, tileLine, |
479 |
I 'dxV','( m - cartesian, degrees - spherical )') |
480 |
CALL WRITE_XY_YLINE_RS( dxV, coordLine, tileLine, |
481 |
I 'dxV','( m - cartesian, degrees - spherical )') |
482 |
CALL WRITE_XY_XLINE_RS( dyU, coordLine, tileLine, |
483 |
I 'dyU','( m - cartesian, degrees - spherical )') |
484 |
CALL WRITE_XY_YLINE_RS( dyU, coordLine, tileLine, |
485 |
I 'dyU','( m - cartesian, degrees - spherical )') |
486 |
CALL WRITE_XY_XLINE_RS( rA, coordLine, tileLine, |
487 |
I 'rA','( m - cartesian, degrees - spherical )') |
488 |
CALL WRITE_XY_YLINE_RS( rA, coordLine, tileLine, |
489 |
I 'rA','( m - cartesian, degrees - spherical )') |
490 |
CALL WRITE_XY_XLINE_RS( rAw, coordLine, tileLine, |
491 |
I 'rAw','( m - cartesian, degrees - spherical )') |
492 |
CALL WRITE_XY_YLINE_RS( rAw, coordLine, tileLine, |
493 |
I 'rAw','( m - cartesian, degrees - spherical )') |
494 |
CALL WRITE_XY_XLINE_RS( rAs, coordLine, tileLine, |
495 |
I 'rAs','( m - cartesian, degrees - spherical )') |
496 |
CALL WRITE_XY_YLINE_RS( rAs, coordLine, tileLine, |
497 |
I 'rAs','( m - cartesian, degrees - spherical )') |
498 |
|
499 |
WRITE(msgBuf,'(A)') ' ' |
500 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
501 |
& SQUEEZE_RIGHT , 1) |
502 |
|
503 |
_END_MASTER(myThid) |
504 |
_BARRIER |
505 |
|
506 |
|
507 |
RETURN |
508 |
100 FORMAT(A, |
509 |
&' ' |
510 |
&) |
511 |
END |
512 |
|