1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/config_summary.F,v 1.35 2002/09/25 19:36:50 mlosch 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 |
CALL WRITE_1D_R8( tRef, Nr, INDEX_K,'tRef =', |
86 |
&' /* Reference temperature profile ( oC or oK ) */') |
87 |
CALL WRITE_1D_R8( sRef, Nr, INDEX_K,'sRef =', |
88 |
&' /* Reference salinity profile ( ppt ) */') |
89 |
CALL WRITE_0D_R8( viscAh, INDEX_NONE,'viscAh =', |
90 |
&' /* Lateral eddy viscosity ( m^2/s ) */') |
91 |
CALL WRITE_0D_R8( viscA4, INDEX_NONE,'viscAh =', |
92 |
&' /* Lateral biharmonic viscosity ( m^4/s ) */') |
93 |
CALL WRITE_0D_L( no_slip_sides, INDEX_NONE, |
94 |
& 'no_slip_sides =', ' /* Viscous BCs: No-slip sides */') |
95 |
IF ( viscAz .NE. UNSET_RL ) THEN |
96 |
CALL WRITE_0D_R8( viscAz, INDEX_NONE,'viscAz =', |
97 |
& ' /* Vertical eddy viscosity ( m^2/s ) */') |
98 |
ENDIF |
99 |
IF ( viscAp .NE. UNSET_RL ) THEN |
100 |
CALL WRITE_0D_R8( viscAp, INDEX_NONE,'viscAp =', |
101 |
& ' /* Vertical eddy viscosity ( Pa^2/s ) */') |
102 |
ENDIF |
103 |
CALL WRITE_0D_R8( viscAr, INDEX_NONE,'viscAr =', |
104 |
&' /* Vertical eddy viscosity ( units of r^2/s ) */') |
105 |
CALL WRITE_0D_R8( diffKhT, INDEX_NONE,'diffKhT =', |
106 |
&' /* Laplacian diffusion of heat laterally ( m^2/s ) */') |
107 |
CALL WRITE_0D_R8( diffK4T, INDEX_NONE,'diffK4T =', |
108 |
&' /* Bihaarmonic diffusion of heat laterally ( m^4/s ) */') |
109 |
CALL WRITE_0D_R8( diffKzT, INDEX_NONE,'diffKzT =', |
110 |
&' /* Laplacian diffusion of heat vertically ( m^2/s ) */') |
111 |
CALL WRITE_0D_R8( diffKrT, INDEX_NONE,'diffKrT =', |
112 |
&' /* Laplacian diffusion of heat vertically ( m^2/s ) */') |
113 |
CALL WRITE_0D_R8( diffKhS, INDEX_NONE,'diffKhS =', |
114 |
&' /* Laplacian diffusion of salt laterally ( m^2/s ) */') |
115 |
CALL WRITE_0D_R8( diffK4S, INDEX_NONE,'diffK4S =', |
116 |
&' /* Bihaarmonic diffusion of salt laterally ( m^4/s ) */') |
117 |
CALL WRITE_0D_R8( diffKzS, INDEX_NONE,'diffKzS =', |
118 |
&' /* Laplacian diffusion of salt vertically ( m^2/s ) */') |
119 |
CALL WRITE_0D_R8( diffKrS, INDEX_NONE,'diffKrS =', |
120 |
&' /* Laplacian diffusion of salt vertically ( m^2/s ) */') |
121 |
CALL WRITE_0D_R8( tAlpha, INDEX_NONE,'tAlpha =', |
122 |
&' /* Linear EOS thermal expansion coefficient ( 1/degree ) */') |
123 |
CALL WRITE_0D_R8( sBeta, INDEX_NONE,'sBeta =', |
124 |
&' /* Linear EOS haline contraction coefficient ( 1/ppt ) */') |
125 |
IF ( eosType .EQ. 'POLY3' ) THEN |
126 |
WRITE(msgBuf,'(A)') |
127 |
& '// Polynomial EQS parameters ( from POLY3.COEFFS ) ' |
128 |
DO K = 1, Nr |
129 |
WRITE(msgBuf,'(I3,13F8.3)') |
130 |
& K,eosRefT(K),eosRefS(K),eosSig0(K), (eosC(I,K),I=1,9) |
131 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
132 |
& SQUEEZE_RIGHT , 1) |
133 |
ENDDO |
134 |
ENDIF |
135 |
CALL WRITE_0D_R8( rhonil, INDEX_NONE,'rhonil =', |
136 |
&' /* Reference density ( kg/m^3 ) */') |
137 |
CALL WRITE_0D_R8( rhoConst, INDEX_NONE,'rhoConst =', |
138 |
&' /* Reference density ( kg/m^3 ) */') |
139 |
CALL WRITE_0D_R8( rhoConstFresh, INDEX_NONE,'rhoConstFresh =', |
140 |
&' /* Reference density ( kg/m^3 ) */') |
141 |
CALL WRITE_0D_R8( gravity, INDEX_NONE,'gravity =', |
142 |
&' /* Gravitational acceleration ( m/s^2 ) */') |
143 |
CALL WRITE_0D_R8( gBaro, INDEX_NONE,'gBaro =', |
144 |
&' /* Barotropic gravity ( m/s^2 ) */') |
145 |
CALL WRITE_0D_R8( f0, INDEX_NONE,'f0 =', |
146 |
&' /* Reference coriolis parameter ( 1/s ) */') |
147 |
CALL WRITE_0D_R8( beta, INDEX_NONE,'beta =', |
148 |
&' /* Beta ( 1/(m.s) ) */') |
149 |
|
150 |
CALL WRITE_0D_R8( freeSurfFac, INDEX_NONE,'freeSurfFac =', |
151 |
&' /* Implicit free surface factor */') |
152 |
CALL WRITE_0D_L( implicitFreeSurface, INDEX_NONE, |
153 |
& 'implicitFreeSurface =', |
154 |
&' /* Implicit free surface on/off flag */') |
155 |
CALL WRITE_0D_L( rigidLid, INDEX_NONE, |
156 |
& 'rigidLid =', |
157 |
&' /* Rigid lid on/off flag */') |
158 |
CALL WRITE_0D_R8( implicSurfPress, INDEX_NONE, |
159 |
&'implicSurfPress =', |
160 |
&' /* Surface Pressure implicit factor (0-1)*/') |
161 |
CALL WRITE_0D_R8( implicDiv2Dflow, INDEX_NONE, |
162 |
&'implicDiv2Dflow =', |
163 |
&' /* Barot. Flow Div. implicit factor (0-1)*/') |
164 |
CALL WRITE_0D_L( exactConserv, INDEX_NONE, |
165 |
&'exactConserv =', |
166 |
&' /* Exact Volume Conservation on/off flag*/') |
167 |
CALL WRITE_0D_L( uniformLin_PhiSurf, INDEX_NONE, |
168 |
&'uniformLin_PhiSurf =', |
169 |
&' /* use uniform Bo_surf on/off flag*/') |
170 |
CALL WRITE_0D_I( nonlinFreeSurf, INDEX_NONE, |
171 |
&'nonlinFreeSurf =', |
172 |
&' /* Non-linear Free Surf. options (-1,0,1,2,3)*/') |
173 |
WRITE(msgBuf,'(2A)') ' -1,0= Off ; 1,2,3= On,', |
174 |
& ' 2=+rescale gU,gV, 3=+update cg2d solv.' |
175 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
176 |
& SQUEEZE_RIGHT , 1) |
177 |
CALL WRITE_0D_R8( hFacInf, INDEX_NONE, |
178 |
&'hFacInf =', |
179 |
&' /* lower threshold for hFac (nonlinFreeSurf only)*/') |
180 |
CALL WRITE_0D_R8( hFacSup, INDEX_NONE, |
181 |
&'hFacSup =', |
182 |
&' /* upper threshold for hFac (nonlinFreeSurf only)*/') |
183 |
CALL WRITE_0D_L( useRealFreshWaterFlux, INDEX_NONE, |
184 |
&'useRealFreshWaterFlux =', |
185 |
&' /* Real Fresh Water Flux on/off flag*/') |
186 |
IF (useRealFreshWaterFlux .AND. nonlinFreeSurf.GT.0) THEN |
187 |
CALL WRITE_0D_R8( temp_EvPrRn, INDEX_NONE, |
188 |
&'temp_EvPrRn =', |
189 |
&' /* Temp. of Evap/Prec/R (UNSET=use local T)(oC)*/') |
190 |
CALL WRITE_0D_R8( salt_EvPrRn, INDEX_NONE, |
191 |
&'salt_EvPrRn =', |
192 |
&' /* Salin. of Evap/Prec/R (UNSET=use local S)(ppt)*/') |
193 |
CALL WRITE_0D_R8( trac_EvPrRn, INDEX_NONE, |
194 |
&'trac_EvPrRn =', |
195 |
&' /* Tracer in Evap/Prec/R (UNSET=use local Tr)*/') |
196 |
ELSE |
197 |
CALL WRITE_0D_R8( convertFW2Salt, INDEX_NONE, |
198 |
&'convertFW2Salt =', |
199 |
&' /* convert F.W. Flux to Salt Flux (-1=use local S)(ppt)*/') |
200 |
ENDIF |
201 |
|
202 |
CALL WRITE_0D_L( multiDimAdvection, INDEX_NONE, |
203 |
& 'multiDimAdvection =', |
204 |
&' /* enable/disable Multi-Dim Advection */') |
205 |
CALL WRITE_0D_L( staggerTimeStep, INDEX_NONE, |
206 |
& 'staggerTimeStep =', |
207 |
&' /* Stagger time stepping on/off flag */') |
208 |
CALL WRITE_0D_L( momStepping, INDEX_NONE, |
209 |
& 'momStepping =', ' /* Momentum equation on/off flag */') |
210 |
CALL WRITE_0D_L( momAdvection, INDEX_NONE, |
211 |
& 'momAdvection =', ' /* Momentum advection on/off flag */') |
212 |
CALL WRITE_0D_L( momViscosity, INDEX_NONE, |
213 |
& 'momViscosity =', ' /* Momentum viscosity on/off flag */') |
214 |
CALL WRITE_0D_L( useCoriolis, INDEX_NONE, |
215 |
& 'useCoriolis =', ' /* Coriolis on/off flag */') |
216 |
CALL WRITE_0D_L( momForcing, INDEX_NONE, |
217 |
& 'momForcing =', ' /* Momentum forcing on/off flag */') |
218 |
CALL WRITE_0D_L( momPressureForcing, INDEX_NONE, |
219 |
& 'momPressureForcing =', |
220 |
& ' /* Momentum pressure term on/off flag */') |
221 |
CALL WRITE_0D_L( tempStepping, INDEX_NONE, |
222 |
& 'tempStepping =', ' /* Temperature equation on/off flag */') |
223 |
CALL WRITE_0D_L( tempAdvection, INDEX_NONE, |
224 |
& 'tempAdvection=', ' /* Temperature advection on/off flag */') |
225 |
CALL WRITE_0D_L( tempForcing, INDEX_NONE, |
226 |
& 'tempForcing =', ' /* Temperature forcing on/off flag */') |
227 |
CALL WRITE_0D_L( saltStepping, INDEX_NONE, |
228 |
& 'saltStepping =', ' /* Salinity equation on/off flag */') |
229 |
CALL WRITE_0D_L( saltAdvection, INDEX_NONE, |
230 |
& 'saltAdvection=', ' /* Salinity advection on/off flag */') |
231 |
CALL WRITE_0D_L( saltForcing, INDEX_NONE, |
232 |
& 'saltForcing =', ' /* Salinity forcing on/off flag */') |
233 |
CALL WRITE_0D_L( nonHydrostatic, INDEX_NONE, |
234 |
& 'nonHydrostatic =', ' /* Non-Hydrostatic on/off flag */') |
235 |
WRITE(msgBuf,'(A)') '// ' |
236 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
237 |
& SQUEEZE_RIGHT , 1) |
238 |
|
239 |
WRITE(msgBuf,'(A)') |
240 |
& '// Elliptic solver(s) paramters ( PARM02 in namelist ) ' |
241 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
242 |
& SQUEEZE_RIGHT , 1) |
243 |
WRITE(msgBuf,'(A)') '// ' |
244 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
245 |
& SQUEEZE_RIGHT , 1) |
246 |
CALL WRITE_0D_I( cg2dMaxIters, INDEX_NONE,'cg2dMaxIters =', |
247 |
&' /* Upper limit on 2d con. grad iterations */') |
248 |
CALL WRITE_0D_I( cg2dChkResFreq, INDEX_NONE,'cg2dChkResFreq =', |
249 |
&' /* 2d con. grad convergence test frequency */') |
250 |
CALL WRITE_0D_R8( cg2dTargetResidual, INDEX_NONE, |
251 |
& 'cg2dTargetResidual =', |
252 |
&' /* 2d con. grad target residual */') |
253 |
|
254 |
WRITE(msgBuf,'(A)') '// ' |
255 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
256 |
& SQUEEZE_RIGHT , 1) |
257 |
WRITE(msgBuf,'(A)') |
258 |
& '// Time stepping paramters ( PARM03 in namelist ) ' |
259 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
260 |
& SQUEEZE_RIGHT , 1) |
261 |
WRITE(msgBuf,'(A)') '// ' |
262 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
263 |
& SQUEEZE_RIGHT , 1) |
264 |
CALL WRITE_0D_I( nIter0, INDEX_NONE,'nIter0 =', |
265 |
&' /* Base timestep number */') |
266 |
CALL WRITE_0D_I( nTimeSteps, INDEX_NONE,'nTimeSteps =', |
267 |
&' /* Number of timesteps */') |
268 |
CALL WRITE_0D_R8( deltaTmom, INDEX_NONE,'deltatTmom =', |
269 |
&' /* Momentum equation timestep ( s ) */') |
270 |
CALL WRITE_0D_R8( deltaTtracer, INDEX_NONE,'deltatTtracer =', |
271 |
&' /* Tracer equation timestep ( s ) */') |
272 |
CALL WRITE_0D_R8( deltaTClock, INDEX_NONE,'deltatTClock =', |
273 |
&' /* Model clock timestep ( s ) */') |
274 |
CALL WRITE_0D_R8( cAdjFreq, INDEX_NONE,'cAdjFreq =', |
275 |
&' /* Convective adjustment interval ( s ) */') |
276 |
CALL WRITE_0D_L( forcing_In_AB,INDEX_NONE,'forcing_In_AB =', |
277 |
&' /* put T,S Forcing in Adams-Bash. stepping */') |
278 |
CALL WRITE_0D_R8( abeps, INDEX_NONE,'abeps =', |
279 |
&' /* Adams-Bashforth stabilizing weight */') |
280 |
CALL WRITE_0D_R8( tauCD, INDEX_NONE,'tauCD =', |
281 |
&' /* CD coupling time-scale ( s ) */') |
282 |
CALL WRITE_0D_R8( rCD, INDEX_NONE,'rCD =', |
283 |
&' /* Normalised CD coupling parameter */') |
284 |
CALL WRITE_0D_R8( startTime, INDEX_NONE,'startTime =', |
285 |
&' /* Run start time ( s ). */') |
286 |
CALL WRITE_0D_R8( endTime, INDEX_NONE,'endTime =', |
287 |
&' /* Integration ending time ( s ). */') |
288 |
CALL WRITE_0D_R8( pChkPtFreq, INDEX_NONE,'pChkPtFreq =', |
289 |
&' /* Permanent restart/checkpoint file interval ( s ). */') |
290 |
CALL WRITE_0D_R8( chkPtFreq, INDEX_NONE,'chkPtFreq =', |
291 |
&' /* Rolling restart/checkpoint file interval ( s ). */') |
292 |
CALL WRITE_0D_R8( dumpFreq, INDEX_NONE,'dumpFreq =', |
293 |
&' /* Model state write out interval ( s ). */') |
294 |
|
295 |
WRITE(msgBuf,'(A)') '// ' |
296 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
297 |
& SQUEEZE_RIGHT , 1) |
298 |
WRITE(msgBuf,'(A)') |
299 |
& '// Gridding paramters ( PARM04 in namelist ) ' |
300 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
301 |
& SQUEEZE_RIGHT , 1) |
302 |
WRITE(msgBuf,'(A)') '// ' |
303 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
304 |
& SQUEEZE_RIGHT , 1) |
305 |
CALL WRITE_0D_L( usingCartesianGrid, INDEX_NONE, |
306 |
& 'usingCartesianGrid =', |
307 |
&' /* Cartesian coordinates flag ( True / False ) */') |
308 |
CALL WRITE_0D_L( usingSphericalPolarGrid, INDEX_NONE, |
309 |
& 'usingSphericalPolarGrid =', |
310 |
&' /* Spherical coordinates flag ( True / False ) */') |
311 |
CALL WRITE_0D_L( groundAtK1, INDEX_NONE, 'groundAtK1 =', |
312 |
&' /* Lower Boundary (ground) at the surface(k=1) ( T / F ) */') |
313 |
CALL WRITE_0D_R8( Ro_SeaLevel, INDEX_NONE,'Ro_SeaLevel =', |
314 |
&' /* r(1) ( units of r ) */') |
315 |
CALL WRITE_0D_R8( rkFac, INDEX_NONE,'rkFac =', |
316 |
&' /* minus Vertical index orientation */') |
317 |
CALL WRITE_0D_R8( horiVertRatio, INDEX_NONE,'horiVertRatio =', |
318 |
&' /* Ratio on units : Horiz - Vertical */') |
319 |
c CALL WRITE_1D_R8( delZ,Nr, INDEX_K,'delZ = ', |
320 |
c &' /* W spacing ( m ) */') |
321 |
c CALL WRITE_1D_R8( delP,Nr, INDEX_K,'delP = ', |
322 |
c &' /* W spacing ( Pa ) */') |
323 |
c CALL WRITE_1D_R8( delR,Nr, INDEX_K,'delR = ', |
324 |
c &' /* W spacing ( units of r ) */') |
325 |
CALL WRITE_1D_R8( drC,Nr, INDEX_K,'drC = ', |
326 |
&' /* C spacing ( units of r ) */') |
327 |
CALL WRITE_1D_R8( drF,Nr, INDEX_K,'drF = ', |
328 |
&' /* W spacing ( units of r ) */') |
329 |
CALL WRITE_1D_R8( delX, Nx, INDEX_I,'delX = ', |
330 |
&' /* U spacing ( m - cartesian, degrees - spherical ) */') |
331 |
CALL WRITE_1D_R8( delY, Ny, INDEX_J,'delY = ', |
332 |
&' /* V spacing ( m - cartesian, degrees - spherical ) */') |
333 |
CALL WRITE_0D_R8( phiMin, INDEX_NONE,'phiMin = ', |
334 |
&' /* South edge (ignored - cartesian, degrees - spherical ) */') |
335 |
CALL WRITE_0D_R8( thetaMin, INDEX_NONE,'thetaMin = ', |
336 |
&' /* West edge ( ignored - cartesian, degrees - spherical ) */') |
337 |
CALL WRITE_0D_R8( rSphere, INDEX_NONE,'rSphere = ', |
338 |
&' /* Radius ( ignored - cartesian, m - spherical ) */') |
339 |
DO bi=1,nSx |
340 |
DO I=1,sNx |
341 |
xcoord((bi-1)*sNx+I) = xC(I,1,bi,1) |
342 |
ENDDO |
343 |
ENDDO |
344 |
CALL WRITE_1D_R8( xcoord, sNx*nSx, INDEX_I,'xcoord = ', |
345 |
&' /* P-point X coord ( m - cartesian, degrees - spherical ) */') |
346 |
DO bj=1,nSy |
347 |
DO J=1,sNy |
348 |
ycoord((bj-1)*sNy+J) = yC(1,J,1,bj) |
349 |
ENDDO |
350 |
ENDDO |
351 |
CALL WRITE_1D_R8( ycoord, sNy*nSy, INDEX_J,'ycoord = ', |
352 |
&' /* P-point Y coord ( m - cartesian, degrees - spherical ) */') |
353 |
DO K=1,Nr |
354 |
rcoord(K) = rC(K) |
355 |
ENDDO |
356 |
CALL WRITE_1D_R8( rcoord, Nr, INDEX_K,'rcoord = ', |
357 |
&' /* P-point R coordinate ( units of r ) */') |
358 |
DO K=1,Nr+1 |
359 |
rcoord(K) = rF(K) |
360 |
ENDDO |
361 |
CALL WRITE_1D_R8( rcoord, Nr+1, INDEX_K,'rF = ', |
362 |
&' /* W-Interf. R coordinate ( units of r ) */') |
363 |
|
364 |
C Grid along selected grid lines |
365 |
coordLine = 1 |
366 |
tileLine = 1 |
367 |
CALL WRITE_XY_XLINE_RS( dxF, coordLine, tileLine, |
368 |
I 'dxF','( m - cartesian, degrees - spherical )') |
369 |
CALL WRITE_XY_YLINE_RS( dxF, coordLine, tileLine, |
370 |
I 'dxF','( m - cartesian, degrees - spherical )') |
371 |
CALL WRITE_XY_XLINE_RS( dyF, coordLine, tileLine, |
372 |
I 'dyF','( m - cartesian, degrees - spherical )') |
373 |
CALL WRITE_XY_YLINE_RS( dyF, coordLine, tileLine, |
374 |
I 'dyF','( m - cartesian, degrees - spherical )') |
375 |
CALL WRITE_XY_XLINE_RS( dxG, coordLine, tileLine, |
376 |
I 'dxG','( m - cartesian, degrees - spherical )') |
377 |
CALL WRITE_XY_YLINE_RS( dxG, coordLine, tileLine, |
378 |
I 'dxG','( m - cartesian, degrees - spherical )') |
379 |
CALL WRITE_XY_XLINE_RS( dyG, coordLine, tileLine, |
380 |
I 'dyG','( m - cartesian, degrees - spherical )') |
381 |
CALL WRITE_XY_YLINE_RS( dyG, coordLine, tileLine, |
382 |
I 'dyG','( m - cartesian, degrees - spherical )') |
383 |
CALL WRITE_XY_XLINE_RS( dxC, coordLine, tileLine, |
384 |
I 'dxC','( m - cartesian, degrees - spherical )') |
385 |
CALL WRITE_XY_YLINE_RS( dxC, coordLine, tileLine, |
386 |
I 'dxC','( m - cartesian, degrees - spherical )') |
387 |
CALL WRITE_XY_XLINE_RS( dyC, coordLine, tileLine, |
388 |
I 'dyC','( m - cartesian, degrees - spherical )') |
389 |
CALL WRITE_XY_YLINE_RS( dyC, coordLine, tileLine, |
390 |
I 'dyC','( m - cartesian, degrees - spherical )') |
391 |
CALL WRITE_XY_XLINE_RS( dxV, coordLine, tileLine, |
392 |
I 'dxV','( m - cartesian, degrees - spherical )') |
393 |
CALL WRITE_XY_YLINE_RS( dxV, coordLine, tileLine, |
394 |
I 'dxV','( m - cartesian, degrees - spherical )') |
395 |
CALL WRITE_XY_XLINE_RS( dyU, coordLine, tileLine, |
396 |
I 'dyU','( m - cartesian, degrees - spherical )') |
397 |
CALL WRITE_XY_YLINE_RS( dyU, coordLine, tileLine, |
398 |
I 'dyU','( m - cartesian, degrees - spherical )') |
399 |
CALL WRITE_XY_XLINE_RS( rA, coordLine, tileLine, |
400 |
I 'rA','( m - cartesian, degrees - spherical )') |
401 |
CALL WRITE_XY_YLINE_RS( rA, coordLine, tileLine, |
402 |
I 'rA','( m - cartesian, degrees - spherical )') |
403 |
CALL WRITE_XY_XLINE_RS( rAw, coordLine, tileLine, |
404 |
I 'rAw','( m - cartesian, degrees - spherical )') |
405 |
CALL WRITE_XY_YLINE_RS( rAw, coordLine, tileLine, |
406 |
I 'rAw','( m - cartesian, degrees - spherical )') |
407 |
CALL WRITE_XY_XLINE_RS( rAs, coordLine, tileLine, |
408 |
I 'rAs','( m - cartesian, degrees - spherical )') |
409 |
CALL WRITE_XY_YLINE_RS( rAs, coordLine, tileLine, |
410 |
I 'rAs','( m - cartesian, degrees - spherical )') |
411 |
|
412 |
WRITE(msgBuf,'(A)') ' ' |
413 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
414 |
& SQUEEZE_RIGHT , 1) |
415 |
|
416 |
_END_MASTER(myThid) |
417 |
_BARRIER |
418 |
|
419 |
|
420 |
RETURN |
421 |
100 FORMAT(A, |
422 |
&' ' |
423 |
&) |
424 |
END |
425 |
|