layerInfoI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 "layerInfo.H"
27 #include "polyMesh.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 inline void Foam::layerInfo::collide() const
32 {
34  << "Layer extrusions collided. Check the patches/zones from which "
35  << "layers are being extruded and ensure that they do not point "
36  << "in opposite directions."
37  << exit(FatalError);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 :
45  layer_(-labelMax),
46  direction_(0),
47  prevFace_(-labelMax)
48 {}
49 
50 
52 :
53  layer_(2*faceLayer - 1),
54  direction_(direction),
55  prevFace_(-labelMax)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
62 {
63  if (layer_ % 2 != 1)
64  {
66  << "Face layer index requested from cell layer info"
67  << exit(FatalError);
68  }
69 
70  return (layer_ + 1)/2;
71 }
72 
73 
75 {
76  if (layer_ % 2 != 0)
77  {
79  << "Cell layer index requested from face layer info"
80  << exit(FatalError);
81  }
82 
83  return layer_/2;
84 }
85 
86 
87 template<class TrackingData>
88 inline bool Foam::layerInfo::valid(TrackingData& td) const
89 {
90  return layer_ != -labelMax;
91 }
92 
93 
94 template<class TrackingData>
96 (
97  const polyMesh&,
98  const layerInfo& l,
99  const scalar tol,
100  TrackingData& td
101 ) const
102 {
103  return layer_ == l.layer_;
104 }
105 
106 
107 template<class TrackingData>
108 inline void Foam::layerInfo::transform
109 (
110  const polyPatch& patch,
111  const label patchFacei,
112  const transformer& transform,
113  TrackingData& td
114 )
115 {}
116 
117 
118 template<class TrackingData>
119 inline bool Foam::layerInfo::updateCell
120 (
121  const polyMesh& mesh,
122  const label thisCelli,
123  const label neighbourFacei,
124  const layerInfo& neighbourLayerInfo,
125  const scalar tol,
126  TrackingData& td
127 )
128 {
129  const bool o = thisCelli == mesh.faceOwner()[neighbourFacei];
130 
131  if (o == (neighbourLayerInfo.direction_ < 0))
132  {
133  if (valid(td) && prevFace_ != neighbourFacei) collide();
134 
135  layer_ = valid(td) ? -labelMax : neighbourLayerInfo.layer_ + 1;
136  direction_ = 0;
137  prevFace_ = neighbourFacei;
138 
139  return true;
140  }
141  else
142  {
143  return false;
144  }
145 }
146 
147 
148 template<class TrackingData>
149 inline bool Foam::layerInfo::updateFace
150 (
151  const polyMesh& mesh,
152  const label thisFacei,
153  const label neighbourCelli,
154  const layerInfo& neighbourLayerInfo,
155  const scalar tol,
156  TrackingData& td
157 )
158 {
159  const cell& c = mesh.cells()[neighbourCelli];
160  const label prevFacei = neighbourLayerInfo.prevFace_;
161  const label nextFacei = c.opposingFaceLabel(prevFacei, mesh.faces());
162 
163  if (nextFacei == thisFacei)
164  {
165  const label direction =
166  mesh.faceOwner()[thisFacei] == neighbourCelli ? +1 : -1;
167 
168  if (valid(td) && direction != direction_) collide();
169 
170  layer_ = valid(td) ? -labelMax : neighbourLayerInfo.layer_ + 1;
171  direction_ = direction;
172  prevFace_ = -labelMax;
173 
174  return true;
175  }
176  else
177  {
178  return false;
179  }
180 }
181 
182 
183 template<class TrackingData>
184 inline bool Foam::layerInfo::updateFace
185 (
186  const polyMesh& mesh,
187  const label thisFacei,
188  const layerInfo& neighbourLayerInfo,
189  const scalar tol,
190  TrackingData& td
191 )
192 {
193  const label direction = -neighbourLayerInfo.direction_;
194 
195  if (valid(td) && direction != direction_) collide();
196 
197  layer_ = valid(td) ? -labelMax : neighbourLayerInfo.layer_;
198  direction_ = direction;
199  prevFace_ = -labelMax;
200 
201  return true;
202 }
203 
204 
205 template<class TrackingData>
206 inline bool Foam::layerInfo::equal
207 (
208  const layerInfo& rhs,
209  TrackingData& td
210 ) const
211 {
212  return operator==(rhs);
213 }
214 
215 
216 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
217 
218 inline bool Foam::layerInfo::operator==(const Foam::layerInfo& rhs) const
219 {
220  return layer_ == rhs.layer_;
221 }
222 
223 
224 inline bool Foam::layerInfo::operator!=(const Foam::layerInfo& rhs) const
225 {
226  return !(*this == rhs);
227 }
228 
229 
230 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
231 
233 {
234  return os << l.layer_ << token::SPACE << l.direction_;
235 }
236 
237 
239 {
240  return is >> l.layer_ >> l.direction_;
241 }
242 
243 
244 // ************************************************************************* //
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool sameGeometry(const polyMesh &, const layerInfo &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for checking.
Definition: layerInfoI.H:96
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label faceLayer() const
Return the face layer index.
Definition: layerInfoI.H:61
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Definition: layerInfoI.H:109
uint8_t direction
Definition: direction.H:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
layerInfo()
Construct null.
Definition: layerInfoI.H:43
label cellLayer() const
Return the face layer index.
Definition: layerInfoI.H:74
const cellList & cells() const
const dimensionedScalar c
Speed of light in a vacuum.
label opposingFaceLabel(const label masterFaceLabel, const faceUList &meshFaces) const
Return index of opposite face.
bool valid(TrackingData &td) const
Check whether the layerInfo has been changed at all or still.
Definition: layerInfoI.H:88
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const layerInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
Definition: layerInfoI.H:120
Istream & operator>>(Istream &, directionInfo &)
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
Class to be used with FaceCellWave which enumerates layers of cells.
Definition: layerInfo.H:59
static const label labelMax
Definition: label.H:62
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
bool equal(const layerInfo &, TrackingData &td) const
Test equality.
Definition: layerInfoI.H:207
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
Ostream & operator<<(Ostream &, const ensightPart &)
bool operator!=(const layerInfo &) const
Definition: layerInfoI.H:224
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const layerInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Definition: layerInfoI.H:150
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
bool operator==(const layerInfo &) const
Definition: layerInfoI.H:218