fieldExtents.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) 2018-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fieldExtents.H"
29 #include "processorPolyPatch.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace functionObjects
37 {
38  defineTypeNameAndDebug(fieldExtents, 0);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
47 {
49  {
50  return;
51  }
52 
53  if (writtenHeader_)
54  {
55  writeBreak(os);
56  }
57  else
58  {
59  writeHeader(os, "Field extents");
60  writeHeaderValue(os, "Reference position", C0_);
61  writeHeaderValue(os, "Threshold", threshold_);
62  }
63 
64  writeCommented(os, "Time");
65 
66  for (const word& fieldName : fieldSet_.selectionNames())
67  {
68  if (internalField_)
69  {
70  writeTabbed(os, fieldName + "_internal");
71  }
72  for (const label patchi : patchIDs_)
73  {
74  const word& patchName = mesh_.boundaryMesh()[patchi].name();
75  writeTabbed(os, fieldName + "_" + patchName);
76  }
77  }
78 
79  os << endl;
80 
81  writtenHeader_ = true;
82 }
83 
84 
85 template<>
87 (
89 ) const
90 {
91  return
92  pos
93  (
94  field
95  - dimensionedScalar("t", field.dimensions(), threshold_)
96  );
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
103 (
104  const word& name,
105  const Time& runTime,
106  const dictionary& dict
107 )
108 :
110  writeFile(mesh_, name, typeName, dict),
111  internalField_(true),
112  threshold_(0),
113  C0_(Zero),
114  fieldSet_(mesh_)
115 {
116  read(dict);
117 
118  // Note: delay creating the output file header to handle field names
119  // specified using regular expressions
120 }
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 {
127  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
128 
130  {
131  dict.readIfPresent<bool>("internalField", internalField_);
132 
133  threshold_ = dict.get<scalar>("threshold");
134 
135  dict.readIfPresent<vector>("referencePosition", C0_);
136 
138  if (dict.readIfPresent("patches", patchNames))
139  {
140  patchIDs_ = pbm.indices(patchNames);
141  }
142  else
143  {
144  labelHashSet patchSet(2*pbm.size());
145 
146  // All non-processor and non-empty patches
147  forAll(pbm, patchi)
148  {
149  const polyPatch& pp = pbm[patchi];
150 
151  if
152  (
153  !isA<processorPolyPatch>(pp)
154  && !isA<emptyPolyPatch>(pp)
155  )
156  {
157  patchSet.insert(patchi);
158  }
159  }
160 
161  patchIDs_ = patchSet.sortedToc();
162  }
163 
164  if (!internalField_ && patchIDs_.empty())
165  {
167  << "No internal field or patches selected - no field extent "
168  << "information will be generated" << endl;
169  }
170 
171  fieldSet_.read(dict);
172 
173  return true;
174  }
175 
176  return false;
177 }
178 
181 {
182  return true;
183 }
184 
185 
187 {
188  writeFileHeader(file());
189 
190  Log << type() << " " << name() << " write:" << nl;
191 
192  for (const word& fieldName : fieldSet_.selectionNames())
193  {
194  calcFieldExtents<scalar>(fieldName, true);
195  calcFieldExtents<vector>(fieldName);
196  calcFieldExtents<sphericalTensor>(fieldName);
197  calcFieldExtents<symmTensor>(fieldName);
198  calcFieldExtents<tensor>(fieldName);
199  }
200 
201  Log << endl;
202 
203  return true;
204 }
205 
206 
207 // ************************************************************************* //
const polyBoundaryMesh & pbm
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
Computes the spatial minimum and maximum extents of an input field.
Definition: fieldExtents.H:179
rDeltaTY field()
labelList patchIDs_
Patches to assess.
Definition: fieldExtents.H:211
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void writeBreak(Ostream &os) const
Write a break marker to the stream.
Definition: writeFile.C:362
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
tmp< volScalarField > calcMask(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return the field mask.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
wordList selectionNames() const
Return the current field selection, in sorted order.
point C0_
Reference position.
Definition: fieldExtents.H:201
bool writtenHeader_
Flag to identify whether the header has been written.
Definition: writeFile.H:157
Abstract base-class for Time/database function objects.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
scalar threshold_
Threshold value.
Definition: fieldExtents.H:196
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dimensionedScalar pos(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
virtual bool write()
Write the fieldExtents.
Definition: fieldExtents.C:179
virtual void writeFileHeader(Ostream &os)
Output file header information.
Definition: fieldExtents.C:39
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
wordList patchNames(nPatches)
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
volFieldSelection fieldSet_
Fields to assess.
Definition: fieldExtents.H:206
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define Log
Definition: PDRblock.C:28
virtual bool updateSelection()
Update the selection using current contents of obr_.
virtual bool read(const dictionary &)
Read the field extents data.
Definition: fieldExtents.C:118
virtual bool read(const dictionary &dict)
Read optional controls.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
bool internalField_
Flag to write the internal field extents.
Definition: fieldExtents.H:191
virtual bool execute()
Execute, currently does nothing.
Definition: fieldExtents.C:173
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Base class for writing single files from the function objects.
Definition: writeFile.H:112
const fvMesh & mesh_
Reference to the fvMesh.
fieldExtents(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: fieldExtents.C:96
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329