interfaceHeight.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-2021 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 \*---------------------------------------------------------------------------*/
25 
26 #include "interfaceHeight.H"
27 #include "fvMesh.H"
28 #include "interpolation.H"
29 #include "IOmanip.H"
30 #include "meshSearch.H"
31 #include "lineCellFace.H"
32 #include "Time.H"
34 #include "volFields.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43  defineTypeNameAndDebug(interfaceHeight, 0);
44  addToRunTimeSelectionTable(functionObject, interfaceHeight, dictionary);
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::functionObjects::interfaceHeight::writePositions()
52 {
54  mesh_.lookupObject<uniformDimensionedVectorField>("g");
55  const vector gHat = g.value()/mag(g.value());
56 
57  const volScalarField& alpha =
58  mesh_.lookupObject<volScalarField>(alphaName_);
59 
60  autoPtr<interpolation<scalar>>
61  interpolator
62  (
63  interpolation<scalar>::New(interpolationScheme_, alpha)
64  );
65 
66  if (Pstream::master())
67  {
68  writeTime(file(fileID::heightFile));
69  writeTime(file(fileID::positionFile));
70  }
71 
72  forAll(locations_, li)
73  {
74  // Create a set along a ray projected in the direction of gravity
75  const sampledSets::lineCellFace set
76  (
77  "",
78  mesh_,
79  meshSearch(mesh_),
80  coordSet::axisTypeNames_[coordSet::axisType::XYZ],
81  locations_[li] + gHat*mesh_.bounds().mag(),
82  locations_[li] - gHat*mesh_.bounds().mag()
83  );
84 
85  // Find the height of the location above the boundary
86  scalar hLB =
87  set.size()
88  ? - gHat & (locations_[li] - set.pointCoord(0))
89  : - vGreat;
90  reduce(hLB, maxOp<scalar>());
91 
92  // Calculate the integrals of length and length*alpha along the sampling
93  // line. The latter is equal to the equivalent length with alpha equal
94  // to one.
95  scalar sumLength = 0, sumLengthAlpha = 0;
96  for(label si = 0; si < set.size() - 1; ++ si)
97  {
98  if (set.segments()[si] != set.segments()[si+1])
99  {
100  continue;
101  }
102 
103  const vector& p0 = set.pointCoord(si), p1 = set.pointCoord(si+1);
104  const label c0 = set.cells()[si], c1 = set.cells()[si+1];
105  const label f0 = set.faces()[si], f1 = set.faces()[si+1];
106  const scalar a0 = interpolator->interpolate(p0, c0, f0);
107  const scalar a1 = interpolator->interpolate(p1, c1, f1);
108 
109  const scalar l = - gHat & (p1 - p0);
110  sumLength += l;
111  sumLengthAlpha += l*(a0 + a1)/2;
112  }
113 
114  reduce(sumLength, sumOp<scalar>());
115  reduce(sumLengthAlpha, sumOp<scalar>());
116 
117  // Write out
118  if (Pstream::master())
119  {
120  // Interface heights above the boundary and location
121  const scalar hIB =
122  liquid_ ? sumLengthAlpha : sumLength - sumLengthAlpha;
123  const scalar hIL = hIB - hLB;
124 
125  // Position of the interface
126  const point p = locations_[li] - gHat*hIL;
127 
128  const Foam::Omanip<int> w = valueWidth(1);
129 
130  file(fileID::heightFile) << w << hIB << w << hIL;
131  file(fileID::positionFile) << '(' << w << p.x() << w << p.y()
132  << valueWidth() << p.z() << ") ";
133  }
134  }
135 
136  if (Pstream::master())
137  {
138  file(fileID::heightFile).endl();
139  file(fileID::positionFile).endl();
140  }
141 }
142 
143 
144 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
145 
147 {
148  forAll(locations_, li)
149  {
150  writeHeaderValue
151  (
152  file(i),
153  "Location " + Foam::name(li),
154  locations_[li]
155  );
156  }
157 
158  switch (fileID(i))
159  {
160  case fileID::heightFile:
161  writeHeaderValue
162  (
163  file(i),
164  "hB",
165  "Interface height above the boundary"
166  );
167  writeHeaderValue
168  (
169  file(i),
170  "hL",
171  "Interface height above the location"
172  );
173  break;
174  case fileID::positionFile:
175  writeHeaderValue(file(i), "p", "Interface position");
176  break;
177  }
178 
179  const Foam::Omanip<int> w = valueWidth(1);
180 
181  writeCommented(file(i), "Location");
182  forAll(locations_, li)
183  {
184  switch (fileID(i))
185  {
186  case fileID::heightFile:
187  file(i) << w << li << w << ' ';
188  break;
189  case fileID::positionFile:
190  file(i) << w << li << w << ' ' << w << ' ' << " ";
191  break;
192  }
193  }
194  file(i).endl();
195 
196  writeCommented(file(i), "Time");
197  forAll(locations_, li)
198  {
199  switch (fileID(i))
200  {
201  case fileID::heightFile:
202  file(i) << w << "hB" << w << "hL";
203  break;
204  case fileID::positionFile:
205  file(i) << w << "p" << w << ' ' << w << ' ' << " ";
206  break;
207  }
208  }
209  file(i).endl();
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
214 
216 (
217  const word& name,
218  const Time& runTime,
219  const dictionary& dict
220 )
221 :
222  fvMeshFunctionObject(name, runTime, dict),
223  logFiles(mesh_, name),
224  alphaName_("alpha"),
225  liquid_(true),
226  locations_(),
227  interpolationScheme_("cellPoint")
228 {
229  read(dict);
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
234 
236 {}
237 
238 
239 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
240 
242 {
244 
245  dict.readIfPresent("alpha", alphaName_);
246  dict.readIfPresent("liquid", liquid_);
247  dict.lookup("locations") >> locations_;
248  dict.readIfPresent("interpolationScheme", interpolationScheme_);
249 
250  resetNames({"height", "position"});
251 
252  return true;
253 }
254 
255 
257 {
258  return wordList(alphaName_);
259 }
260 
261 
263 {
264  return true;
265 }
266 
267 
269 {
270  return true;
271 }
272 
273 
275 {
276  logFiles::write();
277 
278  writePositions();
279 
280  return true;
281 }
282 
283 
284 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Write function.
Definition: logFiles.C:167
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static autoPtr< interpolation< scalar > > New(const word &interpolationType, const GeometricField< scalar, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolation.C:48
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
scalar f1
Definition: createFields.H:15
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual void writeFileHeader(const label i=0)
Output file header information.
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
interfaceHeight(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Istream and Ostream manipulators taking arguments.
virtual bool end()
Execute at the final time-loop.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
List< word > wordList
A List of words.
Definition: fileName.H:54
vector point
Point is a vector.
Definition: point.H:41
virtual bool read(const dictionary &)
Read.
dimensioned< scalar > mag(const dimensioned< Type > &)
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual wordList fields() const
Return the list of fields required.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57