WAVEWATCH III  beta 0.0.1
serv_xnl4v5 Module Reference

Functions/Subroutines

subroutine y_gauleg (x1, x2, x, w, n)
 
subroutine z_cmpcg (sigma, depth, grav, cg)
 
subroutine z_intp1 (x1, y1, x2, y2, n1, n2, ierr)
 
subroutine z_polyarea (xpol, ypol, npol, area)
 
subroutine z_steps (x, dx, nx)
 
real function z_root2 (func, x1, x2, xacc, iprint, ierr)
 
subroutine z_upper (str)
 
real function z_wnumb (w, d, grav)
 

Function/Subroutine Documentation

◆ y_gauleg()

subroutine serv_xnl4v5::y_gauleg ( real, intent(in)  x1,
real, intent(in)  x2,
real, dimension(n), intent(out)  x,
real, dimension(n), intent(out)  w,
integer, intent(in)  n 
)

Definition at line 4 of file serv_xnl4v5.f90.

4 !-------------------------------------------------------------------
5 INTEGER, intent(in) :: n ! Number of intervals
6 real, intent(in) :: x1 ! lower limit of integration interval
7 real, intent(in) :: x2 ! upper limit of integration interval
8 real, intent(out) :: x(n) ! Positions for function evaluations
9 real, intent(out) :: w(n) ! Weights
10 !
11 !-----------------------------------------------------------------------
12 DOUBLE PRECISION EPS
13 parameter(eps=3.d-14)
14 INTEGER i,j,m
15 DOUBLE PRECISION p1,p2,p3,pp,xl,xm,z,z1
16 !-----------------------------------------------------------------------
17  m=(n+1)/2
18  xm=0.5d0*(x2+x1)
19  xl=0.5d0*(x2-x1)
20  do 12 i=1,m
21  z=cos(3.141592654d0*(i-.25d0)/(n+.5d0))
22 1 continue
23  p1=1.d0
24  p2=0.d0
25  do 11 j=1,n
26  p3=p2
27  p2=p1
28  p1=((2.d0*j-1.d0)*z*p2-(j-1.d0)*p3)/j
29 11 continue
30  pp=n*(z*p1-p2)/(z*z-1.d0)
31  z1=z
32 
33  z=z1-p1/pp
34  if(abs(z-z1).gt.eps)goto 1
35  x(i)=xm-xl*z
36  x(n+1-i)=xm+xl*z
37  w(i)=2.d0*xl/((1.d0-z*z)*pp*pp)
38  w(n+1-i)=w(i)
39 12 continue
40 !
41 return

Referenced by m_xnldata::q_modify().

◆ z_cmpcg()

subroutine serv_xnl4v5::z_cmpcg ( real, intent(in)  sigma,
real, intent(in)  depth,
real, intent(in)  grav,
real, intent(out)  cg 
)

Definition at line 46 of file serv_xnl4v5.f90.

46 !-----------------------------------------------------------------------------!
47 !
48 ! +-------+ ALKYON Hydraulic Consultancy & Research
49 ! | | Gerbrant van Vledder
50 ! | +---+
51 ! | | +---+
52 ! +---+ | |
53 ! +---+
54 !
55 implicit none
56 !
57 ! 0. Update history
58 !
59 ! 12/01/2001 Initial version
60 ! 11/04/2001 Check included for the cases tat sigma < 0 or depth <0
61 ! Result is cg = -10
62 !
63 ! 1. Purpose:
64 !
65 ! Compute group velocity for a given radian frequency and depth
66 !
67 ! 2. Method
68 !
69 ! Linear wave theory
70 !
71 ! 3. Parameter list:
72 !
73 !Type I/O Name Description
74 !------------------------------------------------------------------------------
75 real, intent(in) :: sigma ! radian frequency (rad)
76 real, intent(in) :: depth ! water depth (m)
77 real, intent(in) :: grav ! gravitational acceleration (m/s^2)
78 real, intent(out) :: cg ! group velocity (m/s)
79 !
80 real k ! wave number
81 !/A
82 !!real z_wnumb ! compute wave number
83 !/Z
84 !-----------------------------------------------------------------------------
85 k = z_wnumb(sigma,depth,grav)
86 !
87 if(depth <= 0. .or. sigma <= 0.) then
88  cg = -10.
89 else
90  if(depth*k > 30.) then
91  cg = grav/(2.*sigma)
92  else
93  cg = sigma/k*(0.5+depth*k/sinh(2.*depth*k))
94  end if
95 end if
96 !
97 return

References z_wnumb().

Referenced by m_xnldata::q_init(), and m_xnldata::q_xnl4v4().

◆ z_intp1()

subroutine serv_xnl4v5::z_intp1 ( real, dimension(n1), intent(in)  x1,
real, dimension(n1), intent(in)  y1,
real, dimension(n2), intent(in)  x2,
real, dimension(n2), intent(out)  y2,
integer, intent(in)  n1,
integer, intent(in)  n2,
integer, intent(out)  ierr 
)

Definition at line 102 of file serv_xnl4v5.f90.

102 !-----------------------------------------------------------------------------!
103 !
104 ! +-------+ ALKYON Hydraulic Consultancy & Research
105 ! | | Gerbrant van Vledder
106 ! | +---+
107 ! | | +---+
108 ! +---+ | |
109 ! +---+
110 !
111 implicit none
112 !
113 ! 0. Update history
114 !
115 ! 30/03/1999 Initical version
116 ! 9/04/1999 Check included for monotonicity of x-data
117 ! 11/10/1999 Error messages added and updated
118 ! 18/01/2001 Check include if n1==1
119 ! 24/01/2001 Check for equality of y2 data loosened if n2==1
120 ! 13/09/2001 Documentation updated
121 !
122 ! 1. Purpose
123 !
124 ! Interpolate function values
125 !
126 ! 2. Method
127 !
128 ! Linear interpolation
129 
130 ! If a requested point falls outside the input domain, then
131 ! the nearest point is used (viz. begin or end point of x1/y1 array
132 !
133 ! If the input array has only one point. A constant value is assumed
134 !
135 ! 3. Parameter list
136 !
137 ! Name I/O Type Description
138 !
139 integer, intent(in) :: n1 ! number of data points in x1-y1 arrays
140 integer, intent(in) :: n2 ! number of data points in x2-y2 arrays
141 real, intent(in) :: x1(n1) ! x-values of input data
142 real, intent(in) :: y1(n1) ! y-values of input data
143 real, intent(in) :: x2(n2) ! x-values of output data
144 real, intent(out) :: y2(n2) ! y-values of output data
145 integer, intent(out) :: ierr ! Error indicator
146 !
147 ! 4. Subroutines used
148 !
149 ! 5. Error messages
150 !
151 ! ierr = 0 No errors detected
152 ! = 1 x1-data not monotonic increasing
153 ! = 10 x2-data not monotonic increasing
154 ! = 11 x1- and x2 data not monotonic increasing
155 ! = 2 x1-data not monotonic decreasing
156 ! = 20 x1-data not monotonic decreasing
157 ! = 22 x1- and x2 data not monotonic decreasing
158 !
159 ! = 2 No variation in x1-data
160 ! = 3 No variation in x2-data is allowed if n2=1
161 !
162 ! 6. Remarks
163 !
164 ! It is assumed that the x1- and x2-data are either
165 ! monotonic increasing or decreasing
166 !
167 ! If a requested x2-value falls outside the range of x1-values
168 ! it is assumed that the corresponding y2-value is equal to
169 ! the nearest boundary value of the y1-values
170 !
171 ! Example: x1 = [0 1 2 3]
172 ! y1 = [1 2 1 0]
173 !
174 ! x2 = -1, y2 = 1
175 ! x2 = 5, y2 = 0
176 !
177 !------------------------------------------------------------------------------
178 integer i1,i2 ! counters
179 !
180 real ds ! step size
181 real fac ! factor in linear interpolation
182 real s1,s2 ! search values
183 real xmin1,xmax1 ! minimum and maximum of x1-data
184 real xmin2,xmax2 ! minimum and maximum of x2-data
185 !
186 real, parameter :: eps=1.e-20
187 !------------------------------------------------------------------------------
188 ! initialisation
189 !
190 ierr = 0
191 !
192 ! check number of points of input array
193 !
194 if(n1==1) then
195  y2 = y1(1)
196  goto 9999
197 end if
198 !
199 ! check minimum and maximum data values
200 !
201 xmin1 = minval(x1)
202 xmax1 = maxval(x1)
203 xmin2 = minval(x2)
204 xmax2 = maxval(x2)
205 !
206 if (abs(xmin1-xmax1) < eps .or. abs(x1(1)-x1(n1)) < eps) then
207  ierr = 2
208  goto 9999
209 end if
210 !
211 if ((abs(xmin2-xmax2) < eps .or. abs(x2(1)-x2(n2)) < eps) .and. n2 > 1) then
212  ierr = 3
213  goto 9999
214 end if
215 !
216 ! check input data for monotonicity
217 !
218 if(x1(1) < x1(n1)) then ! data increasing
219  do i1=1,n1-1
220  if(x1(i1) > x1(i1+1)) then
221  ierr=1
222  write(*,*) 'z_intp1: i1 x1(i1) x1(i1+1):',i1,x1(i1),x1(i1+1)
223  goto 9999
224  end if
225  end do
226 !
227  do i2=1,n2-1
228  if(x2(i2) > x2(i2+1)) then
229  ierr=ierr+10
230  write(*,*) 'z_intp1: i2 x2(i2) x2(i2+1):',i2,x2(i2),x2(i2+1)
231  goto 9999
232  end if
233  end do
234 !
235 else ! data decreasing
236  do i1=1,n1-1
237  if(x1(i1) < x1(i1+1)) then
238  ierr=2
239  write(*,*) 'z_intp1: i1 x1(i1) x1(i1+1):',i1,x1(i1),x1(i1+1)
240  goto 9999
241  end if
242  end do
243 !
244  do i2=1,n2-1
245  if(x2(i2) < x2(i2+1)) then
246  ierr=ierr + 20
247  write(*,*) 'z_intp1: i2 x2(i2) x2(i2+1):',i2,x2(i2),x2(i2+1)
248  goto 9999
249  end if
250  end do
251 end if
252 !
253 !------------------------------------------------------------------------------
254 ! initialize
255 !------------------------------------------------------------------------------
256 if(ierr==0) then
257  i1 = 1
258  s1 = x1(i1)
259 !
260  do i2 = 1,n2
261  s2 = x2(i2)
262  do while (s1 <= s2 .and. i1 < n1)
263  i1 = i1 + 1
264  s1 = x1(i1)
265  end do
266 !
267 ! special point
268 ! choose lowest s1-value if x2(:) < x1(1)
269 !
270  if(i1 ==1) then
271  y2(i2) = y1(i1)
272  else
273  ds = s2 - x1(i1-1)
274  fac = ds/(x1(i1)-x1(i1-1))
275  y2(i2) = y1(i1-1) + fac*(y1(i1)-y1(i1-1))
276  end if
277 !
278 ! special case at end: choose s2(n2) > s1(n1), choose last value of y1(1)
279 !
280  if(i2==n2 .and. s2>s1) y2(n2) = y1(n1)
281  end do
282 end if
283 !
284 9999 continue
285 !
286 return

Referenced by m_xnldata::q_modify().

◆ z_polyarea()

subroutine serv_xnl4v5::z_polyarea ( real, dimension(npol), intent(in)  xpol,
real, dimension(npol), intent(in)  ypol,
integer, intent(in)  npol,
real, intent(out)  area 
)

Definition at line 291 of file serv_xnl4v5.f90.

291 !-----------------------------------------------------------------------------!
292 !
293 ! +-------+ ALKYON Hydraulic Consultancy & Research
294 ! | | P.O. Box 248
295 ! | +---+ 8300 AE Emmeloord
296 ! | | +---+ Tel: +31 527 620909
297 ! +---+ | | Fax: +31 527 610020
298 ! +---+ http://www.alkyon.nl
299 !
300 ! Gerbrant van Vledder
301 !
302 ! 0. Update history
303 !
304 ! 0.01 12/06/2003 Initial version
305 !
306 ! 1. Purpose
307 !
308 ! Computes area of a closed polygon
309 !
310 ! 2. Method
311 !
312 ! The area of the polygon
313 !
314 ! 3. Parameter list
315 !
316 ! Name I/O Type Description
317 !
318 integer, intent(in) :: npol ! Number of points of polygon
319 real, intent(in) :: xpol(npol) ! x-coodinates of polygon
320 real, intent(in) :: ypol(npol) ! y-coordinates of polygon
321 real, intent(out) :: area ! area of polygon
322 !
323 ! 4. Subroutines used
324 !
325 ! 5. Error messages
326 !
327 ! 6. Remarks
328 !
329 integer ipol,ipol1 ! counters
330 real xmin,xmax,ymin,ymax ! minima and maxima of polygon
331 real xmean,ymean ! mean values
332 real xa,ya,xb,yb ! temporary variables
333 real sumx,sumy ! sums
334 real darea ! piece of area
335 !-------------------------------------------------------------------------------
336 if(npol<=1) then
337  crf = 0.
338  xz = 0.
339  yz = 0.
340  area = 0.
341  return
342 end if
343 !
344 ! compute minimum and maximum coordinates
345 !
346 xmin = minval(xpol)
347 xmax = maxval(xpol)
348 ymin = minval(ypol)
349 ymax = maxval(ypol)
350 !
351 ! compute mean of range of x- and y-coordinates
352 !
353 xmean = 0.5*(xmin + xmax)
354 ymean = 0.5*(ymin + ymax)
355 !
356 ! compute area and center of gravity
357 ! do loop over all line pieces of polygon
358 !
359 area = 0.
360 sumx = 0.
361 sumy = 0.
362 !
363 do ipol=1,npol
364  ipol1 = ipol + 1
365  if(ipol==npol) ipol1 = 1
366  xa = xpol(ipol)
367  ya = ypol(ipol)
368  xb = xpol(ipol1)
369  yb = ypol(ipol1)
370 !
371  darea = 0.5*((xa-xmean)*(yb-ymean) - (xb-xmean)*(ya-ymean))
372  area = area + darea
373  sumx = sumx + darea*(xa+xb+xmean)/3.
374  sumy = sumy + darea*(ya+yb+ymean)/3.
375 end do
376 !
377 return

Referenced by m_xnldata::q_cmplocus().

◆ z_root2()

real function serv_xnl4v5::z_root2 ( external  func,
real, intent(in)  x1,
real, intent(in)  x2,
real, intent(in)  xacc,
integer, intent(in)  iprint,
integer, intent(out)  ierr 
)

Definition at line 432 of file serv_xnl4v5.f90.

432 !-----------------------------------------------------------------------------!
433 !
434 ! +-------+ ALKYON Hydraulic Consultancy & Research
435 ! | |
436 ! | +---+
437 ! | | +---+
438 ! +---+ | |
439 ! +---+
440 !
441 ! 0. Update history
442 !
443 ! Version Date Modification
444 !
445 ! 0.01 29/11/1999 Initial version
446 ! 0.02 07/11/1999 Test added to check boundaries, and reverse if necessary
447 ! Bug fixed in assigning answer
448 ! 0.03 02/09/2002 Maximum number of iterations set to 20, instead of 10
449 !
450 ! 1. Purpose
451 !
452 ! Find zero crossing point of function FUNC between the
453 ! initial values on either side of zero crossing
454 !
455 ! 2. Method
456 !
457 ! Ridders method of root finding
458 !
459 ! adapted from routine zridddr
460 ! Numerical Recipes
461 ! The art if scientific computing, second edition, 1992
462 ! W.H. Press, S.A. Teukolsky, W.T. Vetterling and B.P. Flannery
463 !
464 ! 3. Parameter list
465 !
466 ! Name I/O Type Description
467 !
468 ! func i r real function
469 ! x1 i r initial x-value on left/right side of zero-crossing
470 ! x2 i r initial x-value on right/left side of zero-crossing
471 ! xacc i r accuracy, used as |x1(i)-x2(i)|< xacc
472 ! iprint i i Output channel and test level
473 ! ierr o i Error indicator
474 !
475 ! 4. Subroutines used
476 !
477 ! Func user supplied real function
478 !
479 ! 5. Error messages
480 !
481 ! ierr = 0 No errors occured during iteration process
482 ! 1 Iteration halted in dead end, this combination may NEVER occur
483 ! 2 Maximum number of iterations exceeded
484 ! 3 Solution jumped outside interval
485 !
486 ! 6. Remarks
487 !
488 ! It is assumed that the x1- and x2-coordinate lie
489 ! on different sides of the actual zero crossing
490 !
491 ! The input parameter IPRINT is used to generate test output.
492 ! If IPRINT==0, no test output is created
493 ! > 0, test output is directed to the file connected to unit LUPRINT=IPRINT
494 ! if no file is connected to this unit, no output is written
495 !
496 !
497 implicit none
498 !
499 real func ! external function
500 real, intent (in) :: x1 ! x-value at one side of interval
501 real, intent (in) :: x2 ! x-value at other side of interval
502 real, intent (in) :: xacc ! requested accuracy
503 integer, intent (in) :: iprint ! number of output channel, only used when
504 integer, intent (out) :: ierr ! error indicator
505 !
506 real unused ! default value
507 real zriddr ! intermediate function value
508 real xx1,xx2,xx ! local boundaries during iteration
509 integer maxit ! maximum number of iteration
510 integer luprint ! unit of test output
511 logical lopen ! check if a file is opened
512 
513 parameter(maxit = 20)
514 external func
515 !
516 integer iter ! counter for number of iterations
517 real fh ! function value FUNC(xh)
518 real fl ! function value FUNC(xl)
519 real fm ! function value FUNC(xm)
520 real fnew ! function value FUNC(xnew)
521 real s ! temp. function value, used for inverse quadratic interpolation
522 real xh ! upper (high) boundary of interval
523 real xl ! lower boundary of interval
524 real xm ! middle point of interval
525 real xnew ! new estimate according to Ridders method
526 !
527 ierr = 0 ! set error level
528 unused =-1.11e30 ! set start value
529 !
530 xx1 = x1 ! copy boundaries of interval to local variables
531 xx2 = x2
532 !
533 luprint = iprint
534 !
535 if(luprint > 0) then
536  inquire(unit=luprint,opened=lopen)
537  if(.not.lopen) then
538  luprint = 0
539  write(*,'(a,i4)') 'Z_ROOT2: invalid unit number:',iprint
540  end if
541 end if
542 !
543 ! check boundaries on requirement x2 > x1
544 !
545 if(xx1 > xx2) then
546  xx = xx1
547  xx1 = xx2
548  xx2 = xx
549 end if
550 !
551 fl = func(xx1)
552 fh = func(xx2)
553 !
554 !if(luprint > 0) write(luprint,'(a,4e13.5)') &
555 !& 'Z_ROOT2: xx1 xx2 fl fh:',xx1,xx2,fl,fh
556 !
557 if((fl > 0. .and. fh < 0.) .or. (fl < 0. .and. fh > 0.))then
558  xl = xx1
559  xh = xx2
560  zriddr = unused
561 !
562  do iter=1,maxit
563  xm = 0.5*(xl+xh)
564  fm = func(xm)
565  s = sqrt(fm**2-fl*fh)
566  if(s == 0.) goto 9000
567  xnew = xm+(xm-xl)*(sign(1.,fl-fh)*fm/s)
568 !
569 ! if(luprint>0) write(luprint,'(a,4e13.5)') &
570 !& 'Z_ROOT2: xm,fm,s,xnew:',xm,fm,s,xnew
571 !
572  if (abs(xnew-zriddr) <= xacc) then
573 ! if(luprint>0) write(luprint,'(a)') 'Z_ROOT2: xnew=zriddr'
574  goto 9000
575  end if
576 !
577  zriddr = xnew
578  fnew = func(zriddr)
579  if (fnew == 0.) goto 9000
580 !
581  if(sign(fm,fnew) /= fm) then
582  xl = xm
583  fl = fm
584  xh = zriddr
585  fh = fnew
586  elseif(sign(fl,fnew) /= fl) then
587  xh = zriddr
588  fh = fnew
589  elseif(sign(fh,fnew) /= fh) then
590  xl = zriddr
591  fl = fnew
592  else
593  ierr = 1
594  goto 9000
595  endif
596 !
597  if(abs(xh-xl) <= xacc) goto 9000
598 !
599  if(luprint > 0) write(luprint,'(a,i4,5e14.6)') &
600 & 'Z_ROOT2: iter,x1,x2,|x1-x2|,xacc,z:', iter,xl,xh,abs(xl-xh),xacc,fnew
601 !
602  end do
603  ierr = 2
604  if(luprint > 0) write(luprint,'(a)') 'Z_ROOT2: -> ierr=2'
605  goto 9000
606 else if (fl == 0.) then
607  zriddr = xx1
608 else if (fh == 0.) then
609  zriddr = xx2
610 else
611  ierr = 3
612  goto 9999
613 ! 'root must be bracketed in zriddr'
614 endif
615 !
616 9000 continue
617 !
618 z_root2 = zriddr
619 !
620 if(luprint > 0) write(luprint,'(a,2i3,5e13.5)') &
621 & 'Z_ROOT2: ierr,iter,xl,xh,acc,x0,z0:', ierr,iter,xl,xh,xacc,z_root2,func(z_root2)
622 !
623 9999 continue
624 !
625 return

Referenced by m_xnldata::q_locpos().

◆ z_steps()

subroutine serv_xnl4v5::z_steps ( real, dimension(nx), intent(in)  x,
real, dimension(nx), intent(out)  dx,
integer, intent(in)  nx 
)

Definition at line 382 of file serv_xnl4v5.f90.

382 !-----------------------------------------------------------------------------!
383 !
384 !
385 ! +-------+ ALKYON Hydraulic Consultancy & Research
386 ! | | Gerbrant van Vledder
387 ! | +---+
388 ! | | +---+ Creation date: September 28, 1998
389 ! +---+ | | Last Update: march 19, 2003
390 ! +---+
391 !
392 ! 0. Update history
393 !
394 ! 19/03/2003 Input argument defined using intent option
395 ! check included nx > 0
396 !
397 ! 1. Purpose
398 !
399 ! Compute bandwidth of spectral discretization
400 !
401 implicit none
402 !
403 integer, intent(in) :: nx ! Number of elements in array
404 real, intent(in) :: x(nx) ! Input data array with elements
405 real, intent(out) :: dx(nx) ! Output array with step sizes
406 !
407 integer ix ! counter
408 !------------------------------------------------------------------------------
409 if (nx<1) then
410  return
411 !
412 elseif (nx==1) then
413  dx = 0
414 else
415  do ix=2,nx-1
416  dx(ix) = 0.5 * (x(ix+1) - x(ix-1))
417  end do
418 !
419  if (nx >= 4) then
420  dx(1) = dx(2)*dx(2)/dx(3)
421  dx(nx) = dx(nx-1)*dx(nx-1)/dx(nx-2)
422  else
423  dx(1) = dx(2)
424  dx(nx) = dx(nx-1)
425  end if
426 end if
427 !
428 return

Referenced by m_xnldata::q_dscale(), and m_xnldata::q_init().

◆ z_upper()

subroutine serv_xnl4v5::z_upper ( character(len=*), intent(inout)  str)

Definition at line 630 of file serv_xnl4v5.f90.

630 !-----------------------------------------------------------------------------!
631 !
632 ! +-------+ ALKYON Hydraulic Consultancy & Research
633 ! | | Gerbrant van Vledder
634 ! | +---+
635 ! | | +---+ Creation date: July 3,1998
636 ! +---+ | | Last Update:
637 ! +---+
638 !
639 ! 0. Update history
640 !
641 ! 1. Purpose
642 !
643 ! Transform all lower capitals to UPPER CAPITALS in string STR
644 !
645 ! 2. Method
646 !
647 ! 3. Parameter list
648 !
649 ! Name I/O Type Description
650 !
651 implicit none
652 character(len=*), intent(inout) :: str ! Character string to be converted
653 !
654 ! 4. Subroutines used
655 !
656 ! 5. Error messages
657 !
658 ! 6. Remarks
659 !
660 integer nlen
661 integer i,ial,iau,izl
662 !
663 nlen = len(str)
664 !
665 ial = ichar('a')
666 iau = ichar('A')
667 izl = ichar('z')
668 !
669 do i=1,nlen
670  if(ichar(str(i:i)) >= ial.and. ichar(str(i:i)) <= izl) then
671  str(i:i) = char(ichar(str(i:i))-ial+iau)
672  end if
673 end do
674 !
675 return

Referenced by m_xnldata::q_setconfig().

◆ z_wnumb()

real function serv_xnl4v5::z_wnumb ( real, intent(in)  w,
real, intent(in)  d,
real, intent(in)  grav 
)

Definition at line 680 of file serv_xnl4v5.f90.

680 !-----------------------------------------------------------------------------!
681 !
682 ! +-------+ ALKYON Hydraulic Consultancy & Research
683 ! | | Gerbrant van Vledder
684 ! | +---+
685 ! | | +---+
686 ! +---+ | |
687 ! +---+
688 !
689 implicit none
690 !
691 ! 0. Update history
692 !
693 ! 01/04/1999 Initial version
694 ! 12/01/2001 grav added as input parameter
695 ! 11/04/2001 Check included for the case w < 0 or d < 0
696 !
697 ! 1. Purpose:
698 !
699 ! Compute wave number k for a given radian frequency and water depth
700 !
701 ! 2. Method
702 !
703 ! finite depth linear dispersion relation, using a Pade approximation
704 !
705 ! 3. Parameter list:
706 !
707 !Type I/O Name Description
708 !------------------------------------------------------------------------------
709 real, intent(in) :: w ! radian frequency (rad)
710 real, intent(in) :: d ! water depth (m)
711 real, intent(in) :: grav ! graviational acceleration (m/s^2)
712 !
713 ! 4. Subroutines used
714 !
715 ! 5. Error messages
716 !
717 ! 6. Remarks
718 !
719 ! The Pade approximation has been described in Hunt, 198.
720 !
721 !
722 real x,xx,y,omega
723 !
724 if(d<=0 .or. w<= 0.) then
725  z_wnumb = -10.
726 else
727  omega = w**2/grav
728  y = omega*d
729  xx = y*(y+1./(1.+y*(0.66667+y*(0.35550+y*(0.16084+y*(0.06320+y* &
730 & (0.02174+y*(0.00654+y*(0.00171+y*(0.00039+y*0.00011))))))))))
731  x = sqrt(xx)
732  z_wnumb = x/d
733 end if
734 !
735 return

Referenced by m_xnldata::q_dscale(), m_xnldata::q_init(), m_xnldata::q_polar2(), m_xnldata::x_cosk(), and z_cmpcg().

serv_xnl4v5::z_root2
real function z_root2(func, x1, x2, xacc, iprint, ierr)
Definition: serv_xnl4v5.f90:432
w3snl4md::ds
real, dimension(:), allocatable ds
Definition: w3snl4md.F90:173
serv_xnl4v5::z_wnumb
real function z_wnumb(w, d, grav)
Definition: serv_xnl4v5.f90:680
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61