tetrahedron.C
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) 2011-2017 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 "tetrahedron.H"
29 #include "scalarField.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Point, class PointRef>
35 (
36  const scalar tol
37 ) const
38 {
39  // (Probably very inefficient) minimum containment sphere calculation.
40  // From http://www.imr.sandia.gov/papers/imr11/shewchuk2.pdf:
41  // Sphere ctr is smallest one of
42  // - tet circumcentre
43  // - triangle circumcentre
44  // - edge mids
45 
46  const scalar fac = 1 + tol;
47 
48  // Halve order of tolerance for comparisons of sqr.
49  const scalar facSqr = Foam::sqrt(fac);
50 
51 
52  // 1. Circumcentre itself.
53 
54  pointHit pHit(circumCentre());
55  pHit.setHit();
56  scalar minRadiusSqr = magSqr(pHit.point() - a_);
57 
58 
59  // 2. Try circumcentre of tet triangles. Create circumcircle for triFace and
60  // check if 4th point is inside.
61 
62  {
63  point ctr = triPointRef(a_, b_, c_).circumCentre();
64  scalar radiusSqr = magSqr(ctr - a_);
65 
66  if
67  (
68  radiusSqr < minRadiusSqr
69  && Foam::magSqr(d_-ctr) <= facSqr*radiusSqr
70  )
71  {
72  pHit.setMiss(false);
73  pHit.setPoint(ctr);
74  minRadiusSqr = radiusSqr;
75  }
76  }
77  {
78  point ctr = triPointRef(a_, b_, d_).circumCentre();
79  scalar radiusSqr = magSqr(ctr - a_);
80 
81  if
82  (
83  radiusSqr < minRadiusSqr
84  && Foam::magSqr(c_-ctr) <= facSqr*radiusSqr
85  )
86  {
87  pHit.setMiss(false);
88  pHit.setPoint(ctr);
89  minRadiusSqr = radiusSqr;
90  }
91  }
92  {
93  point ctr = triPointRef(a_, c_, d_).circumCentre();
94  scalar radiusSqr = magSqr(ctr - a_);
95 
96  if
97  (
98  radiusSqr < minRadiusSqr
99  && Foam::magSqr(b_-ctr) <= facSqr*radiusSqr
100  )
101  {
102  pHit.setMiss(false);
103  pHit.setPoint(ctr);
104  minRadiusSqr = radiusSqr;
105  }
106  }
107  {
108  point ctr = triPointRef(b_, c_, d_).circumCentre();
109  scalar radiusSqr = magSqr(ctr - b_);
110 
111  if
112  (
113  radiusSqr < minRadiusSqr
114  && Foam::magSqr(a_-ctr) <= facSqr*radiusSqr
115  )
116  {
117  pHit.setMiss(false);
118  pHit.setPoint(ctr);
119  minRadiusSqr = radiusSqr;
120  }
121  }
122 
123 
124  // 3. Try midpoints of edges
125 
126  // mid of edge A-B
127  {
128  point ctr = 0.5*(a_ + b_);
129  scalar radiusSqr = magSqr(a_ - ctr);
130  scalar testRadSrq = facSqr*radiusSqr;
131 
132  if
133  (
134  radiusSqr < minRadiusSqr
135  && magSqr(c_-ctr) <= testRadSrq
136  && magSqr(d_-ctr) <= testRadSrq)
137  {
138  pHit.setMiss(false);
139  pHit.setPoint(ctr);
140  minRadiusSqr = radiusSqr;
141  }
142  }
143 
144  // mid of edge A-C
145  {
146  point ctr = 0.5*(a_ + c_);
147  scalar radiusSqr = magSqr(a_ - ctr);
148  scalar testRadSrq = facSqr*radiusSqr;
149 
150  if
151  (
152  radiusSqr < minRadiusSqr
153  && magSqr(b_-ctr) <= testRadSrq
154  && magSqr(d_-ctr) <= testRadSrq
155  )
156  {
157  pHit.setMiss(false);
158  pHit.setPoint(ctr);
159  minRadiusSqr = radiusSqr;
160  }
161  }
162 
163  // mid of edge A-D
164  {
165  point ctr = 0.5*(a_ + d_);
166  scalar radiusSqr = magSqr(a_ - ctr);
167  scalar testRadSrq = facSqr*radiusSqr;
168 
169  if
170  (
171  radiusSqr < minRadiusSqr
172  && magSqr(b_-ctr) <= testRadSrq
173  && magSqr(c_-ctr) <= testRadSrq
174  )
175  {
176  pHit.setMiss(false);
177  pHit.setPoint(ctr);
178  minRadiusSqr = radiusSqr;
179  }
180  }
181 
182  // mid of edge B-C
183  {
184  point ctr = 0.5*(b_ + c_);
185  scalar radiusSqr = magSqr(b_ - ctr);
186  scalar testRadSrq = facSqr*radiusSqr;
187 
188  if
189  (
190  radiusSqr < minRadiusSqr
191  && magSqr(a_-ctr) <= testRadSrq
192  && magSqr(d_-ctr) <= testRadSrq
193  )
194  {
195  pHit.setMiss(false);
196  pHit.setPoint(ctr);
197  minRadiusSqr = radiusSqr;
198  }
199  }
200 
201  // mid of edge B-D
202  {
203  point ctr = 0.5*(b_ + d_);
204  scalar radiusSqr = magSqr(b_ - ctr);
205  scalar testRadSrq = facSqr*radiusSqr;
206 
207  if
208  (
209  radiusSqr < minRadiusSqr
210  && magSqr(a_-ctr) <= testRadSrq
211  && magSqr(c_-ctr) <= testRadSrq)
212  {
213  pHit.setMiss(false);
214  pHit.setPoint(ctr);
215  minRadiusSqr = radiusSqr;
216  }
217  }
218 
219  // mid of edge C-D
220  {
221  point ctr = 0.5*(c_ + d_);
222  scalar radiusSqr = magSqr(c_ - ctr);
223  scalar testRadSrq = facSqr*radiusSqr;
224 
225  if
226  (
227  radiusSqr < minRadiusSqr
228  && magSqr(a_-ctr) <= testRadSrq
229  && magSqr(b_-ctr) <= testRadSrq
230  )
231  {
232  pHit.setMiss(false);
233  pHit.setPoint(ctr);
234  minRadiusSqr = radiusSqr;
235  }
236  }
237 
238 
239  pHit.setDistance(sqrt(minRadiusSqr));
241  return pHit;
242 }
243 
244 
245 template<class Point, class PointRef>
247 (
248  scalarField& buffer
249 ) const
250 {
251  // Change of sign between face area vector and gradient
252  // does not matter because of square
253 
254  // Warning: Added a mag to produce positive coefficients even if
255  // the tetrahedron is twisted inside out. This is pretty
256  // dangerous, but essential for mesh motion.
257  scalar magVol = Foam::mag(mag());
258 
259  buffer[0] = (1.0/9.0)*magSqr(Sa())/magVol;
260  buffer[1] = (1.0/9.0)*magSqr(Sb())/magVol;
261  buffer[2] = (1.0/9.0)*magSqr(Sc())/magVol;
262  buffer[3] = (1.0/9.0)*magSqr(Sd())/magVol;
263 }
264 
265 
266 template<class Point, class PointRef>
268 (
269  scalarField& buffer
270 ) const
271 {
272  // Warning. Ordering of edges needs to be the same for a tetrahedron
273  // class, a tetrahedron cell shape model and a tetCell
274 
275  // Warning: Added a mag to produce positive coefficients even if
276  // the tetrahedron is twisted inside out. This is pretty
277  // dangerous, but essential for mesh motion.
278 
279  // Double change of sign between face area vector and gradient
280 
281  scalar magVol = Foam::mag(mag());
282  vector sa = Sa();
283  vector sb = Sb();
284  vector sc = Sc();
285  vector sd = Sd();
286 
287  buffer[0] = (1.0/9.0)*(sa & sb)/magVol;
288  buffer[1] = (1.0/9.0)*(sa & sc)/magVol;
289  buffer[2] = (1.0/9.0)*(sa & sd)/magVol;
290  buffer[3] = (1.0/9.0)*(sd & sb)/magVol;
291  buffer[4] = (1.0/9.0)*(sb & sc)/magVol;
292  buffer[5] = (1.0/9.0)*(sd & sc)/magVol;
293 }
294 
295 
296 template<class Point, class PointRef>
298 (
299  tensorField& buffer
300 ) const
301 {
302  // Change of sign between face area vector and gradient
303  // does not matter because of square
304 
305  scalar magVol = Foam::mag(mag());
306 
307  buffer[0] = (1.0/9.0)*sqr(Sa())/magVol;
308  buffer[1] = (1.0/9.0)*sqr(Sb())/magVol;
309  buffer[2] = (1.0/9.0)*sqr(Sc())/magVol;
310  buffer[3] = (1.0/9.0)*sqr(Sd())/magVol;
311 }
312 
313 
314 template<class Point, class PointRef>
316 (
317  tensorField& buffer
318 ) const
319 {
320  // Warning. Ordering of edges needs to be the same for a tetrahedron
321  // class, a tetrahedron cell shape model and a tetCell
322 
323  // Double change of sign between face area vector and gradient
324 
325  scalar magVol = Foam::mag(mag());
326  vector sa = Sa();
327  vector sb = Sb();
328  vector sc = Sc();
329  vector sd = Sd();
330 
331  buffer[0] = (1.0/9.0)*(sa * sb)/magVol;
332  buffer[1] = (1.0/9.0)*(sa * sc)/magVol;
333  buffer[2] = (1.0/9.0)*(sa * sd)/magVol;
334  buffer[3] = (1.0/9.0)*(sd * sb)/magVol;
335  buffer[4] = (1.0/9.0)*(sb * sc)/magVol;
336  buffer[5] = (1.0/9.0)*(sd * sc)/magVol;
337 }
338 
339 
340 // ************************************************************************* //
Calculate the second temporal derivative.
pointHit containmentSphere(const scalar tol) const
Return (min)containment sphere, i.e. the smallest sphere with.
Definition: tetrahedron.C:28
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void setPoint(const point_type &p)
Set the point.
Definition: pointHit.H:235
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void gradNiSquared(scalarField &buffer) const
Fill buffer with shape function products.
Definition: tetrahedron.C:240
dimensionedScalar sqrt(const dimensionedScalar &ds)
void setHit() noexcept
Set the hit status on.
Definition: pointHit.H:217
void gradNiGradNi(tensorField &buffer) const
Definition: tetrahedron.C:291
void setDistance(const scalar d) noexcept
Set the distance.
Definition: pointHit.H:243
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
void setMiss(const bool eligible) noexcept
Set the hit status off and set the eligible miss status.
Definition: pointHit.H:226
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
void gradNiGradNj(tensorField &buffer) const
Definition: tetrahedron.C:309
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
const point_type & point() const noexcept
Return the point, no checks.
Definition: pointHit.H:161
void gradNiDotGradNj(scalarField &buffer) const
Definition: tetrahedron.C:261
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)