ROMS
Loading...
Searching...
No Matches
set_scoord.F File Reference
#include "cppdefs.h"
Include dependency graph for set_scoord.F:

Go to the source code of this file.

Functions/Subroutines

subroutine set_scoord (ng)
 

Function/Subroutine Documentation

◆ set_scoord()

subroutine set_scoord ( integer, intent(in) ng)

Definition at line 3 of file set_scoord.F.

4!
5!git $Id$
6!=======================================================================
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md Hernan G. Arango !
10!========================================== Alexander F. Shchepetkin ===
11! !
12! This routine sets and initializes relevant variables associated !
13! with the vertical terrain-following coordinates transformation. !
14! !
15! Definitions: !
16! !
17! N(ng) : Number of vertical levels for each nested grid. !
18! !
19! zeta : time-varying free-surface, zeta(x,y,t), (m) !
20! !
21! h : bathymetry, h(x,y), (m, positive, maybe time-varying) !
22! !
23! hc : critical (thermocline, pycnocline) depth (m, positive) !
24! !
25! z : vertical depths, z(x,y,s,t), meters, negative !
26! z_w(x,y,0:N(ng)) at W-points (top/bottom cell) !
27! z_r(z,y,1:N(ng)) at RHO-points (cell center) !
28! !
29! z_w(x,y,0 ) = -h(x,y) !
30! z_w(x,y,N(ng)) = zeta(x,y,t) !
31! !
32! s : nondimensional stretched vertical coordinate, !
33! -1 <= s <= 0 !
34! !
35! s = 0 at the free-surface, z(x,y, 0,t) = zeta(x,y,t) !
36! s = -1 at the bottom, z(x,y,-1,t) = - h(x,y,t) !
37! !
38! sc_w(k) = (k-N(ng))/N(ng) k=0:N, W-points !
39! sc_r(k) = (k-N(ng)-0.5)/N(ng) k=1:N, RHO-points !
40! !
41! C : nondimensional vertical stretching function, C(s), !
42! -1 <= C(s) <= 0 !
43! !
44! C(s) = 0 for s = 0, at the free-surface !
45! C(s) = -1 for s = -1, at the bottom !
46! !
47! Cs_w(k) = F(s,theta_s,theta_b) k=0:N, W-points !
48! Cs_r(k) = C(s,theta_s,theta_b) k=1:N, RHO-points !
49! !
50! Zo : vertical transformation functional, Zo(x,y,s): !
51! !
52! Zo(x,y,s) = H(x,y)C(s) separable functions !
53! !
54! !
55! Two vertical transformations are supported, z => z(x,y,s,t): !
56! !
57! (1) Original transformation (Shchepetkin and McWilliams, 2005): In !
58! ROMS since 1999 (version 1.8): !
59! !
60! z(x,y,s,t) = Zo(x,y,s) + zeta(x,y,t) * [1 + Zo(x,y,s)/h(x,y)] !
61! !
62! where !
63! !
64! Zo(x,y,s) = hc * s + [h(x,y) - hc] * C(s) !
65! !
66! Zo(x,y,s) = 0 for s = 0, C(s) = 0, at the surface !
67! Zo(x,y,s) = -h(x,y) for s = -1, C(s) = -1, at the bottom !
68! !
69! (2) New transformation: In UCLA-ROMS since 2005: !
70! !
71! z(x,y,s,t) = zeta(x,y,t) + [zeta(x,y,t) + h(x,y)] * Zo(x,y,s) !
72! !
73! where !
74! !
75! Zo(x,y,s) = [hc * s(k) + h(x,y) * C(k)] / [hc + h(x,y)] !
76! !
77! Zo(x,y,s) = 0 for s = 0, C(s) = 0, at the surface !
78! Zo(x,y,s) = -1 for s = -1, C(s) = -1, at the bottom !
79! !
80! At the rest state, corresponding to zero free-surface, this !
81! transformation yields the following unperturbed depths, zhat: !
82! !
83! zhat = z(x,y,s,0) = h(x,y) * Zo(x,y,s) !
84! !
85! = h(x,y) * [hc * s(k) + h(x,y) * C(k)] / [hc + h(x,y)] !
86! !
87! and !
88! !
89! d(zhat) = ds * h(x,y) * hc / [hc + h(x,y)] !
90! !
91! As a consequence, the uppermost grid box retains very little !
92! dependency from bathymetry in the areas where hc << h(x,y), !
93! that is deep areas. For example, if hc=250 m, and h(x,y) !
94! changes from 2000 to 6000 meters, the uppermost grid box !
95! changes only by a factor of 1.08 (less than 10%). !
96! !
97! Notice that: !
98! !
99! * Regardless of the design of C(s), transformation (2) behaves !
100! like equally-spaced sigma-coordinates in shallow areas, where !
101! h(x,y) << hc. This is advantageous because high vertical !
102! resolution and associated CFL limitation is avoided in these !
103! areas. !
104! !
105! * Near-surface refinement is close to geopotential coordinates !
106! in deep areas (level thickness do not depend or weakly-depend !
107! on the bathymetry). Contrarily, near-bottom refinement is !
108! like sigma-coordinates with thicknesses roughly proportional !
109! to depth reducing high r-factors in these areas. !
110! !
111! !
112! This generic transformation design facilitates numerous vertical !
113! stretching functions, C(s). These functions are set-up in this !
114! routine in terms of several stretching parameters specified in !
115! the standard input file. !
116! !
117! C(s) vertical stretching function properties: !
118! !
119! * a nonlinear, monotonic function !
120! * a continuous differentiable function, or !
121! * a piecewise function with smooth transition and differentiable !
122! * must be constrained by -1 <= C(s) <= 0, with C(0)=0 at the !
123! free-surface and C(-1)=-1 at the bottom (bathymetry). !
124! !
125! References: !
126! !
127! Shchepetkin, A.F. and J.C. McWilliams, 2005: The regional oceanic !
128! modeling system (ROMS): a split-explicit, free-surface, !
129! topography-following-coordinate oceanic model, Ocean !
130! Modelling, 9, 347-404. !
131! !
132! Song, Y. and D. Haidvogel, 1994: A semi-implicit ocean !
133! circulation model using a generalized topography- !
134! following coordinate system, J. Comp. Physics, !
135! 115, 228-244. !
136! !
137!=======================================================================
138!
139 USE mod_param
140 USE mod_parallel
141 USE mod_grid
142 USE mod_iounits
143 USE mod_scalars
144!
145 implicit none
146!
147! Imported variable declarations.
148!
149 integer, intent(in) :: ng
150!
151! Local variable declarations.
152!
153 integer :: k
154
155 real(dp) :: Aweight, Bweight, Cweight, Cbot, Csur, Hscale
156 real(dp) :: ds, exp_bot, exp_sur, rk, rN, sc_r, sc_w
157 real(dp) :: cff, cff1, cff2
158 real(dp) :: zhc, z1, z2, z3
159!
160!-----------------------------------------------------------------------
161! Set thickness controlling vertical coordinate stretching.
162!-----------------------------------------------------------------------
163!
164! Set hc <= hmin, in the original formulation (Vtransform=1) to avoid
165! [h(x,y)-hc] to be negative which results in dz/ds to be negative.
166! Notice that this restriction is REMOVED in the new transformation
167! (Vtransform=2): hc can be any value. It works for both hc < hmin
168! and hc > hmin.
169!
170 IF (vtransform(ng).eq.1) THEN
171# if defined WET_DRY
172 hc(ng)=min(max(hmin(ng),0.0_dp),tcline(ng))
173# else
174 hc(ng)=min(hmin(ng),tcline(ng))
175# endif
176 ELSE IF (vtransform(ng).eq.2) THEN
177 hc(ng)=tcline(ng)
178 END IF
179!
180!-----------------------------------------------------------------------
181! Original vertical strectching function, Song and Haidvogel (1994).
182!-----------------------------------------------------------------------
183!
184 IF (vstretching(ng).eq.1) THEN
185!
186! This vertical stretching function is defined as:
187!
188! C(s) = (1 - b) * [SINH(s * a) / SINH(a)] +
189!
190! b * [-0.5 + 0.5 * TANH(a * (s + 0.5)) / TANH(0.5 * a)]
191!
192! where the stretching parameters (a, b) are specify at input:
193!
194! a = theta_s 0 < theta_s <= 8
195! b = theta_b 0 <= theta_b <= 1
196!
197! If theta_b=0, the refinement is surface intensified as theta_s is
198! increased.
199! If theta_b=1, the refinement is both bottom ans surface intensified
200! as theta_s is increased.
201!
202 IF (theta_s(ng).ne.0.0_dp) THEN
203 cff1=1.0_dp/sinh(theta_s(ng))
204 cff2=0.5_dp/tanh(0.5_dp*theta_s(ng))
205 END IF
206 scalars(ng)%sc_w(0)=-1.0_dp
207 scalars(ng)%Cs_w(0)=-1.0_dp
208 ds=1.0_dp/real(n(ng),dp)
209 DO k=1,n(ng)
210 scalars(ng)%sc_w(k)=ds*real(k-n(ng),dp)
211 scalars(ng)%sc_r(k)=ds*(real(k-n(ng),dp)-0.5_dp)
212 IF (theta_s(ng).ne.0.0_dp) THEN
213 scalars(ng)%Cs_w(k)=(1.0_dp-theta_b(ng))* &
214 & cff1*sinh(theta_s(ng)* &
215 & scalars(ng)%sc_w(k))+ &
216 & theta_b(ng)* &
217 & (cff2*tanh(theta_s(ng)* &
218 & (scalars(ng)%sc_w(k)+ &
219 & 0.5_dp))- &
220 & 0.5_dp)
221 scalars(ng)%Cs_r(k)=(1.0_dp-theta_b(ng))* &
222 & cff1*sinh(theta_s(ng)* &
223 & scalars(ng)%sc_r(k))+ &
224 & theta_b(ng)* &
225 & (cff2*tanh(theta_s(ng)* &
226 & (scalars(ng)%sc_r(k)+ &
227 & 0.5_dp))- &
228 & 0.5_dp)
229 ELSE
230 scalars(ng)%Cs_w(k)=scalars(ng)%sc_w(k)
231 scalars(ng)%Cs_r(k)=scalars(ng)%sc_r(k)
232 END IF
233 END DO
234!
235!-----------------------------------------------------------------------
236! A. Shchepetkin vertical stretching function. This function was
237! improved further to allow bottom refiment (see Vstretching=4).
238!-----------------------------------------------------------------------
239!
240 ELSE IF (vstretching(ng).eq.2) THEN
241!
242! This vertical stretching function is defined, in the simplest form,
243! as:
244!
245! C(s) = [1.0 - COSH(theta_s * s)] / [COSH(theta_s) - 1.0]
246!
247! it is similar in meaning to the original vertical stretcing function
248! (Song and Haidvogel, 1994), but note that hyperbolic functions are
249! COSH, and not SINH.
250!
251! Note that the above definition results in
252!
253! -1 <= C(s) <= 0
254!
255! as long as
256!
257! -1 <= s <= 0
258!
259! and, unlike in any previous definition
260!
261! d[C(s)]/ds --> 0 if s --> 0
262!
263! For the purpose of bottom boundary layer C(s) is further modified
264! to allow near-bottom refinement. This is done by blending it with
265! another function.
266!
267 aweight=1.0_dp
268 bweight=1.0_dp
269 ds=1.0_dp/real(n(ng),dp)
270!
271 scalars(ng)%sc_w(n(ng))=0.0_dp
272 scalars(ng)%Cs_w(n(ng))=0.0_dp
273 DO k=n(ng)-1,1,-1
274 sc_w=ds*real(k-n(ng),dp)
275 scalars(ng)%sc_w(k)=sc_w
276 IF (theta_s(ng).gt.0.0_dp) THEN
277 csur=(1.0_dp-cosh(theta_s(ng)*sc_w))/ &
278 & (cosh(theta_s(ng))-1.0_dp)
279 IF (theta_b(ng).gt.0.0_dp) THEN
280 cbot=sinh(theta_b(ng)*(sc_w+1.0_dp))/ &
281 & sinh(theta_b(ng))-1.0_dp
282 cweight=(sc_w+1.0_dp)**aweight* &
283 & (1.0_dp+(aweight/bweight)* &
284 & (1.0_dp-(sc_w+1.0_dp)**bweight))
285 scalars(ng)%Cs_w(k)=cweight*csur+(1.0_dp-cweight)*cbot
286 ELSE
287 scalars(ng)%Cs_w(k)=csur
288 END IF
289 ELSE
290 scalars(ng)%Cs_w(k)=sc_w
291 END IF
292 END DO
293 scalars(ng)%sc_w(0)=-1.0_dp
294 scalars(ng)%Cs_w(0)=-1.0_dp
295!
296 DO k=1,n(ng)
297 sc_r=ds*(real(k-n(ng),dp)-0.5_dp)
298 scalars(ng)%sc_r(k)=sc_r
299 IF (theta_s(ng).gt.0.0_dp) THEN
300 csur=(1.0_dp-cosh(theta_s(ng)*sc_r))/ &
301 & (cosh(theta_s(ng))-1.0_dp)
302 IF (theta_b(ng).gt.0.0_dp) THEN
303 cbot=sinh(theta_b(ng)*(sc_r+1.0_dp))/ &
304 & sinh(theta_b(ng))-1.0_dp
305 cweight=(sc_r+1.0_dp)**aweight* &
306 & (1.0_dp+(aweight/bweight)* &
307 & (1.0_dp-(sc_r+1.0_dp)**bweight))
308 scalars(ng)%Cs_r(k)=cweight*csur+(1.0_dp-cweight)*cbot
309 ELSE
310 scalars(ng)%Cs_r(k)=csur
311 END IF
312 ELSE
313 scalars(ng)%Cs_r(k)=sc_r
314 END IF
315 END DO
316!
317!-----------------------------------------------------------------------
318! R. Geyer stretching function for high bottom boundary layer
319! resolution.
320!-----------------------------------------------------------------------
321!
322 ELSE IF (vstretching(ng).eq.3) THEN
323!
324! This stretching function is intended for very shallow coastal
325! applications, like gravity sediment flows.
326!
327! At the surface, C(s=0)=0
328!
329! C(s) = - LOG(COSH(Hscale * ABS(s) ** alpha)) /
330! LOG(COSH(Hscale))
331!
332! At the bottom, C(s=-1)=-1
333!
334! C(s) = LOG(COSH(Hscale * (s + 1) ** beta)) /
335! LOG(COSH(Hscale)) - 1
336!
337! where
338!
339! Hscale : scale value for all hyperbolic functions
340! Hscale = 3.0 set internally here
341! alpha : surface stretching exponent
342! alpha = 0.65 minimal increase of surface resolution
343! 1.0 significant amplification
344! beta : bottoom stretching exponent
345! beta = 0.58 no amplification
346! 1.0 significant amplification
347! 3.0 super-high bottom resolution
348! s : stretched vertical coordinate, -1 <= s <= 0
349! s(k) = (k-N)/N k=0:N, W-points (s_w)
350! s(k) = (k-N-0.5)/N k=1:N, RHO-points (s_rho)
351!
352! The stretching exponents (alpha, beta) are specify at input:
353!
354! alpha = theta_s
355! beta = theta_b
356!
357 exp_sur=theta_s(ng)
358 exp_bot=theta_b(ng)
359 hscale=3.0_dp
360 ds=1.0_dp/real(n(ng),dp)
361!
362 scalars(ng)%sc_w(n(ng))=0.0_dp
363 scalars(ng)%Cs_w(n(ng))=0.0_dp
364 DO k=n(ng)-1,1,-1
365 sc_w=ds*real(k-n(ng),dp)
366 scalars(ng)%sc_w(k)=sc_w
367 cbot= log(cosh(hscale*(sc_w+1.0_dp)**exp_bot))/ &
368 & log(cosh(hscale))-1.0_dp
369 csur=-log(cosh(hscale*abs(sc_w)**exp_sur))/ &
370 & log(cosh(hscale))
371 cweight=0.5_dp*(1.0_dp-tanh(hscale*(sc_w+0.5_dp)))
372 scalars(ng)%Cs_w(k)=cweight*cbot+(1.0_dp-cweight)*csur
373 END DO
374 scalars(ng)%sc_w(0)=-1.0_dp
375 scalars(ng)%Cs_w(0)=-1.0_dp
376!
377 DO k=1,n(ng)
378 sc_r=ds*(real(k-n(ng),dp)-0.5_dp)
379 scalars(ng)%sc_r(k)=sc_r
380 cbot= log(cosh(hscale*(sc_r+1.0_dp)**exp_bot))/ &
381 & log(cosh(hscale))-1.0_dp
382 csur=-log(cosh(hscale*abs(sc_r)**exp_sur))/ &
383 & log(cosh(hscale))
384 cweight=0.5_dp*(1.0_dp-tanh(hscale*(sc_r+0.5_dp)))
385 scalars(ng)%Cs_r(k)=cweight*cbot+(1.0_dp-cweight)*csur
386 END DO
387!
388!-----------------------------------------------------------------------
389! A. Shchepetkin improved double vertical stretching functions with
390! bottom refiment.
391!-----------------------------------------------------------------------
392!
393 ELSE IF (vstretching(ng).eq.4) THEN
394!
395! The range of meaningful values for the control parameters are:
396!
397! 0 < theta_s <= 10.0
398! 0 <= theta_b <= 3.0
399!
400! Users need to pay attention to extreme r-factor (rx1) values near
401! the bottom.
402!
403! This vertical stretching function is defined, in the simplest form,
404! as:
405!
406! C(s) = [1.0 - COSH(theta_s * s)] / [COSH(theta_s) - 1.0]
407!
408! it is similar in meaning to the original vertical stretcing function
409! (Song and Haidvogel, 1994), but note that hyperbolic functions are
410! COSH, and not SINH.
411!
412! Note that the above definition results in
413!
414! -1 <= C(s) <= 0
415!
416! as long as
417!
418! -1 <= s <= 0
419!
420! and
421!
422! d[C(s)]/ds --> 0 if s --> 0
423!
424! For the purpose of bottom boundary layer C(s) is further modified
425! to allow near-bottom refinement by using a continuous, second
426! stretching function
427!
428! C(s) = [EXP(theta_b * C(s)) - 1.0] / [1.0 - EXP(-theta_b)]
429!
430! This double transformation is continuous with respect to "theta_s"
431! and "theta_b", as both values approach to zero.
432!
433 ds=1.0_dp/real(n(ng),dp)
434!
435 scalars(ng)%sc_w(n(ng))=0.0_dp
436 scalars(ng)%Cs_w(n(ng))=0.0_dp
437 DO k=n(ng)-1,1,-1
438 sc_w=ds*real(k-n(ng),dp)
439 scalars(ng)%sc_w(k)=sc_w
440 IF (theta_s(ng).gt.0.0_dp) THEN
441 csur=(1.0_dp-cosh(theta_s(ng)*sc_w))/ &
442 & (cosh(theta_s(ng))-1.0_dp)
443 ELSE
444 csur=-sc_w**2
445 END IF
446 IF (theta_b(ng).gt.0.0_dp) THEN
447 cbot=(exp(theta_b(ng)*csur)-1.0_dp)/ &
448 & (1.0_dp-exp(-theta_b(ng)))
449 scalars(ng)%Cs_w(k)=cbot
450 ELSE
451 scalars(ng)%Cs_w(k)=csur
452 END IF
453 END DO
454 scalars(ng)%sc_w(0)=-1.0_dp
455 scalars(ng)%Cs_w(0)=-1.0_dp
456!
457 DO k=1,n(ng)
458 sc_r=ds*(real(k-n(ng),dp)-0.5_dp)
459 scalars(ng)%sc_r(k)=sc_r
460 IF (theta_s(ng).gt.0.0_dp) THEN
461 csur=(1.0_dp-cosh(theta_s(ng)*sc_r))/ &
462 & (cosh(theta_s(ng))-1.0_dp)
463 ELSE
464 csur=-sc_r**2
465 END IF
466 IF (theta_b(ng).gt.0.0_dp) THEN
467 cbot=(exp(theta_b(ng)*csur)-1.0_dp)/ &
468 & (1.0_dp-exp(-theta_b(ng)))
469 scalars(ng)%Cs_r(k)=cbot
470 ELSE
471 scalars(ng)%Cs_r(k)=csur
472 END IF
473 END DO
474!
475!----------------------------------------------------------------------
476! Stretching 5 case using a quadratic Legendre polynomial function
477! aproach for the s-coordinate to enhance the surface exchange layer.
478!
479! J. Souza, B.S. Powell, A.C. Castillo-Trujillo, and P. Flament, 2015:
480! The Vorticity Balance of the Ocean Surface in Hawaii from a
481! Regional Reanalysis.'' J. Phys. Oceanogr., 45, 424-440.
482!
483! Added by Joao Marcos Souza - SOEST - 05/07/2012.
484!-----------------------------------------------------------------------
485!
486 ELSE IF (vstretching(ng).eq.5) THEN
487 scalars(ng)%sc_w(n(ng))=0.0_dp
488 scalars(ng)%Cs_w(n(ng))=0.0_dp
489 DO k=n(ng)-1,1,-1
490 rk=real(k,dp)
491 rn=real(n(ng),dp)
492 sc_w=-(rk*rk - 2.0_dp*rk*rn + rk + rn*rn - rn)/(rn*rn - rn)- &
493 0.01_dp*(rk*rk - rk*rn)/(1.0_dp - rn)
494 scalars(ng)%sc_w(k)=sc_w
495 IF (theta_s(ng).gt.0.0_dp) THEN
496 csur=(1.0_dp-cosh(theta_s(ng)*sc_w))/ &
497 & (cosh(theta_s(ng))-1.0_dp)
498 ELSE
499 csur=-sc_w**2
500 END IF
501 IF (theta_b(ng).gt.0.0_dp) THEN
502 cbot=(exp(theta_b(ng)*csur)-1.0_dp)/ &
503 & (1.0_dp-exp(-theta_b(ng)))
504 scalars(ng)%Cs_w(k)=cbot
505 ELSE
506 scalars(ng)%Cs_w(k)=csur
507 END IF
508 END DO
509 scalars(ng)%sc_w(0)=-1.0_dp
510 scalars(ng)%Cs_w(0)=-1.0_dp
511!
512 DO k=1,n(ng)
513 rk=real(k,dp)-0.5_dp
514 rn=real(n(ng),dp)
515 sc_r=-(rk*rk - 2.0_dp*rk*rn + rk + rn*rn - rn)/(rn*rn - rn)- &
516 0.01_dp*(rk*rk - rk*rn)/(1.0_dp - rn)
517 scalars(ng)%sc_r(k)=sc_r
518 IF (theta_s(ng).gt.0.0_dp) THEN
519 csur=(1.0_dp-cosh(theta_s(ng)*sc_r))/ &
520 & (cosh(theta_s(ng))-1.0_dp)
521 ELSE
522 csur=-sc_r**2
523 END IF
524 IF (theta_b(ng).gt.0.0_dp) THEN
525 cbot=(exp(theta_b(ng)*csur)-1.0_dp)/ &
526 & (1.0_dp-exp(-theta_b(ng)))
527 scalars(ng)%Cs_r(k)=cbot
528 ELSE
529 scalars(ng)%Cs_r(k)=csur
530 END IF
531 END DO
532 END IF
533!
534!-----------------------------------------------------------------------
535! Report information about vertical transformation.
536!-----------------------------------------------------------------------
537!
538 IF (master.and.lwrtinfo(ng)) THEN
539 WRITE (stdout,10) ng
540 cff=0.5_dp*(hmax(ng)+hmin(ng))
541 DO k=n(ng),0,-1
542 IF (vtransform(ng).eq.1) THEN
543 zhc=hc(ng)*scalars(ng)%sc_w(k)
544 z1=zhc+(hmin(ng)-hc(ng))*scalars(ng)%Cs_w(k)
545 z2=zhc+(cff -hc(ng))*scalars(ng)%Cs_w(k)
546 z3=zhc+(hmax(ng)-hc(ng))*scalars(ng)%Cs_w(k)
547 ELSE IF (vtransform(ng).eq.2) THEN
548 z1=hmin(ng)*(hc(ng) *scalars(ng)%sc_w(k)+ &
549 & hmin(ng)*scalars(ng)%Cs_w(k))/(hc(ng)+hmin(ng))
550 z2=cff *(hc(ng) *scalars(ng)%sc_w(k)+ &
551 & cff *scalars(ng)%Cs_w(k))/(hc(ng)+cff)
552 z3=hmax(ng)*(hc(ng) *scalars(ng)%sc_w(k)+ &
553 & hmax(ng)*scalars(ng)%Cs_w(k))/(hc(ng)+hmax(ng))
554 IF (hc(ng).gt.hmax(ng)) THEN
555 zhc=z3 ! same as hmax, other values do not make sense
556 ELSE
557 zhc=0.5_dp*hc(ng)*(scalars(ng)%sc_w(k)+ &
558 & scalars(ng)%Cs_w(k))
559 END IF
560 END IF
561 WRITE (stdout,20) k, scalars(ng)%sc_w(k), &
562 & scalars(ng)%Cs_w(k), z1, zhc, z2, z3
563 END DO
564 END IF
565
566 10 FORMAT (/,' Vertical S-coordinate System, Grid ',i2.2,':'/,/, &
567 & ' level S-coord Cs-curve Z',3x, &
568 & 'at hmin at hc half way at hmax',/)
569 20 FORMAT (i6,2f12.7,1x,4f12.3)
570
571 RETURN
integer stdout
logical master
integer, dimension(:), allocatable n
Definition mod_param.F:479
real(dp), dimension(:), allocatable hmin
real(dp), dimension(:), allocatable theta_s
real(dp), dimension(:), allocatable tcline
real(dp), dimension(:), allocatable theta_b
logical, dimension(:), allocatable lwrtinfo
real(dp), dimension(:), allocatable hc
type(t_scalars), dimension(:), allocatable scalars
Definition mod_scalars.F:65
real(dp), dimension(:), allocatable hmax
integer, dimension(:), allocatable vstretching
integer, dimension(:), allocatable vtransform

References mod_scalars::hc, mod_scalars::hmax, mod_scalars::hmin, mod_scalars::lwrtinfo, mod_parallel::master, mod_param::n, mod_scalars::scalars, mod_iounits::stdout, mod_scalars::tcline, mod_scalars::theta_b, mod_scalars::theta_s, mod_scalars::vstretching, and mod_scalars::vtransform.

Referenced by set_grid().

Here is the caller graph for this function: