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-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
34namespace Foam
35{
36namespace functionObjects
37{
40}
41}
42
43
44// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45
47{
49 {
50 return;
51 }
52
54 {
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 {
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
85template<>
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 patchIDs_()
116{
117 read(dict);
118
119 // Note: delay creating the output file header to handle field names
120 // specified using regular expressions
121}
122
123
124// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125
127{
129 {
130 dict.readIfPresent<bool>("internalField", internalField_);
131
132 threshold_ = dict.get<scalar>("threshold");
133
134 dict.readIfPresent<vector>("referencePosition", C0_);
135
136 patchIDs_.clear();
137 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
138
140 if (dict.readIfPresent("patches", patchNames))
141 {
142 for (const wordRe& name : patchNames)
143 {
144 patchIDs_.insert(pbm.indices(name));
145 }
146 }
147 else
148 {
149 // Add all non-processor and non-empty patches
150 forAll(pbm, patchi)
151 {
152 const polyPatch& pp = pbm[patchi];
153 if (!isA<processorPolyPatch>(pp) && !isA<emptyPolyPatch>(pp))
154 {
155 patchIDs_.insert(patchi);
156 }
157 }
158 }
159
160 if (!internalField_ && patchIDs_.empty())
161 {
163 << "No internal field or patches selected - no field extent "
164 << "information will be generated" << endl;
165 }
166
167 fieldSet_.read(dict);
168
169 return true;
170 }
171
172 return false;
173}
174
175
177{
178 return true;
179}
180
181
183{
184 writeFileHeader(file());
185
186 Log << type() << " " << name() << " write:" << nl;
187
188 for (const word& fieldName : fieldSet_.selectionNames())
189 {
190 calcFieldExtents<scalar>(fieldName, true);
191 calcFieldExtents<vector>(fieldName);
192 calcFieldExtents<sphericalTensor>(fieldName);
193 calcFieldExtents<symmTensor>(fieldName);
194 calcFieldExtents<tensor>(fieldName);
195 }
196
197 Log << endl;
198
199 return true;
200}
201
202
203// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Generic GeometricField class.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
Computes the spatial minimum and maximum extents of an input field.
Definition: fieldExtents.H:187
volFieldSelection fieldSet_
Fields to assess.
Definition: fieldExtents.H:203
scalar threshold_
Threshold value.
Definition: fieldExtents.H:197
labelHashSet patchIDs_
Patches to assess.
Definition: fieldExtents.H:206
bool internalField_
Flag to write the internal field extents.
Definition: fieldExtents.H:194
tmp< volScalarField > calcMask(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return the field mask.
point C0_
Reference position.
Definition: fieldExtents.H:200
virtual void writeFileHeader(Ostream &os)
Output file header information.
Definition: fieldExtents.C:46
virtual bool execute()
Execute, currently does nothing.
Definition: fieldExtents.C:176
virtual bool write()
Write the fieldExtents.
Definition: fieldExtents.C:182
virtual bool read(const dictionary &)
Read the field extents data.
Definition: fieldExtents.C:126
wordHashSet selectionNames() const
Return the current field selection.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
virtual bool updateSelection()
Update the selection using current contents of obr_.
Base class for writing single files from the function objects.
Definition: writeFile.H:120
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:285
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:295
virtual void writeBreak(Ostream &os) const
Write a break marker to the stream.
Definition: writeFile.C:318
bool writtenHeader_
Flag to identify whether the header has been written.
Definition: writeFile.H:148
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:269
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) patch indices for all matches.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
A class for managing temporary objects.
Definition: tmp.H:65
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:83
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
rDeltaTY field()
engineTime & runTime
OBJstream os(runTime.globalPath()/outputName)
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
wordList patchNames(nPatches)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333