patchEdgeFaceRegionsI.H
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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2020 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 {}
35 
36 
38 (
39  const labelList& regions
40 )
41 :
42  regions_(regions)
43 {}
44 
45 
47 (
48  const labelPair& regions
49 )
50 :
51  regions_(2)
52 {
53  regions_[0] = regions[0];
54  regions_[1] = regions[1];
55 }
56 
57 
58 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 
60 template<class TrackingData>
61 inline bool Foam::patchEdgeFaceRegions::valid(TrackingData& td) const
62 {
63  return regions_.size() && !regions_.found(labelMax);
64 }
65 
66 
67 template<class Patch, class TrackingData>
69 (
70  const polyMesh& mesh,
71  const Patch& patch,
72  const tensor& rotTensor,
73  const scalar tol,
74  TrackingData& td
75 )
76 {}
77 
78 
79 template<class Patch, class TrackingData>
81 (
82  const polyMesh& mesh,
83  const Patch& patch,
84  const label edgeI,
85  const label facei,
86  const patchEdgeFaceRegions& faceInfo,
87  const scalar tol,
88  TrackingData& td
89 )
90 {
91  const face& f = patch.localFaces()[facei];
92  const edge& e = patch.edges()[edgeI];
93 
94  label index = patch.faceEdges()[facei].find(edgeI);
95  bool sameOrientation = (f[index] == e.start());
96 
97  // Get information in edge-order
98  edge orientedInfo
99  (
100  faceInfo.regions_[index],
101  faceInfo.regions_[f.fcIndex(index)]
102  );
103  if (!sameOrientation)
104  {
105  orientedInfo.flip();
106  }
107 
108  if (!faceInfo.valid(td))
109  {
111  << "problem." << abort(FatalError);
112  }
113 
114  if (orientedInfo.found(-1) || regions_.found(-1))
115  {
116  // Blocked edge/face
117  return false;
118  }
119 
120 
121  bool changed = false;
122 
123  regions_.setSize(orientedInfo.size(), labelMax);
124  forAll(orientedInfo, i)
125  {
126  if (orientedInfo[i] != -1 && orientedInfo[i] < regions_[i])
127  {
128  regions_[i] = orientedInfo[i];
129  changed = true;
130  }
131  }
132  return changed;
133 }
134 
135 
136 template<class Patch, class TrackingData>
138 (
139  const polyMesh& mesh,
140  const Patch& patch,
141  const patchEdgeFaceRegions& edgeInfo,
142  const bool sameOrientation,
143  const scalar tol,
144  TrackingData& td
145 )
146 {
147  // Get information in edge-order
148  edge orientedInfo(edgeInfo.regions_[0], edgeInfo.regions_[1]);
149  if (!sameOrientation)
150  {
151  orientedInfo.flip();
152  }
153 
154 
155  if (!edgeInfo.valid(td))
156  {
158  << "problem." << abort(FatalError);
159  }
160 
161  if (orientedInfo.found(-1) || regions_.found(-1))
162  {
163  // Blocked edge/face
164  return false;
165  }
166 
167 
168  bool changed = false;
169 
170  regions_.setSize(orientedInfo.size(), labelMax);
171  forAll(orientedInfo, i)
172  {
173  if (orientedInfo[i] != -1 && orientedInfo[i] < regions_[i])
174  {
175  regions_[i] = orientedInfo[i];
176  changed = true;
177  }
178  }
179  return changed;
180 }
181 
182 
183 template<class Patch, class TrackingData>
185 (
186  const polyMesh& mesh,
187  const Patch& patch,
188  const label facei,
189  const label edgeI,
190  const patchEdgeFaceRegions& edgeInfo,
191  const scalar tol,
192  TrackingData& td
193 )
194 {
195  const face& f = patch.localFaces()[facei];
196  const edge& e = patch.edges()[edgeI];
197 
198  // Find starting point of edge on face.
199  label index0 = patch.faceEdges()[facei].find(edgeI);
200  label index1 = f.fcIndex(index0);
201  bool sameOrientation = (f[index0] == e.start());
202 
203 
204  // Get information in face-order
205  edge orientedInfo
206  (
207  edgeInfo.regions_[0],
208  edgeInfo.regions_[1]
209  );
210  if (!sameOrientation)
211  {
212  orientedInfo.flip();
213  }
214 
215  if (!edgeInfo.valid(td))
216  {
218  << "problem." << abort(FatalError);
219  }
220 
221  if (orientedInfo.found(-1) || regions_.found(-1))
222  {
223  // Blocked edge/face
224  return false;
225  }
226 
227 
228  bool changed = false;
229 
230  // Scale if needed
231  regions_.setSize(f.size(), labelMax);
232 
233  if (orientedInfo[0] < regions_[index0])
234  {
235  regions_[index0] = orientedInfo[0];
236  changed = true;
237  }
238  if (orientedInfo[1] < regions_[index1])
239  {
240  regions_[index1] = orientedInfo[1];
241  changed = true;
242  }
243 
244  return changed;
245 }
246 
247 
248 template<class TrackingData>
250 (
251  const patchEdgeFaceRegions& rhs,
252  TrackingData& td
253 ) const
254 {
255  return operator==(rhs);
256 }
257 
258 
259 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
260 
261 inline bool Foam::patchEdgeFaceRegions::operator==
262 (
263  const patchEdgeFaceRegions& rhs
264 ) const
265 {
266  return regions_ == rhs.regions_;
267 }
268 
269 
270 inline bool Foam::patchEdgeFaceRegions::operator!=
271 (
272  const patchEdgeFaceRegions& rhs
273 ) const
274 {
275  return !(*this == rhs);
276 }
277 
278 
279 // ************************************************************************* //
Foam::Tensor< scalar >
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
Foam::patchEdgeFaceRegions::updateFace
bool updateFace(const polyMesh &mesh, const Patch &patch, const label facei, const label edgeI, const patchEdgeFaceRegions &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on face.
Definition: patchEdgeFaceRegionsI.H:185
Foam::patchEdgeFaceRegions::transform
void transform(const polyMesh &mesh, const Patch &patch, const tensor &rotTensor, const scalar tol, TrackingData &td)
Apply rotation matrix.
Definition: patchEdgeFaceRegionsI.H:69
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::patchEdgeFaceRegions::updateEdge
bool updateEdge(const polyMesh &mesh, const Patch &patch, const label edgeI, const label facei, const patchEdgeFaceRegions &faceInfo, const scalar tol, TrackingData &td)
Influence of face on edge.
Definition: patchEdgeFaceRegionsI.H:81
Foam::patchEdgeFaceRegions::valid
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
Definition: patchEdgeFaceRegionsI.H:61
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::patchEdgeFaceRegions
Transport of regions for use in PatchEdgeFaceWave.
Definition: patchEdgeFaceRegions.H:63
Foam::edge::found
bool found(const label pointLabel) const
Return true if point label is found in edge.
Definition: edgeI.H:136
Foam::patchEdgeFaceRegions::equal
bool equal(const patchEdgeFaceRegions &, TrackingData &) const
Test for equality, with TrackingData.
Definition: patchEdgeFaceRegionsI.H:250
Foam::FatalError
error FatalError
Foam::Pair::flip
void flip()
Flip the Pair in-place.
Definition: PairI.H:160
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pair< label >
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::patchEdgeFaceRegions::patchEdgeFaceRegions
patchEdgeFaceRegions()
Default construct.
Definition: patchEdgeFaceRegionsI.H:33
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72