WAVEWATCH III  beta 0.0.1
constants Module Reference

Define some much-used constants for global use (all defined as PARAMETER). More...

Functions/Subroutines

subroutine tabu_fw
 Estimate friction coefficients in oscillatory boundary layers using tabulation on Kelvin functions. More...
 
subroutine kzeone (X, Y, RE0, IM0, RE1, IM1)
 June 1999 adaptation to CRESTb, all tests on range of (x,y) have been bypassed, we implicitly expect X to be positive or |x,y| non zero. More...
 
subroutine kerkei (X, KER, KEI)
 Computes the values of the zeroth order Kelvin function Ker and Kei. More...
 

Variables

logical, parameter tstout = .FALSE.
 TSTOUT Flag for generation of test files. More...
 
real, parameter grav = 9.806
 GRAV Acc. More...
 
real, parameter dwat = 1000.
 DWAT Density of water (kg/m3). More...
 
real, parameter dair = 1.225
 DAIR Density of air (kg/m3). More...
 
real, parameter nu_air = 1.4E-5
 NU_AIR Kinematic viscosity of air (m2/s). More...
 
real, parameter nu_water = 1.31E-6
 NU_WATER Kinematic viscosity of water (m2/s). More...
 
real, parameter sed_sg = 2.65
 SED_SG Specific gravity of sediments (N.D.). More...
 
real, parameter kappa = 0.40
 KAPPA von Karman's constant (N.D.). More...
 
real, parameter pi = 3.141592653589793
 PI Value of Pi. More...
 
real, parameter tpi = 2.0 * PI
 TPI 2*Pi. More...
 
real, parameter hpi = 0.5 * PI
 HPI 1/2*Pi. More...
 
real, parameter tpiinv = 1. / TPI
 TPIINV Inverse of 2*Pi. More...
 
real, parameter hpiinv = 1. / HPI
 HPIINV Inverse of 1/2*Pi. More...
 
real, parameter rade = 180. / PI
 RADE Conversion factor from radians to degrees. More...
 
real, parameter dera = PI / 180.
 DERA Conversion factor from degrees to radians. More...
 
real, parameter radius = 4.E7 * TPIINV
 RADIUS Radius of the earth (m). More...
 
real, parameter g2pi3i = 1. / ( GRAV**2 * TPI**3 )
 G2PI3I Inverse of gravity^2 * (2*Pi)^3. More...
 
real, parameter g1pi1i = 1. / ( GRAV * TPI )
 G1PI1I Inverse of gravity * 2 * Pi. More...
 
real undef = -999.9
 UNDEF Value for undefined variable in output. More...
 
character(*), parameter file_endian = ENDIANNESS
 FILE_ENDIAN Filled by preprocessor with 'big_endian', 'little_endian', or 'native'. More...
 
integer, parameter sizefwtable =300
 SIZEFWTABLE. More...
 
real, dimension(0:sizefwtablefwtable
 FWTABLE. More...
 
real delab
 DELAB. More...
 
real, parameter abmin = -1.
 ABMIN. More...
 
integer, parameter srce_direct = 0
 srce_direct More...
 
integer, parameter srce_imp_post = 1
 srce_imp_post More...
 
integer, parameter srce_imp_pre = 2
 srce_imp_pre More...
 
integer, parameter debug_node = 1014
 DEBUG_NODE Node number used for debugging. More...
 
integer, parameter debug_element = 50
 DEBUG_ELEMENT Element number used for debug. More...
 
logical lpdlib = .FALSE.
 LPDLIB Logical for using the PDLIB library. More...
 
logical lsetup = .FALSE.
 LSETUP Logical LSETUP is not used. More...
 
logical is_esmf_component = .FALSE.
 IS_ESMF_COMPONENT Flag for model invoked via ESMF. More...
 

Detailed Description

Define some much-used constants for global use (all defined as PARAMETER).

Author
H. L. Tolman
Date
05-Jun-2018

Function/Subroutine Documentation

◆ kerkei()

subroutine constants::kerkei ( real  X,
real  KER,
real  KEI 
)

Computes the values of the zeroth order Kelvin function Ker and Kei.

These functions are used to determine the friction factor fw as a function of the bottom roughness length assuming a linear profile of eddy viscosity (See Grant and Madsen, 1979).

Parameters
X
KER
KEI
Author
N/A
Date
N/A

Definition at line 424 of file constants.F90.

424  !**********************************************************************
425  ! Computes the values of the zeroth order Kelvin function Ker and Kei
426  ! These functions are used to determine the friction factor fw as a
427  ! function of the bottom roughness length assuming a linear profile
428  ! of eddy viscosity (See Grant and Madsen, 1979)
429  !**********************************************************************
430  IMPLICIT NONE
431 
432  DOUBLE PRECISION ZR,ZI,CYR,CYI,CYR1,CYI1
433  REAL X,KER,KEI
434 
435  zr=x*.50d0*sqrt(2.0d0)
436  zi=zr
437  CALL kzeone(zr, zi, cyr, cyi,cyr1,cyi1)
438  ker=cyr/exp(zr)
439  kei=cyi/exp(zr)

References kzeone().

Referenced by tabu_fw().

◆ kzeone()

subroutine constants::kzeone ( double precision  X,
double precision  Y,
double precision  RE0,
double precision  IM0,
double precision  RE1,
double precision  IM1 
)

June 1999 adaptation to CRESTb, all tests on range of (x,y) have been bypassed, we implicitly expect X to be positive or |x,y| non zero.

The variables X and Y are the real and imaginary parts of the argument of the first two modified bessel functions of the second kind,k0 and k1. Re0,im0,re1 and im1 give the real and imaginary parts of exp(x)*k0 and exp(x)*k1, respectively. Although the real notation used in this subroutine may seem inelegant when compared with the complex notation that fortran allows, this version runs about 30 percent faster than one written using complex variables.

Parameters
XReal part of argument to modified Bessel functions.
YImaginary part of argument to modified Bessel functions.
RE0Real part of exp(x)*k0.
IM0Imaginary part of exp(x)*k0.
RE1Real part of exp(x)*k1.
IM1Imaginary part of exp(x)*k1.
Author
N/A
Date
N/A

Definition at line 244 of file constants.F90.

244  ! June 1999 adaptation to CRESTb, all tests on range of (x,y) have been
245  ! bypassed, we implicitly expect X to be positive or |x,y| non zero
246  !
247  ! This subroutine is copyright by ACM
248  ! see http://www.acm.org/pubs/copyright_policy/softwareCRnotice.html
249  ! ACM declines any responsibility of any kind
250  !
251  ! THE VARIABLES X AND Y ARE THE REAL AND IMAGINARY PARTS OF
252  ! THE ARGUMENT OF THE FIRST TWO MODIFIED BESSEL FUNCTIONS
253  ! OF THE SECOND KIND,K0 AND K1. RE0,IM0,RE1 AND IM1 GIVE
254  ! THE REAL AND IMAGINARY PARTS OF EXP(X)*K0 AND EXP(X)*K1,
255  ! RESPECTIVELY. ALTHOUGH THE REAL NOTATION USED IN THIS
256  ! SUBROUTINE MAY SEEM INELEGANT WHEN COMPARED WITH THE
257  ! COMPLEX NOTATION THAT FORTRAN ALLOWS, THIS VERSION RUNS
258  ! ABOUT 30 PERCENT FASTER THAN ONE WRITTEN USING COMPLEX
259  ! VARIABLES.
260  ! ACM Libraries
261  ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
262  IMPLICIT NONE
263  DOUBLE PRECISION X, Y, X2, Y2, RE0, IM0, RE1, IM1, &
264  R1, R2, T1, T2, P1, P2, RTERM, ITERM, L
265  DOUBLE PRECISION , PARAMETER, DIMENSION(8) :: EXSQ = &
266  (/ 0.5641003087264d0,0.4120286874989d0,0.1584889157959d0, &
267  0.3078003387255d-1,0.2778068842913d-2,0.1000044412325d-3, &
268  0.1059115547711d-5,0.1522475804254d-8 /)
269  DOUBLE PRECISION , PARAMETER, DIMENSION(8) :: TSQ = &
270  (/ 0.0d0,3.19303633920635d-1,1.29075862295915d0, &
271  2.95837445869665d0,5.40903159724444d0,8.80407957805676d0, &
272  1.34685357432515d1,2.02499163658709d1 /)
273  INTEGER N,M,K
274  ! THE ARRAYS TSQ AND EXSQ CONTAIN THE SQUARE OF THE
275  ! ABSCISSAS AND THE WEIGHT FACTORS USED IN THE GAUSS-
276  ! HERMITE QUADRATURE.
277  r2 = x*x + y*y
278  IF (r2.GE.1.96d2) GO TO 50
279  IF (r2.GE.1.849d1) GO TO 30
280  ! THIS SECTION CALCULATES THE FUNCTIONS USING THE SERIES
281  ! EXPANSIONS
282  x2 = x/2.0d0
283  y2 = y/2.0d0
284  p1 = x2*x2
285  p2 = y2*y2
286  t1 = -(dlog(p1+p2)/2.0d0+0.5772156649015329d0)
287  ! THE CONSTANT IN THE PRECEDING STATEMENT IS EULER*S
288  ! CONSTANT
289  t2 = -datan2(y,x)
290  x2 = p1 - p2
291  y2 = x*y2
292  rterm = 1.0d0
293  iterm = 0.0d0
294  re0 = t1
295  im0 = t2
296  t1 = t1 + 0.5d0
297  re1 = t1
298  im1 = t2
299  p2 = dsqrt(r2)
300  l = 2.106d0*p2 + 4.4d0
301  IF (p2.LT.8.0d-1) l = 2.129d0*p2 + 4.0d0
302  DO n=1,int(l)
303  p1 = n
304  p2 = n*n
305  r1 = rterm
306  rterm = (r1*x2-iterm*y2)/p2
307  iterm = (r1*y2+iterm*x2)/p2
308  t1 = t1 + 0.5d0/p1
309  re0 = re0 + t1*rterm - t2*iterm
310  im0 = im0 + t1*iterm + t2*rterm
311  p1 = p1 + 1.0d0
312  t1 = t1 + 0.5d0/p1
313  re1 = re1 + (t1*rterm-t2*iterm)/p1
314  im1 = im1 + (t1*iterm+t2*rterm)/p1
315  END DO
316  r1 = x/r2 - 0.5d0*(x*re1-y*im1)
317  r2 = -y/r2 - 0.5d0*(x*im1+y*re1)
318  p1 = dexp(x)
319  re0 = p1*re0
320  im0 = p1*im0
321  re1 = p1*r1
322  im1 = p1*r2
323  RETURN
324  ! THIS SECTION CALCULATES THE FUNCTIONS USING THE INTEGRAL
325  ! REPRESENTATION, EQN 3, EVALUATED WITH 15 POINT GAUSS-
326  ! HERMITE QUADRATURE
327 30 x2 = 2.0d0*x
328  y2 = 2.0d0*y
329  r1 = y2*y2
330  p1 = dsqrt(x2*x2+r1)
331  p2 = dsqrt(p1+x2)
332  t1 = exsq(1)/(2.0d0*p1)
333  re0 = t1*p2
334  im0 = t1/p2
335  re1 = 0.0d0
336  im1 = 0.0d0
337  DO n=2,8
338  t2 = x2 + tsq(n)
339  p1 = dsqrt(t2*t2+r1)
340  p2 = dsqrt(p1+t2)
341  t1 = exsq(n)/p1
342  re0 = re0 + t1*p2
343  im0 = im0 + t1/p2
344  t1 = exsq(n)*tsq(n)
345  re1 = re1 + t1*p2
346  im1 = im1 + t1/p2
347  END DO
348  t2 = -y2*im0
349  re1 = re1/r2
350  r2 = y2*im1/r2
351  rterm = 1.41421356237309d0*dcos(y)
352  iterm = -1.41421356237309d0*dsin(y)
353  ! THE CONSTANT IN THE PREVIOUS STATEMENTS IS,OF COURSE,
354  ! SQRT(2.0).
355  im0 = re0*iterm + t2*rterm
356  re0 = re0*rterm - t2*iterm
357  t1 = re1*rterm - r2*iterm
358  t2 = re1*iterm + r2*rterm
359  re1 = t1*x + t2*y
360  im1 = -t1*y + t2*x
361  RETURN
362  ! THIS SECTION CALCULATES THE FUNCTIONS USING THE
363  ! ASYMPTOTIC EXPANSIONS
364 50 rterm = 1.0d0
365  iterm = 0.0d0
366  re0 = 1.0d0
367  im0 = 0.0d0
368  re1 = 1.0d0
369  im1 = 0.0d0
370  p1 = 8.0d0*r2
371  p2 = dsqrt(r2)
372  l = 3.91d0+8.12d1/p2
373  r1 = 1.0d0
374  r2 = 1.0d0
375  m = -8
376  k = 3
377  DO n=1,int(l)
378  m = m + 8
379  k = k - m
380  r1 = float(k-4)*r1
381  r2 = float(k)*r2
382  t1 = float(n)*p1
383  t2 = rterm
384  rterm = (t2*x+iterm*y)/t1
385  iterm = (-t2*y+iterm*x)/t1
386  re0 = re0 + r1*rterm
387  im0 = im0 + r1*iterm
388  re1 = re1 + r2*rterm
389  im1 = im1 + r2*iterm
390  END DO
391  t1 = dsqrt(p2+x)
392  t2 = -y/t1
393  p1 = 8.86226925452758d-1/p2
394  ! THIS CONSTANT IS SQRT(PI)/2.0, WITH PI=3.14159...
395  rterm = p1*dcos(y)
396  iterm = -p1*dsin(y)
397  r1 = re0*rterm - im0*iterm
398  r2 = re0*iterm + im0*rterm
399  re0 = t1*r1 - t2*r2
400  im0 = t1*r2 + t2*r1
401  r1 = re1*rterm - im1*iterm
402  r2 = re1*iterm + im1*rterm
403  re1 = t1*r1 - t2*r2
404  im1 = t1*r2 + t2*r1
405  RETURN

Referenced by kerkei().

◆ tabu_fw()

subroutine constants::tabu_fw

Estimate friction coefficients in oscillatory boundary layers using tabulation on Kelvin functions.

Author
F. Ardhuin
Date
28-Feb-2013

Definition at line 120 of file constants.F90.

120  !/
121  !/ +-----------------------------------+
122  !/ | WAVEWATCH III NOAA/NCEP |
123  !/ | F. Ardhuin |
124  !/ | FORTRAN 90 |
125  !/ | Last update : 28-Feb-2013 |
126  !/ +-----------------------------------+
127  !/
128  !/ 19-Oct-2007 : Origination. ( version 3.13 )
129  !/ 28-Feb-2013 : Caps the friction factor to 0.5 ( version 4.08 )
130  !/
131  ! 1. Purpose :
132  ! TO estimate friction coefficients in oscillatory boundary layers
133  ! METHOD.
134  ! tabulation on Kelvin functions
135  !
136  ! 2. Method :
137  !
138  ! 3. Parameters :
139  !
140  ! Parameter list
141  ! ----------------------------------------------------------------
142  ! ----------------------------------------------------------------
143  !
144  ! 4. Subroutines used :
145  !
146  ! Name Type Module Description
147  ! ----------------------------------------------------------------
148  ! STRACE Subr. W3SERVMD Subroutine tracing.
149  ! ----------------------------------------------------------------
150  !
151  ! 5. Called by :
152  !
153  ! Name Type Module Description
154  ! ----------------------------------------------------------------
155  ! WW3_GRID Prog. WW3_GRID Model grid initialization
156  ! ----------------------------------------------------------------
157  !
158  ! 6. Error messages :
159  !
160  ! None.
161  !
162  ! 7. Remarks :
163  !
164  ! 8. Structure :
165  !
166  ! See source code.
167  !
168  ! 9. Switches :
169  !
170  ! !/S Enable subroutine tracing.
171  !
172  ! 10. Source code :
173  !
174  !/ ------------------------------------------------------------------- /
175  IMPLICIT NONE
176  INTEGER, PARAMETER :: NITER=100
177  REAL , PARAMETER :: XM=0.50
178  REAL , PARAMETER :: EPS1=0.00001
179  ! ----------------------------------------------------------------------
180  INTEGER I,ITER
181  REAL KER, KEI
182  REAL ABR,ABRLOG,L10,FACT,FSUBW,FSUBWMEMO,dzeta0,dzeta0memo
183  !
184  delab = (abmax-abmin)/real(sizefwtable)
185  l10=alog(10.)
186  DO i=0,sizefwtable
187  !
188  ! index I in this table corresponds to a normalized roughness z0/ABR = 10^ABMIN+REAL(I)*DELAB
189  !
190  abrlog=abmin+real(i)*delab
191  abr=exp(abrlog*l10)
192  fact=1/abr/(21.2*kappa)
193  fsubw=0.05
194  dzeta0=0.
195  DO iter=1,niter
196  fsubwmemo=fsubw
197  dzeta0memo=dzeta0
198  dzeta0=fact*fsubw**(-.5)
199  CALL kerkei(2.*sqrt(dzeta0),ker,kei)
200  fsubw=.08/(ker**2+kei**2)
201  fsubw=.5*(fsubwmemo+fsubw)
202  dzeta0=.5*(dzeta0memo+dzeta0)
203  END DO
204  !
205  ! Maximum value of 0.5 for fe is based on field
206  ! and lab experiment by Lowe et al. JGR 2005, 2007
207  !
208  fwtable(i) = min(fsubw,0.5)
209  ! WRITE(994,*) 'Friction factor:',I,ABR,FWTABLE(I)
210  END DO
211  RETURN

References abmin, delab, fwtable, kappa, kerkei(), and sizefwtable.

Referenced by w3gridmd::w3grid().

Variable Documentation

◆ abmin

real, parameter constants::abmin = -1.

ABMIN.

Definition at line 94 of file constants.F90.

94  REAL, PARAMETER :: ABMIN = -1.

Referenced by tabu_fw(), w3sbt4md::w3sbt4(), w3sic2md::w3sic2(), and w3sic3md::w3sic3().

◆ dair

real, parameter constants::dair = 1.225

DAIR Density of air (kg/m3).

Definition at line 63 of file constants.F90.

63  REAL, PARAMETER :: DAIR = 1.225

Referenced by w3wdatmd::w3dimw(), w3exnc(), w3expo(), w3src6md::w3sds6(), w3wavemd::w3wave(), and wmupdtmd::wmupd2().

◆ debug_element

integer, parameter constants::debug_element = 50

DEBUG_ELEMENT Element number used for debug.

Definition at line 100 of file constants.F90.

100  INTEGER, PARAMETER :: DEBUG_ELEMENT = 50

◆ debug_node

integer, parameter constants::debug_node = 1014

DEBUG_NODE Node number used for debugging.

Definition at line 99 of file constants.F90.

99  INTEGER, PARAMETER :: DEBUG_NODE = 1014

Referenced by w3sdb1md::w3sdb1(), w3src4md::w3sds4(), and w3wavemd::w3wave().

◆ delab

real constants::delab

DELAB.

Definition at line 93 of file constants.F90.

93  REAL :: DELAB

Referenced by tabu_fw(), w3sbt4md::w3sbt4(), w3sic2md::w3sic2(), and w3sic3md::w3sic3().

◆ dera

◆ dwat

◆ file_endian

character(*), parameter constants::file_endian = ENDIANNESS

◆ fwtable

real, dimension(0:sizefwtable) constants::fwtable

FWTABLE.

Definition at line 92 of file constants.F90.

92  REAL :: FWTABLE(0:SIZEFWTABLE)

Referenced by tabu_fw(), w3sbt4md::w3sbt4(), w3sic2md::w3sic2(), w3sic3md::w3sic3(), and w3src4md::w3sin4().

◆ g1pi1i

real, parameter constants::g1pi1i = 1. / ( GRAV * TPI )

G1PI1I Inverse of gravity * 2 * Pi.

Definition at line 82 of file constants.F90.

82  REAL, PARAMETER :: G1PI1I = 1. / ( grav * tpi )

Referenced by w3src2md::w3sds2().

◆ g2pi3i

real, parameter constants::g2pi3i = 1. / ( GRAV**2 * TPI**3 )

G2PI3I Inverse of gravity^2 * (2*Pi)^3.

Definition at line 81 of file constants.F90.

81  REAL, PARAMETER :: G2PI3I = 1. / ( grav**2 * tpi**3 )

Referenced by w3src2md::w3sds2().

◆ grav

real, parameter constants::grav = 9.806

GRAV Acc.

of gravity (m/s2).

Definition at line 61 of file constants.F90.

61  REAL, PARAMETER :: GRAV = 9.806

Referenced by pdlib_w3profsmd::action_limiter_local(), w3iogomd::calc_u3stokes(), w3src3md::calc_ustar(), w3src4md::calc_ustar(), pdlib_w3profsmd::calcarray_jacobi_source_1(), pdlib_w3profsmd::calcarray_jacobi_source_2(), w3snl1md::couple(), w3gig1md::df1f2theta(), w3dispmd::distab(), w3sic3md::drfun_dble_cheng(), w3sic3md::drfun_quad_cheng(), wmesmfmd::fieldindex(), w3sic5md::fsdisp(), w3sic3md::ic3table_cheng(), w3src4md::insin4(), w3snl2md::insnl2(), w3snl3md::insnl3(), w3snl1md::insnlgqm(), w3snlsmd::insnls(), w3dispmd::liu_reverse_dispersion(), pdlib_w3profsmd::pdlib_init(), pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block(), w3iogomd::secondhh(), w3fld1md::sig2wn(), w3iogomd::skewness(), w3src3md::tabu_stress(), w3src4md::tabu_stress(), w3src3md::tabu_tauhf(), w3src4md::tabu_tauhf(), w3src4md::tabu_tauhf2(), w3src6md::tau_wave_atmos(), w3wavset::trig_compute_lh_stress(), w3canomd::w3add2ndorder(), w3exnc(), w3expo(), w3fld1md::w3fld1(), w3fld2md::w3fld2(), w3flx2md::w3flx2(), w3flx3md::w3flx3(), w3flx5md::w3flx5(), w3sic3md::w3ic3wncg_v1(), w3iogomd::w3outg(), w3psmcmd::w3psmc(), w3sis2md::w3rpwnice(), w3sbs1md::w3sbs1(), w3sbt4md::w3sbt4(), w3sbt8md::w3sbt8(), w3sbt9md::w3sbt9(), w3src2md::w3sds2(), w3src3md::w3sds3(), w3src4md::w3sds4(), w3src6md::w3sds6(), w3sic2md::w3sic2(), w3sic3md::w3sic3(), w3src2md::w3sin2(), w3src3md::w3sin3(), w3src4md::w3sin4(), w3src6md::w3sin6(), w3sis2md::w3sis2(), w3sln1md::w3sln1(), w3snl3md::w3snl3(), w3snl5md::w3snl5(), w3snlsmd::w3snls(), w3src2md::w3spr2(), w3src4md::w3spr4(), w3srcemd::w3srce(), w3ref1md::w3sref(), w3str1md::w3str1(), w3swldmd::w3swl4(), w3swldmd::w3swl6(), w3updtmd::w3uini(), w3updtmd::w3uwnd(), w3pro1md::w3xyp1(), w3pro2md::w3xyp2(), w3pro3md::w3xyp3(), w3dispmd::wavnu1(), w3dispmd::wavnu2(), w3dispmd::wavnu3(), wmgridmd::wmgeql(), and wmgridmd::wmghgh().

◆ hpi

real, parameter constants::hpi = 0.5 * PI

HPI 1/2*Pi.

Definition at line 73 of file constants.F90.

73  REAL, PARAMETER :: HPI = 0.5 * pi

◆ hpiinv

real, parameter constants::hpiinv = 1. / HPI

HPIINV Inverse of 1/2*Pi.

Definition at line 75 of file constants.F90.

75  REAL, PARAMETER :: HPIINV = 1. / hpi

◆ is_esmf_component

logical constants::is_esmf_component = .FALSE.

IS_ESMF_COMPONENT Flag for model invoked via ESMF.

Definition at line 109 of file constants.F90.

109  LOGICAL :: IS_ESMF_COMPONENT = .false.

Referenced by wmesmfmd::setservices(), wminitmd::wminit(), and wminitmd::wminitnml().

◆ kappa

real, parameter constants::kappa = 0.40

KAPPA von Karman's constant (N.D.).

Definition at line 69 of file constants.F90.

69  REAL, PARAMETER :: KAPPA = 0.40

Referenced by w3src4md::calc_ustar(), tabu_fw(), w3src4md::tabu_stress(), w3src4md::tabu_tauhf(), w3src4md::tabu_tauhf2(), w3fld1md::w3fld1(), w3fld2md::w3fld2(), w3flx5md::w3flx5(), and w3src4md::w3sin4().

◆ lpdlib

◆ lsetup

logical constants::lsetup = .FALSE.

LSETUP Logical LSETUP is not used.

Definition at line 102 of file constants.F90.

102  LOGICAL :: LSETUP = .false.

◆ nu_air

real, parameter constants::nu_air = 1.4E-5

NU_AIR Kinematic viscosity of air (m2/s).

Definition at line 64 of file constants.F90.

64  REAL, PARAMETER :: NU_AIR = 1.4e-5

Referenced by w3src4md::calc_ustar(), w3flx5md::w3flx5(), w3src4md::w3sin4(), and w3src4md::w3spr4().

◆ nu_water

real, parameter constants::nu_water = 1.31E-6

NU_WATER Kinematic viscosity of water (m2/s).

Definition at line 67 of file constants.F90.

67  REAL, PARAMETER :: NU_WATER = 1.31e-6

Referenced by w3sbt8md::w3sbt8(), and w3sbt9md::w3sbt9().

◆ pi

◆ rade

◆ radius

◆ sed_sg

real, parameter constants::sed_sg = 2.65

SED_SG Specific gravity of sediments (N.D.).

Definition at line 68 of file constants.F90.

68  REAL, PARAMETER :: SED_SG = 2.65

Referenced by w3sbt4md::w3sbt4().

◆ sizefwtable

integer, parameter constants::sizefwtable =300

SIZEFWTABLE.

Definition at line 91 of file constants.F90.

91  INTEGER, PARAMETER :: SIZEFWTABLE=300

Referenced by tabu_fw(), w3sbt4md::w3sbt4(), w3sic2md::w3sic2(), w3sic3md::w3sic3(), and w3src4md::w3sin4().

◆ srce_direct

integer, parameter constants::srce_direct = 0

srce_direct

Definition at line 96 of file constants.F90.

96  INTEGER, PARAMETER :: srce_direct = 0

Referenced by w3srcemd::w3srce(), and w3wavemd::w3wave().

◆ srce_imp_post

integer, parameter constants::srce_imp_post = 1

srce_imp_post

Definition at line 97 of file constants.F90.

97  INTEGER, PARAMETER :: srce_imp_post = 1

Referenced by w3srcemd::w3srce(), and w3wavemd::w3wave().

◆ srce_imp_pre

integer, parameter constants::srce_imp_pre = 2

srce_imp_pre

Definition at line 98 of file constants.F90.

98  INTEGER, PARAMETER :: srce_imp_pre = 2

Referenced by w3srcemd::w3srce(), and w3wavemd::w3wave().

◆ tpi

real, parameter constants::tpi = 2.0 * PI

TPI 2*Pi.

Definition at line 72 of file constants.F90.

72  REAL, PARAMETER :: TPI = 2.0 * pi

Referenced by pdlib_w3profsmd::action_limiter_local(), w3fld1md::appendtail(), pdlib_w3profsmd::block_solver_init(), w3iogomd::calc_u3stokes(), pdlib_w3profsmd::calcarray_jacobi_source_1(), pdlib_w3profsmd::calcarray_jacobi_source_2(), expand(), wmesmfmd::fieldindex(), w3sic5md::fsdisp(), w3sis2md::insis2(), w3snl5md::insnl5(), w3dispmd::liu_forward_dispersion(), w3dispmd::liu_reverse_dispersion(), pdlib_w3profsmd::pdlib_init(), pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block(), w3partmd::ptmean(), w3iogomd::secondhh(), w3iogomd::skewness(), w3src3md::tabu_tauhf(), w3src6md::tau_wave_atmos(), uvtocart(), w3gig1md::w3addig(), w3bounc(), w3bound(), w3bullmd::w3bull(), w3iosfmd::w3cprt(), w3cspcmd::w3cspc(), w3exnc(), w3expo(), w3fld1md::w3fld1(), w3fld2md::w3fld2(), w3flx2md::w3flx2(), w3flx3md::w3flx3(), w3sic3md::w3ic3wncg_v1(), w3sic5md::w3ic5wncg(), w3initmd::w3init(), w3iopomd::w3iope(), w3pro1md::w3ktp1(), w3pro2md::w3ktp2(), w3pro3md::w3ktp3(), w3iogomd::w3outg(), w3outp(), w3partmd::w3part(), w3sbt4md::w3sbt4(), w3src3md::w3sds3(), w3src4md::w3sds4(), w3src6md::w3sds6(), w3gdatmd::w3setref(), w3sic1md::w3sic1(), w3sic2md::w3sic2(), w3sic3md::w3sic3(), w3sic4md::w3sic4(), w3sic5md::w3sic5(), w3src2md::w3sin2(), w3src3md::w3sin3(), w3src4md::w3sin4(), w3src6md::w3sin6(), w3sis2md::w3sis2(), w3snl1md::w3snl1(), w3snl2md::w3snl2(), w3snl3md::w3snl3(), w3snl4md::w3snl4(), w3snl5md::w3snl5(), w3snl1md::w3snlgqm(), w3snlsmd::w3snls(), w3srcemd::w3srce(), w3ref1md::w3sref(), w3str1md::w3str1(), w3strt(), w3servmd::w3thrtn(), w3updtmd::w3ucur(), w3updtmd::w3ulev(), w3updtmd::w3utau(), w3updtmd::w3uwnd(), w3servmd::w3xyrtn(), and w3sic3md::wn_precalc_cheng().

◆ tpiinv

◆ tstout

logical, parameter constants::tstout = .FALSE.

TSTOUT Flag for generation of test files.

Definition at line 56 of file constants.F90.

56  LOGICAL, PARAMETER :: TSTOUT = .false.

Referenced by w3initmd::w3init(), wminitmd::wminit(), and wminitmd::wminitnml().

◆ undef

w3snlsmd::abmax
real, parameter abmax
Definition: w3snlsmd.F90:101
constants::tpi
real, parameter tpi
TPI 2*Pi.
Definition: constants.F90:72
w3adatmd::iter
integer, dimension(:,:), pointer iter
Definition: w3adatmd.F90:654