polyMeshGenModifierReplaceBoundary.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 "polyMeshGenModifier.H"
29 #include "demandDrivenData.H"
30 #include "DynList.H"
31 
32 // #define DEBUGSearch
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
42 (
43  const wordList& patchNames,
44  const VRWGraph& boundaryFaces,
45  const labelLongList& faceOwners,
46  const labelLongList& facePatches
47 )
48 {
49  const label nIntFaces = mesh_.nInternalFaces();
50 
51  faceListPMG& faces = this->facesAccess();
52  cellListPMG& cells = this->cellsAccess();
53 
54  labelLongList newFaceLabel(faces.size(), -1);
55  for(label faceI=0;faceI<nIntFaces;++faceI)
56  newFaceLabel[faceI] = faceI;
57 
58  if( Pstream::parRun() )
59  {
60  //- shift processor faces
61  PtrList<processorBoundaryPatch>& procBoundaries =
62  mesh_.procBoundaries_;
63 
64  label nProcFaces(0);
65  forAll(procBoundaries, patchI)
66  nProcFaces += procBoundaries[patchI].patchSize();
67 
68  const label procStart = nIntFaces + boundaryFaces.size();
69  const label shift = procStart - procBoundaries[0].patchStart();
70 
71  if( shift > 0 )
72  {
73  faces.setSize(procStart+nProcFaces);
74  forAllReverse(procBoundaries, patchI)
75  {
76  const label start = procBoundaries[patchI].patchStart();
77  const label end = procBoundaries[patchI].patchSize() + start;
78  for(label faceI=end-1;faceI>=start;--faceI)
79  {
80  faces[faceI+shift].transfer(faces[faceI]);
81  newFaceLabel[faceI] = faceI+shift;
82  }
83 
84  procBoundaries[patchI].patchStart() += shift;
85  }
86  }
87  else if( shift < 0 )
88  {
89  forAll(procBoundaries, patchI)
90  {
91  const label start = procBoundaries[patchI].patchStart();
92  const label end = procBoundaries[patchI].patchSize() + start;
93  for(label faceI=start;faceI<end;++faceI)
94  {
95  faces[faceI+shift].transfer(faces[faceI]);
96  newFaceLabel[faceI] = faceI+shift;
97  }
98 
99  procBoundaries[patchI].patchStart() += shift;
100  }
101 
102  faces.setSize(procStart+nProcFaces);
103  }
104  else
105  {
106  //- processor faces are not moved
107  for(label fI=0;fI<nProcFaces;++fI)
108  newFaceLabel[procStart+fI] = procStart+fI;
109  }
110  }
111  else
112  {
113  faces.setSize(nIntFaces + boundaryFaces.size());
114  }
115 
116  //- change cells according to the new face ordering
117  List<direction> nFacesInCell(cells.size(), direction(0));
118  forAll(cells, cellI)
119  {
120  const cell& c = cells[cellI];
121 
122  cell newC(c.size());
123  forAll(c, fI)
124  if( newFaceLabel[c[fI]] != -1 )
125  newC[nFacesInCell[cellI]++] = newFaceLabel[c[fI]];
126 
127  cells[cellI].transfer(newC);
128  }
129 
130  mesh_.updateFaceSubsets(newFaceLabel);
131  newFaceLabel.setSize(0);
132 
133  //- store boundary faces
134  labelList newPatchStart(patchNames.size());
135  labelList newPatchSize(patchNames.size(), 0);
136  forAll(facePatches, bfI)
137  ++newPatchSize[facePatches[bfI]];
138 
139  newPatchStart[0] = nIntFaces;
140  for(label i=1;i<newPatchSize.size();++i)
141  newPatchStart[i] = newPatchStart[i-1] + newPatchSize[i-1];
142 
143  //- store boundary faces
144  newPatchSize = 0;
145  forAll(boundaryFaces, faceI)
146  {
147  const label fPatch = facePatches[faceI];
148  const label fOwn = faceOwners[faceI];
149 
150  const label fLabel = newPatchStart[fPatch] + newPatchSize[fPatch]++;
151 
152  cells[fOwn].newElmt(nFacesInCell[fOwn]++) = fLabel;
153  faces[fLabel].setSize(boundaryFaces.sizeOfRow(faceI));
154  forAllRow(boundaryFaces, faceI, pI)
155  faces[fLabel][pI] = boundaryFaces(faceI, pI);
156  }
157 
158  forAll(cells, cellI)
159  cells[cellI].setSize(nFacesInCell[cellI]);
160 
161  PtrList<boundaryPatch>& boundaries = mesh_.boundaries_;
162  if( boundaries.size() == patchNames.size() )
163  {
164  forAll(boundaries, patchI)
165  {
166  boundaries[patchI].patchName() = patchNames[patchI];
167  boundaries[patchI].patchStart() = newPatchStart[patchI];
168  boundaries[patchI].patchSize() = newPatchSize[patchI];
169  }
170  }
171  else
172  {
173  boundaries.clear();
174  boundaries.setSize(patchNames.size());
175  forAll(boundaries, patchI)
176  boundaries.set
177  (
178  patchI,
179  new boundaryPatch
180  (
181  patchNames[patchI],
182  "patch",
183  newPatchSize[patchI],
184  newPatchStart[patchI]
185  )
186  );
187  }
188 
189  mesh_.clearOut();
190  this->clearOut();
191 }
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 } // End namespace Foam
196 
197 // ************************************************************************* //
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::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::List::newElmt
T & newElmt(const label)
Return subscript-checked element of UList.
Definition: ListI.H:64
polyMeshGenModifier.H
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Foam::polyMeshGenModifier::replaceBoundary
void replaceBoundary(const wordList &patchNames, const VRWGraph &boundaryFaces, const labelLongList &faceOwners, const labelLongList &facePatches)
replace the boundary with new boundary faces
Definition: polyMeshGenModifierReplaceBoundary.C:42
Foam::cellListPMG
Definition: cellListPMG.H:49
Foam::LongList::setSize
void setSize(const label)
Reset size of List.
Definition: LongListI.H:223
Foam::boundaryPatch
Like polyPatch but without reference to mesh. patchIdentifier::index is not used. Used in boundaryMes...
Definition: boundaryPatch.H:50
Foam::LongList< label >
Foam::PtrList::set
bool set(const label) const
Is element set.
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::VRWGraph::size
label size() const
Returns the number of rows.
Definition: VRWGraphI.H:122
Foam::faceListPMG::setSize
void setSize(const label nElmts)
set the number of used elements
Definition: faceListPMGI.H:78
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
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
patchNames
wordList patchNames(nPatches)
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:418
Foam::PtrList::clear
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:185
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::List::setSize
void setSize(const label)
Reset size of List.
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
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::faceListPMG::transfer
void transfer(faceList &)
DynList.H