-
Notifications
You must be signed in to change notification settings - Fork 54
/
WaveTools.py
3323 lines (2844 loc) · 110 KB
/
WaveTools.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""Tools for working with water waves.
The primary objective of this module is to provide solutions (exact and
approximate) for the free surface deformation and subsurface velocity
components of water waves. These can be used as boundary conditions, wave
generation sources, and validation solutions for numerical wave codes.
"""
from __future__ import print_function
from __future__ import absolute_import
from __future__ import division
from builtins import str
from builtins import zip
from builtins import range
#from builtins import object
from past.utils import old_div
import cython
import numpy as np
import cmath as cmat
from .Profiling import logEvent, logFile
from proteus import Comm
import time as tt
import sys as sys
__all__ = ['SteadyCurrent',
'SolitaryWave',
'MonochromaticWaves',
'NewWave',
'RandomWaves',
'MultiSpectraRandomWaves',
'DirectionalWaves',
'TimeSeries',
'RandomWavesFast',
'RandomNLWaves',
'RandomNLWavesFast',
'CombineWaves',
'fastcos_test',
'fastcosh_test',
'fastsinh_test',
'coshkzd_test',
'sinhkzd_test',
'loadExistingFunction',
'setVertDir',
'loadExistingFunction',
'setVertDir',
'setDirVector',
'dirCheck',
'reduceToIntervals',
'returnRectangles',
'returnRectangles3D',
'normIntegral',
'eta_mode',
'Udrift',
'vel_mode',
'sigma',
'JONSWAP',
'PM_mod',
'cos2s',
'mitsuyasu',
'dispersion',
'tophat',
'costap',
'decompose_tseries']
def fastcos_test(phase,sinus=False):
"""Fast cosine function with Taylor approximation - TO BE USED FOR TESTING"
Parameters
----------
phase : double
Phase
sinus : bool
Switch for cosine or sine calculation
Returns
--------
cos(phi) or sin(phi)
"""
if(sinus):
phase = old_div(np.pi,2.) - phase
return fastcos(phase,True)
def fastcosh_test(k,Z,fast=True):
"""Fast hyperbolic cosine function with Taylor approximation - TO BE USED FOR TESTING"
Parameters
----------
k : double
Wavenumber
Z : double
Z coordinate
Returns
--------
cosh(k*z)
"""
cython.declare(xx=cython.double[2])
fastcosh(xx,k,Z,fast)
return xx[0]
def fastsinh_test(k,Z,fast=True):
"""Fast hyperbolic sine function with Taylor approximation - TO BE USED FOR TESTING"
Parameters
----------
k : double
Wavenumber
Z : double
Z coordinate
Returns
--------
sinh(k*z)
"""
cython.declare(xx=cython.double[2])
fastcosh(xx,k,Z,fast)
return xx[1]
def coshkzd_test(k,Z,d, fast=True):
"""Calculation of u horizontal profile cosh(k(d+Z))/sinh(kd) using fast appoximaitons
and hyp trig relation cosh(a+b) = cosha*coshb+sinha*sinhb
Parameters
----------
----------
k : double
Wavenumber
Z : double
Z coordinate
d : double
depth
Returns
--------
cosh(k*(z+d))/sinh(kd) for Z>-d/2, 0 otherwise
"""
if (Z > old_div(-d,2.)):
return old_div(fastcosh_test(k,Z,fast), np.tanh(k*d)) + fastsinh_test(k,Z,fast)
else:
return 0.
def sinhkzd_test(k,Z,d,fast=True):
"""Calculation of v vertical profile cosh(k(d+Z))/sinh(kd) using fast appoximaitons
and hyp trig relation sinh(a+b) = sinha*coshb+cosha*sinhb
Parameters
----------
----------
k : double
Wavenumber
Z : double
Z coordinate
d : double
depth
Returns
--------
sinh(k*(z+d))/sinh(kd) for Z>-d/2, 0 otherwise
"""
if (Z> old_div(-d,2.)):
return fastcosh_test(k,Z,fast) + old_div(fastsinh_test(k,Z,fast), np.tanh(k*d))
else:
return 0.
def loadExistingFunction(funcName, validFunctions):
"""Checks if a function name is known function and returns it
Checks if a function name is present in a list of functions.
If True, the function is returned. If False, raises SystemExit.
Parameters
----------
funcName : string
Function name
validFunctions : List[function]
List of valid functions (list of objects)
Returns
--------
function
Raises
---------
SystemExit
"""
funcNames = []
for func in validFunctions:
funcNames.append(func.__name__)
if func.__name__ == funcName:
func_ret = func
if funcName not in funcNames:
logEvent("ERROR! Wavetools.py: Wrong function type (%s) given: Valid wavetypes are %s" %(funcName,funcNames), level=0)
sys.exit(1)
return func_ret
def setVertDir(g):
""" Returns the unit vector for the vertical direction
The vertical direction is opposite to the gravity direction
Parameters
----------
g : numpy.ndarray
Gravitational acceleration vector (must have 3 components)
Returns
--------
numpy.ndarray
"""
return -np.array(old_div(g,(sqrt(g[0]**2 + g[1]**2 + g[2]**2))))
def setDirVector(vector):
""" Returns the direction of a vector in the form of a unit vector
Parameters
----------
vector : numpy.ndarray
1D numpy array with three components
Returns
--------
numpy.ndarray
"""
return old_div(vector,(sqrt(vector[0]**2 + vector[1]**2 + vector[2]**2)))
def dirCheck(v1, v2):
""" Checks if two vectors are vertical raises SystemError if True
Parameters
----------
v1 : numpy.ndarray
1st vector with three components
v2 : numpy.ndarray
2nd vector with three components
Returns
--------
None
Raises
---------
SystemExit
"""
dircheck = abs(v1[0]*v2[0]+v1[1]*v2[1]+v1[2]*v2[2])
#print self.dircheck
if dircheck > 1e-10:
logEvent("Wave direction is not perpendicular to gravity vector. Check input",level=0)
return sys.exit(1)
else:
return None
def reduceToIntervals(fi,df):
""" Prepares the x-axis array with N elements for numerical integration
Integration is performed along he x- axis. If
fi = [a1, a2, a3,...,a_N-1 a_N] then it returns the array
[a1, 0.5(a1+a2), 0.5(a2+a3),...0.5(a_N-1+a_N), a_N].
Input array must have constant step
Parameters
----------
fi : numpy.ndarray
x- array with N elements
df : float
Constant step of array
Returns
--------
numpy.ndarray
"""
fim_tmp = (0.5*(fi[1:]+fi[:-1])).tolist()
return np.array([fim_tmp[0]-0.5*df]+fim_tmp+[fim_tmp[-1]+0.5*df])
def returnRectangles(a,x):
""" Returns 2D discrete integral array using the rectangle method
The calculation for each array element is
:math:`(\Delta y_i = 0.5(a_{n-1}+a_{n})*(x_{n-1}-x_{n})`
Parameters
----------
a : numpy.ndarray
Description: Array of y(x) function with N+1 elements
x : numpy.ndarray
Description: x- coordinate array with N elements
Returns
--------
numpy.ndarray
"""
return 0.5*(a[1:]+a[:-1])*(x[1:]-x[:-1])
def returnRectangles3D(a,x,y):
""" Returns 3D discrete integrals using the rectangle method
The calculation for each array element is
:math: `(\Delta y = 0.25*(a_(n-1,m-1)+a_(n,m-1)+a_(n-1,m) ...
...+a_(n,m))*(x_n-1-x_n) *(z_m-1-z_m))`
Parameters
----------
a : numpy.ndarray
2D Array of y(x,y) function with (N+1)x(M+1)elements
x : numpy.ndarray
Description: x- coordinate array with N+1 elements
y : numpy.ndarray
Description: x- coordinate array with N+1 elements
Returns
--------
numpy.ndarray
"""
ai = 0.5*(a[1:,:]+a[:-1,:])
ai = 0.5*(ai[:,1:]+ai[:,:-1])
for ii in range(len(x)-1):
ai[ii,:] *= (y[1:]-y[:-1])
for jj in range(len(y) - 1):
ai[:,jj] *= (x[1:]-x[:-1])
return ai
def normIntegral(f,dom):
"""Returns a normalised 2D function
The calculation is :math: `(\int_\Omega f d\Omega =1)`
Parameters
----------
f : numpy.ndarray
Discrete 2D function
Numpy array or list
dom : float
Discrete function step
Returns
--------
numpy.ndarray
"""
G0 = old_div(1.,sum(returnRectangles(f,dom)))
return G0*f
def eta_mode(x, t, kDir, omega, phi, amplitude):
"""Calculates the free surface elevation for a single frequency mode
Parameters
----------
x : numpy.ndarray
Position vector
t : float
Time variable
kDir : numpy.ndarray
Wave number vector
omega : float
Angular frequency
phi : float
Description: Wave phase
amp : float
Description: Wave amplitude
Returns
--------
float
The free surface elevation at x,t
"""
phase = x[0]*kDir[0]+x[1]*kDir[1]+x[2]*kDir[2] - omega*t + phi
return amplitude*cos(phase)
def Udrift(amp,gAbs,c,d):
"""Calculates the 2nd order Stokes drift for a linear mode
Parameters
----------
amp : float
Description: Wave amplitude
gAbs : float
Magnitude of gravitational acceleration
c : float
Wave celerity
d : float
Water depth
Returns
--------
float
Magnitude of the mean velocity drift
"""
return 0.5*gAbs*amp*amp/c/d
def vel_mode(x, t, kDir, kAbs, omega, phi, amplitude, mwl, depth, vDir, gAbs):
"""Calculates the wave velocity components for a single frequency mode
Parameters
----------
x : numpy.ndarray
Position vector
t : float
Time variable
kDir : numpy.ndarray
Wave number vector
kAbs : floatkAbs
Wave number magnitude
omega : float
Angular frequency
phi : float
Description: Wave phase
amplidute : float
Description: Wave amplitude
mwl : float
Mean water level
depth : float
Water depth
vDir : numpy.ndarray
Unit vector aligned with vertical direction
Returns
--------
numpy.ndarray
1D Numpy array of the velocity vector at x,t
"""
phase = x[0]*kDir[0]+x[1]*kDir[1]+x[2]*kDir[2] - omega*t + phi
Z = (vDir[0]*x[0] + vDir[1]*x[1]+ vDir[2]*x[2]) - mwl
UH = 0.
UV=0.
ii=0.
UH=amplitude*omega*cosh(kAbs*(Z + depth))*cos( phase )/sinh(kAbs*depth)
UV=amplitude*omega*sinh(kAbs*(Z + depth))*sin( phase )/sinh(kAbs*depth)
waveDir = old_div(kDir,kAbs)
UH = UH - Udrift(amplitude,gAbs,old_div(omega,kAbs),depth)
#waves(period = 1./self.fi[ii], waveHeight = 2.*self.ai[ii],mwl = self.mwl, depth = self.d,g = self.g,waveDir = self.waveDir,wavelength=self.wi[ii], phi0 = self.phi[ii]).u(x,y,z,t)
V = np.array([UH*waveDir[0]+UV*vDir[0],
UH*waveDir[1]+UV*vDir[1],
UH*waveDir[2]+UV*vDir[2]])
return V
def sigma(omega,omega0):
"""Calculates sigma function for JONSWAP spectrum
See http://www.wikiwaves.org/Ocean-Wave_Sectra
Parameters
----------
omega : numpy.ndarray
Angular frequency array
omega0 : numpy.ndarray
Peak angular frequency
Returns
--------
numpy.ndarray
1D Numpy array of simga function with respect to f
"""
sigmaReturn = np.where(omega > omega0,0.09,0.07)
return sigmaReturn
def JONSWAP(f,f0,Hs,gamma=3.3,TMA=False, depth = None):
"""Calculates the JONSWAP frequency spectrum (Goda 2009)
The calculation includes the TMA modification, if TMA =True
Parameters
----------
f : numpy.ndarray
Frequency array
f0 : float
Peak frequency
Hs : float
Significant wave height
gamma : Optional[float]
Peak enhancement factor
TMA : bool
Description: TMA switch
depth : Optional[float]
Water depth
Returns
--------
numpy.ndarray
1D Numpy array of the spectrum in frequency domain
"""
Tp = old_div(1.,f0)
bj = 0.0624*(1.094-0.01915*log(gamma))/(0.23+0.0336*gamma-old_div(0.185,(1.9+gamma)))
r = np.exp(old_div(-(Tp*f-1.)**2,(2.*sigma(f,f0)**2)))
tma = 1.
if TMA:
if (depth is None):
logEvent("Wavetools:py. Provide valid depth definition definition for TMA spectrum")
logEvent("Wavetools:py. Stopping simulation")
sys.exit(1)
k = dispersion(2*M_PI*f,depth)
tma = np.tanh(k*depth)*np.tanh(k*depth)/(1.+ 2.*k*depth/np.sinh(2.*k*depth))
return tma * bj*(Hs**2)*(old_div(1.,((Tp**4) *(f**5))))*np.exp(-1.25*(old_div(1.,(Tp*f)**(4.))))*(gamma**r)
def PM_mod(f,f0,Hs):
"""Calculates the Pierson-Moskovitz spectrum (or Bretschneider or ISSC)
Reference:
http://www.orcina.com/SoftwareProducts/OrcaFlex/Documentation/Help/Content/html/Waves,WaveSpectra.htm
And then to Tucker M J, 1991. Waves in Ocean Engineering. Ellis Horwood Ltd. (Chichester).
Parameters
--------
f : numpy.ndarray
Frequency array
f0 : float
Peak frequency
Hs : float
Significant wave height
Returns
--------
numpy.ndarray
1D Numpy array of the spectrum in frequency domain
"""
return (old_div(5.0,16.0))*Hs**2*(old_div(f0**4,f**5))*np.exp((old_div(-5.0,4.0))*(old_div(f0,f))**4)
def cos2s(theta,f,s=10):
"""Calculates the cos-2s directional spreading function
see USACE - CETN-I-28 http://chl.erdc.usace.army.mil/library/publications/chetn/pdf/cetn-i-28.pdf
Parameters
----------
theta : numpy.ndarray
Wave angle array
f : numpy.ndarray
Frequency array
s : Optional[float]
Spreading parameter
Returns
--------
numpy.ndarray
2D Numpy array of cos-2s spectrum
"""
fun = np.zeros((len(theta),len(f)),)
for ii in range(len(fun[0,:])):
fun[:,ii] = np.cos(old_div(theta,2))**(2*s)
return fun
def mitsuyasu(theta,fi,f0,smax=10):
"""The cos2s wave directional spread with wave frequency dependency
Equation from "Random Seas and Design of Maritime Structures" - Y. Goda - 2010
(3rd ed) eq. 2.22 - 2.25
Parameters
----------
theta : numpy.ndarray
Wave angle array
fi : numpy.ndarray
Frequency array
f0 : float
Peak frequency
smax : Optional[float]
Spreading parameter
Returns
--------
numpy.ndarray
2D Numpy array of Mitsuyashu-type spectrum
"""
s = smax * (old_div(fi,f0))**(5)
ii = np.where(fi>f0)[0][0]
s[ii:] = smax * (old_div(fi[ii:],f0))**(-2.5)
fun = np.zeros((len(theta),len(fi)),)
for ii in range(len(fun[0,:])):
fun[:,ii] = np.cos(old_div(theta,2))**(2.*s[ii])
return fun
def dispersion(w,d, g = 9.81,niter = 1000):
"""Calculates the wave number for single or multiple frequencies using linear dispersion relation.
Parameters
----------
w : float or np.ndarray
Angular frequency
d : float
Water depth
g : Optional[float]
Gravitational acceleration
niter : Optional[int]
Solution iterations
Returns
--------
float or numpy.ndarray
Wavenumber as a float or 1D array for multiple frequencies
"""
w_aux = np.array(w)
K = old_div(w_aux**2,g)
for jj in range(niter):
K = old_div(w_aux**2,(g*np.tanh(K*d)))
if type(K) is float:
return K[0]
else:
return K
def tophat(l,cutoff):
""" Calculates and returns a top hat filter array
Parameters
----------
l : int
Length of array
cutoff : float
Cut off fraction at both the leading and tailing part of the array
Returns
--------
numpy.ndarray
"""
a = np.zeros(l,)
cut = int(cutoff*l)
a[cut:-cut] = 1.
return a
def costap(l,cutoff=0.1):
""" Calculates and returns a top hat filter array
Parameters
----------
l : int
Length of array
cutoff : float
Cut off fraction at both the leading and tailing part of the array
Returns
--------
numpy.ndarray
"""
npoints = int(cutoff*l)
wind = np.ones(l)
for k in range(l): # (k,np) = (n,N) normally used
if k < npoints:
wind[k] = 0.5*(1.-cos(M_PI*float(k)/float(npoints)))
if k > l - npoints -1:
wind[k] = 0.5*(1.-cos(M_PI*float(l-k-1)/float(npoints)))
return wind
def decompose_tseries(time,eta,dt):
""" Performs spectral analysis and calculates angular frequency components, amplitude, phase and mean level power
of a time series with constant sampling.
Parameters
----------
time : numpy.ndarray
Time array
eta :numpy.ndarray
Signal array
dt : float
Sampling interval
Returns
--------
List
A list with results with four entries:
0 -> numpy array with angular frequency components
1 -> numpy array with amplitude of each component aa
2 -> numpy array with phase of each component pp
3 -> float of the 0th fourier mode (wave setup)
"""
nfft = len(time)
results = []
fft_x = np.fft.fft(eta,nfft)
freq = np.fft.fftfreq(nfft,dt) #%complex spectrum
iend = np.where(freq<0)[0][0]
setup = old_div(np.real(fft_x[0]),nfft)
fft_x = fft_x[1:iend]
freq = freq[1:iend]
#%retaining only first half of the spectrum
aa = 2.*abs(fft_x)/nfft #%amplitudes (only the ones related to positive frequencies)
ww = 2*M_PI*freq
pp = np.zeros(len(aa),complex)
for k in range(len(aa)):
pp[k]=cmat.phase(fft_x[k]) #% Calculating phases phases
pp = np.real(pp) # Append results to list
results.append(ww)
results.append(aa)
results.append(pp)
results.append(setup)
return results
class SteadyCurrent(object):
"""
This class is used for generating a steady current
Parameters
----------
U: numpy.ndarray
Current velocity in vector form
mwl : float
Still water level
rampTime : float
Ramp time for current
"""
def __init__(self,
U,
mwl,
rampTime = 0.):
self.mwl = mwl
try:
if len(U)!=3:
logEvent("ERROR! Wavetools.py: SteadyCurrent velocity argument needs to be a vector with length = 3")
sys.exit(1)
except:
logEvent("ERROR! Wavetools.py: SteadyCurrent velocity argument needs to be a vector with length = 3")
sys.exit(1)
self.U = np.array(U)
self.ramp = rampTime
def eta(self,x,t):
"""Calculates free surface elevation (SolitaryWave class)
Parameters
----------
x : numpy.ndarray
Position vector
t : float
Time variable
Returns
--------
float
Free-surface elevation as a float
"""
return 0.
def u(self,x,t):
"""Calculates wave velocity vector (SolitaryWave class).
Parameters
----------
x : numpy.ndarray
Position vector
t : float
Time variable
Returns
--------
numpy.ndarray
Velocity vector as 1D array
"""
U = np.zeros(3,)
if(t<self.ramp):
U=self.U*t/self.ramp
else:
U=self.U
return U
class SolitaryWave(object):
"""
This class is used for generating 1st order solitary wave
Parameters
----------
waveHeight: float
Regular wave height
mwl : float
Still water level
depth : float
Water depth
g : numpy.ndarray
Gravitational acceleration vector
waveDir : numpy.ndarray
Wave direction in vector form
trans : numpy.ndarray
Position vector of the peak
fast : bool
Switch for optimised functions
"""
def __init__(self,
waveHeight,
mwl,
depth,
g,
waveDir,
trans = np.zeros(3,"d"),
fast = True):
self.H = waveHeight
self.fast = fast
self.g = np.array(g)
self.waveDir = setDirVector(np.array(waveDir))
self.vDir = setVertDir(g)
self.gAbs = sqrt(self.g[0]*self.g[0]+self.g[1]*self.g[1]+self.g[2]*self.g[2])
self.trans = trans
self.c = np.sqrt(self.gAbs * (depth+self.H))
self.mwl = mwl
self.depth = depth
self.K = old_div(np.sqrt(3. *self.H/ (4. * self.depth)),self.depth)
self.d2 = depth*depth
self.d3 = self.d2 * depth
#Checking if g and waveDir are perpendicular
dirCheck(self.waveDir,self.vDir)
def eta(self,x,t):
"""Calculates free surface elevation (SolitaryWave class)
Parameters
----------
x : numpy.ndarray
Position vector
t : float
Time variable
Returns
--------
float
Free-surface elevation as a float
"""
phase = sum( (x[:]-self.trans[:])*self.waveDir[:]) - self.c * t
a1 = self.K*phase
return self.H*1.0/ cosh(a1)**2
def u(self,x,t):
"""Calculates wave velocity vector (SolitaryWave class).
Parameters
----------
x : numpy.ndarray
Position vector
t : float
Time variable
Returns
--------
numpy.ndarray
Velocity vector as 1D array
"""
phase = sum( (x[:]-self.trans[:])*self.waveDir[:]) - self.c * t
a1 = cosh(self.K*phase*2.)
a2 = cosh(self.K*phase)
Z = (self.vDir[0]*x[0] + self.vDir[1]*x[1]+ self.vDir[2]*x[2]) - self.mwl
Uhorz = 1.0 /(4.0 * self.depth**4 ) * np.sqrt(self.gAbs * self.depth) * self.H *(
2.0 * self.d3 + self.d2 * self.H + 12.0 * self.depth * self.H * Z + 6.0 * self.H * Z**2.0 +
(2.0 * self.d3 - self.d2 * self.H - 6.0 * self.depth * self.H * Z - 3.0 * self.H * Z**2 ) * a1)/(a2)**4
Uvert = 1.0 / ( 4.0 * np.sqrt(self.gAbs* self.depth) ) * np.sqrt(3.0) * self.gAbs * (old_div(self.H, self.depth**3.0))** 1.5 * (self.depth + Z)*(
2.0 * self.depth**3 - 7.0 * self.depth**2.0 * self.H + 10.0 * self.depth * self.H * Z + 5.0 * self.H * Z**2.0 +
(2.0 * self.depth**3.0 + self.depth**2.0 * self.H - 2.0 * self.depth * self.H * Z - self.H * Z**2.0)*
cosh(np.sqrt( 3.0 * self.H / self.depth**3.0) * phase ))/(
cosh(np.sqrt( 3.0 * self.H / ( 4.0 * self.depth**3.0))*
phase ) )** 4.0*( tanh( np.sqrt( 3.0 * self.H / ( 4.0 * self.depth**3.0))*phase ))
"""
phase = sum( (x[:]-self.trans[:])*self.waveDir[:]) - self.c * t
a1 = cosh(self.K * phase)
a2 = tanh( self.K * phase)
Z = (self.vDir[0]*x[0] + self.vDir[1]*x[1]+ self.vDir[2]*x[2]) - self.mwl
Uhorz = np.sqrt( self.gAbs * self.depth) * ( self.H / self.depth) * ( 1 / ( a1**2)) * ( 1 - ( self.H / ( 4 * self.depth)) * ( 1 / ( a1**2)))
Uvert = -np.sqrt( self.gAbs * self.depth) * ( Z / self.depth) * ( 1 - ( self.H / ( 4 * self.depth)) * ( 1 / ( a1**2))) * ( ( 2 * self.H / self.depth) * self.K * ( a2 / ( a1**2)))
"""
return self.waveDir*Uhorz + self.vDir*Uvert
class MonochromaticWaves(object):
"""
This class is used for generating regular waves in both linear and nonlinear regimes. See Dean and Dalrymple 1994 for equations.
Parameters
----------
period : float
Regular wave period
waveHeight: float
Regular wave height
mwl : float
Still water level
depth : float
Water depth
g : numpy.ndarray
Gravitational acceleration vector
waveDir : numpy.ndarray
Wave direction in vector form
wavelength : float
Regular wave length, calculated from linear dispersion if set to None
waveType : string
Defines regular wave theory ("Linear", "Fenton")
Fenton: uses BCoeffs/YCoeffs provided by user
autoFenton: bool
autoFenton=True: uses waveheight, period, depth, and g to
calculate coeffs
autoFenton=False: uses BCoeffs/YCoeffs provided by user
autoFentonOpts: dict
options for autoFenton. The dictionary must contain the following
entries (here the default values if autoFentonOpts is None):
autoFentonOpts = {'mode': 'Period',
'current_criterion': 1,
'height_steps': 1,
'niter': 40,
'conv_crit': 1e-05,
'points_freesurface': 50,
'points_velocity': 16,
'points_vertical': 20}
Ycoeff : numpy.ndarray
Fenton Fourier coefficients for free-surface elevation
Bcoeff : numpy.ndarray
Fenton Fourier coefficients for velocity (set to None for linear wave theory)
Nf : integer
Fenton Fourier components for reconstruction (set to 1000, needs to be equal to the size of Bcoeff and Ycoeff)
meanVelocity : numpy.ndarray
Mean velocity for Fenton Fourier approximation
phi0 : float
Regular wave phase (0 by default)
fast : bool
Switch for optimised functions
"""
def __init__(self,
period,
waveHeight,
mwl,
depth,
g,
waveDir,
wavelength=None,
waveType="Linear",
autoFenton=True,
autoFentonOpts=None,
Ycoeff = np.zeros(1000,),
Bcoeff =np.zeros(1000,),
Nf = 1000,
meanVelocity = np.array([0.,0,0.]),
phi0 = 0.,
fast = True):
self.fast = fast
knownWaveTypes = ["Linear","Fenton"]
self.waveType = waveType
if waveType not in knownWaveTypes:
logEvent("ERROR!!: Wrong wavetype given: Valid wavetypes are %s" %(knownWaveTypes), level=0)
sys.exit(1)
self.g = np.array(g)
self.waveDir = setDirVector(np.array(waveDir))
self.vDir = setVertDir(g)
self.gAbs = sqrt(self.g[0]*self.g[0]+self.g[1]*self.g[1]+self.g[2]*self.g[2])
#Checking if g and waveDir are perpendicular
dirCheck(self.waveDir,self.vDir)
self.phi0=phi0
self.mwl = mwl
self.depth = depth
self.omega = 2.0*M_PI/period
self.Nf = Nf
self.Ycoeff = Ycoeff
self.Bcoeff = Bcoeff
self.tanhF = np.zeros(Nf,"d")
#Calculating / checking wavelength data
if waveType== "Linear":
self.k = dispersion(w=self.omega,d=self.depth,g=self.gAbs)
self.wavelength = 2.0*M_PI/self.k
elif waveType == "Fenton":
if autoFenton is False:
try:
self.k = 2.0*M_PI/wavelength
self.wavelength=wavelength
except:
logEvent("ERROR! Wavetools.py: Wavelenght is not defined for nonlinear waves. Enter wavelength in class arguments",level=0)
sys.exit(1)
if ( (len(self.Ycoeff)!=self.Nf) or (len(self.Bcoeff)!=self.Nf) or (Ycoeff[0]==0.) or (Bcoeff[0]==0.) ):
logEvent("ERROR! Wavetools.py: Ycoeff and Bcoeff must have the same length and equal to Nf and the 1st order harmonic must not be zero",level=0)
sys.exit(1)
else:
for ii in range(len(self.tanhF)):
kk = (ii+1)*self.k
self.tanhF[ii] = float(np.tanh(kk*self.depth) )
elif autoFenton is True:
from proteus.fenton import Fenton