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-2018 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 
43 using namespace Foam;
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  timeSelector::addOptions(false, false);
50 
51  #include "addDictOption.H"
52  #include "addRegionOption.H"
53 
55  (
56  "alpha",
57  "name",
58  "name of the volume fraction field, default is \"alpha\""
59  );
60 
62  (
63  "U",
64  "name",
65  "name of the velocity field, default is \"U\""
66  );
67 
69  (
70  "gas",
71  "the volume fraction field is that that of the gas above the wave"
72  );
73 
74  #include "setRootCase.H"
75  #include "createTime.H"
76 
78 
79  #include "createNamedMesh.H"
80 
81  const word dictName("setWavesDict");
83  Info<< "Reading " << dictName << "\n" << endl;
84  IOdictionary setWavesDict(dictIO);
85 
86  #include "readGravitationalAcceleration.H"
87 
88  const pointMesh& pMesh = pointMesh::New(mesh);
89 
90  // Parse options
91  const word alphaName = setWavesDict.lookupOrDefault<word>("alpha", "alpha");
92  const word UName = setWavesDict.lookupOrDefault<word>("U", "U");
93  const bool liquid = setWavesDict.lookupOrDefault<bool>("liquid", true);
94 
95  // Get the wave models
97 
98  forAll(timeDirs, timeI)
99  {
100  runTime.setTime(timeDirs[timeI], timeI);
101  const scalar t = runTime.value();
102 
103  Info<< "Time = " << runTime.timeName() << nl << endl;
104 
105  mesh.readUpdate();
106 
107  // Read the fields which are to be set
109  (
110  IOobject
111  (
112  alphaName,
113  runTime.timeName(),
114  mesh,
116  ),
117  mesh
118  );
120  (
121  IOobject
122  (
123  UName,
124  runTime.timeName(),
125  mesh,
127  ),
128  mesh
129  );
130 
131  // Create modelled fields on both cells and points
133  (
134  IOobject("h", runTime.timeName(), mesh),
135  mesh,
137  );
139  (
140  IOobject("hp", runTime.timeName(), mesh),
141  pMesh,
143  );
144  volVectorField uGas
145  (
146  IOobject("uGas", runTime.timeName(), mesh),
147  mesh,
149  );
150  pointVectorField uGasp
151  (
152  IOobject("uGasp", runTime.timeName(), mesh),
153  pMesh,
155  );
156  volVectorField uLiq
157  (
158  IOobject("uLiq", runTime.timeName(), mesh),
159  mesh,
161  );
162  pointVectorField uLiqp
163  (
164  IOobject("uLiqp", runTime.timeName(), mesh),
165  pMesh,
167  );
168 
169  // Cell centres and points
170  const pointField& ccs = mesh.cellCentres();
171  const pointField& pts = mesh.points();
172 
173  // Internal field
174  h.primitiveFieldRef() = waves.height(t, ccs);
175  hp.primitiveFieldRef() = waves.height(t, pts);
176  uGas.primitiveFieldRef() = waves.UGas(t, ccs);
177  uGasp.primitiveFieldRef() = waves.UGas(t, pts);
178  uLiq.primitiveFieldRef() = waves.ULiquid(t, ccs);
179  uLiqp.primitiveFieldRef() = waves.ULiquid(t, pts);
180 
181  // Boundary fields
182  forAll(mesh.boundary(), patchj)
183  {
184  const pointField& fcs = mesh.boundary()[patchj].Cf();
185  h.boundaryFieldRef()[patchj] = waves.height(t, fcs);
186  uGas.boundaryFieldRef()[patchj] = waves.UGas(t, fcs);
187  uLiq.boundaryFieldRef()[patchj] = waves.ULiquid(t, fcs);
188  }
189 
190  // Calculate the fields
191  volScalarField alphaNoBCs(levelSetFraction(h, hp, !liquid));
192  volVectorField UNoBCs(levelSetAverage(h, hp, uGas, uGasp, uLiq, uLiqp));
193 
194  // Set the wave and non-wall fixed-value patch fields
196  {
197  const polyPatch& patch = mesh.boundaryMesh()[patchi];
198 
200  fvPatchVectorField& Up = U.boundaryFieldRef()[patchi];
201  if
202  (
203  !isA<wallPolyPatch>(patch)
204  || isA<waveAlphaFvPatchScalarField>(alphap)
205  || isA<waveVelocityFvPatchVectorField>(Up)
206  )
207  {
208  alphap == alphaNoBCs.boundaryField()[patchi];
209  Up == UNoBCs.boundaryField()[patchi];
210  }
211  }
212 
213  // Set the internal fields and all non-fixed value patch fields
214  alpha = alphaNoBCs;
215  U = UNoBCs;
216 
217  // Output
218  Info<< "Writing " << alpha.name() << nl;
219  alpha.write();
220  Info<< "Writing " << U.name() << nl << endl;
221  U.write();
222  }
223 
224  Info<< "End\n" << endl;
225 
226  return 0;
227 }
228 
229 
230 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const word & name() const
Return name.
Definition: IOobject.H:303
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
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
dynamicFvMesh & mesh
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:489
A class for handling words, derived from string.
Definition: word.H:59
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
const Type & value() const
Return const reference to value.
const word dictName("particleTrackDict")
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0)
Definition: timeSelector.C:283
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:879
const vectorField & cellCentres() const
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
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.
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
label patchi
U
Definition: pEqn.H:72
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
const dimensionedScalar & h
Planck constant.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual bool write(const bool write=true) const
Write using setting from DB.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540
const dimensionSet dimVelocity