pointEdgeCollapseI.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) 2012-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 "polyMesh.H"
27 #include "transformer.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 // Update this with w2.
32 template<class TrackingData>
33 inline bool Foam::pointEdgeCollapse::update
34 (
35  const pointEdgeCollapse& w2,
36  const scalar tol,
37  TrackingData& td
38 )
39 {
40  if (!w2.valid(td))
41  {
43  << "problem." << abort(FatalError);
44  }
45 
46  if (!valid(td))
47  {
48  operator=(w2);
49  return true;
50  }
51 
52  if (w2.collapseIndex_ == -1 || collapseIndex_ == -1)
53  {
54  // Not marked for collapse; only happens on edges.
55  return false;
56  }
57 
58  if (w2.collapsePriority_ < collapsePriority_)
59  {
60  return false;
61  }
62  else if (w2.collapsePriority_ > collapsePriority_)
63  {
64  operator=(w2);
65  return true;
66  }
67 
68  // Get overwritten by w2 if it has a higher priority
69  if (w2.collapseIndex_ < collapseIndex_)
70  {
71  operator=(w2);
72  return true;
73  }
74  else if (w2.collapseIndex_ == collapseIndex_)
75  {
76  bool identicalPoint = samePoint(w2.collapsePoint_);
77 
78  bool nearer = (magSqr(w2.collapsePoint_) < magSqr(collapsePoint_));
79 
80  if (nearer)
81  {
82  operator=(w2);
83  }
84 
85  if (identicalPoint)
86  {
87  return false;
88  }
89  else
90  {
91  return nearer;
92  }
93  }
94  else
95  {
96  return false;
97  }
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
103 // Null constructor
105 :
106  collapsePoint_(great, great, great),
107  collapseIndex_(-2),
108  collapsePriority_(-2)
109 {}
110 
111 
112 // Construct from origin, distance
114 (
115  const point& collapsePoint,
116  const label collapseIndex,
117  const label collapsePriority
118 )
119 :
120  collapsePoint_(collapsePoint),
121  collapseIndex_(collapseIndex),
122  collapsePriority_(collapsePriority)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 {
130  return collapsePoint_;
131 }
132 
133 
135 {
136  return collapseIndex_;
137 }
138 
139 
141 {
142  return collapsePriority_;
143 }
144 
145 
146 inline bool Foam::pointEdgeCollapse::samePoint(const point& pt) const
147 {
148  bool isLegal1 = (cmptMin(collapsePoint_) < 0.5*great);
149  bool isLegal2 = (cmptMin(pt) < 0.5*great);
150 
151  if (isLegal1 && isLegal2)
152  {
153  return mag(collapsePoint_ - pt) < 1e-9;//small;
154  }
155  else
156  {
157  return isLegal1 == isLegal2;
158  }
159 }
160 
161 
162 template<class TrackingData>
163 inline bool Foam::pointEdgeCollapse::valid(TrackingData& td) const
164 {
165  return collapseIndex_ != -2;
166 }
167 
168 
169 template<class TrackingData>
171 (
172  const polyPatch& patch,
173  const label patchFacei,
174  const transformer& transform,
175  TrackingData& td
176 )
177 {
178  collapsePoint_ = transform.transformPosition(collapsePoint_);
179 }
180 
181 
182 // Update this with information from connected edge
183 template<class TrackingData>
185 (
186  const polyMesh& mesh,
187  const label pointi,
188  const label edgeI,
189  const pointEdgeCollapse& edgeInfo,
190  const scalar tol,
191  TrackingData& td
192 )
193 {
194  return update(edgeInfo, tol, td);
195 }
196 
197 
198 // Update this with new information on same point
199 template<class TrackingData>
201 (
202  const polyMesh& mesh,
203  const label pointi,
204  const pointEdgeCollapse& newPointInfo,
205  const scalar tol,
206  TrackingData& td
207 )
208 {
209  return update(newPointInfo, tol, td);
210 }
211 
212 
213 // Update this with new information on same point. No extra information.
214 template<class TrackingData>
216 (
217  const pointEdgeCollapse& newPointInfo,
218  const scalar tol,
219  TrackingData& td
220 )
221 {
222  return update(newPointInfo, tol, td);
223 }
224 
225 
226 // Update this with information from connected point
227 template<class TrackingData>
229 (
230  const polyMesh& mesh,
231  const label edgeI,
232  const label pointi,
233  const pointEdgeCollapse& pointInfo,
234  const scalar tol,
235  TrackingData& td
236 )
237 {
238  return update(pointInfo, tol, td);
239 }
240 
241 
242 template<class TrackingData>
244 (
245  const pointEdgeCollapse& rhs,
246  TrackingData& td
247 ) const
248 {
249  return operator==(rhs);
250 }
251 
252 
253 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
254 
255 inline bool Foam::pointEdgeCollapse::operator==
256 (
257  const Foam::pointEdgeCollapse& rhs
258 ) const
259 {
260  return
261  collapseIndex_ == rhs.collapseIndex_
262  && collapsePriority_ == rhs.collapsePriority_
263  && samePoint(rhs.collapsePoint_);
264 }
265 
266 
267 inline bool Foam::pointEdgeCollapse::operator!=
268 (
269  const Foam::pointEdgeCollapse& rhs
270 ) const
271 {
272  return !(*this == rhs);
273 }
274 
275 
276 // ************************************************************************* //
#define w2
Definition: blockCreate.C:32
Determines length of string of edges walked to point.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgeCollapse &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
pointEdgeCollapse()
Construct null.
bool equal(const pointEdgeCollapse &, TrackingData &) const
Same (like operator==)
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgeCollapse &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
const point & collapsePoint() const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const doubleScalar e
Definition: doubleScalar.H:105
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
dimensioned< scalar > magSqr(const dimensioned< Type > &)