Loading...
Searching...
No Matches
interfaceHeight.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2019 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "interfaceHeight.H"
30#include "fvMesh.H"
31#include "interpolation.H"
32#include "IOmanip.H"
33#include "meshSearch.H"
34#include "midPointAndFaceSet.H"
35#include "Time.H"
37#include "volFields.H"
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace functionObjects
45{
48}
49}
50
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
54void Foam::functionObjects::interfaceHeight::writePositions()
55{
58 vector gHat = vector::zero;
59
60 if (mag(direction_) > 0.0)
61 {
62 gHat = direction_/mag(direction_);
63 }
64 else
65 {
66 gHat = g.value()/mag(g.value());
67 }
68
69 const volScalarField& alpha =
70 mesh_.lookupObject<volScalarField>(alphaName_);
71
72 autoPtr<interpolation<scalar>>
73 interpolator
74 (
75 interpolation<scalar>::New(interpolationScheme_, alpha)
76 );
77
78 if (Pstream::master())
79 {
80 files(fileID::heightFile) << mesh_.time().timeName() << tab;
81 files(fileID::positionFile) << mesh_.time().timeName() << tab;
82 }
83
84 forAll(locations_, li)
85 {
86 // Create a set along a ray projected in the direction of gravity
87 const midPointAndFaceSet set
88 (
89 "",
90 mesh_,
91 meshSearch(mesh_),
92 "xyz",
93 locations_[li] + gHat*mesh_.bounds().mag(),
94 locations_[li] - gHat*mesh_.bounds().mag()
95 );
96
97 // Find the height of the location above the boundary
98 scalar hLB = set.size() ? - gHat & (locations_[li] - set[0]) : - GREAT;
99 reduce(hLB, maxOp<scalar>());
100
101 // Calculate the integrals of length and length*alpha along the sampling
102 // line. The latter is equal to the equivalent length with alpha equal
103 // to one.
104 scalar sumLength = 0, sumLengthAlpha = 0;
105 for(label si = 0; si < set.size() - 1; ++ si)
106 {
107 if (set.segments()[si] != set.segments()[si+1])
108 {
109 continue;
110 }
111
112 const vector& p0 = set[si], p1 = set[si+1];
113 const label c0 = set.cells()[si], c1 = set.cells()[si+1];
114 const label f0 = set.faces()[si], f1 = set.faces()[si+1];
115 const scalar a0 = interpolator->interpolate(p0, c0, f0);
116 const scalar a1 = interpolator->interpolate(p1, c1, f1);
117
118 const scalar l = - gHat & (p1 - p0);
119 sumLength += l;
120 sumLengthAlpha += l*(a0 + a1)/2;
121 }
122
123 reduce(sumLength, sumOp<scalar>());
124 reduce(sumLengthAlpha, sumOp<scalar>());
125
126 // Write out
127 if (Pstream::master())
128 {
129 // Interface heights above the boundary and location
130 const scalar hIB =
131 liquid_ ? sumLengthAlpha : sumLength - sumLengthAlpha;
132 const scalar hIL = hIB - hLB;
133
134 // Position of the interface
135 const point p = locations_[li] - gHat*hIL;
136
137 const Foam::Omanip<int> w = valueWidth(1);
138
139 files(fileID::heightFile) << w << hIB << w << hIL;
140 files(fileID::positionFile) << '(' << w << p.x() << w << p.y()
141 << valueWidth() << p.z() << ") ";
142 }
143 }
144
145 if (Pstream::master())
146 {
147 files(fileID::heightFile).endl();
148 files(fileID::positionFile).endl();
149 }
150}
151
152
153// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
154
156{
157 forAll(locations_, li)
158 {
159 writeHeaderValue
160 (
161 files(fid),
162 "Location " + Foam::name(li),
163 locations_[li]
164 );
165 }
166
167 switch (fileID(fid))
168 {
169 case fileID::heightFile:
170 {
171 writeHeaderValue
172 (
173 files(fid),
174 "hB",
175 "Interface height above the boundary"
176 );
177 writeHeaderValue
178 (
179 files(fid),
180 "hL",
181 "Interface height above the location"
182 );
183 break;
184 }
185 case fileID::positionFile:
186 {
187 writeHeaderValue(files(fid), "p", "Interface position");
188 break;
189 }
190 }
191
192 const Foam::Omanip<int> w = valueWidth(1);
193
194 writeCommented(files(fid), "Location");
195 forAll(locations_, li)
196 {
197 switch (fid)
198 {
199 case fileID::heightFile:
200 files(fid) << w << li << w << ' ';
201 break;
202 case fileID::positionFile:
203 files(fid) << w << li << w << ' ' << w << ' ' << " ";
204 break;
205 }
206 }
207 files(fid).endl();
208
209 writeCommented(files(fid), "Time");
210 forAll(locations_, li)
211 {
212 switch (fid)
213 {
214 case fileID::heightFile:
215 files(fid) << w << "hB" << w << "hL";
216 break;
217 case fileID::positionFile:
218 files(fid) << w << "p" << w << ' ' << w << ' ' << " ";
219 break;
220 }
222 files(fid).endl();
223}
224
225
226// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
227
229(
230 const word& name,
231 const Time& runTime,
232 const dictionary& dict
233)
234:
235 fvMeshFunctionObject(name, runTime, dict),
236 logFiles(obr_, name),
237 liquid_(true),
238 alphaName_("alpha"),
239 interpolationScheme_("cellPoint"),
240 direction_(vector::zero),
241 locations_()
242{
243 read(dict);
244 resetNames({"height", "position"});
245
246 if (Pstream::master())
247 {
248 writeFileHeader(fileID::heightFile);
249 writeFileHeader(fileID::positionFile);
250 }
251}
252
253
254// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
255
257{
258 dict.readIfPresent("alpha", alphaName_);
259 dict.readIfPresent("liquid", liquid_);
260 dict.readEntry("locations", locations_);
261 dict.readIfPresent("interpolationScheme", interpolationScheme_);
262 dict.readIfPresent("direction", direction_);
263
264 return true;
265}
266
269{
270 return true;
271}
272
281{
283
284 writePositions();
285
286 return true;
287}
288
289
290// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
An Ostream manipulator taking arguments.
Definition IOmanip.H:143
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static const Form zero
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
This function object reports the height of the interface above a set of locations.
interfaceHeight(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
PtrList< OFstream > & files()
Inherit logFiles methods.
Definition logFiles.C:109
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual void writeFileHeader(const fileID fid)
Output file header information.
virtual bool execute()
Execute the function-object operations.
OFstream & files(const fileID fid)
Return file corresponding to enumeration.
virtual bool write()
Write the function-object results.
virtual bool end()
Execute at the final time-loop.
virtual void resetNames(const wordList &names)
Reset the list of names from a wordList.
Definition logFiles.C:48
virtual bool write()
Write function.
Definition logFiles.C:142
const objectRegistry & obr_
Reference to the region objectRegistry.
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
Omanip< int > valueWidth(const label offset=0) const
Return the value width when writing to stream with optional offset.
Definition writeFile.C:173
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
engineTime & runTime
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
const dimensionedScalar a0
Bohr radius: default SI units: [m].
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
vector point
Point is a vector.
Definition point.H:37
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Various UniformDimensionedField types.