patchSummary.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  patchSummary
26 
27 Description
28  Writes fields and boundary condition info for each patch at each requested
29  time instance.
30 
31  Default action is to write a single entry for patches/patchGroups with the
32  same boundary conditions. Use the -expand option to print every patch
33  separately. In case of multiple groups matching it will print only the
34  first one.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "volFields.H"
40 #include "pointFields.H"
41 #include "IOobjectList.H"
42 #include "timeSelector.H"
43 #include "patchSummaryTemplates.H"
44 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
52 
53  #include "addRegionOption.H"
55  (
56  "expand",
57  "Do not combine patches"
58  );
59  #include "setRootCase.H"
60  #include "createTime.H"
61 
63 
64  const bool expand = args.optionFound("expand");
65 
66 
67  #include "createNamedMesh.H"
68  const polyBoundaryMesh& bm = mesh.boundaryMesh();
69 
70 
71  forAll(timeDirs, timeI)
72  {
73  runTime.setTime(timeDirs[timeI], timeI);
74 
75  Info<< "Time = " << runTime.userTimeName() << nl << endl;
76 
77  // Update the mesh if changed
78  if (mesh.readUpdate() == polyMesh::TOPO_PATCH_CHANGE)
79  {
80  Info<< "Detected changed patches. Recreating patch group table."
81  << endl;
82  }
83 
84 
85  const IOobjectList fieldObjs(mesh, runTime.name());
86  const wordList objNames = fieldObjs.names();
87 
88  PtrList<volScalarField> vsf(objNames.size());
89  PtrList<volVectorField> vvf(objNames.size());
90  PtrList<volSphericalTensorField> vsptf(objNames.size());
91  PtrList<volSymmTensorField> vsytf(objNames.size());
92  PtrList<volTensorField> vtf(objNames.size());
93 
94  PtrList<pointScalarField> psf(objNames.size());
95  PtrList<pointVectorField> pvf(objNames.size());
96  PtrList<pointSphericalTensorField> psptf(objNames.size());
97  PtrList<pointSymmTensorField> psytf(objNames.size());
98  PtrList<pointTensorField> ptf(objNames.size());
99 
100  Info<< "Valid fields:" << endl;
101 
102  forAll(objNames, objI)
103  {
104  IOobject obj
105  (
106  objNames[objI],
107  runTime.name(),
108  mesh,
110  );
111 
112  if (obj.headerOk())
113  {
114  addToFieldList(vsf, obj, objI, mesh);
115  addToFieldList(vvf, obj, objI, mesh);
116  addToFieldList(vsptf, obj, objI, mesh);
117  addToFieldList(vsytf, obj, objI, mesh);
118  addToFieldList(vtf, obj, objI, mesh);
119 
120  addToFieldList(psf, obj, objI, pointMesh::New(mesh));
121  addToFieldList(pvf, obj, objI, pointMesh::New(mesh));
122  addToFieldList(psptf, obj, objI, pointMesh::New(mesh));
123  addToFieldList(psytf, obj, objI, pointMesh::New(mesh));
124  addToFieldList(ptf, obj, objI, pointMesh::New(mesh));
125  }
126  }
127 
128  Info<< endl;
129 
130 
131  if (expand)
132  {
133  // Print each patch separately
134 
135  forAll(bm, patchi)
136  {
137  Info<< bm[patchi].type() << "\t: " << bm[patchi].name() << nl;
138  outputFieldList(vsf, patchi);
139  outputFieldList(vvf, patchi);
140  outputFieldList(vsptf, patchi);
141  outputFieldList(vsytf, patchi);
142  outputFieldList(vtf, patchi);
143 
144  outputFieldList(psf, patchi);
145  outputFieldList(pvf, patchi);
146  outputFieldList(psptf, patchi);
147  outputFieldList(psytf, patchi);
148  outputFieldList(ptf, patchi);
149  Info<< endl;
150  }
151  }
152  else
153  {
154  // Collect for each patch the bc type per field. Merge similar
155  // patches.
156 
157  // Per 'group', the map from fieldname to patchfield type
158  DynamicList<HashTable<word>> fieldToTypes(bm.size());
159  // Per 'group' the patches
160  DynamicList<DynamicList<label>> groupToPatches(bm.size());
161  forAll(bm, patchi)
162  {
163  HashTable<word> fieldToType;
164  collectFieldList(vsf, patchi, fieldToType);
165  collectFieldList(vvf, patchi, fieldToType);
166  collectFieldList(vsptf, patchi, fieldToType);
167  collectFieldList(vsytf, patchi, fieldToType);
168  collectFieldList(vtf, patchi, fieldToType);
169 
170  collectFieldList(psf, patchi, fieldToType);
171  collectFieldList(pvf, patchi, fieldToType);
172  collectFieldList(psptf, patchi, fieldToType);
173  collectFieldList(psytf, patchi, fieldToType);
174  collectFieldList(ptf, patchi, fieldToType);
175 
176  label groupI = findIndex(fieldToTypes, fieldToType);
177  if (groupI == -1)
178  {
180  group.append(patchi);
181  groupToPatches.append(group);
182  fieldToTypes.append(fieldToType);
183  }
184  else
185  {
186  groupToPatches[groupI].append(patchi);
187  }
188  }
189 
190 
191  forAll(groupToPatches, groupI)
192  {
193  const DynamicList<label>& patchIDs = groupToPatches[groupI];
194 
195  if (patchIDs.size() > 1)
196  {
197  // Check if part of a group
198  wordList groups;
199  labelHashSet nonGroupPatches;
200  bm.matchGroups(patchIDs, groups, nonGroupPatches);
201 
202  const labelList sortedPatches(nonGroupPatches.sortedToc());
203  forAll(sortedPatches, i)
204  {
205  Info<< bm[sortedPatches[i]].type()
206  << "\t: " << bm[sortedPatches[i]].name() << nl;
207  }
208  if (groups.size())
209  {
210  forAll(groups, i)
211  {
212  Info<< "group\t: " << groups[i] << nl;
213  }
214  }
215  outputFieldList(vsf, patchIDs[0]);
216  outputFieldList(vvf, patchIDs[0]);
217  outputFieldList(vsptf, patchIDs[0]);
218  outputFieldList(vsytf, patchIDs[0]);
219  outputFieldList(vtf, patchIDs[0]);
220 
221  outputFieldList(psf, patchIDs[0]);
222  outputFieldList(pvf, patchIDs[0]);
223  outputFieldList(psptf, patchIDs[0]);
224  outputFieldList(psytf, patchIDs[0]);
225  outputFieldList(ptf, patchIDs[0]);
226  Info<< endl;
227  }
228  else
229  {
230  // No group.
231  forAll(patchIDs, i)
232  {
233  label patchi = patchIDs[i];
234  Info<< bm[patchi].type()
235  << "\t: " << bm[patchi].name() << nl;
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 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
An STL-conforming hash table.
Definition: HashTable.H:127
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:217
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
Foam::polyBoundaryMesh.
void matchGroups(const labelUList &patchIDs, wordList &groups, labelHashSet &nonGroupPatches) const
Match the patches to groups. Returns all the (fully matched) groups.
@ TOPO_PATCH_CHANGE
Definition: polyMesh.H:95
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static instantList timeDirs
Definition: globalFoam.H:44
const char *const group
Group name for atomic constants.
Namespace for OpenFOAM.
void addToFieldList(PtrList< GeoField > &fieldList, const IOobject &obj, const label fieldi, const typename GeoField::Mesh &mesh)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void collectFieldList(const PtrList< GeoField > &fieldList, const label patchi, HashTable< word > &fieldToType)
messageStream Info
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void outputFieldList(const PtrList< GeoField > &fieldList, const label patchi)
static const char nl
Definition: Ostream.H:260
string expand(const string &s, string::size_type &index, const dictionary &dict, const bool allowEnvVars, const bool allowEmpty)
Definition: stringOps.C:146
Foam::argList args(argc, argv)