polyMeshInitMesh.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "polyMesh.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::polyMesh::initMesh()
34 {
35  if (debug)
36  {
37  InfoInFunction << "initialising primitiveMesh" << endl;
38  }
39 
40  // For backward compatibility check if the neighbour array is the same
41  // length as the owner and shrink to remove the -1s padding
42  if (neighbour_.size() == owner_.size())
43  {
45 
46  forAll(neighbour_, facei)
47  {
48  if (neighbour_[facei] == -1)
49  {
50  break;
51  }
52  else
53  {
55  }
56  }
57 
58  neighbour_.setSize(nInternalFaces);
59  }
60 
61  label nCells = -1;
62 
63  forAll(owner_, facei)
64  {
65  if (owner_[facei] < 0)
66  {
68  << "Illegal cell label " << owner_[facei]
69  << " in neighbour addressing for face " << facei
70  << exit(FatalError);
71  }
72  nCells = max(nCells, owner_[facei]);
73  }
74 
75  // The neighbour array may or may not be the same length as the owner
76  forAll(neighbour_, facei)
77  {
78  if (neighbour_[facei] < 0)
79  {
81  << "Illegal cell label " << neighbour_[facei]
82  << " in neighbour addressing for face " << facei
83  << exit(FatalError);
84  }
85  nCells = max(nCells, neighbour_[facei]);
86  }
87 
88  nCells++;
89 
90  // Reset the primitiveMesh with the sizes of the primitive arrays
92  (
93  points_.size(),
94  neighbour_.size(),
95  owner_.size(),
96  nCells
97  );
98 
99  const string meshInfo
100  (
101  "nPoints:" + Foam::name(nPoints())
102  + " nCells:" + Foam::name(this->nCells())
103  + " nFaces:" + Foam::name(nFaces())
104  + " nInternalFaces:" + Foam::name(nInternalFaces())
105  );
106 
107  owner_.note() = meshInfo;
108  neighbour_.note() = meshInfo;
109 }
110 
111 
112 void Foam::polyMesh::initMesh(cellList& c)
113 {
114  if (debug)
115  {
116  InfoInFunction << "Calculating owner-neighbour arrays" << endl;
117  }
118 
119  owner_.setSize(faces_.size(), -1);
120  neighbour_.setSize(faces_.size(), -1);
121 
122  boolList markedFaces(faces_.size(), false);
123 
124  label nInternalFaces = 0;
125 
126  forAll(c, celli)
127  {
128  // get reference to face labels for current cell
129  const labelList& cellfaces = c[celli];
130 
131  forAll(cellfaces, facei)
132  {
133  if (cellfaces[facei] < 0)
134  {
136  << "Illegal face label " << cellfaces[facei]
137  << " in cell " << celli
138  << exit(FatalError);
139  }
140 
141  if (!markedFaces[cellfaces[facei]])
142  {
143  // First visit: owner
144  owner_[cellfaces[facei]] = celli;
145  markedFaces[cellfaces[facei]] = true;
146  }
147  else
148  {
149  // Second visit: neighbour
150  neighbour_[cellfaces[facei]] = celli;
151  nInternalFaces++;
152  }
153  }
154  }
155 
156  // The neighbour array is initialised with the same length as the owner
157  // padded with -1s and here it is truncated to the correct size of the
158  // number of internal faces.
159  neighbour_.setSize(nInternalFaces);
160 
161  // Reset the primitiveMesh
163  (
164  points_.size(),
165  neighbour_.size(),
166  owner_.size(),
167  c.size(),
168  c
169  );
170 
171  const string meshInfo
172  (
173  "nPoints:" + Foam::name(nPoints())
174  + " nCells:" + Foam::name(nCells())
175  + " nFaces:" + Foam::name(nFaces())
176  + " nInternalFaces:" + Foam::name(this->nInternalFaces())
177  );
178 
179  owner_.note() = meshInfo;
180  neighbour_.note() = meshInfo;
181 }
182 
183 
184 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:72
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
polyMesh.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
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:62
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::FatalError
error FatalError
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::IOobject::note
const string & note() const
Return the optional note.
Definition: IOobjectI.H:76
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
Foam::primitiveMesh::reset
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
Definition: primitiveMesh.C:208
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146