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-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 \*---------------------------------------------------------------------------*/
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  "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 = set.size() ? - gHat & (locations_[li] - set[0]) : - vGreat;
87  reduce(hLB, maxOp<scalar>());
88 
89  // Calculate the integrals of length and length*alpha along the sampling
90  // line. The latter is equal to the equivalent length with alpha equal
91  // to one.
92  scalar sumLength = 0, sumLengthAlpha = 0;
93  for(label si = 0; si < set.size() - 1; ++ si)
94  {
95  if (set.segments()[si] != set.segments()[si+1])
96  {
97  continue;
98  }
99 
100  const vector& p0 = set[si], p1 = set[si+1];
101  const label c0 = set.cells()[si], c1 = set.cells()[si+1];
102  const label f0 = set.faces()[si], f1 = set.faces()[si+1];
103  const scalar a0 = interpolator->interpolate(p0, c0, f0);
104  const scalar a1 = interpolator->interpolate(p1, c1, f1);
105 
106  const scalar l = - gHat & (p1 - p0);
107  sumLength += l;
108  sumLengthAlpha += l*(a0 + a1)/2;
109  }
110 
111  reduce(sumLength, sumOp<scalar>());
112  reduce(sumLengthAlpha, sumOp<scalar>());
113 
114  // Write out
115  if (Pstream::master())
116  {
117  // Interface heights above the boundary and location
118  const scalar hIB =
119  liquid_ ? sumLengthAlpha : sumLength - sumLengthAlpha;
120  const scalar hIL = hIB - hLB;
121 
122  // Position of the interface
123  const point p = locations_[li] - gHat*hIL;
124 
125  const Foam::Omanip<int> w = valueWidth(1);
126 
127  file(fileID::heightFile) << w << hIB << w << hIL;
128  file(fileID::positionFile) << '(' << w << p.x() << w << p.y()
129  << valueWidth() << p.z() << ") ";
130  }
131  }
132 
133  if (Pstream::master())
134  {
135  file(fileID::heightFile).endl();
136  file(fileID::positionFile).endl();
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
142 
144 {
145  forAll(locations_, li)
146  {
147  writeHeaderValue
148  (
149  file(i),
150  "Location " + Foam::name(li),
151  locations_[li]
152  );
153  }
154 
155  switch (fileID(i))
156  {
157  case fileID::heightFile:
158  writeHeaderValue
159  (
160  file(i),
161  "hB",
162  "Interface height above the boundary"
163  );
164  writeHeaderValue
165  (
166  file(i),
167  "hL",
168  "Interface height above the location"
169  );
170  break;
171  case fileID::positionFile:
172  writeHeaderValue(file(i), "p", "Interface position");
173  break;
174  }
175 
176  const Foam::Omanip<int> w = valueWidth(1);
177 
178  writeCommented(file(i), "Location");
179  forAll(locations_, li)
180  {
181  switch (fileID(i))
182  {
183  case fileID::heightFile:
184  file(i) << w << li << w << ' ';
185  break;
186  case fileID::positionFile:
187  file(i) << w << li << w << ' ' << w << ' ' << " ";
188  break;
189  }
190  }
191  file(i).endl();
192 
193  writeCommented(file(i), "Time");
194  forAll(locations_, li)
195  {
196  switch (fileID(i))
197  {
198  case fileID::heightFile:
199  file(i) << w << "hB" << w << "hL";
200  break;
201  case fileID::positionFile:
202  file(i) << w << "p" << w << ' ' << w << ' ' << " ";
203  break;
204  }
205  }
206  file(i).endl();
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
211 
213 (
214  const word& name,
215  const Time& runTime,
216  const dictionary& dict
217 )
218 :
219  fvMeshFunctionObject(name, runTime, dict),
220  logFiles(mesh_, name),
221  alphaName_("alpha"),
222  liquid_(true),
223  locations_(),
224  interpolationScheme_("cellPoint")
225 {
226  read(dict);
227  resetNames({"height", "position"});
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
232 
234 {}
235 
236 
237 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
238 
240 {
241  dict.readIfPresent("alpha", alphaName_);
242  dict.readIfPresent("liquid", liquid_);
243  dict.lookup("locations") >> locations_;
244  dict.readIfPresent("interpolationScheme", interpolationScheme_);
245 
246  return true;
247 }
248 
249 
251 {
252  return true;
253 }
254 
255 
257 {
258  return true;
259 }
260 
261 
263 {
264  logFiles::write();
265 
266  writePositions();
267 
268  return true;
269 }
270 
271 
272 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Write function.
Definition: logFiles.C:180
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
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
static autoPtr< interpolation< scalar > > New(const word &interpolationType, const GeometricField< scalar, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
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:68
scalar f1
Definition: createFields.H:28
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
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)
vector point
Point is a vector.
Definition: point.H:41
virtual bool read(const dictionary &)
Read.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
dimensioned< scalar > mag(const dimensioned< Type > &)
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57