CV2DI.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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
28 inline Foam::label Foam::CV2D::insertPoint
29 (
30  const point2D& p,
31  const label type
32 )
33 {
34  unsigned int nVert = number_of_vertices();
35 
36  return insertPoint(toPoint(p), nVert, type);
37 }
38 
39 
40 inline Foam::label Foam::CV2D::insertPoint
41 (
42  const point2D& p,
43  const label index,
44  const label type
45 )
46 {
47  return insertPoint(toPoint(p), index, type);
48 }
49 
50 
51 inline Foam::label Foam::CV2D::insertPoint
52 (
53  const Point& p,
54  const label index,
55  const label type
56 )
57 {
58  unsigned int nVert = number_of_vertices();
59 
60  Vertex_handle vh = insert(p);
61 
62  if (nVert == number_of_vertices())
63  {
65  << "Failed to insert point " << toPoint2D(p) << endl;
66  }
67  else
68  {
69  vh->index() = index;
70  vh->type() = type;
71  }
72 
73  return vh->index();
74 }
75 
76 
77 inline bool Foam::CV2D::insertMirrorPoint
78 (
79  const point2D& nearSurfPt,
80  const point2D& surfPt
81 )
82 {
83  point2D mirrorPoint(2*surfPt - nearSurfPt);
84 
85  if (qSurf_.outside(toPoint3D(mirrorPoint)))
86  {
87  insertPoint(mirrorPoint, Vb::MIRROR_POINT);
88  return true;
89  }
90  else
91  {
92  return false;
93  }
94 }
95 
96 
97 inline void Foam::CV2D::insertPointPair
98 (
99  const scalar ppDist,
100  const point2D& surfPt,
101  const vector2D& n
102 )
103 {
104  vector2D ppDistn = ppDist*n;
105 
106  label master = insertPoint
107  (
108  surfPt - ppDistn,
109  number_of_vertices() + 1
110  );
111 
112  insertPoint(surfPt + ppDistn, master);
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  return controls_;
121 }
122 
123 
124 inline const Foam::point2D& Foam::CV2D::toPoint2D(const point& p) const
125 {
126  return reinterpret_cast<const point2D&>(p);
127 }
128 
129 
131 {
132  point2DField temp(p.size());
133  forAll(temp, pointi)
134  {
135  temp[pointi] = point2D(p[pointi].x(), p[pointi].y());
136  }
137  return temp;
138 }
139 
140 
142 {
143  return point(p.x(), p.y(), z_);
144 }
145 
146 
147 #ifdef CGAL_INEXACT
148 
149 inline Foam::CV2D::point2DFromPoint Foam::CV2D::toPoint2D(const Point& P) const
150 {
151  return reinterpret_cast<point2DFromPoint>(P);
152 }
153 
154 
156 {
157  return reinterpret_cast<PointFromPoint2D>(p);
158 }
159 
160 #else
161 
163 {
164  return point2D(CGAL::to_double(P.x()), CGAL::to_double(P.y()));
165 }
166 
167 
169 {
170  return Point(p.x(), p.y());
171 }
172 
173 #endif
174 
175 
176 inline Foam::point Foam::CV2D::toPoint3D(const Point& P) const
177 {
178  return point(CGAL::to_double(P.x()), CGAL::to_double(P.y()), z_);
179 }
180 
181 
182 inline void Foam::CV2D::movePoint(const Vertex_handle& vh, const Point& P)
183 {
184  int i = vh->index();
185  int t = vh->type();
186 
187  remove(vh);
188 
189  Vertex_handle newVh = insert(P);
190 
191  newVh->index() = i;
192  newVh->type() = t;
193 
194  // label i = vh->index();
195  // move(vh, P);
196  // vh->index() = i;
197 
198  // vh->set_point(P);
199  // fast_restore_Delaunay(vh);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
204 
205 inline bool Foam::boundaryTriangle(const CV2D::Face_handle fc)
206 {
207  return boundaryTriangle
208  (
209  *fc->vertex(0),
210  *fc->vertex(1),
211  *fc->vertex(2)
212  );
213 }
214 
215 
216 inline bool Foam::outsideTriangle(const CV2D::Face_handle fc)
217 {
218  return outsideTriangle
219  (
220  *fc->vertex(0),
221  *fc->vertex(1),
222  *fc->vertex(2)
223  );
224 }
225 
226 
227 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field< bool > outside(const pointField &samplePts) const
Check if points are outside surfaces to conform to.
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
vector2D point2D
Point2D is a vector.
Definition: point2D.H:41
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Vector2D< scalar > vector2D
vector2D obtained from generic Vector2D
Definition: vector2D.H:49
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
point toPoint3D(const point2D &) const
Definition: CV2DI.H:141
bool outsideTriangle(const CV2D::Face_handle fc)
Definition: CV2DI.H:216
const cv2DControls & meshControls() const
Definition: CV2DI.H:118
const Cmpt & x() const
Definition: Vector2DI.H:68
const point2D & toPoint2D(const point &) const
Definition: CV2DI.H:124
scalar y
PointFromPoint2D toPoint(const point2D &) const
Definition: CV2DI.H:168
const Cmpt & y() const
Definition: Vector2DI.H:74
Pre-declare SubField and related Field type.
Definition: Field.H:56
timeIndices insert(timeIndex, timeDirs[timeI].value())
Controls for the 2D CV mesh generator.
Definition: cv2DControls.H:58
#define WarningInFunction
Report a warning using Foam::Warning.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void movePoint(const Vertex_handle &vh, const Point &P)
Definition: CV2DI.H:182
label n
bool boundaryTriangle(const CV2D::Face_handle fc)
Definition: CV2DI.H:205
volScalarField & p
const Point & PointFromPoint2D
Definition: CV2D.H:351