Test-patchRegion.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 | Copyright (C) 2013-2015 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 Description
25  Detect point pinches
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "argList.H"
30 #include "PatchTools.H"
31 #include "Time.H"
32 #include "polyMesh.H"
33 #include "patchEdgeFaceRegions.H"
34 #include "PatchEdgeFaceWave.H"
35 #include "globalIndex.H"
36 #include "syncTools.H"
37 
38 using namespace Foam;
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 
43 // Main program:
44 
45 int main(int argc, char *argv[])
46 {
47  argList::validArgs.append("patch");
48 
49  #include "setRootCase.H"
50  #include "createTime.H"
51 
52  const word patchName = args[1];
53 
54  #include "createPolyMesh.H"
55 
56  Info<< "Mesh read in = "
57  << runTime.cpuTimeIncrement()
58  << " s\n" << endl << endl;
59 
60 
61  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
62  label patchI = pbm.findPatchID(patchName);
63  const polyPatch& patch = pbm[patchI];
64 
65  Info<< "Patch:" << patch.name() << endl;
66 
67 
68 
69  // Data on all edges and faces
70  List<patchEdgeFaceRegions> allEdgeInfo(patch.nEdges());
71  List<patchEdgeFaceRegions> allFaceInfo(patch.size());
72 
73  // Determine parallel global indexing
74  const globalIndex globalNumbering(patch.size());
75 
76  DynamicList<label> changedEdges(4*patch.size());
77  DynamicList<patchEdgeFaceRegions> changedInfo(changedEdges.size());
78 
79  forAll(patch, faceI)
80  {
81  const labelList& fEdges = patch.faceEdges()[faceI];
82 
83  label globalFaceI = globalNumbering.toGlobal(faceI);
84 
85  forAll(fEdges, i)
86  {
87  changedEdges.append(fEdges[i]);
88  changedInfo.append
89  (
90  patchEdgeFaceRegions(labelPair(globalFaceI, globalFaceI))
91  );
92  }
93  }
94 
95  // Walk
97  <
100  > calc
101  (
102  mesh,
103  patch,
104  changedEdges,
105  changedInfo,
106  allEdgeInfo,
107  allFaceInfo,
108  returnReduce(patch.nEdges(), sumOp<label>())
109  );
110 
111 
112  Info<< "Time now = " << runTime.timeName() << endl;
113 
114 
115  // Detect points with multiple regions
116  labelList duplicateRegion(patch.nPoints(), -1);
117  {
118  labelList currentRegion(patch.nPoints(), -1);
119 
120  forAll(patch.localFaces(), faceI)
121  {
122  const face& f = patch.localFaces()[faceI];
123 
124  forAll(f, fp)
125  {
126  label faceRegion = allFaceInfo[faceI].regions()[fp];
127 
128  label pointI = f[fp];
129 
130  if (currentRegion[pointI] == -1)
131  {
132  currentRegion[pointI] = faceRegion;
133  }
134  else if (currentRegion[pointI] != faceRegion)
135  {
136  if (duplicateRegion[pointI] == -1)
137  {
138  Pout<< "Multi region point:"
139  << patch.localPoints()[pointI]
140  << " with region:" << currentRegion[pointI]
141  << " and region:" << faceRegion
142  << endl;
143  duplicateRegion[pointI] = currentRegion[pointI];
144  }
145  }
146  }
147  }
148  }
149 
150 
151  Info<< "End\n" << endl;
152 
153  return 0;
154 }
155 
156 
157 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
patchEdgeFaceRegions.H
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
main
int main(int argc, char *argv[])
Definition: Test-patchRegion.C:45
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList< label >
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
globalIndex.H
Foam::PatchEdgeFaceWave
Wave propagation of information along patch. Every iteration information goes through one layer of fa...
Definition: PatchEdgeFaceWave.H:69
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::calc
void calc(const argList &args, const Time &runTime, const fvMesh &mesh)
Definition: Lambda2.C:38
Foam::primitivePatch
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Foam::primitivePatch.
Definition: primitivePatch.H:45
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::PrimitivePatch::faceEdges
const labelListList & faceEdges() const
Return face-edge addressing.
Definition: PrimitivePatchTemplate.C:312
polyMesh.H
syncTools.H
PatchEdgeFaceWave.H
Foam::label
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
Foam::patchEdgeFaceRegions
Transport of regions for use in PatchEdgeFaceWave.
Definition: patchEdgeFaceRegions.H:57
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
argList.H
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
setRootCase.H
Foam::sumOp
Definition: ops.H:162
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
createTime.H
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
args
Foam::argList args(argc, argv)
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
createPolyMesh.H
Foam::labelPair
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82