setWaves.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 | Copyright (C) 2017 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 "fvCFD.H"
34 #include "levelSet.H"
35 #include "pointFields.H"
36 #include "timeSelector.H"
37 #include "wallDist.H"
40 #include "waveSuperposition.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 int main(int argc, char *argv[])
45 {
46  timeSelector::addOptions(false, false);
47 
49  (
50  "U",
51  "name",
52  "name of the velocity field, default is \"U\""
53  );
54 
56  (
57  "alpha",
58  "name",
59  "name of the volume fraction field, default is \"alpha\""
60  );
61 
62  #include "setRootCase.H"
63  #include "createTime.H"
64 
65  instantList timeDirs = timeSelector::selectIfPresent(runTime, args);
66 
67  #include "createMesh.H"
68  #include "readGravitationalAcceleration.H"
69 
70  const pointMesh& pMesh = pointMesh::New(mesh);
71 
72  forAll(timeDirs, timeI)
73  {
74  runTime.setTime(timeDirs[timeI], timeI);
75 
76  Info<< "Time = " << runTime.timeName() << nl << endl;
77 
78  mesh.readUpdate();
79 
80  // Read the fields which are to be set
82  (
83  IOobject
84  (
85  args.optionFound("alpha") ? args["alpha"] : "alpha",
86  runTime.timeName(),
87  mesh,
88  IOobject::MUST_READ
89  ),
90  mesh
91  );
93  (
94  IOobject
95  (
96  args.optionFound("U") ? args["U"] : "U",
97  runTime.timeName(),
98  mesh,
99  IOobject::MUST_READ
100  ),
101  mesh
102  );
103 
104  // Create modelled fields on both cells and points
105  volScalarField height
106  (
107  IOobject("height", runTime.timeName(), mesh),
108  mesh,
110  );
111  pointScalarField heightp
112  (
113  IOobject("heightp", runTime.timeName(), mesh),
114  pMesh,
116  );
117  volVectorField uGas
118  (
119  IOobject("uGas", runTime.timeName(), mesh),
120  mesh,
121  dimensionedVector("0", dimVelocity, vector::zero)
122  );
123  pointVectorField uGasp
124  (
125  IOobject("uGasp", runTime.timeName(), mesh),
126  pMesh,
127  dimensionedVector("0", dimLength, vector::zero)
128  );
129  volVectorField uLiquid
130  (
131  IOobject("uLiquid", runTime.timeName(), mesh),
132  mesh,
133  dimensionedVector("0", dimVelocity, vector::zero)
134  );
135  pointVectorField uLiquidp
136  (
137  IOobject("uLiquidp", runTime.timeName(), mesh),
138  pMesh,
139  dimensionedVector("0", dimLength, vector::zero)
140  );
141 
142  // The number of wave patches
143  label nWaves = 0;
144 
145  // Whether the alpha conditions refer to the liquid phase
146  bool liquid = false;
147 
148  // The mean velocity of one of the wave patches
149  vector UMeanp = vector::zero;
150 
151  // Loop the patches, averaging and superimposing wave model data
152  forAll(mesh.boundary(), patchi)
153  {
154  fvPatchScalarField& alphap = alpha.boundaryFieldRef()[patchi];
155  fvPatchVectorField& Up = U.boundaryFieldRef()[patchi];
156 
157  const bool isWave = isA<waveAlphaFvPatchScalarField>(alphap);
158 
159  if (isA<waveVelocityFvPatchVectorField>(Up) != isWave)
160  {
162  << "The alpha condition on patch " << Up.patch().name()
163  << " is " << alphap.type() << " and the velocity condition"
164  << " is " << Up.type() << ". Wave boundary conditions must"
165  << " be set in pairs. If the alpha condition is "
166  << waveAlphaFvPatchScalarField::typeName
167  << " then the velocity condition must be "
168  << waveVelocityFvPatchVectorField::typeName
169  << " and vice-versa." << exit(FatalError);
170  }
171 
172  if (!isWave)
173  {
174  continue;
175  }
176 
177  Info<< "Adding waves from patch " << Up.patch().name() << endl;
178 
179  const waveSuperposition& waves =
180  refCast<waveVelocityFvPatchVectorField>(Up).waves();
181 
182  UMeanp = waves.UMean();
183 
184  const bool liquidp =
185  refCast<waveAlphaFvPatchScalarField>(alphap).liquid();
186  if (nWaves > 0 && liquidp != liquid)
187  {
189  << "All " << waveAlphaFvPatchScalarField::typeName
190  << "patch fields must be configured for the same phase,"
191  << " i.e., the liquid switch must have the same value."
192  << exit(FatalError);
193  }
194  liquid = liquidp;
195 
196  const scalar t = runTime.value();
197  const pointField& ccs = mesh.cellCentres();
198  const pointField& pts = mesh.points();
199 
200  // Internal field superposition
201  height.primitiveFieldRef() += waves.height(t, ccs);
202  heightp.primitiveFieldRef() += waves.height(t, pts);
203  uGas.primitiveFieldRef() += waves.UGas(t, ccs) - UMeanp;
204  uGasp.primitiveFieldRef() += waves.UGas(t, pts) - UMeanp;
205  uLiquid.primitiveFieldRef() += waves.ULiquid(t, ccs) - UMeanp;
206  uLiquidp.primitiveFieldRef() += waves.ULiquid(t, pts) - UMeanp;
207 
208  // Boundary field superposition
209  forAll(mesh.boundary(), patchj)
210  {
211  const pointField& fcs = mesh.boundary()[patchj].Cf();
212  height.boundaryFieldRef()[patchj] += waves.height(t, fcs);
213  uGas.boundaryFieldRef()[patchj] += waves.UGas(t, fcs) - UMeanp;
214  uLiquid.boundaryFieldRef()[patchj] +=
215  waves.ULiquid(t, fcs) - UMeanp;
216  }
217 
218  ++ nWaves;
219  }
220 
221  // Warn and skip to the next time if no wave patches were found
222  if (nWaves == 0)
223  {
225  << "No " << waveAlphaFvPatchScalarField::typeName << " or "
226  << waveVelocityFvPatchVectorField::typeName << " patch fields "
227  << "were found. No waves have been set." << endl;
228 
229  continue;
230  }
231 
232  // Create the mean velocity field
234  (
235  IOobject("UMean", runTime.timeName(), mesh),
236  mesh,
237  dimensionedVector("UMean", dimVelocity, UMeanp)
238  );
239 
240  if (nWaves > 1)
241  {
242  // Create weighted average fields for the mean velocity
243  volScalarField weight
244  (
245  IOobject("weight", runTime.timeName(), mesh),
246  mesh,
248  );
249  volVectorField weightUMean
250  (
251  IOobject("weightUMean", runTime.timeName(), mesh),
252  mesh,
253  dimensionedVector("0", dimVelocity/dimLength, vector::zero)
254  );
255 
256  // Loop the patches, inverse-distance weighting the mean velocities
257  forAll(mesh.boundary(), patchi)
258  {
259  fvPatchScalarField& alphap = alpha.boundaryFieldRef()[patchi];
260  fvPatchVectorField& Up = U.boundaryFieldRef()[patchi];
261 
262  const bool isWave = isA<waveAlphaFvPatchScalarField>(alphap);
263 
264  if (!isWave)
265  {
266  continue;
267  }
268 
269  const waveSuperposition& waves =
270  refCast<waveVelocityFvPatchVectorField>(Up).waves();
271 
272  UMeanp = waves.UMean();
273 
274  const volScalarField w
275  (
276  1
277  /(
278  wallDist(mesh, labelList(1, patchi)).y()
279  + dimensionedScalar("ySmall", dimLength, SMALL)
280  )
281  );
282  weight += w;
283  weightUMean +=
284  w*dimensionedVector("UMeanp", dimVelocity, UMeanp);
285  }
286 
287  // Complete the average for the mean velocity
288  UMean = weightUMean/weight;
289  }
290 
291  // Set the internal fields
292  alpha.ref() = levelSetFraction(mesh, height, heightp, !liquid);
293  U.ref() =
294  UMean
296  (
297  mesh,
298  height,
299  heightp,
300  uGas,
301  uGasp,
302  uLiquid,
303  uLiquidp
304  );
305 
306  // Set the boundary fields
307  forAll(mesh.boundary(), patchi)
308  {
309  fvPatchScalarField& alphap = alpha.boundaryFieldRef()[patchi];
310  fvPatchVectorField& Up = U.boundaryFieldRef()[patchi];
311 
312  const bool isWave = isA<waveAlphaFvPatchScalarField>(alphap);
313 
314  if (!isWave)
315  {
316  alphap ==
318  (
319  mesh.boundary()[patchi],
320  height.boundaryField()[patchi],
321  heightp.boundaryField()[patchi].patchInternalField(),
322  !liquid
323  );
324  Up ==
325  UMean.boundaryField()[patchi]
327  (
328  mesh.boundary()[patchi],
329  height.boundaryField()[patchi],
330  heightp.boundaryField()[patchi].patchInternalField()(),
331  uGas.boundaryField()[patchi],
332  uGasp.boundaryField()[patchi].patchInternalField()(),
333  uLiquid.boundaryField()[patchi],
334  uLiquidp.boundaryField()[patchi].patchInternalField()()
335  );
336  }
337  else
338  {
339  alphap == refCast<waveAlphaFvPatchScalarField>(alphap).alpha();
340  Up == refCast<waveVelocityFvPatchVectorField>(Up).U();
341  }
342  }
343 
344  // Output
345  Info<< "Writing " << alpha.name() << nl;
346  alpha.write();
347  Info<< "Writing " << U.name() << nl << endl;
348  U.write();
349  }
350 
351  Info<< "End\n" << endl;
352 
353  return 0;
354 }
355 
356 
357 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
List< instant > instantList
List of instants.
Definition: instantList.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
volVectorField UMean(UMeanHeader, mesh)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
tmp< DimensionedField< Type, volMesh > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const DimensionedField< Type, volMesh > &positiveC, const DimensionedField< Type, pointMesh > &positiveP, const DimensionedField< Type, volMesh > &negativeC, const DimensionedField< Type, pointMesh > &negativeP)
Calculate the average value of two fields, one on each side of a level set.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
scalar y
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:95
List< label > labelList
A List of labels.
Definition: labelList.H:56
tmp< DimensionedField< scalar, volMesh > > 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 the.
Definition: levelSet.C:35
static const char nl
Definition: Ostream.H:262
const word & name() const
Return const reference to name.
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Definition: pointFields.H:50
messageStream Info
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Foam::argList args(argc, argv)
const dimensionSet dimVelocity