StokesVWaveModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015 IH-Cantabria
9  Copyright (C) 2016-2017 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "StokesVWaveModel.H"
30 #include "mathematicalConstants.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace waveModels
40 {
41  defineTypeNameAndDebug(StokesV, 0);
43  (
44  waveModel,
45  StokesV,
46  patch
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53 
54 Foam::scalar Foam::waveModels::StokesV::A11
55 (
56  const scalar h,
57  const scalar k
58 ) const
59 {
60  const scalar s = sinh(k*h);
61  return 1.0/s;
62 }
63 
64 
65 Foam::scalar Foam::waveModels::StokesV::A13
66 (
67  const scalar h,
68  const scalar k
69 ) const
70 {
71  const scalar s = sinh(k*h);
72  const scalar c = cosh(k*h);
73  return -sqr(c)*(5*sqr(c) + 1)/(8*pow5(s));
74 }
75 
76 
77 Foam::scalar Foam::waveModels::StokesV::A15
78 (
79  const scalar h,
80  const scalar k
81 ) const
82 {
83  const scalar s = sinh(k*h);
84  const scalar c = cosh(k*h);
85 
86  return
87  -(
88  1184*pow(c, 10)
89  - 1440*pow(c, 8)
90  - 1992*pow6(c)
91  + 2641*pow4(c)
92  - 249*sqr(c) + 18
93  )
94  /(1536*pow(s, 11));
95 }
96 
97 
98 Foam::scalar Foam::waveModels::StokesV::A22
99 (
100  const scalar h,
101  const scalar k
102 ) const
103 {
104  const scalar s = sinh(k*h);
105  return 3/(8*pow4(s));
106 }
107 
108 
109 Foam::scalar Foam::waveModels::StokesV::A24
110 (
111  const scalar h,
112  const scalar k
113 ) const
114 {
115  const scalar s = sinh(k*h);
116  const scalar c = cosh(k*h);
117 
118  return
119  (192*pow(c, 8) - 424*pow(c, 6) - 312*pow4(c) + 480*sqr(c) - 17)
120  /(768*pow(s, 10));
121 }
122 
123 
124 Foam::scalar Foam::waveModels::StokesV::A33
125 (
126  const scalar h,
127  const scalar k
128 ) const
129 {
130  const scalar s = sinh(k*h);
131  const scalar c = cosh(k*h);
132 
133  return (13 - 4*sqr(c))/(64*pow(s, 7));
134 }
135 
136 
137 Foam::scalar Foam::waveModels::StokesV::A35
138 (
139  const scalar h,
140  const scalar k
141 ) const
142 {
143  const scalar s = sinh(k*h);
144  const scalar c = cosh(k*h);
145 
146  return
147  (
148  512*pow(c, 12)
149  + 4224*pow(c, 10)
150  - 6800*pow(c, 8)
151  - 12808*pow(c, 6)
152  + 16704.0*pow4(c)
153  - 3154*sqr(c)
154  + 107
155  )
156  /(4096*pow(s, 13)*(6*sqr(c) - 1));
157 }
158 
159 
160 Foam::scalar Foam::waveModels::StokesV::A44
161 (
162  const scalar h,
163  const scalar k
164 ) const
165 {
166  const scalar s = sinh(k*h);
167  const scalar c = cosh(k*h);
168 
169  return
170  (80*pow(c, 6) - 816*pow4(c) + 1338*sqr(c) - 197)
171  /(1536*pow(s, 10)*(6*sqr(c) - 1));
172 }
173 
174 
175 Foam::scalar Foam::waveModels::StokesV::A55
176 (
177  const scalar h,
178  const scalar k
179 ) const
180 {
181  const scalar s = sinh(k*h);
182  const scalar c = cosh(k*h);
183 
184  return
185  -(
186  2880*pow(c, 10)
187  - 72480*pow(c, 8)
188  + 324000*pow(c, 6)
189  - 432000*pow4(c)
190  + 163470*sqr(c)
191  - 16245
192  )
193  /(61440*pow(s, 11)*(6*sqr(c) - 1)*(8*pow4(c) - 11*sqr(c) + 3));
194 }
195 
196 
197 Foam::scalar Foam::waveModels::StokesV::B22
198 (
199  const scalar h,
200  const scalar k
201 ) const
202 {
203  const scalar s = sinh(k*h);
204  const scalar c = cosh(k*h);
205 
206  return (2*sqr(c) + 1)*c/(4*pow3(s));
207 }
208 
209 
210 Foam::scalar Foam::waveModels::StokesV::B24
211 (
212  const scalar h,
213  const scalar k
214 ) const
215 {
216  const scalar s = sinh(k*h);
217  const scalar c = cosh(k*h);
218 
219  return
220  (272*pow(c, 8) - 504*pow(c, 6) - 192*pow4(c) + 322*sqr(c) + 21)*c
221  /(384*pow(s, 9));
222 }
223 
224 
225 Foam::scalar Foam::waveModels::StokesV::B33
226 (
227  const scalar h,
228  const scalar k
229 ) const
230 {
231  const scalar s = sinh(k*h);
232  const scalar c = cosh(k*h);
233 
234  return (8*pow6(c) + 1)*3/(64*pow6(s));
235 }
236 
237 
238 Foam::scalar Foam::waveModels::StokesV::B33k
239 (
240  const scalar h,
241  const scalar k
242 ) const // d B33 / d k
243 {
244  const scalar s = sinh(k*h);
245  const scalar c = cosh(k*h);
246 
247  const scalar sk = h*s;
248  const scalar ck = h*c;
249 
250  return 9.*pow5(c)*ck/(4*pow6(s)) - (9*(8*pow6(c) + 1))/(32*pow(s, 7))*sk;
251 }
252 
253 
254 Foam::scalar Foam::waveModels::StokesV::B35
255 (
256  const scalar h,
257  const scalar k
258 ) const
259 {
260  const scalar s = sinh(k*h);
261  const scalar c = cosh(k*h);
262 
263  return
264  (
265  88128*pow(c, 14)
266  - 208224*pow(c, 12)
267  + 70848*pow(c, 10)
268  + 54000*pow(c, 8)
269  - 21816*pow6(c)
270  + 6264*pow4(c)
271  - 54*sqr(c)
272  - 81
273  )
274  /(12288*pow(s, 12)*(6*sqr(c) - 1));
275 }
276 
277 
278 Foam::scalar Foam::waveModels::StokesV::B35k
279 (
280  const scalar h,
281  const scalar k
282 ) const // d B35 / d k
283 {
284  const scalar s = sinh(k*h);
285  const scalar c = cosh(k*h);
286 
287  const scalar sk = h*s;
288  const scalar ck = h*c;
289 
290  return
291  (
292  14*88128*pow(c, 13)*ck
293  - 12*208224*pow(c, 11)*ck
294  + 10*70848*pow(c, 9)*ck
295  + 8*54000.0*pow(c, 7)*ck
296  - 6*21816*pow5(c)*ck
297  + 4*6264*pow3(c)*ck
298  - 2*54*c*ck
299  )
300  /(12288*pow(s, 12)*(6*sqr(c) - 1))
301  - (
302  88128*pow(c, 14)
303  - 208224*pow(c, 12)
304  + 70848*pow(c, 10)
305  + 54000*pow(c, 8)
306  - 21816*pow6(c)
307  + 6264*pow4(c)
308  - 54*sqr(c)
309  - 81
310  )*12
311  /(12288*pow(s, 13)*(6*sqr(c) - 1))*sk
312  - (
313  88128*pow(c,14)
314  - 208224*pow(c, 12)
315  + 70848*pow(c, 10)
316  + 54000*pow(c, 8)
317  - 21816*pow6(c)
318  + 6264*pow4(c)
319  - 54*sqr(c)
320  - 81
321  )*12*c*ck
322  /(12288*pow(s, 12)*sqr(6*sqr(c) - 1));
323 }
324 
325 
326 Foam::scalar Foam::waveModels::StokesV::B44
327 (
328  const scalar h,
329  const scalar k
330 ) const
331 {
332  const scalar s = sinh(k*h);
333  const scalar c = cosh(k*h);
334 
335  return
336  (
337  768*pow(c, 10)
338  - 448*pow(c, 8)
339  - 48*pow6(c)
340  + 48*pow4(c)
341  + 106*sqr(c)
342  - 21
343  )*c
344  /(384*pow(s, 9)*(6*sqr(c) - 1));
345 }
346 
347 
348 Foam::scalar Foam::waveModels::StokesV::B55
349 (
350  const scalar h,
351  const scalar k
352 ) const
353 {
354  const scalar s = sinh(k*h);
355  const scalar c = cosh(k*h);
356 
357  return
358  (
359  192000*pow(c, 16)
360  - 262720*pow(c, 14)
361  + 83680*pow(c, 12)
362  + 20160*pow(c, 10)
363  - 7280*pow(c, 8)
364  + 7160*pow(c, 6)
365  - 1800*pow(c, 4)
366  - 1050*sqr(c)
367  + 225
368  )
369  /(12288*pow(s, 10)*(6*sqr(c) - 1)*(8*pow4(c) - 11*sqr(c) + 3));
370 }
371 
372 
373 Foam::scalar Foam::waveModels::StokesV::B55k
374 (
375  const scalar h,
376  const scalar k
377 ) const // d B55 / d k
378 {
379  const scalar s = sinh(k*h);
380  const scalar c = cosh(k*h);
381 
382  const scalar sk = h*s;
383  const scalar ck = h*c;
384 
385  return
386  (
387  16*192000*pow(c, 15)*ck
388  - 14*262720*pow(c, 13)*ck
389  + 12*83680*pow(c, 11)*ck
390  + 10*20160*pow(c, 9)*ck
391  - 8*7280*pow(c, 7)*ck
392  + 6*7160*pow(c, 5)*ck
393  - 4*1800*pow(c, 3)*ck
394  - 2*1050*pow(c, 1)*ck
395  )
396  /(12288*pow(s, 10)*(6*sqr(c) - 1)*(8*pow(c, 4) - 11*sqr(c) + 3))
397  - (
398  192000*pow(c, 16)
399  - 262720*pow(c, 14)
400  + 83680*pow(c, 12)
401  + 20160*pow(c, 10)
402  - 7280*pow(c, 8)
403  + 7160*pow(c, 6)
404  - 1800*pow(c, 4)
405  - 1050*pow(c, 2)
406  + 225
407  )*10.0
408  /(12288*pow(s, 11)*(6*sqr(c) - 1)*(8*pow4(c) - 11*sqr(c) + 3))*sk
409  - (
410  192000*pow(c, 16)
411  - 262720*pow(c, 14)
412  + 83680*pow(c, 12)
413  + 20160*pow(c,10)
414  - 7280*pow(c, 8)
415  + 7160*pow(c, 6)
416  - 1800*pow(c,4)
417  - 1050*pow(c, 2)
418  + 225
419  )*12*c*ck
420  /(12288*pow(s, 10)*sqr(6*sqr(c) - 1)*(8*pow4(c) - 11*sqr(c) + 3))
421  - (
422  192000*pow(c, 16)
423  - 262720*pow(c, 14)
424  + 83680*pow(c, 12)
425  + 20160*pow(c, 10)
426  - 7280*pow(c, 8)
427  + 7160*pow(c, 6)
428  - 1800*pow(c, 4)
429  - 1050*pow(c, 2)
430  + 225
431  )*(32*pow3(c) - 22*c)*ck
432  /(12288*pow(s, 10)*(6*sqr(c) - 1)*sqr(8*pow4(c) - 11*sqr(c) + 3));
433 }
434 
435 
436 Foam::scalar Foam::waveModels::StokesV::C1
437 (
438  const scalar h,
439  const scalar k
440 ) const
441 {
442  const scalar s = sinh(k*h);
443  const scalar c = cosh(k*h);
444 
445  return (8*pow4(c) - 8*sqr(c) + 9)/(8*pow4(s));
446 }
447 
448 
449 Foam::scalar Foam::waveModels::StokesV::C1k
450 (
451  const scalar h,
452  const scalar k
453 ) const
454 {
455  const scalar s = sinh(k*h);
456  const scalar c = cosh(k*h);
457 
458  const scalar sk = h*s;
459  const scalar ck = h*c;
460 
461  return
462  (4*8*pow3(c)*ck - 2*8*c*ck)/(8*pow4(s))
463  - (8*pow4(c) - 8*sqr(c) + 9)*4*sk/(8*pow5(s));
464 }
465 
466 
467 Foam::scalar Foam::waveModels::StokesV::C2
468 (
469  const scalar h,
470  const scalar k
471 ) const
472 {
473  const scalar s = sinh(k*h);
474  const scalar c = cosh(k*h);
475 
476  return
477  (
478  3840*pow(c, 12)
479  - 4096*pow(c, 10)
480  + 2592*pow(c, 8)
481  - 1008*pow(c, 6)
482  + 5944*pow(c, 4)
483  - 1830*pow(c, 2)
484  + 147
485  ) // - 2592
486  /(512*pow(s, 10)*(6*sqr(c) - 1));
487 }
488 
489 
490 Foam::scalar Foam::waveModels::StokesV::C2k
491 (
492  const scalar h,
493  const scalar k
494 ) const
495 {
496  const scalar s = sinh(k*h);
497  const scalar c = cosh(k*h);
498 
499  const scalar sk = h*s;
500  const scalar ck = h*c;
501 
502  return
503  (
504  12*3840*pow(c, 11)*ck
505  - 10*4096*pow(c,9)*ck
506  + 8*2592*pow(c, 7)*ck
507  - 6*1008*pow(c, 5)*ck
508  + 4*5944*pow(c, 3)*ck
509  - 2*1830*c*ck
510  )
511  /(512*pow(s, 10)*(6*sqr(c) - 1))
512  - (
513  3840*pow(c, 12)
514  - 4096*pow(c, 10)
515  + 2592*pow(c, 8)
516  - 1008*pow(c, 6)
517  + 5944*pow(c, 4)
518  - 1830*pow(c, 2)
519  + 147
520  )*10*sk
521  /(512*pow(s, 11)*(6*sqr(c) - 1))
522  - (
523  3840*pow(c, 12)
524  - 4096*pow(c, 10)
525  + 2592*pow(c, 8)
526  - 1008*pow(c, 6)
527  + 5944*pow(c, 4)
528  - 1830*pow(c, 2)
529  + 147
530  )*12*c*ck
531  /(512*pow(s, 10)*sqr(6*sqr(c) - 1));
532 }
533 
534 
535 Foam::scalar Foam::waveModels::StokesV::C3
536 (
537  const scalar h,
538  const scalar k
539 ) const
540 {
541  const scalar s = sinh(k*h);
542  const scalar c = cosh(k*h);
543 
544  return -1/(4*s*c);
545 }
546 
547 
548 Foam::scalar Foam::waveModels::StokesV::C4
549 (
550  const scalar h,
551  const scalar k
552 ) const
553 {
554  const scalar s = sinh(k*h);
555  const scalar c = cosh(k*h);
556 
557  return
558  (12*pow(c, 8) + 36*pow(c, 6) - 162*pow(c, 4) + 141*sqr(c) - 27)
559  /(192*c*pow(s, 9));
560 }
561 
562 
563 void Foam::waveModels::StokesV::initialise
564 (
565  const scalar H,
566  const scalar d,
567  const scalar T,
568  scalar& kOut,
569  scalar& LambdaOut,
570  scalar& f1Out,
571  scalar& f2Out
572 ) const
573 {
574  scalar f1 = 1;
575  scalar f2 = 1;
576 
577  const scalar pi = mathematical::pi;
578  scalar k = 2.0*pi/(sqrt(mag(g_)*d)*T);
579  scalar lambda = H/2.0*k;
580 
581  label n = 0;
582 
583  static const scalar tolerance = 1e-12;
584  static const label iterMax = 10000;
585 
586  while ((mag(f1) > tolerance || mag(f2) > tolerance) && (n < iterMax))
587  {
588  const scalar b33 = B33(d, k);
589  const scalar b35 = B35(d, k);
590  const scalar b55 = B55(d, k);
591  const scalar c1 = C1(d, k);
592  const scalar c2 = C2(d, k);
593 
594  const scalar b33k = B33k(d, k);
595  const scalar b35k = B35k(d, k);
596  const scalar b55k = B55k(d, k);
597  const scalar c1k = C1k(d, k);
598  const scalar c2k = C2k(d, k);
599 
600  const scalar l2 = sqr(lambda);
601  const scalar l3 = l2*lambda;
602  const scalar l4 = l3*lambda;
603  const scalar l5 = l4*lambda;
604 
605  const scalar Bmat11 =
606  2*pi/(sqr(k)*d)*(lambda + l3*b33 + l5*(b35 + b55))
607  - 2*pi/(k*d)*(l3*b33k + l5*(b35k + b55k));
608 
609  const scalar Bmat12 =
610  - 2*pi/(k*d)*(1 + 3*l2*b33 + 5*l4*(b35 + b55));
611 
612  const scalar Bmat21 =
613  - d/(2*pi)*tanh(k*d)*(1 + l2*c1 + l4*c2)
614  - k*d/(2*pi)*(1 - sqr(tanh(k*d)))*d*(1 + l2*c1 + l4*c2)
615  - k*d/(2*pi)*tanh(k*d)*(l2*c1k + l4*c2k);
616 
617  const scalar Bmat22 = - k*d/(2.0*pi)*tanh(k*d)*(2*lambda*c1 + 4*l3*c2);
618 
619  f1 = pi*H/d - 2*pi/(k*d)*(lambda + l3*b33 + l5*(b35 + b55));
620 
621  f2 =
622  (
623  (2*pi*d)/(mag(g_)*sqr(T))
624  - k*d/(2*pi)*tanh(k*d)*(1 + l2*c1 + l4*c2)
625  );
626 
627  const scalar lambdaPr =
628  (f1*Bmat21 - f2*Bmat11)/(Bmat11*Bmat22 - Bmat12*Bmat21);
629  const scalar kPr =
630  (f2*Bmat12 - f1*Bmat22)/(Bmat11*Bmat22 - Bmat12*Bmat21);
631 
632  lambda += lambdaPr;
633  k += kPr;
634 
635  n++;
636  }
637 
638  kOut = k;
639  LambdaOut = lambda;
640 
641  f1Out = mag(f1);
642  f2Out = mag(f2);
643 }
644 
645 
646 Foam::scalar Foam::waveModels::StokesV::eta
647 (
648  const scalar h,
649  const scalar kx,
650  const scalar ky,
651  const scalar lambda,
652  const scalar T,
653  const scalar x,
654  const scalar y,
655  const scalar t,
656  const scalar phase
657 ) const
658 {
659  const scalar k = sqrt(kx*kx + ky*ky);
660 
661  const scalar b22 = B22(h, k);
662  const scalar b24 = B24(h, k);
663  const scalar b33 = B33(h, k);
664  const scalar b35 = B35(h, k);
665  const scalar b44 = B44(h, k);
666  const scalar b55 = B55(h, k);
667 
668  const scalar l2 = sqr(lambda);
669  const scalar l3 = l2*lambda;
670  const scalar l4 = l3*lambda;
671  const scalar l5 = l4*lambda;
672 
673  const scalar amp1 = lambda/k;
674  const scalar amp2 = (b22*l2 + b24*l4)/k;
675  const scalar amp3 = (b33*l3 + b35*l5)/k;
676  const scalar amp4 = b44*l4/k;
677  const scalar amp5 = b55*l5/k;
678 
679  const scalar theta = kx*x + ky*y - 2.0*mathematical::pi/T*t + phase;
680 
681  return
682  (
683  amp1*cos(theta)
684  + amp2*cos(2*theta)
685  + amp3*cos(3*theta)
686  + amp4*cos(4*theta)
687  + amp5*cos(5*theta)
688  );
689 }
690 
691 
692 Foam::vector Foam::waveModels::StokesV::Uf
693 (
694  const scalar d,
695  const scalar kx,
696  const scalar ky,
697  const scalar lambda,
698  const scalar T,
699  const scalar x,
700  const scalar y,
701  const scalar t,
702  const scalar phase,
703  const scalar z
704 ) const
705 {
706  const scalar k = sqrt(kx*kx + ky*ky);
707 
708  const scalar a11 = A11(d, k);
709  const scalar a13 = A13(d, k);
710  const scalar a15 = A15(d, k);
711  const scalar a22 = A22(d, k);
712  const scalar a24 = A24(d, k);
713  const scalar a33 = A33(d, k);
714  const scalar a35 = A35(d, k);
715  const scalar a44 = A44(d, k);
716  const scalar a55 = A55(d, k);
717 
718  const scalar pi = mathematical::pi;
719  const scalar l2 = sqr(lambda);
720  const scalar l3 = l2*lambda;
721  const scalar l4 = l3*lambda;
722  const scalar l5 = l4*lambda;
723 
724  const scalar a1u = 2*pi/T/k*(lambda*a11 + l3*a13 + l5*a15);
725  const scalar a2u = 2*2*pi/T/k*(l2*a22 + l4*a24);
726  const scalar a3u = 3*2*pi/T/k*(l3*a33 + l5*a35);
727  const scalar a4u = 4*2*pi/T/k*(l4*a44);
728  const scalar a5u = 5*2*pi/T/k*(l5*a55);
729 
730  const scalar theta = kx*x + ky*y - 2*pi/T*t + phase;
731 
732  scalar u =
733  a1u*cosh(k*z)*cos(theta)
734  + a2u*cosh(2*k*z)*cos(2*theta)
735  + a3u*cosh(3*k*z)*cos(3*theta)
736  + a4u*cosh(4*k*z)*cos(4*theta)
737  + a5u*cosh(5*k*z)*cos(5*theta);
738 
739  scalar w =
740  a1u*sinh(k*z)*sin(theta)
741  + a2u*sinh(2*k*z)*sin(2*theta)
742  + a3u*sinh(3*k*z)*sin(3*theta)
743  + a4u*sinh(4*k*z)*sin(4*theta)
744  + a5u*sinh(5*k*z)*sin(5*theta);
745 
746  scalar v = u*sin(waveAngle_);
747  u *= cos(waveAngle_);
748 
749  return vector(u, v, w);
750 }
751 
752 
753 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
754 
756 (
757  const scalar t,
758  const scalar tCoeff,
759  scalarField& level
760 ) const
761 {
762  const scalar waveK = mathematical::twoPi/waveLength_;
763  const scalar waveKx = waveK*cos(waveAngle_);
764  const scalar waveKy = waveK*sin(waveAngle_);
765 
766  forAll(level, paddlei)
767  {
768  const scalar eta =
769  this->eta
770  (
771  waterDepthRef_,
772  waveKx,
773  waveKy,
774  lambda_,
775  wavePeriod_,
776  xPaddle_[paddlei],
777  yPaddle_[paddlei],
778  t,
779  wavePhase_
780  );
781 
782  level[paddlei] = waterDepthRef_ + tCoeff*eta;
783  }
784 }
785 
786 
788 (
789  const scalar t,
790  const scalar tCoeff,
791  const scalarField& level
792 )
793 {
794  const scalar waveK = mathematical::twoPi/waveLength_;
795  const scalar waveKx = waveK*cos(waveAngle_);
796  const scalar waveKy = waveK*sin(waveAngle_);
797 
798  forAll(U_, facei)
799  {
800  // Fraction of geometry represented by paddle - to be set
801  scalar fraction = 1;
802 
803  // Height - to be set
804  scalar z = 0;
805 
806  setPaddlePropeties(level, facei, fraction, z);
807 
808  if (fraction > 0)
809  {
810  const label paddlei = faceToPaddle_[facei];
811 
812  const vector Uf = this->Uf
813  (
814  waterDepthRef_,
815  waveKx,
816  waveKy,
817  lambda_,
818  wavePeriod_,
819  xPaddle_[paddlei],
820  yPaddle_[paddlei],
821  t,
822  wavePhase_,
823  z
824  );
825 
826  U_[facei] = fraction*Uf*tCoeff;
827  }
828  }
829 }
830 
831 
832 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
833 
835 (
836  const dictionary& dict,
837  const fvMesh& mesh,
838  const polyPatch& patch,
839  const bool readFields
840 )
841 :
842  StokesI(dict, mesh, patch, false),
843  lambda_(0)
844 {
845  if (readFields)
846  {
847  readDict(dict);
848  }
849 }
850 
851 
852 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
853 
855 {
856  if (StokesI::readDict(overrideDict))
857  {
858  scalar f1;
859  scalar f2;
860  scalar waveK;
861 
862  initialise
863  (
864  waveHeight_,
865  waterDepthRef_,
866  wavePeriod_,
867  waveK,
868  lambda_,
869  f1,
870  f2
871  );
872 
873  if (f1 > 0.001 || f2 > 0.001)
874  {
876  << "No convergence for Stokes V wave theory" << nl
877  << " f1: " << f1 << nl
878  << " f2: " << f2 << nl
879  << exit(FatalError);
880  }
881 
882  return true;
883  }
884 
885  return false;
886 }
887 
888 
890 {
891  StokesI::info(os);
892 
893  os << " Lambda : " << lambda_ << nl
894  << " Wave type : " << waveType() << nl;
895 }
896 
897 
898 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::waveModels::StokesV::setVelocity
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
Definition: StokesVWaveModel.C:788
mathematicalConstants.H
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::cosh
dimensionedScalar cosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:271
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:689
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:122
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::Field< scalar >
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::waveModels::StokesV::setLevel
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
Definition: StokesVWaveModel.C:756
Foam::waveModels::StokesV::StokesV
StokesV(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition: StokesVWaveModel.C:835
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:111
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::constant::physicoChemical::c2
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::waveModels::StokesV::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: StokesVWaveModel.C:854
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
StokesVWaveModel.H
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::waveModels::StokesI
Stokes I wave model.
Definition: StokesIWaveModel.H:50
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::sinh
dimensionedScalar sinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:270