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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "polyMesh.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 {}
34 
35 
37 (
38  const labelList& regions
39 )
40 :
41  regions_(regions)
42 {}
43 
44 
46 (
47  const labelPair& regions
48 )
49 :
50  regions_(2)
51 {
52  regions_[0] = regions[0];
53  regions_[1] = regions[1];
54 }
55 
56 
57 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58 
60 {
61  return regions_;
62 }
63 
64 
65 template<class TrackingData>
66 inline bool Foam::patchEdgeFaceRegions::valid(TrackingData& td) const
67 {
68  return regions_.size() && !regions_.found(labelMax);
69 }
70 
71 
72 template<class Patch, class TrackingData>
74 (
75  const polyMesh& mesh,
76  const Patch& patch,
77  const tensor& rotTensor,
78  const scalar tol,
79  TrackingData& td
80 )
81 {}
82 
83 
84 template<class Patch, class TrackingData>
86 (
87  const polyMesh& mesh,
88  const Patch& patch,
89  const label edgeI,
90  const label facei,
91  const patchEdgeFaceRegions& faceInfo,
92  const scalar tol,
93  TrackingData& td
94 )
95 {
96  const face& f = patch.localFaces()[facei];
97  const edge& e = patch.edges()[edgeI];
98 
99  label index = patch.faceEdges()[facei].find(edgeI);
100  bool sameOrientation = (f[index] == e.start());
101 
102  // Get information in edge-order
103  edge orientedInfo
104  (
105  faceInfo.regions_[index],
106  faceInfo.regions_[f.fcIndex(index)]
107  );
108  if (!sameOrientation)
109  {
110  orientedInfo.flip();
111  }
112 
113  if (!faceInfo.valid(td))
114  {
116  << "problem." << abort(FatalError);
117  }
118 
119  if (orientedInfo.found(-1) || regions_.found(-1))
120  {
121  // Blocked edge/face
122  return false;
123  }
124 
125 
126  bool changed = false;
127 
128  regions_.setSize(orientedInfo.size(), labelMax);
129  forAll(orientedInfo, i)
130  {
131  if (orientedInfo[i] != -1 && orientedInfo[i] < regions_[i])
132  {
133  regions_[i] = orientedInfo[i];
134  changed = true;
135  }
136  }
137  return changed;
138 }
139 
140 
141 template<class Patch, class TrackingData>
143 (
144  const polyMesh& mesh,
145  const Patch& patch,
146  const patchEdgeFaceRegions& edgeInfo,
147  const bool sameOrientation,
148  const scalar tol,
149  TrackingData& td
150 )
151 {
152  // Get information in edge-order
153  edge orientedInfo(edgeInfo.regions_[0], edgeInfo.regions_[1]);
154  if (!sameOrientation)
155  {
156  orientedInfo.flip();
157  }
158 
159 
160  if (!edgeInfo.valid(td))
161  {
163  << "problem." << abort(FatalError);
164  }
165 
166  if (orientedInfo.found(-1) || regions_.found(-1))
167  {
168  // Blocked edge/face
169  return false;
170  }
171 
172 
173  bool changed = false;
174 
175  regions_.setSize(orientedInfo.size(), labelMax);
176  forAll(orientedInfo, i)
177  {
178  if (orientedInfo[i] != -1 && orientedInfo[i] < regions_[i])
179  {
180  regions_[i] = orientedInfo[i];
181  changed = true;
182  }
183  }
184  return changed;
185 }
186 
187 
188 template<class Patch, class TrackingData>
190 (
191  const polyMesh& mesh,
192  const Patch& patch,
193  const label facei,
194  const label edgeI,
195  const patchEdgeFaceRegions& edgeInfo,
196  const scalar tol,
197  TrackingData& td
198 )
199 {
200  const face& f = patch.localFaces()[facei];
201  const edge& e = patch.edges()[edgeI];
202 
203  // Find starting point of edge on face.
204  label index0 = patch.faceEdges()[facei].find(edgeI);
205  label index1 = f.fcIndex(index0);
206  bool sameOrientation = (f[index0] == e.start());
207 
208 
209  // Get information in face-order
210  edge orientedInfo
211  (
212  edgeInfo.regions_[0],
213  edgeInfo.regions_[1]
214  );
215  if (!sameOrientation)
216  {
217  orientedInfo.flip();
218  }
219 
220  if (!edgeInfo.valid(td))
221  {
223  << "problem." << abort(FatalError);
224  }
225 
226  if (orientedInfo.found(-1) || regions_.found(-1))
227  {
228  // Blocked edge/face
229  return false;
230  }
231 
232 
233  bool changed = false;
234 
235  // Scale if needed
236  regions_.setSize(f.size(), labelMax);
237 
238  if (orientedInfo[0] < regions_[index0])
239  {
240  regions_[index0] = orientedInfo[0];
241  changed = true;
242  }
243  if (orientedInfo[1] < regions_[index1])
244  {
245  regions_[index1] = orientedInfo[1];
246  changed = true;
247  }
248 
249  return changed;
250 }
251 
252 
253 template<class TrackingData>
255 (
256  const patchEdgeFaceRegions& rhs,
257  TrackingData& td
258 ) const
259 {
260  return operator==(rhs);
261 }
262 
263 
264 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
265 
266 inline bool Foam::patchEdgeFaceRegions::operator==
267 (
268  const Foam::patchEdgeFaceRegions& rhs
269 ) const
270 {
271  return regions() == rhs.regions();
272 }
273 
274 
275 inline bool Foam::patchEdgeFaceRegions::operator!=
276 (
277  const Foam::patchEdgeFaceRegions& rhs
278 ) const
279 {
280  return !(*this == rhs);
281 }
282 
283 
284 // ************************************************************************* //
Foam::Tensor< scalar >
Foam::labelMax
constexpr label labelMax
Definition: label.H:65
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:190
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:74
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:290
Foam::patchEdgeFaceRegions::regions
const labelList & regions() const
Definition: patchEdgeFaceRegionsI.H:59
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:86
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::patchEdgeFaceRegions::valid
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Definition: patchEdgeFaceRegionsI.H:66
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::patchEdgeFaceRegions
Transport of regions for use in PatchEdgeFaceWave.
Definition: patchEdgeFaceRegions.H:68
Foam::edge::found
bool found(const label pointLabel) const
Return true if point label is found in edge.
Definition: edgeI.H:126
Foam::patchEdgeFaceRegions::equal
bool equal(const patchEdgeFaceRegions &, TrackingData &) const
Same (like operator==)
Definition: patchEdgeFaceRegionsI.H:255
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:137
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
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()
Construct null.
Definition: patchEdgeFaceRegionsI.H:32
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:74