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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM 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 OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMesh.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
31 {
32  if (debug)
33  {
34  Info<< "void polyMesh::initMesh() : "
35  << "initialising primitiveMesh" << endl;
36  }
37 
38  // For backward compatibility check if the neighbour array is the same
39  // length as the owner and shrink to remove the -1s padding
40  if (neighbour_.size() == owner_.size())
41  {
43 
44  forAll(neighbour_, faceI)
45  {
46  if (neighbour_[faceI] == -1)
47  {
48  break;
49  }
50  else
51  {
53  }
54  }
55 
57  }
58 
59  label nCells = -1;
60 
61  forAll(owner_, facei)
62  {
63  if (owner_[facei] < 0)
64  {
66  << "Illegal cell label " << owner_[facei]
67  << " in neighbour addressing for face " << facei
68  << exit(FatalError);
69  }
70  nCells = max(nCells, owner_[facei]);
71  }
72 
73  // The neighbour array may or may not be the same length as the owner
74  forAll(neighbour_, facei)
75  {
76  if (neighbour_[facei] < 0)
77  {
79  << "Illegal cell label " << neighbour_[facei]
80  << " in neighbour addressing for face " << facei
81  << exit(FatalError);
82  }
83  nCells = max(nCells, neighbour_[facei]);
84  }
85 
86  nCells++;
87 
88  // Reset the primitiveMesh with the sizes of the primitive arrays
90  (
91  points_.size(),
92  neighbour_.size(),
93  owner_.size(),
94  nCells
95  );
96 
97  string meshInfo =
98  "nPoints:" + Foam::name(nPoints())
99  + " nCells:" + Foam::name(this->nCells())
100  + " nFaces:" + Foam::name(nFaces())
101  + " nInternalFaces:" + Foam::name(nInternalFaces());
102 
103  owner_.note() = meshInfo;
104  neighbour_.note() = meshInfo;
105 }
106 
107 
109 {
110  if (debug)
111  {
112  Info<< "void polyMesh::initMesh(cellList& c) : "
113  << "calculating owner-neighbour arrays" << endl;
114  }
115 
116  owner_.setSize(faces_.size(), -1);
117  neighbour_.setSize(faces_.size(), -1);
118 
119  boolList markedFaces(faces_.size(), false);
120 
121  label nInternalFaces = 0;
122 
123  forAll(c, cellI)
124  {
125  // get reference to face labels for current cell
126  const labelList& cellfaces = c[cellI];
127 
128  forAll(cellfaces, faceI)
129  {
130  if (cellfaces[faceI] < 0)
131  {
133  << "Illegal face label " << cellfaces[faceI]
134  << " in cell " << cellI
135  << exit(FatalError);
136  }
137 
138  if (!markedFaces[cellfaces[faceI]])
139  {
140  // First visit: owner
141  owner_[cellfaces[faceI]] = cellI;
142  markedFaces[cellfaces[faceI]] = true;
143  }
144  else
145  {
146  // Second visit: neighbour
147  neighbour_[cellfaces[faceI]] = cellI;
148  nInternalFaces++;
149  }
150  }
151  }
152 
153  // The neighbour array is initialised with the same length as the owner
154  // padded with -1s and here it is truncated to the correct size of the
155  // number of internal faces.
156  neighbour_.setSize(nInternalFaces);
157 
158  // Reset the primitiveMesh
160  (
161  points_.size(),
162  neighbour_.size(),
163  owner_.size(),
164  c.size(),
165  c
166  );
167 
168  string meshInfo =
169  "nPoints: " + Foam::name(nPoints())
170  + " nCells: " + Foam::name(nCells())
171  + " nFaces: " + Foam::name(nFaces())
172  + " nInternalFaces: " + Foam::name(this->nInternalFaces());
173 
174  owner_.note() = meshInfo;
175  neighbour_.note() = meshInfo;
176 }
177 
178 // ************************************************************************* //
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::polyMesh::neighbour_
labelIOList neighbour_
Face neighbour.
Definition: polyMesh.H:127
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyMesh.H
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::polyMesh::owner_
labelIOList owner_
Face owner.
Definition: polyMesh.H:124
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
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::FatalError
error FatalError
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
Foam::polyMesh::points_
pointIOField points_
Points.
Definition: polyMesh.H:118
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
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::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:121
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::IOobject::note
string & note()
Return non-constant access to the optional note.
Definition: IOobject.H:285
Foam::polyMesh::initMesh
void initMesh()
Initialise the polyMesh from the primitive data.
Definition: polyMeshInitMesh.C:30