patchEdgeFaceRegionsI.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) 2013-2018 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 {}
32 
33 
35 (
36  const labelList& regions
37 )
38 :
39  regions_(regions)
40 {}
41 
42 
44 (
45  const labelPair& regions
46 )
47 :
48  regions_(2)
49 {
50  regions_[0] = regions[0];
51  regions_[1] = regions[1];
52 }
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
58 {
59  return regions_;
60 }
61 
62 
63 template<class TrackingData>
64 inline bool Foam::patchEdgeFaceRegions::valid(TrackingData& td) const
65 {
66  return regions_.size() && (findIndex(regions_, labelMax) == -1);
67 }
68 
69 
70 template<class Patch, class TrackingData>
72 (
73  const polyMesh& mesh,
74  const Patch& patch,
75  const tensor& rotTensor,
76  const scalar tol,
77  TrackingData& td
78 )
79 {}
80 
81 
82 template<class Patch, class TrackingData>
84 (
85  const polyMesh& mesh,
86  const Patch& patch,
87  const label edgeI,
88  const label facei,
89  const patchEdgeFaceRegions& faceInfo,
90  const scalar tol,
91  TrackingData& td
92 )
93 {
94  const face& f = patch.localFaces()[facei];
95  const edge& e = patch.edges()[edgeI];
96 
97  label index = findIndex(patch.faceEdges()[facei], edgeI);
98  bool sameOrientation = (f[index] == e.start());
99 
100  // Get information in edge-order
101  edge orientedInfo
102  (
103  faceInfo.regions_[index],
104  faceInfo.regions_[f.fcIndex(index)]
105  );
106  if (!sameOrientation)
107  {
108  orientedInfo.flip();
109  }
110 
111  if (!faceInfo.valid(td))
112  {
114  << "problem." << abort(FatalError);
115  }
116 
117  if ((findIndex(orientedInfo, -1) != -1) || (findIndex(regions_, -1) != -1))
118  {
119  // Blocked edge/face
120  return false;
121  }
122 
123 
124  bool changed = false;
125 
126  regions_.setSize(orientedInfo.size(), labelMax);
127  forAll(orientedInfo, i)
128  {
129  if (orientedInfo[i] != -1 && orientedInfo[i] < regions_[i])
130  {
131  regions_[i] = orientedInfo[i];
132  changed = true;
133  }
134  }
135  return changed;
136 }
137 
138 
139 template<class Patch, class TrackingData>
141 (
142  const polyMesh& mesh,
143  const Patch& patch,
144  const patchEdgeFaceRegions& edgeInfo,
145  const bool sameOrientation,
146  const scalar tol,
147  TrackingData& td
148 )
149 {
150  // Get information in edge-order
151  edge orientedInfo(edgeInfo.regions_[0], edgeInfo.regions_[1]);
152  if (!sameOrientation)
153  {
154  orientedInfo.flip();
155  }
156 
157 
158  if (!edgeInfo.valid(td))
159  {
161  << "problem." << abort(FatalError);
162  }
163 
164  if ((findIndex(orientedInfo, -1) != -1) || (findIndex(regions_, -1) != -1))
165  {
166  // Blocked edge/face
167  return false;
168  }
169 
170 
171  bool changed = false;
172 
173  regions_.setSize(orientedInfo.size(), labelMax);
174  forAll(orientedInfo, i)
175  {
176  if (orientedInfo[i] != -1 && orientedInfo[i] < regions_[i])
177  {
178  regions_[i] = orientedInfo[i];
179  changed = true;
180  }
181  }
182  return changed;
183 }
184 
185 
186 template<class Patch, class TrackingData>
188 (
189  const polyMesh& mesh,
190  const Patch& patch,
191  const label facei,
192  const label edgeI,
193  const patchEdgeFaceRegions& edgeInfo,
194  const scalar tol,
195  TrackingData& td
196 )
197 {
198  const face& f = patch.localFaces()[facei];
199  const edge& e = patch.edges()[edgeI];
200 
201  // Find starting point of edge on face.
202  label index0 = findIndex(patch.faceEdges()[facei], edgeI);
203  label index1 = f.fcIndex(index0);
204  bool sameOrientation = (f[index0] == e.start());
205 
206 
207  // Get information in face-order
208  edge orientedInfo
209  (
210  edgeInfo.regions_[0],
211  edgeInfo.regions_[1]
212  );
213  if (!sameOrientation)
214  {
215  orientedInfo.flip();
216  }
217 
218  if (!edgeInfo.valid(td))
219  {
221  << "problem." << abort(FatalError);
222  }
223 
224  if ((findIndex(orientedInfo, -1) != -1) || (findIndex(regions_, -1) != -1))
225  {
226  // Blocked edge/face
227  return false;
228  }
229 
230 
231  bool changed = false;
232 
233  // Scale if needed
234  regions_.setSize(f.size(), labelMax);
235 
236  if (orientedInfo[0] < regions_[index0])
237  {
238  regions_[index0] = orientedInfo[0];
239  changed = true;
240  }
241  if (orientedInfo[1] < regions_[index1])
242  {
243  regions_[index1] = orientedInfo[1];
244  changed = true;
245  }
246 
247  return changed;
248 }
249 
250 
251 template<class TrackingData>
253 (
254  const patchEdgeFaceRegions& rhs,
255  TrackingData& td
256 ) const
257 {
258  return operator==(rhs);
259 }
260 
261 
262 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
263 
264 inline bool Foam::patchEdgeFaceRegions::operator==
265 (
266  const Foam::patchEdgeFaceRegions& rhs
267 ) const
268 {
269  return regions() == rhs.regions();
270 }
271 
272 
273 inline bool Foam::patchEdgeFaceRegions::operator!=
274 (
275  const Foam::patchEdgeFaceRegions& rhs
276 ) const
277 {
278  return !(*this == rhs);
279 }
280 
281 
282 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
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.
Transport of regions for use in PatchEdgeFaceWave.
patchEdgeFaceRegions()
Construct null.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
void transform(const polyMesh &mesh, const Patch &patch, const tensor &rotTensor, const scalar tol, TrackingData &td)
Apply rotation matrix.
static const label labelMax
Definition: label.H:62
errorManip< error > abort(error &err)
Definition: errorManip.H:131
bool operator==(const patchEdgeFaceRegions &) const
label size() const
Return the number of elements in the FixedList.
Definition: FixedListI.H:429
const labelList & regions() const
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
bool equal(const patchEdgeFaceRegions &, TrackingData &) const
Same (like operator==)
void setSize(const label)
Reset size of List.
Definition: List.C:281
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.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void flip()
Flip the edge in-place.
Definition: edgeI.H:157
label start() const
Return start vertex label.
Definition: edgeI.H:81