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