Loading...
Searching...
No Matches
patchSummary.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2023 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
27Application
28 patchSummary
29
30Group
31 grpMiscUtilities
32
33Description
34 Write field and boundary condition info for each patch at each requested
35 time instance.
36
37 Default action is to write a single entry for patches/patchGroups with the
38 same boundary conditions. Use the -expand option to print every patch
39 separately. In case of multiple groups matching it will print only the
40 first one.
41
42\*---------------------------------------------------------------------------*/
43
44#include "fvCFD.H"
45#include "volFields.H"
46#include "pointFields.H"
47#include "IOobjectList.H"
49
50// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52int main(int argc, char *argv[])
53{
54 argList::addNote
55 (
56 "Write field and boundary condition info for each patch"
57 " at each requested time instance"
58 );
59
60 timeSelector::addOptions();
61
62 #include "addRegionOption.H"
63 argList::addBoolOption
64 (
65 "expand",
66 "Do not combine patches"
67 );
68 #include "setRootCase.H"
69 #include "createTime.H"
70
71 instantList timeDirs = timeSelector::select0(runTime, args);
72
73 const bool expand = args.found("expand");
74
75
76 #include "createNamedMesh.H"
77 const polyBoundaryMesh& bm = mesh.boundaryMesh();
78
79
80 forAll(timeDirs, timeI)
81 {
82 runTime.setTime(timeDirs[timeI], timeI);
83
84 Info<< "Time = " << runTime.timeName() << nl << endl;
85
86 // Update the mesh if changed
87 if (mesh.readUpdate() == polyMesh::TOPO_PATCH_CHANGE)
88 {
89 Info<< "Detected changed patches. Recreating patch group table."
90 << endl;
91 }
92
93 const IOobjectList objects(mesh, runTime.timeName());
94
95 Info<< "Reading fields:" << endl;
96
97 // Read fields
98 #undef createFields
99 #define createFields(FieldType, Variable) \
100 PtrList<FieldType> Variable \
101 ( \
102 readFields<FieldType>(objects, mesh) \
103 );
104
105 createFields(volScalarField, vsf);
106 createFields(volVectorField, vvf);
107 createFields(volSphericalTensorField, vsptf);
108 createFields(volSymmTensorField, vsytf);
109 createFields(volTensorField, vtf);
110
111 // Point fields
112 const pointMesh& pMesh = pointMesh::New(mesh);
113
114 #undef createFields
115 #define createFields(FieldType, Variable) \
116 PtrList<FieldType> Variable \
117 ( \
118 readFields<FieldType>(objects, pMesh) \
119 );
120
121 createFields(pointScalarField, psf);
122 createFields(pointVectorField, pvf);
123 createFields(pointSphericalTensorField, psptf);
124 createFields(pointSymmTensorField, psytf);
125 createFields(pointTensorField, ptf);
126
127 #undef createFields
128
129 Info<< endl;
130
131
132 if (expand)
133 {
134 // Print each patch separately
135
136 forAll(bm, patchi)
137 {
138 Info<< bm[patchi].type() << "\t: " << bm[patchi].name() << nl;
139
140 outputFieldList(vsf, patchi);
141 outputFieldList(vvf, patchi);
142 outputFieldList(vsptf, patchi);
143 outputFieldList(vsytf, patchi);
144 outputFieldList(vtf, patchi);
145
146 outputFieldList(psf, patchi);
147 outputFieldList(pvf, patchi);
148 outputFieldList(psptf, patchi);
149 outputFieldList(psytf, patchi);
150 outputFieldList(ptf, patchi);
151 Info<< endl;
152 }
153 }
154 else
155 {
156 // Collect for each patch the bc type per field. Merge similar
157 // patches.
158
159 // Per 'group', the map from fieldname to patchfield type
160 DynamicList<HashTable<word>> fieldToTypes(bm.size());
161 // Per 'group' the patches
162 DynamicList<DynamicList<label>> groupToPatches(bm.size());
163
164 forAll(bm, patchi)
165 {
166 HashTable<word> fieldToType;
167 collectFieldList(vsf, patchi, fieldToType);
168 collectFieldList(vvf, patchi, fieldToType);
169 collectFieldList(vsptf, patchi, fieldToType);
170 collectFieldList(vsytf, patchi, fieldToType);
171 collectFieldList(vtf, patchi, fieldToType);
172
173 collectFieldList(psf, patchi, fieldToType);
174 collectFieldList(pvf, patchi, fieldToType);
175 collectFieldList(psptf, patchi, fieldToType);
176 collectFieldList(psytf, patchi, fieldToType);
177 collectFieldList(ptf, patchi, fieldToType);
178
179 label groupI = fieldToTypes.find(fieldToType);
180 if (groupI == -1)
181 {
182 DynamicList<label> group(1);
183 group.append(patchi);
184 groupToPatches.append(group);
185 fieldToTypes.append(fieldToType);
186 }
187 else
188 {
189 groupToPatches[groupI].append(patchi);
190 }
191 }
192
193
194 for (const auto& patchIDs : groupToPatches)
195 {
196 if (patchIDs.size() > 1)
197 {
198 // Check if part of a group
199 wordList groups;
200 labelHashSet nonGroupPatches;
201 bm.matchGroups(patchIDs, groups, nonGroupPatches);
202
203 for (const label patchi : nonGroupPatches.sortedToc())
204 {
205 Info<< bm[patchi].type()
206 << "\t: " << bm[patchi].name() << nl;
207 }
208 for (const word& groupName : groups)
209 {
210 Info<< "group\t: " << groupName << nl;
211 }
212
213 const label patchi = patchIDs[0];
214
215 outputFieldList(vsf, patchi);
216 outputFieldList(vvf, patchi);
217 outputFieldList(vsptf, patchi);
218 outputFieldList(vsytf, patchi);
219 outputFieldList(vtf, patchi);
220
221 outputFieldList(psf, patchi);
222 outputFieldList(pvf, patchi);
223 outputFieldList(psptf, patchi);
224 outputFieldList(psytf, patchi);
225 outputFieldList(ptf, patchi);
226 Info<< endl;
227 }
228 else
229 {
230 // No group.
231 for (const label patchi : patchIDs)
232 {
233 Info<< bm[patchi].type()
234 << "\t: " << bm[patchi].name() << nl;
235
236 outputFieldList(vsf, patchi);
237 outputFieldList(vvf, patchi);
238 outputFieldList(vsptf, patchi);
239 outputFieldList(vsytf, patchi);
240 outputFieldList(vtf, patchi);
241
242 outputFieldList(psf, patchi);
243 outputFieldList(pvf, patchi);
244 outputFieldList(psptf, patchi);
245 outputFieldList(psytf, patchi);
246 outputFieldList(ptf, patchi);
247 Info<< endl;
248 }
249 }
250 }
251 }
252 }
253
254 Info<< "End\n" << endl;
255
256 return 0;
257}
258
259
260// ************************************************************************* //
Required Classes.
labelList patchIDs
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
constexpr const char *const group
Group name for atomic constants.
string expand(const std::string &s, const HashTable< string > &mapping, const char sigil='$')
Expand occurrences of variables according to the mapping and return the expanded string.
List< word > wordList
List of word.
Definition fileName.H:60
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
void collectFieldList(const UPtrList< GeoField > &fieldList, const label patchi, HashTable< word > &fieldToType)
List< instant > instantList
List of instants.
Definition instantList.H:41
void outputFieldList(const UPtrList< GeoField > &fieldList, const label patchi)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299