Test-PrimitivePatch.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) 2011-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  Test new primitive patches.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "PrimitivePatch.H"
30 #include "argList.H"
31 #include "Time.H"
32 #include "polyMesh.H"
33 #include "primitivePatch.H"
34 #include "IFstream.H"
35 #include "OFstream.H"
36 
37 using namespace Foam;
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
42 
43 
45 {
46  forAll(points, pointI)
47  {
48  const point& pt = points[pointI];
49 
50  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
51  }
52 }
53 
54 
55 void checkFaceEdges
56 (
57  const faceList& localFaces,
58  const edgeList& edges,
59  const labelListList& faceEdges
60 )
61 {
62  forAll(faceEdges, faceI)
63  {
64  const face& f = localFaces[faceI];
65 
66  const labelList& myEdges = faceEdges[faceI];
67 
68  forAll(f, fp)
69  {
70  label fp1 = f.fcIndex(fp);
71 
72  if (edges[myEdges[fp]] != edge(f[fp], f[fp1]))
73  {
75  << "Edges of face not in face point order:"
76  << "face:" << faceI << " localF:" << f
77  << " faceEdges:" << myEdges
78  << abort(FatalError);
79  }
80  }
81  }
82 }
83 
84 
85 void writeEdges
86 (
87  const pointField& localPoints,
88  const edgeList& edges,
89  const label nInternalEdges
90 )
91 {
92  Info<< "Writing internal edges to internalEdges.obj" << nl << endl;
93 
94  OFstream intStream("internalEdges.obj");
95 
96  writeObj(intStream, localPoints);
97 
98  for (label edgeI = 0; edgeI < nInternalEdges; edgeI++)
99  {
100  const edge& e = edges[edgeI];
101 
102  intStream << "l " << e.start()+1 << ' ' << e.end()+1 << endl;
103  }
104 
105  Info<< "Writing boundary edges to boundaryEdges.obj" << nl << endl;
106 
107  OFstream bndStream("boundaryEdges.obj");
108 
109  writeObj(bndStream, localPoints);
110 
111  for (label edgeI = nInternalEdges; edgeI < edges.size(); edgeI++)
112  {
113  const edge& e = edges[edgeI];
114 
115  bndStream << "l " << e.start()+1 << ' ' << e.end()+1 << endl;
116  }
117 }
118 
119 
120 void writeFaceEdges
121 (
122  const pointField& localPoints,
123  const edgeList& edges,
124  const labelListList& faceEdges
125 )
126 {
127  Info<< "Writing faceEdges per face to faceEdges.obj" << nl << endl;
128 
129  OFstream feStream("faceEdges.obj");
130 
131  writeObj(feStream, localPoints);
132 
133  forAll(faceEdges, faceI)
134  {
135  const labelList& myEdges = faceEdges[faceI];
136 
137  forAll(myEdges, i)
138  {
139  const edge& e = edges[myEdges[i]];
140 
141  feStream<< "l " << e.start()+1 << ' ' << e.end()+1 << endl;
142  }
143  }
144 }
145 
146 
147 void writeEdgeFaces
148 (
149  const pointField& localPoints,
150  const faceList& localFaces,
151  const edgeList& edges,
152  const labelListList& edgeFaces
153 )
154 {
155  Info<< "Writing edgeFaces per face to edgeFaces.obj" << nl << endl;
156 
157  OFstream efStream("edgeFaces.obj");
158 
159  pointField ctrs(localFaces.size(), vector::zero);
160 
161  forAll(localFaces, faceI)
162  {
163  ctrs[faceI] = localFaces[faceI].centre(localPoints);
164  }
165  writeObj(efStream, ctrs);
166 
167 
168  forAll(edgeFaces, edgeI)
169  {
170  const labelList& myFaces = edgeFaces[edgeI];
171 
172  forAll(myFaces, i)
173  {
174  efStream<< "l " << myFaces[0]+1 << ' ' << myFaces[i]+1 << endl;
175  }
176  }
177 }
178 
179 
180 void writeFaceFaces
181 (
182  const pointField& localPoints,
183  const faceList& localFaces,
184  const labelListList& faceFaces
185 )
186 {
187  Info<< "Writing faceFaces per face to faceFaces.obj" << nl << endl;
188 
189  OFstream ffStream("faceFaces.obj");
190 
191  pointField ctrs(localFaces.size(), vector::zero);
192 
193  forAll(localFaces, faceI)
194  {
195  ctrs[faceI] = localFaces[faceI].centre(localPoints);
196  }
197  writeObj(ffStream, ctrs);
198 
199  forAll(faceFaces, faceI)
200  {
201  const labelList& nbrs = faceFaces[faceI];
202 
203  forAll(nbrs, nbI)
204  {
205  ffStream << "l " << faceI+1 << ' ' << nbrs[nbI]+1 << endl;
206  }
207  }
208 }
209 
210 
211 // Main program:
212 
213 int main(int argc, char *argv[])
214 {
216  argList::validArgs.append("patch");
217 
218  #include "setRootCase.H"
219  #include "createTime.H"
220  #include "createPolyMesh.H"
221 
222  const word patchName = args[1];
223  const polyPatch& patch = mesh.boundaryMesh()[patchName];
224 
225  Info<< "Patch:" << patch.name() << endl;
226 
227  // Test addressing
228  {
229  myPrimitivePatch pp(patch, patch.points());
230 
231  const pointField& localPoints = pp.localPoints();
232  const faceList& localFaces = pp.localFaces();
233  const labelListList& faceFaces = pp.faceFaces();
234  const edgeList& edges = pp.edges();
235  const labelListList& edgeFaces = pp.edgeFaces();
236  const labelListList& faceEdges = pp.faceEdges();
237 
238 
239  checkFaceEdges(localFaces, edges, faceEdges);
240 
241  writeEdges(localPoints, edges, pp.nInternalEdges());
242 
243  writeFaceEdges(localPoints, edges, faceEdges);
244 
245  writeEdgeFaces(localPoints, localFaces, edges, edgeFaces);
246 
247  writeFaceFaces(localPoints, localFaces, faceFaces);
248  }
249 
250  // Test construction from Xfer
251  {
252  faceList patchFaces = patch;
253  pointField allPoints = patch.points();
254 
256  (
257  patchFaces.xfer(),
258  allPoints.xfer()
259  );
260  }
261 
262  Info<< "End\n" << endl;
263 
264  return 0;
265 }
266 
267 
268 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatchTemplate.C:292
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::List::xfer
Xfer< List< T > > xfer()
Transfer contents to the Xfer container.
Definition: ListI.H:90
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
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
writeFaceEdges
void writeFaceEdges(const pointField &localPoints, const edgeList &edges, const labelListList &faceEdges)
Definition: Test-PrimitivePatch.C:121
polyMesh.H
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
myPrimitivePatch
PrimitivePatch< face, List, const pointField & > myPrimitivePatch
Definition: Test-PrimitivePatch.C:41
OFstream.H
main
int main(int argc, char *argv[])
Definition: Test-PrimitivePatch.C:213
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
PrimitivePatch.H
Foam::nl
static const char nl
Definition: Ostream.H:260
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
writeEdgeFaces
void writeEdgeFaces(const pointField &localPoints, const faceList &localFaces, const edgeList &edges, const labelListList &edgeFaces)
Definition: Test-PrimitivePatch.C:148
argList.H
IFstream.H
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
writeObj
void writeObj(Ostream &os, const pointField &points)
Definition: Test-PrimitivePatch.C:44
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
checkFaceEdges
void checkFaceEdges(const faceList &localFaces, const edgeList &edges, const labelListList &faceEdges)
Definition: Test-PrimitivePatch.C:56
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatchTemplate.C:232
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
writeEdges
void writeEdges(const pointField &localPoints, const edgeList &edges, const label nInternalEdges)
Definition: Test-PrimitivePatch.C:86
f
labelList f(nPoints)
Foam::Vector< scalar >
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
createTime.H
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::argList::noParallel
static void noParallel()
Remove the parallel options.
Definition: argList.C:161
writeFaceFaces
void writeFaceFaces(const pointField &localPoints, const faceList &localFaces, const labelListList &faceFaces)
Definition: Test-PrimitivePatch.C:181
Foam::PrimitivePatch::faceFaces
const labelListList & faceFaces() const
Return face-face addressing.
Definition: PrimitivePatchTemplate.C:272
args
Foam::argList args(argc, argv)
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
createPolyMesh.H
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88