waveForcing.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2022-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::fv::waveForcing
26 
27 Description
28  This fvModel applies forcing to the liquid phase-fraction field and all
29  components of the vector field to relax the fields towards those
30  calculated from the current wave distribution.
31 
32  The force coefficient \f$\lambda\f$ should be set based on the desired level
33  of forcing and the residence time the waves through the forcing zone. For
34  example, if waves moving at 2 [m/s] are travelling through a forcing zone 8
35  [m] in length, then the residence time is 4 [s]. If it is deemed necessary
36  to force for 5 time-scales, then \f$\lambda\f$ should be set to equal 5/(4
37  [s]) = 1.2 [1/s]. If more aggressive forcing is required adjacent to the
38  boundaries, which is often the case if wave boundary conditions are
39  specified at outflow boundaries, the optional \c lambdaBoundary
40  coefficient can be specified higher than the value of \f$\lambda\f$.
41 
42  Alternatively the forcing force coefficient \c lambdaCoeff can be
43  specified in which case \f$\lambda\f$ is evaluated by multiplying the
44  maximum wave speed by \c lambdaCoeff and dividing by the forcing region
45  length. \c lambdaBoundary is similarly evaluated from
46  \c lambdaBoundaryCoeff if specified.
47 
48 Usage
49  Example usage:
50  \verbatim
51  waveForcing1
52  {
53  type waveForcing;
54 
55  libs ("libwaves.so");
56 
57  liquidPhase water;
58 
59  // Define the line along which to apply the graduation
60  origin (600 0 0);
61  direction (-1 0 0);
62 
63  // // Or, define multiple lines
64  // origins ((600 0 0) (600 -300 0) (600 300 0));
65  // directions ((-1 0 0) (0 1 0) (0 -1 0));
66 
67  scale
68  {
69  type halfCosineRamp;
70  start 0;
71  duration 300;
72  }
73 
74  lambdaCoeff 2; // Forcing coefficient
75 
76  // lambdaBoundaryCoeff 10; // Optional boundary forcing coefficient
77  // Useful when wave BCs are specified at outlets
78 
79  // Write the forcing fields: forcing:scale, forcing:forceCoeff
80  writeForceFields true;
81  }
82  \endverbatim
83 
84 See also
85  Foam::fv::forcing
86 
87 SourceFiles
88  waveForcing.C
89 
90 \*---------------------------------------------------------------------------*/
91 
92 #ifndef waveForcing_H
93 #define waveForcing_H
94 
95 #include "forcing.H"
96 #include "waveSuperposition.H"
97 
98 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100 namespace Foam
101 {
102 namespace fv
103 {
104 
105 /*---------------------------------------------------------------------------*\
106  Class waveForcing Declaration
107 \*---------------------------------------------------------------------------*/
108 
109 class waveForcing
110 :
111  public forcing
112 {
113  // Private Data
114 
115  scalar lambdaCoeff_;
116 
117  scalar lambdaBoundaryCoeff_;
118 
119  dimensionedScalar regionLength_;
120 
121  //- Reference to the waves
122  const waveSuperposition& waves_;
123 
124  //- Name of the liquid phase
125  const word liquidPhaseName_;
126 
127  //- Name of the liquid phase-fraction field
128  const word alphaName_;
129 
130  //- Name of the velocity field
131  const word UName_;
132 
133  //- Phase-fraction field calculated from the current wave form
134  tmp<volScalarField::Internal> alphaWaves_;
135 
136  //- Velocity field calculated from the current wave form
138 
139  //- Forcing coefficient scale field
141 
142 
143  // Private Member Functions
144 
145  //- Non-virtual read
146  void readCoeffs();
147 
148  //- Return the scale distribution
149  virtual tmp<volScalarField::Internal> scale() const;
150 
151  //- Return the force coefficient
152  virtual tmp<volScalarField::Internal> forceCoeff() const;
153 
154 
155 public:
156 
157  //- Runtime type information
158  TypeName("waveForcing");
159 
160 
161  // Constructors
162 
163  //- Construct from components
165  (
166  const word& name,
167  const word& modelType,
168  const fvMesh& mesh,
169  const dictionary& dict
170  );
171 
172 
173  //- Destructor
174  virtual ~waveForcing()
175  {}
176 
177 
178  // Member Functions
179 
180  //- Read dictionary
181  virtual bool read(const dictionary& dict);
182 
183 
184  // Checks
185 
186  //- Return the list of fields for which the fvModel adds source term
187  // to the transport equation
188  virtual wordList addSupFields() const;
189 
190 
191  // Add explicit and implicit contributions
192 
193  //- Source term to VoF phase-fraction equation
194  virtual void addSup
195  (
196  const volScalarField& alpha,
197  fvMatrix<scalar>& eqn
198  ) const;
199 
200  //- Source term to momentum equation
201  virtual void addSup
202  (
203  const volScalarField& rho,
204  const volVectorField& U,
205  fvMatrix<vector>& eqn
206  ) const;
207 
208 
209  // Mesh changes
210 
211  //- Update for mesh motion
212  virtual bool movePoints();
213 
214  //- Update topology using the given map
215  virtual void topoChange(const polyTopoChangeMap&);
216 
217  //- Update from another mesh using the given map
218  virtual void mapMesh(const polyMeshMap&);
219 
220  //- Redistribute or update using the given distribution map
221  virtual void distribute(const polyDistributionMap&);
222 
223 
224  //- Correct the wave forcing coefficients
225  virtual void correct();
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace fv
232 } // End namespace Foam
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 #endif
237 
238 // ************************************************************************* //
Generic GeometricField class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:47
Base fvModel for forcing functions.
Definition: forcing.H:59
This fvModel applies forcing to the liquid phase-fraction field and all components of the vector fiel...
Definition: waveForcing.H:111
virtual bool movePoints()
Update for mesh motion.
Definition: waveForcing.C:170
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: waveForcing.C:122
virtual void correct()
Correct the wave forcing coefficients.
Definition: waveForcing.C:197
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: waveForcing.C:179
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: waveForcing.C:191
virtual ~waveForcing()
Destructor.
Definition: waveForcing.H:173
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: waveForcing.C:232
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: waveForcing.C:185
TypeName("waveForcing")
Runtime type information.
virtual void addSup(const volScalarField &alpha, fvMatrix< scalar > &eqn) const
Source term to VoF phase-fraction equation.
Definition: waveForcing.C:129
waveForcing(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: waveForcing.C:98
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity....
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
labelList fv(nPoints)
dictionary dict