columnAverage.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-2021 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 "columnAverage.H"
29#include "volFields.H"
31#include "meshStructure.H"
32#include "globalIndex.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace functionObjects
39{
42}
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
49Foam::functionObjects::columnAverage::meshAddressing(const polyMesh& mesh) const
50{
51 if (!meshStructurePtr_)
52 {
53 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
54
55 // Count
56 label sz = 0;
57 for (const label patchi : patchIDs_)
58 {
59 sz += pbm[patchi].size();
60 }
61
62 // Fill
63 labelList meshFaces(sz);
64 sz = 0;
65 for (const label patchi : patchIDs_)
66 {
67 label start = pbm[patchi].start();
68 label size = pbm[patchi].size();
69 for (label i = 0; i < size; ++i)
70 {
71 meshFaces[sz++] = start+i;
72 }
73 }
74
75 if (sz == 0)
76 {
77 // TODO: If source patch is a cyclic it may have have been
78 // converted to a processorCyclic for parallel runs
79
81 << "Requested patches have zero faces"
82 << endl;
83 }
84
86 (
87 UIndirectList<face>(mesh.faces(), meshFaces),
88 mesh.points()
89 );
90
91 globalFaces_.reset(new globalIndex(uip.size()));
92 globalEdges_.reset(new globalIndex(uip.nEdges()));
93 globalPoints_.reset(new globalIndex(uip.nPoints()));
94 meshStructurePtr_.reset
95 (
96 new meshStructure
97 (
98 mesh,
99 uip,
100 globalFaces_(),
101 globalEdges_(),
102 globalPoints_()
103 )
104 );
105 }
106
107 return *meshStructurePtr_;
108}
109
110
111const Foam::word Foam::functionObjects::columnAverage::averageName
112(
113 const word& fieldName
114) const
115{
116 return scopedName("columnAverage(" + fieldName + ")");
117}
118
119
120// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121
123(
124 const word& name,
125 const Time& runTime,
126 const dictionary& dict
127)
128:
130 patchIDs_(),
131 fieldSet_(mesh_)
132{
133 read(dict);
134}
135
136
137// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138
140{
142
143 patchIDs_ =
144 mesh_.boundaryMesh().patchSet
145 (
146 dict.get<wordRes>("patches")
147 ).sortedToc();
148
149 fieldSet_.read(dict);
150
151 return true;
152}
153
154
156{
157 // Make fields up to date with current selection
158 fieldSet_.updateSelection();
159
160 for (const word& fieldName : fieldSet_.selectionNames())
161 {
162 columnAverageField<scalar>(fieldName);
163 columnAverageField<vector>(fieldName);
164 columnAverageField<sphericalTensor>(fieldName);
165 columnAverageField<symmTensor>(fieldName);
166 columnAverageField<tensor>(fieldName);
167 }
168
169 return true;
170}
171
172
174{
175 for (const word& fieldName : fieldSet_.selectionNames())
176 {
177 const regIOobject* obj =
178 obr_.cfindObject<regIOobject>(averageName(fieldName));
179
180 if (obj)
181 {
182 obj->write();
183 }
184 }
185
186 return true;
187}
188
189
190// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
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
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
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 arithmetic average of given quantities along columns of cells in a given direction for s...
virtual bool read(const dictionary &dict)
Read the settings.
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Write the results.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Detect extruded mesh structure given a set of faces (uindirectPrimitivePatch).
Definition: meshStructure.H:62
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
virtual bool write(const bool valid=true) const
Write using setting from DB.
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
dynamicFvMesh & mesh
engineTime & runTime
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict