setWaves.C
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) 2017-2023 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 Application
25  setWaves
26 
27 Description
28  Applies wave models to the entire domain for case initialisation using
29  level sets for second-order accuracy.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "levelSet.H"
35 #include "pointFields.H"
36 #include "timeSelector.H"
38 #include "volFields.H"
39 #include "wallPolyPatch.H"
42 #include "systemDict.H"
43 
44 using namespace Foam;
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 int main(int argc, char *argv[])
49 {
50  timeSelector::addOptions(false, false);
51 
52  #include "addDictOption.H"
53  #include "addRegionOption.H"
54 
56  (
57  "alpha",
58  "name",
59  "name of the volume fraction field, default is \"alpha\""
60  );
61 
63  (
64  "U",
65  "name",
66  "name of the velocity field, default is \"U\""
67  );
68 
70  (
71  "gas",
72  "the volume fraction field is that of the gas above the wave"
73  );
74 
75  #include "setRootCase.H"
76  #include "createTime.H"
77 
79 
80  #include "createNamedMesh.H"
81 
82  const dictionary setWavesDict(systemDict("setWavesDict", args, mesh));
83 
84  #include "readGravitationalAcceleration.H"
85 
86  const pointMesh& pMesh = pointMesh::New(mesh);
87 
88  // Parse options
89  const word alphaName = setWavesDict.lookupOrDefault<word>("alpha", "alpha");
90  const word UName = setWavesDict.lookupOrDefault<word>("U", "U");
91  const bool liquid = setWavesDict.lookupOrDefault<bool>("liquid", true);
92 
93  // Get the wave models
94  const waveSuperposition& waves = waveSuperposition::New(mesh);
95 
96  forAll(timeDirs, timeI)
97  {
98  runTime.setTime(timeDirs[timeI], timeI);
99  const scalar t = runTime.value();
100 
101  Info<< "Time = " << runTime.userTimeName() << nl << endl;
102 
103  mesh.readUpdate();
104 
105  // Read the fields which are to be set
107  (
108  IOobject
109  (
110  alphaName,
111  runTime.name(),
112  mesh,
114  ),
115  mesh
116  );
118  (
119  IOobject
120  (
121  UName,
122  runTime.name(),
123  mesh,
125  ),
126  mesh
127  );
128 
129  // Create modelled fields on both cells and points
131  (
132  IOobject("h", runTime.name(), mesh),
133  mesh,
135  );
137  (
138  IOobject("hp", runTime.name(), mesh),
139  pMesh,
141  );
142  volVectorField uGas
143  (
144  IOobject("uGas", runTime.name(), mesh),
145  mesh,
147  );
148  pointVectorField uGasp
149  (
150  IOobject("uGasp", runTime.name(), mesh),
151  pMesh,
153  );
154  volVectorField uLiq
155  (
156  IOobject("uLiq", runTime.name(), mesh),
157  mesh,
159  );
160  pointVectorField uLiqp
161  (
162  IOobject("uLiqp", runTime.name(), mesh),
163  pMesh,
165  );
166 
167  // Cell centres and points
168  const pointField& ccs = mesh.cellCentres();
169  const pointField& pts = mesh.points();
170 
171  // Internal field
172  h.primitiveFieldRef() = waves.height(t, ccs);
173  hp.primitiveFieldRef() = waves.height(t, pts);
174  uGas.primitiveFieldRef() = waves.UGas(t, ccs);
175  uGasp.primitiveFieldRef() = waves.UGas(t, pts);
176  uLiq.primitiveFieldRef() = waves.ULiquid(t, ccs);
177  uLiqp.primitiveFieldRef() = waves.ULiquid(t, pts);
178 
179  // Boundary fields
180  forAll(mesh.boundary(), patchj)
181  {
182  const pointField& fcs = mesh.boundary()[patchj].Cf();
183  h.boundaryFieldRef()[patchj] = waves.height(t, fcs);
184  uGas.boundaryFieldRef()[patchj] = waves.UGas(t, fcs);
185  uLiq.boundaryFieldRef()[patchj] = waves.ULiquid(t, fcs);
186  }
187 
188  // Calculate the fields
189  volScalarField alphaNoBCs(levelSetFraction(h, hp, !liquid));
190  volVectorField UNoBCs(levelSetAverage(h, hp, uGas, uGasp, uLiq, uLiqp));
191 
192  // Set the wave and non-wall fixed-value patch fields
193  forAll(mesh.boundary(), patchi)
194  {
195  const polyPatch& patch = mesh.boundaryMesh()[patchi];
196 
198  fvPatchVectorField& Up = U.boundaryFieldRef()[patchi];
199  if
200  (
201  !isA<wallPolyPatch>(patch)
202  || isA<waveAlphaFvPatchScalarField>(alphap)
203  || isA<waveVelocityFvPatchVectorField>(Up)
204  )
205  {
206  alphap == alphaNoBCs.boundaryField()[patchi];
207  Up == UNoBCs.boundaryField()[patchi];
208  }
209  }
210 
211  // Set the internal fields and all non-fixed value patch fields
212  alpha = alphaNoBCs;
213  U = UNoBCs;
214 
215  // Output
216  Info<< "Writing " << alpha.name() << nl;
217  alpha.write();
218  Info<< "Writing " << U.name() << nl << endl;
219  U.write();
220  }
221 
222  Info<< "End\n" << endl;
223 
224  return 0;
225 }
226 
227 
228 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
static const Form zero
Definition: VectorSpace.H:113
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:270
virtual bool write(const bool write=true) const
Write using setting from DB.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0)
Definition: timeSelector.C:283
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity....
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static instantList timeDirs
Definition: globalFoam.H:44
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar h
Planck constant.
Namespace for OpenFOAM.
tmp< scalarField > levelSetFraction(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const bool above)
Calculate the volume-fraction that a level set occupies. This gives the.
Definition: levelSet.C:34
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
const dimensionSet dimLength
const dimensionSet dimVelocity
static const char nl
Definition: Ostream.H:260
tmp< Field< Type > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const Field< Type > &positiveC, const Field< Type > &positiveP, const Field< Type > &negativeC, const Field< Type > &negativeP)
Calculate the average value of two fields, one on each side of a level set.
Foam::argList args(argc, argv)