Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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"
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
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 (
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// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Generic GeometricField class.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Computes the spatial minimum and maximum extents of an input field.
volFieldSelection fieldSet_
Fields to assess.
scalar threshold_
Threshold value.
void calcFieldExtents(const word &fieldName, const bool calcMag=false)
Main calculation.
fieldExtents(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
bool internalField_
Flag to write the internal field extents.
tmp< volScalarField > calcMask(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return the field mask.
point C0_
Reference position.
virtual void writeFileHeader(Ostream &os)
Output file header information.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
labelList patchIDs_
Patches to assess.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
virtual bool updateSelection()
Update the selection using current contents of obr_.
Base class for writing single files from the function objects.
Definition writeFile.H:113
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition writeFile.C:344
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
virtual void writeBreak(Ostream &os) const
Write a break marker to the stream.
Definition writeFile.C:367
bool writtenHeader_
Flag to identify whether the header has been written.
Definition writeFile.H:157
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
rDeltaTY field()
engineTime & runTime
OBJstream os(runTime.globalPath()/outputName)
auto & name
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299