meshOptimizerI.H
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 
29 #include "meshOptimizer.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 template<class labelListType>
39 inline void meshOptimizer::lockCells(const labelListType& l)
40 {
41  boolList lockedFace(mesh_.faces().size(), false);
42  const cellListPMG& cells = mesh_.cells();
43  forAll(l, lcI)
44  {
45  const cell& c = cells[l[lcI]];
46 
47  forAll(c, fI)
48  lockedFace[c[fI]] = true;
49  }
50 
51  if( Pstream::parRun() )
52  {
53  const PtrList<processorBoundaryPatch>& procBoundaries =
55 
56  forAll(procBoundaries, patchI)
57  {
58  labelLongList dataToSend;
59 
60  const label start = procBoundaries[patchI].patchStart();
61  const label end = start+procBoundaries[patchI].patchSize();
62 
63  for(label faceI=start;faceI<end;++faceI)
64  if( lockedFace[faceI] )
65  dataToSend.append(faceI-start);
66 
67  OPstream toOtherProc
68  (
70  procBoundaries[patchI].neiProcNo(),
71  dataToSend.byteSize()
72  );
73 
74  toOtherProc << dataToSend;
75  }
76 
77  forAll(procBoundaries, patchI)
78  {
79  const label start = procBoundaries[patchI].patchStart();
80 
81  IPstream fromOtherProc
82  (
84  procBoundaries[patchI].neiProcNo()
85  );
86 
87  labelList receivedData;
88  fromOtherProc >> receivedData;
89 
90  forAll(receivedData, i)
91  lockedFace[start+receivedData[i]];
92  }
93  }
94 
95  //- Finally, mark locked points and faces
96  const faceListPMG& faces = mesh_.faces();
97  forAll(lockedFace, faceI)
98  {
99  if( lockedFace[faceI] )
100  {
101  lockedFaces_.append(faceI);
102 
103  const face& f = faces[faceI];
104 
105  forAll(f, pI)
106  vertexLocation_[f[pI]] |= LOCKED;
107  }
108  }
109 
110  # ifdef DEBUGSmoothing
111  const label lockedFacesI = mesh_.addFaceSubset("lockedFaces");
112  forAll(lockedFaces_, lfI)
113  mesh_.addFaceToSubset(lockedFacesI, lockedFaces_[lfI]);
114 
115  const label lockPointsI = mesh_.addPointSubset("lockedPoints");
116  forAll(vertexLocation_, pointI)
117  if( vertexLocation_[pointI] & LOCKED )
118  mesh_.addPointToSubset(lockPointsI, pointI);
119  # endif
120 }
121 
122 template<class labelListType>
123 inline void meshOptimizer::lockFaces(const labelListType& lf)
124 {
125  boolList lockedFace(mesh_.faces().size(), false);
126  forAll(lf, lfI)
127  {
128  lockedFace[lf[lfI]] = true;
129  }
130 
131  if( Pstream::parRun() )
132  {
133  const PtrList<processorBoundaryPatch>& procBoundaries =
135 
136  forAll(procBoundaries, patchI)
137  {
138  labelLongList dataToSend;
139 
140  const label start = procBoundaries[patchI].patchStart();
141  const label end = start+procBoundaries[patchI].patchSize();
142 
143  for(label faceI=start;faceI<end;++faceI)
144  if( lockedFace[faceI] )
145  dataToSend.append(faceI-start);
146 
147  OPstream toOtherProc
148  (
150  procBoundaries[patchI].neiProcNo(),
151  dataToSend.byteSize()
152  );
153 
154  toOtherProc << dataToSend;
155  }
156 
157  forAll(procBoundaries, patchI)
158  {
159  const label start = procBoundaries[patchI].patchStart();
160 
161  IPstream fromOtherProc
162  (
164  procBoundaries[patchI].neiProcNo()
165  );
166 
167  labelList receivedData;
168  fromOtherProc >> receivedData;
169 
170  forAll(receivedData, i)
171  lockedFace[start+receivedData[i]];
172  }
173  }
174 
175  //- Finally, mark locked points and faces
176  const faceListPMG& faces = mesh_.faces();
177  forAll(lockedFace, faceI)
178  {
179  if( lockedFace[faceI] )
180  {
181  lockedFaces_.append(faceI);
182 
183  const face& f = faces[faceI];
184 
185  forAll(f, pI)
186  vertexLocation_[f[pI]] |= LOCKED;
187  }
188  }
189 
190  # ifdef DEBUGSmoothing
191  const label lockedFacesI = mesh_.addFaceSubset("lockedFaces");
192  forAll(lockedFaces_, lfI)
193  mesh_.addFaceToSubset(lockedFacesI, lockedFaces_[lfI]);
194 
195  const label lockPointsI = mesh_.addPointSubset("lockedPoints");
196  forAll(vertexLocation_, pointI)
197  if( vertexLocation_[pointI] & LOCKED )
198  mesh_.addPointToSubset(lockPointsI, pointI);
199  # endif
200 }
201 
202 template<class labelListType>
203 inline void meshOptimizer::lockPoints(const labelListType& lp)
204 {
205  forAll(lp, lpI)
206  vertexLocation_[lp[lpI]] |= LOCKED;
207 
208  if( Pstream::parRun() )
209  {
210  const PtrList<processorBoundaryPatch>& procBoundaries =
212  const faceListPMG& faces = mesh_.faces();
213 
214  forAll(procBoundaries, patchI)
215  {
216  labelLongList dataToSend;
217 
218  const label start = procBoundaries[patchI].patchStart();
219  const label end = start+procBoundaries[patchI].patchSize();
220 
221  for(label faceI=start;faceI<end;++faceI)
222  {
223  const face& f = faces[faceI];
224 
225  forAll(f, pI)
226  {
227  if( vertexLocation_[pI] & LOCKED )
228  {
229  // send the face lable and the location in the face
230  dataToSend.append(faceI-start);
231  dataToSend.append((f.size()-pI)%f.size());
232  }
233  }
234  }
235 
236  OPstream toOtherProc
237  (
239  procBoundaries[patchI].neiProcNo(),
240  dataToSend.byteSize()
241  );
242 
243  toOtherProc << dataToSend;
244  }
245 
246  forAll(procBoundaries, patchI)
247  {
248  const label start = procBoundaries[patchI].patchStart();
249 
250  IPstream fromOtherProc
251  (
253  procBoundaries[patchI].neiProcNo()
254  );
255 
256  labelList receivedData;
257  fromOtherProc >> receivedData;
258 
259  label counter(0);
260  while( counter < receivedData.size() )
261  {
262  const face& f = faces[start+receivedData[counter++]];
263  vertexLocation_[f[receivedData[counter++]]] |= LOCKED;
264  }
265  }
266  }
267 
268  # ifdef DEBUGSmoothing
269  const label lockPointsI = mesh_.addPointSubset("lockedPoints");
270  forAll(vertexLocation_, pointI)
271  if( vertexLocation_[pointI] & LOCKED )
272  mesh_.addPointToSubset(lockPointsI, pointI);
273  # endif
274 }
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Foam::meshOptimizer::lockCells
void lockCells(const labelListType &)
lock the cells which shall not be modified
Definition: meshOptimizerI.H:39
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:50
Foam::meshOptimizer::vertexLocation_
List< direction > vertexLocation_
location of vertex (internal, boundary, edge, corner)
Definition: meshOptimizer.H:67
Foam::meshOptimizer::lockPoints
void lockPoints(const labelListType &)
lock points which shall not be moved
Definition: meshOptimizerI.H:203
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::polyMeshGenFaces::addFaceSubset
label addFaceSubset(const word &)
Definition: polyMeshGenFaces.C:262
Foam::cellListPMG
Definition: cellListPMG.H:49
Foam::polyMeshGenFaces::faces
const faceListPMG & faces() const
access to faces
Definition: polyMeshGenFacesI.H:43
Foam::meshOptimizer::lockFaces
void lockFaces(const labelListType &)
lock the faces whih shall not be modified
Definition: meshOptimizerI.H:123
Foam::LongList< label >
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
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::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::meshOptimizer::lockedFaces_
labelLongList lockedFaces_
locked faces which shall not be changed
Definition: meshOptimizer.H:70
Foam::polyMeshGenPoints::addPointToSubset
void addPointToSubset(const label, const label)
Definition: polyMeshGenPointsI.H:60
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::meshOptimizer::mesh_
polyMeshGen & mesh_
reference to the mesh
Definition: meshOptimizer.H:64
Foam::meshOptimizer::LOCKED
@ LOCKED
Definition: meshOptimizer.H:198
Foam::LongList::byteSize
label byteSize() const
Return the binary size in number of characters of the UList.
Definition: LongListI.H:209
f
labelList f(nPoints)
Foam::faceListPMG::size
label size() const
return the number of used elements
Definition: faceListPMGI.H:73
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
meshOptimizer.H
Foam::polyMeshGenFaces::procBoundaries
const PtrList< processorBoundaryPatch > & procBoundaries() const
inter-processor boundaries
Definition: polyMeshGenFacesI.H:106
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::polyMeshGenFaces::addFaceToSubset
void addFaceToSubset(const label, const label)
Definition: polyMeshGenFacesI.H:117
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::polyMeshGenPoints::addPointSubset
label addPointSubset(const word &)
point subsets
Definition: polyMeshGenPoints.C:88
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56