waveGenerationModel.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 "waveGenerationModel.H"
30 #include "unitConversion.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace waveModels
37 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  const scalar h(get<scalar>("waveHeight"));
48  if (h < 0)
49  {
51  << "Wave height must be greater than zero. Supplied"
52  << " value waveHeight = " << h
53  << exit(FatalIOError);
54  }
55 
56  return h;
57 }
58 
59 
61 {
62  return degToRad(get<scalar>("waveAngle"));
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
69 (
70  const dictionary& dict,
71  const fvMesh& mesh,
72  const polyPatch& patch,
73  const bool readFields
74 )
75 :
76  waveModel(dict, mesh, patch, false)
77 {
78  if (readFields)
79  {
80  readDict(dict);
81  }
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
88 (
89  const dictionary& overrideDict
90 )
91 {
92  if (waveModel::readDict(overrideDict))
93  {
94  readEntry("activeAbsorption", activeAbsorption_);
95 
96  return true;
97  }
98 
99  return false;
100 }
101 
102 
104 {
106 }
107 
108 
109 // ************************************************************************* //
Foam::waveModel
Base class for waveModels.
Definition: waveModel.H:57
Foam::waveModels::waveGenerationModel::readWaveHeight
scalar readWaveHeight() const
Helper function to read the wave height from the coeff dictionary.
Definition: waveGenerationModel.C:45
Foam::waveModels::waveGenerationModel::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: waveGenerationModel.C:88
Foam::waveModels::waveGenerationModel::waveGenerationModel
waveGenerationModel(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition: waveGenerationModel.C:69
unitConversion.H
Unit conversion functions.
Foam::FatalIOError
IOerror FatalIOError
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:689
waveGenerationModel
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::waveModels::waveGenerationModel::readWaveAngle
scalar readWaveAngle() const
Helper function to read the wave angle from the coeff dictionary.
Definition: waveGenerationModel.C:60
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
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::waveModel::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: waveModel.C:295
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
waveGenerationModel.H
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::waveModels::defineTypeNameAndDebug
defineTypeNameAndDebug(waveAbsorptionModel, 0)