meshSurfaceEdgeCreateEdgeVertices.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "demandDrivenData.H"
29 #include "polyMeshGenAddressing.H"
31 #include "meshOctree.H"
32 #include "Map.H"
33 
34 #define DEBUGMapping
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
44 {
46  const faceListPMG& faces = mesh_.faces();
47 
48  const edgeList& edges = mesh_.addressingData().edges();
49  const VRWGraph& faceEdges = mesh_.addressingData().faceEdges();
50 
51  Map<label> newEdgePoint;
52 
53  const label nIntFaces = mesh_.nInternalFaces();
54  const label nFaces = faces.size();
55 
56  nPoints_ = points.size();
57 
58  for(label faceI=nIntFaces;faceI<nFaces;++faceI)
59  {
60  const face& f = faces[faceI];
61 
62  forAll(f, pI)
63  {
64  const label edgeI = faceEdges(faceI, pI);
65 
66  if( newEdgePoint.found(edgeI) ) continue;
67 
68  const label s = f[pI];
69  const label e = f.nextLabel(pI);
70 
72  {
73  Warning << "Boundary vertices " << s << " and " << e
74  << " are not mapped to the boundary!" << endl;
75 
76  continue;
77  }
78 
79  if( pointRegions_(s, 0) != pointRegions_(e, 0) )
80  {
81  point newP;
82  scalar distSq;
83  label nse;
84 
85  FixedList<point, 2> edgePoints;
87 
88  edgePoints[0] = points[s];
89  edgePoints[1] = points[e];
90  patches[0] = pointRegions_(s, 0);
91  patches[1] = pointRegions_(e, 0);
92 
93  const bool found =
95  (
96  newP,
97  distSq,
98  nse,
99  edgePoints,
100  patches
101  );
102 
103  if( found )
104  {
105  points.append(newP);
106  }
107  else
108  {
109  points.append
110  (
111  edges[faceEdges(faceI, pI)].centre(points)
112  );
113  }
114 
116 
117  newEdgePoint.insert(edgeI, nPoints_);
118  ++nPoints_;
119  }
120  }
121  }
122 
123  points.setSize(nPoints_);
124 
125  //- create new faces
126  DynList<label> newF;
127  forAll(faces, faceI)
128  {
129  const face& f = faces[faceI];
130 
131  newF.clear();
132 
133  forAll(f, eI)
134  {
135  newF.append(f[eI]);
136  if( newEdgePoint.found(faceEdges(faceI, eI)) )
137  newF.append(newEdgePoint[faceEdges(faceI, eI)]);
138  }
139 
140  if( newF.size() > f.size() )
141  {
142  //- face must be changed
143  face& mf = const_cast<face&>(f);
144  mf.setSize(newF.size());
145  forAll(mf, pI)
146  mf[pI] = newF[pI];
147  }
148  }
149 
151 
152  Info << "Finished creating mesh edges" << endl;
153 }
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 } // End namespace Foam
158 
159 // ************************************************************************* //
Foam::polyMeshGenCells::addressingData
const polyMeshGenAddressing & addressingData() const
addressing which may be needed
Definition: polyMeshGenCells.C:327
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::meshSurfaceEdgeExtractor::createEdgeVertices
void createEdgeVertices()
create vertices on surface edges
Definition: meshSurfaceEdgeCreateEdgeVertices.C:43
Foam::Warning
messageStream Warning
Foam::polyMeshGenAddressing::faceEdges
const VRWGraph & faceEdges() const
Definition: polyMeshGenAddressingFaceEdges.C:112
Foam::Map< label >
Foam::meshOctree::findNearestPointToEdge
bool findNearestPointToEdge(point &nearest, scalar &distSq, label &nearestEdge, const FixedList< point, 2 > &edgePoints, const FixedList< label, 2 > &edgePointRegions) const
Definition: meshOctreeFindNearestSurfacePoint.C:321
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMeshGenPoints::points
const pointFieldPMG & points() const
access to points
Definition: polyMeshGenPointsI.H:44
Foam::meshSurfaceEdgeExtractor::meshOctree_
const meshOctree & meshOctree_
reference to the octree
Definition: meshSurfaceEdgeExtractor.H:66
meshOctree.H
Foam::polyMeshGenFaces::faces
const faceListPMG & faces() const
access to faces
Definition: polyMeshGenFacesI.H:43
Foam::meshSurfaceEdgeExtractor::mesh_
polyMeshGen & mesh_
mesh
Definition: meshSurfaceEdgeExtractor.H:57
Foam::polyMeshGenFaces::nInternalFaces
label nInternalFaces() const
return number of internal faces
Definition: polyMeshGenFacesI.H:48
Foam::meshSurfaceEdgeExtractor::nPoints_
label nPoints_
Definition: meshSurfaceEdgeExtractor.H:59
Map.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::Info
messageStream Info
Foam::polyMeshGenCells::clearAddressingData
void clearAddressingData() const
clear addressing data
Definition: polyMeshGenCells.C:346
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::meshSurfaceEdgeExtractor::pointRegions_
VRWGraph pointRegions_
regions for boundary vertices
Definition: meshSurfaceEdgeExtractor.H:69
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::polyMeshGenAddressing::edges
const edgeList & edges() const
Return mesh edges.
Definition: polyMeshGenAddressingEdges.C:157
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::DynList< label >
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::List::setSize
void setSize(const label)
Reset size of List.
f
labelList f(nPoints)
Foam::faceListPMG::size
label size() const
return the number of used elements
Definition: faceListPMGI.H:73
Foam::VRWGraph::appendList
void appendList(const ListType &l)
Append a list as a row at the end of the graph.
Definition: VRWGraphI.H:286
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
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
points
const pointField & points
Definition: gmvOutputHeader.H:1
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::DynList::size
label size() const
Definition: DynListI.H:235
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::pointFieldPMG
Definition: pointFieldPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
polyMeshGenAddressing.H
Foam::DynList::clear
void clear()
Clear the list, i.e. set next free to zero.
Definition: DynListI.H:279
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304
meshSurfaceEdgeExtractor.H