meshUntanglerCutRegion.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 "meshUntangler.H"
30 #include "plane.H"
31 #include "primitiveMesh.H"
32 #include "polyMeshGenModifier.H"
33 #include "sortEdgesIntoChains.H"
34 
35 //#define DEBUGSmooth
36 
37 #ifdef DEBUGSmooth
38 #include "Time.H"
39 #include "objectRegistry.H"
40 #endif
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
50 (
51  const boundBox& bb
52 )
53 {
54  pointsPtr_ = new DynList<point, 64>();
55  DynList<point, 64>& bVertices = *pointsPtr_;
56  edgesPtr_ = new DynList<edge, 128>();
57  DynList<edge, 128>& bEdges = *edgesPtr_;
58  facesPtr_ = new DynList<DynList<label, 8>, 64>();
59  DynList<DynList<label, 8>, 64>& bFaces = *facesPtr_;
60 
61  //- set vertices
62  const point c = (bb.max() + bb.min()) / 2.0;
63  const point vec = (bb.max() - bb.min()) / 2.0;
64  bVertices.append
65  (
66  point(c.x() - vec.x(), c.y() - vec.y(), c.z() - vec.z())
67  );
68  bVertices.append
69  (
70  point(c.x() + vec.x(), c.y() - vec.y(), c.z() - vec.z())
71  );
72  bVertices.append
73  (
74  point(c.x() + vec.x(), c.y() + vec.y(), c.z() - vec.z())
75  );
76  bVertices.append
77  (
78  point(c.x() - vec.x(), c.y() + vec.y(), c.z() - vec.z())
79  );
80  bVertices.append
81  (
82  point(c.x() - vec.x(), c.y() - vec.y(), c.z() + vec.z())
83  );
84  bVertices.append
85  (
86  point(c.x() + vec.x(), c.y() - vec.y(), c.z() + vec.z())
87  );
88  bVertices.append
89  (
90  point(c.x() + vec.x(), c.y() + vec.y(), c.z() + vec.z())
91  );
92  bVertices.append
93  (
94  point(c.x() - vec.x(), c.y() + vec.y(), c.z() + vec.z())
95  );
96 
97  //- set edges
98 
99  //- edges in x direction
100  bEdges.append(edge(0, 1));
101  bEdges.append(edge(3, 2));
102  bEdges.append(edge(7, 6));
103  bEdges.append(edge(4, 5));
104 
105  //- edges in y direction
106  bEdges.append(edge(1, 2));
107  bEdges.append(edge(0, 3));
108  bEdges.append(edge(4, 7));
109  bEdges.append(edge(5, 6));
110 
111  //- edges in z direction
112  bEdges.append(edge(0, 4));
113  bEdges.append(edge(1, 5));
114  bEdges.append(edge(2, 6));
115  bEdges.append(edge(3, 7));
116 
117  //- set faces
119  f.setSize(4);
120 
121  //- faces in x direction
122  f[0] = 5;
123  f[1] = 11;
124  f[2] = 6;
125  f[3] = 8;
126  bFaces.append(f);
127  f[0] = 4;
128  f[1] = 10;
129  f[2] = 7;
130  f[3] = 9;
131  bFaces.append(f);
132  //- faces in y direction
133  f[0] = 0;
134  f[1] = 8;
135  f[2] = 3;
136  f[3] = 9;
137  bFaces.append(f);
138  f[0] = 1;
139  f[1] = 11;
140  f[2] = 2;
141  f[3] = 10;
142  bFaces.append(f);
143  //- faces in z direction
144  f[0] = 0;
145  f[1] = 4;
146  f[2] = 1;
147  f[3] = 5;
148  bFaces.append(f);
149  f[0] = 3;
150  f[1] = 7;
151  f[2] = 2;
152  f[3] = 6;
153  bFaces.append(f);
154 
155  # ifdef DEBUGSmooth
156  Info << "Original vertices " << *pointsPtr_ << endl;
157  Info << "Original edges " << *edgesPtr_ << endl;
158  Info << "Original faces " << *facesPtr_ << endl;
159  # endif
160 }
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
165 :
166  pointsPtr_(NULL),
167  edgesPtr_(NULL),
168  facesPtr_(NULL),
169  cPtsPtr_(NULL),
170  cEdgesPtr_(NULL),
171  cFacesPtr_(NULL),
172  newVertexLabel_(),
173  vertexDistance_(),
174  vertexTypes_(),
175  newEdgeLabel_(),
176  origNumVertices_(),
177  tol_(SMALL * bb.mag()),
178  valid_(true)
179 {
181 }
182 
184 {
185  deleteDemandDrivenData(pointsPtr_);
186  deleteDemandDrivenData(edgesPtr_);
187  deleteDemandDrivenData(facesPtr_);
188  deleteDemandDrivenData(cPtsPtr_);
189  deleteDemandDrivenData(cEdgesPtr_);
190  deleteDemandDrivenData(cFacesPtr_);
191 }
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
196 {
197  if( !valid_ )
198  return;
199 
200  # ifdef DEBUGSmooth
201  if( cFacesPtr_ || cPtsPtr_ || cEdgesPtr_ )
202  {
204  (
205  "void meshUntangler::"
206  "cutRegion::planeCut(const plane& plane)"
207  ) << "Pointers should not be allocated!" << abort(FatalError);
208  }
209 
210  Foam::Time runTime
211  (
213  "../.",
214  "testSmoothing"
215  );
216 
217  objectRegistry oR(runTime);
218 
219  polyMeshGen pmg
220  (
221  oR
222  );
223  this->createPolyMeshFromRegion(pmg);
224  # endif
225 
226  if( findNewVertices(plane) )
227  {
228  findNewEdges();
229 
230  findNewFaces();
231 
232  if( !valid_ ) return;
233 
234  deleteDemandDrivenData(pointsPtr_);
235  pointsPtr_ = cPtsPtr_;
236  cPtsPtr_ = NULL;
237 
238  deleteDemandDrivenData(edgesPtr_);
239  edgesPtr_ = cEdgesPtr_;
240  cEdgesPtr_ = NULL;
241 
242  deleteDemandDrivenData(facesPtr_);
243  facesPtr_ = cFacesPtr_;
244  cFacesPtr_ = NULL;
245  }
246 }
247 
249 (
251 ) const
252 {
253  polyMeshGenModifier meshModifier(mesh);
254  pointFieldPMG& points = meshModifier.pointsAccess();
255  points.setSize(pointsPtr_->size());
256  forAll(points, pI)
257  points[pI] = (*pointsPtr_)[pI];
258 
259  faceListPMG& faces = meshModifier.facesAccess();
260  cellListPMG& cells = meshModifier.cellsAccess();
261  cells.setSize(1);
262  cells[0].setSize(facesPtr_->size());
263  faces.setSize(facesPtr_->size());
264 
265  const DynList<edge, 128>& edges = *edgesPtr_;
266  const DynList<DynList<label, 8>, 64>& fcs = *facesPtr_;
267  forAll(faces, fI)
268  {
269  DynList<edge> fEdges;
270  const DynList<label, 8>& f = fcs[fI];
271  forAll(f, eI)
272  fEdges.append(edges[f[eI]]);
273 
274  Info << "Edges forming face " << fI << " are " << fEdges << endl;
276  if( sf.size() != 1 )
278  (
279  "void meshOptimizer::meshUntangler::"
280  "cutRegion::createPolyMeshFromRegion(polyMesgGen&)"
281  ) << "More than one face created!" << abort(FatalError);
282 
283  faces[fI].setSize(sf[0].size());
284  forAll(sf[0], pI)
285  faces[fI][pI] = sf[0][pI];
286 
287  cells[0][fI] = fI;
288  }
289 }
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 } // End namespace Foam
294 
295 // ************************************************************************* //
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::meshUntangler::cutRegion::planeCut
void planeCut(const plane &plane)
cut the region woth the plane
Definition: meshUntanglerCutRegion.C:195
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
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::sortEdgesIntoChains::sortedChains
const DynList< DynList< label > > & sortedChains() const
a list of points which have not yet been resolved
Definition: sortEdgesIntoChains.C:264
polyMeshGenModifier.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
objectRegistry.H
Foam::polyMeshGen
Definition: polyMeshGen.H:46
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::cellListPMG
Definition: cellListPMG.H:49
Foam::plane
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::polyMeshGenModifier::cellsAccess
cellListPMG & cellsAccess()
access to cells
Definition: polyMeshGenModifier.H:119
plane.H
Foam::meshUntangler::cutRegion::createInitialConfiguration
void createInitialConfiguration(const boundBox &)
Definition: meshUntanglerCutRegion.C:50
Foam::Info
messageStream Info
sortEdgesIntoChains.H
Foam::faceListPMG::setSize
void setSize(const label nElmts)
set the number of used elements
Definition: faceListPMGI.H:78
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::FatalError
error FatalError
meshUntangler.H
Foam::polyMeshGenModifier::pointsAccess
pointFieldPMG & pointsAccess()
access to mesh points
Definition: polyMeshGenModifier.H:107
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
sf
volScalarField sf(fieldObject, mesh)
Foam::DynList
Definition: DynList.H:53
Foam::polyMeshGenModifier
Definition: polyMeshGenModifier.H:52
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::sortEdgesIntoChains
Definition: sortEdgesIntoChains.H:49
f
labelList f(nPoints)
Foam::Vector< scalar >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::polyMeshGenModifier::facesAccess
faceListPMG & facesAccess()
access to mesh faces
Definition: polyMeshGenModifier.H:113
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::Time::controlDictName
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:213
Foam::pointFieldPMG
Definition: pointFieldPMG.H:50
Foam::meshUntangler::cutRegion::cutRegion
cutRegion(const boundBox &)
Construct from boundBox.
Definition: meshUntanglerCutRegion.C:164
Foam::meshUntangler::cutRegion::~cutRegion
~cutRegion()
Definition: meshUntanglerCutRegion.C:183
Foam::meshUntangler::cutRegion::createPolyMeshFromRegion
void createPolyMeshFromRegion(polyMeshGen &) const
export the feasible region as polyMeshGen
Definition: meshUntanglerCutRegion.C:249
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304