faMeshToolsChecks.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) 2021-2022 OpenCFD Ltd.
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 \*---------------------------------------------------------------------------*/
27 
28 #include "faMeshTools.H"
29 #include "areaFields.H"
30 #include "edgeFields.H"
31 #include "processorFaPatch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
36 (
37  const faMesh& mesh,
38  const int verbose
39 )
40 {
41  const faBoundaryMesh& patches = mesh.boundary();
42  const label nNonProcessor = patches.nNonProcessor();
43  const label nPatches = patches.size();
44 
45  label nLocalProcEdges = 0;
46  if (Pstream::parRun())
47  {
48  for (const faPatch& fap : patches)
49  {
50  const auto* cpp = isA<processorFaPatch>(fap);
51 
52  if (cpp)
53  {
54  nLocalProcEdges += fap.nEdges();
55  }
56  }
57  }
58 
59  const labelList nFaces
60  (
61  UPstream::listGatherValues<label>(mesh.nFaces())
62  );
63  const labelList nPoints
64  (
65  UPstream::listGatherValues<label>(mesh.nPoints())
66  );
67  const labelList nEdges
68  (
69  UPstream::listGatherValues<label>(mesh.nEdges())
70  );
71  const labelList nIntEdges
72  (
73  UPstream::listGatherValues<label>(mesh.nInternalEdges())
74  );
75 
76  // The "real" (non-processor) boundary edges
77  const labelList nBndEdges
78  (
79  UPstream::listGatherValues<label>
80  (
81  mesh.nBoundaryEdges() - nLocalProcEdges
82  )
83  );
84  const labelList nProcEdges
85  (
86  UPstream::listGatherValues<label>(nLocalProcEdges)
87  );
88 
89 
90  // Format output as
91  // Number of faces: ...
92  // per-proc: (...)
93 
94  const auto reporter =
95  [&,verbose](const char* tag, const labelList& list)
96  {
97  Info<< " Number of " << tag << ": " << sum(list) << nl;
98  if (Pstream::parRun() && verbose)
99  {
100  int padding = static_cast<int>
101  (
102  // strlen(" Number of ") - strlen("per-proc")
103  (12 - 8)
104  + strlen(tag)
105  );
106 
107  do { Info<< ' '; } while (--padding > 0);
108 
109  Info<< "per-proc: " << flatOutput(list) << nl;
110  }
111  };
112 
113 
114  Info<< "----------------" << nl
115  << "Mesh Information" << nl
116  << "----------------" << nl
117  << " " << "boundingBox: " << boundBox(mesh.points()) << nl;
118 
119  if (Pstream::master())
120  {
121  reporter("faces", nFaces);
122  reporter("points", nPoints);
123  reporter("edges", nEdges);
124  reporter("internal edges", nIntEdges);
125  reporter("boundary edges", nBndEdges);
126 
127  if (Pstream::parRun())
128  {
129  reporter("processor edges", nProcEdges);
130  }
131  }
132 
133  Info<< "----------------" << nl
134  << "Patches" << nl
135  << "----------------" << nl;
136 
137  for (label patchi = 0; patchi < nNonProcessor; ++patchi)
138  {
139  const faPatch& p = patches[patchi];
140 
141  // Report physical size (nEdges) not virtual size
142  Info<< " " << "patch " << p.index()
143  << " (size: " << returnReduce(p.nEdges(), sumOp<label>())
144  << ") name: " << p.name()
145  << nl;
146  }
147 
148  Info<< "----------------" << nl
149  << "Used polyPatches: " << flatOutput(mesh.whichPolyPatches()) << nl;
150 
151 
152  // Geometry information
153  Info<< nl;
154  {
155  scalarMinMax limit(gMinMax(mesh.S().field()));
156  Info<< "Face area:" << nl
157  << " min = " << limit.min() << " max = " << limit.max() << nl;
158  }
159 
160  {
161  scalarMinMax limit(minMax(mesh.magLe().primitiveField()));
162 
163  // Include processor boundaries into 'internal' edges
164  if (Pstream::parRun())
165  {
166  for (label patchi = nNonProcessor; patchi < nPatches; ++patchi)
167  {
168  limit.add(minMax(mesh.magLe().boundaryField()[patchi]));
169  }
170 
172  }
173 
174  Info<< "Edge length (internal):" << nl
175  << " min = " << limit.min() << " max = " << limit.max() << nl;
176 
177 
178  // Include (non-processor) boundaries
179  for (label patchi = 0; patchi < nNonProcessor; ++patchi)
180  {
181  limit.add(minMax(mesh.magLe().boundaryField()[patchi]));
182  }
183 
184  if (Pstream::parRun())
185  {
187  }
188 
189  Info<< "Edge length:" << nl
190  << " min = " << limit.min() << " max = " << limit.max() << nl;
191  }
192 
193  // Not particularly meaningful
194  #if 0
195  {
196  MinMax<vector> limit(gMinMax(mesh.faceAreaNormals().field()));
197 
198  Info<< "Face area normals:" << nl
199  << " min = " << limit.min() << " max = " << limit.max() << nl;
200  }
201  #endif
202 }
203 
204 
205 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:54
dynamicFvMesh & mesh
label nPoints
complex limit(const complex &c1, const complex &c2)
Definition: complexI.H:235
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:72
static void printMeshChecks(const faMesh &mesh, const int verbose=1)
Report mesh information.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
Finite area boundary mesh.
messageStream Info
Information stream (stdout output on master, null elsewhere)
volScalarField & p
Combine values and/or MinMax ranges.
Definition: MinMaxOps.H:153
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225