sweptFaceAreaWeightAMI.C
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) 2017-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 "sweptFaceAreaWeightAMI.H"
27 #include "cut.H"
28 #include "linearEqn.H"
29 #include "quadraticEqn.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 template<class SourcePatch, class TargetPatch>
35 minCutRatio_ = 10*small;
36 
37 template<class SourcePatch, class TargetPatch>
39 maxDot_ = - cos(degToRad(89.9));
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 template<class SourcePatch, class TargetPatch>
45 template<unsigned Size>
47 (
48  const cutTriList<Size>& tris,
49  const word& name
50 ) const
51 {
52  OFstream obj(typeName + "_" + name + ".vtk");
53 
54  obj << "# vtk DataFile Version 2.0" << endl
55  << "triSurface" << endl
56  << "ASCII" << endl
57  << "DATASET POLYDATA" << endl
58  << "POINTS " << 3*tris.size() << " float" << endl;
59 
60  for (label i = 0; i < tris.size(); ++ i)
61  {
62  for (label j = 0; j < 3; ++ j)
63  {
64  const vector& x = tris[i][j];
65  obj << x.x() << ' ' << x.y() << ' ' << x.z() << endl;
66  }
67  }
68 
69  obj << "POLYGONS " << tris.size() << ' ' << 4*tris.size() << endl;
70 
71  for (label i = 0; i < tris.size(); ++ i)
72  {
73  obj << 3 << ' ' << 3*i << ' ' << 3*i + 1 << ' ' << 3*i + 2 << endl;
74  }
75 }
76 
77 
78 template<class SourcePatch, class TargetPatch>
80 (
81  const face& f,
82  const pointField& ps,
83  const string& name
84 ) const
85 {
86  OFstream obj(typeName + "_" + name + ".obj");
87 
88  for (label i = 0; i < f.size(); ++ i)
89  {
90  const vector& x = ps[f[i]];
91  obj << "v " << x.x() << ' ' << x.y() << ' ' << x.z() << endl;
92  }
93 
94  obj << 'f';
95  for (label i = 0; i < f.size(); ++ i)
96  {
97  obj << ' ' << i + 1;
98  }
99  obj << endl;
100 }
101 
102 
103 template<class SourcePatch, class TargetPatch>
105 (
106  const label srcN,
107  const FixedList<point, 4>& srcTri,
108  const FixedList<point, 4>& srcNrm
109 ) const
110 {
111  scalar l = 0;
112  for (label i = 0; i < srcN - 1; ++ i)
113  {
114  l = max(l, mag(srcTri[i] - srcTri[i + 1]));
115  }
116 
117  const label nu = 20, nv = 20;
118  const scalar u0 = 0, u1 = 1, v0 = -l, v1 = l;
119 
120  OFstream obj(typeName + "_projection.obj");
121 
122  for (label i = 0; i < srcN; ++ i)
123  {
124  const point& p0 = srcTri[i], p1 = srcTri[(i + 1) % srcN];
125  const vector& n0 = srcNrm[i], n1 = srcNrm[(i + 1) % srcN];
126 
127  for (label iu = 0; iu <= nu; ++ iu)
128  {
129  const scalar u = u0 + (u1 - u0)*scalar(iu)/nu;
130  for (label iv = 0; iv <= nv; ++ iv)
131  {
132  const scalar v = v0 + (v1 - v0)*scalar(iv)/nv;
133  const vector x = p0 + (p1 - p0)*u + (n0 + (n1 - n0)*u)*v;
134  obj << "v " << x.x() << ' ' << x.y() << ' ' << x.z() << endl;
135  }
136  }
137  }
138 
139  for (label i = 0; i < srcN; ++ i)
140  {
141  for (label iu = 0; iu < nu; ++ iu)
142  {
143  for (label iv = 0; iv < nv; ++ iv)
144  {
145  obj << "f "
146  << i*(nu + 1)*(nv + 1) + (nv + 1)*iu + iv + 1 << ' '
147  << i*(nu + 1)*(nv + 1) + (nv + 1)*iu + iv + 2 << ' '
148  << i*(nu + 1)*(nv + 1) + (nv + 1)*(iu + 1) + iv + 2 << ' '
149  << i*(nu + 1)*(nv + 1) + (nv + 1)*(iu + 1) + iv + 1<< endl;
150  }
151  }
152  }
153 }
154 
155 
156 template<class SourcePatch, class TargetPatch>
159 (
160  FixedList<point, 4>& srcTri,
161  FixedList<point, 4>& srcNrm,
162  const FixedList<point, 3>& tgtTri
163 ) const
164 {
165  // The target normal
166  const vector np = plane(tgtTri[0], tgtTri[1], tgtTri[2]).normal();
167 
168  // Dot products between the projection normals and the target plane
169  FixedList<scalar, 4> dots;
170  for (label i = 0; i < 3; ++ i)
171  {
172  dots[i] = srcNrm[i] & np;
173  }
174  dots[3] = 0;
175 
176  // Search though for the number of reversed directions and an index for
177  // both a reverse and a forward direction
178  label nR = 0, iR = -1, iF = -1;
179  for (label i = 0; i < 3; ++ i)
180  {
181  if (dots[i] > maxDot_)
182  {
183  ++ nR;
184  iR = i;
185  }
186  else
187  {
188  iF = i;
189  }
190  }
191 
192  // If all the normals hit in the forward or reverse direction then return
193  // the existing projection or no projection respectively
194  if (nR == 0 || nR == 3)
195  {
196  return 3 - nR;
197  }
198 
199  // One normal hits in the reverse direction
200  if (nR == 1)
201  {
202  /*
203  o j1
204  / \
205  / \
206  / \
207  / o jR
208  / / .
209  / / .
210  / / .
211  o - > - o . . . . (old iR)
212  i0 iR
213  */
214 
215  // Make a gap by duplicating the reverse the point
216  for (label j = 2; j >= iR; -- j)
217  {
218  srcTri[j + 1] = srcTri[j];
219  srcNrm[j + 1] = srcNrm[j];
220  dots[j + 1] = dots[j];
221  }
222 
223  const label jR = iR + 1;
224  const label i0 = (iR + 3) % 4, j1 = (jR + 1) % 4;
225 
226  const scalar wi = (dots[iR] - maxDot_)/(dots[iR] - dots[i0]);
227  const scalar wj = (dots[jR] - maxDot_)/(dots[jR] - dots[j1]);
228 
229  srcTri[iR] = srcTri[iR] + (srcTri[i0] - srcTri[iR])*wi;
230  srcTri[jR] = srcTri[jR] + (srcTri[j1] - srcTri[jR])*wj;
231 
232  srcNrm[iR] = srcNrm[iR] + (srcNrm[i0] - srcNrm[iR])*wi;
233  srcNrm[jR] = srcNrm[jR] + (srcNrm[j1] - srcNrm[jR])*wj;
234 
235  return 4;
236  }
237 
238  // Two normals hit in the reverse direction
239  if (nR == 2)
240  {
241  /*
242  . (old jR)
243  . .
244  jR o .
245  / \ .
246  / \ .
247  / \ .
248  / \ .
249  / \ .
250  o - - > - - o . . (old iR)
251  iF iR
252  */
253 
254  const label iR = (iF + 1) % 3, jR = (iF + 2) % 3;
255  const label i0 = (iR + 2) % 3, j1 = (jR + 1) % 3;
256 
257  const scalar wi = (dots[iR] - maxDot_)/(dots[iR] - dots[i0]);
258  const scalar wj = (dots[jR] - maxDot_)/(dots[jR] - dots[j1]);
259 
260  srcTri[iR] = srcTri[iR] + (srcTri[i0] - srcTri[iR])*wi;
261  srcTri[jR] = srcTri[jR] + (srcTri[j1] - srcTri[jR])*wj;
262 
263  srcNrm[iR] = srcNrm[iR] + (srcNrm[i0] - srcNrm[iR])*wi;
264  srcNrm[jR] = srcNrm[jR] + (srcNrm[j1] - srcNrm[jR])*wj;
265 
266  return 3;
267  }
268 
269  // This cannot happen
270  return -1;
271 }
272 
273 
274 template<class SourcePatch, class TargetPatch>
276 (
277  const point& p0,
278  const point& p1,
279  const vector& n0,
280  const vector& n1,
281  const FixedList<point, 3>& tgtTri
282 ) const
283 {
284  // The target plane
285  const plane tgtPln(tgtTri[0], tgtTri[1], tgtTri[2]);
286  const point& pp = tgtPln.refPoint();
287  const vector& np = tgtPln.normal();
288 
289  // Calculate the bounding intersection points. These are the locations at
290  // which the bounding lines of the projected surface intersect with the
291  // target plane.
292  const scalar v0 = ((pp - p0) & np)/(n0 & np);
293  const scalar v1 = ((pp - p1) & np)/(n1 & np);
294  Pair<point> cutPnts(p0 + v0*n0, p1 + v1*n1);
295  Pair<vector> cutNrms((p1 - p0 + n1*v0) ^ n0, (p1 - p0 - n0*v1) ^ n1);
296 
297  // Calculate edge intersections with the surface
298  for (label i = 0; i < 3; ++ i)
299  {
300  // The current target edge
301  const point& l0 = tgtTri[i], l1 = tgtTri[(i + 1)%3];
302 
303  // Form the quadratic in the surface parameter u
304  const vector k = l0 - p0;
305  const vector a = l0 - l1, b = p1 - p0, c = n0, d = n1 - n0;
306  const vector ka = k ^ a, ba = b ^ a, ca = c ^ a, da = d ^ a;
307  const scalar A = d & ba;
308  const scalar B = (c & ba) - (d & ka);
309  const scalar C = - c & ka;
310 
311  // Solve and extract the other parameters
312  const Roots<2> us = quadraticEqn(A, B, C).roots();
313  Pair<scalar> vs, ws;
314  Pair<vector> ns;
315  Pair<bool> cuts(false, false);
316  for (label j = 0; j < 2; ++ j)
317  {
318  if (us.type(j) == roots::real)
319  {
320  const vector den = ca + da*us[j];
321  if (magSqr(den) > vSmall)
322  {
323  const vector vNum = ka - ba*us[j];
324  const vector wNum = (- k + b*us[j]) ^ (c + d*us[j]);
325  vs[j] = (vNum & den)/magSqr(den);
326  ws[j] = (wNum & den)/magSqr(den);
327  ns[j] = (b ^ (c + d*us[j])) + (n1 ^ n0)*vs[j];
328  const bool cutu = 0 < us[j] && us[j] < 1;
329  const bool cutw = 0 < ws[j] && ws[j] < 1;
330  cuts[j] = cutu && cutw;
331  }
332  }
333  }
334 
335  // If we have just one intersection in bounds then store it in the
336  // result list based on its direction
337  if (cuts[0] != cuts[1])
338  {
339  const label j = cuts[0] ? 0 : 1;
340  const scalar na = ns[j] & a;
341  cutPnts[na < 0] = l0 + (l1 - l0)*ws[j];
342  cutNrms[na < 0] = ns[j];
343  }
344  }
345 
346  // Generate and return the cut plane. If the cut points are not coincident
347  // then form a plane normal by crossing the displacement between the points
348  // by the target plane normal. If the points are coincident then use the
349  // projected surface normal evaluated at the first cut point.
350  const vector cutDelta = cutPnts[1] - cutPnts[0];
351  const bool coincident = mag(cutDelta) < minCutRatio_*mag(p1 - p0);
352  return plane(cutPnts[0], coincident ? cutNrms[0] : np ^ cutDelta);
353 };
354 
355 
356 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
357 
358 template<class SourcePatch, class TargetPatch>
360 (
361  const label srcFacei,
362  const label tgtFacei
363 ) const
364 {
365  const label debugTgtFace =
366  (- 1 - debug) % this->tgtPatch_.size() == tgtFacei
367  ? (- 1 - debug) / this->tgtPatch_.size() + 1 : 0;
368  const bool debugPrint = debug > 0 || debugTgtFace > 0;
369  const bool debugWrite = debug > 1 || debugTgtFace > 1;
370 
371  if (debugPrint)
372  {
373  Info<< "Inter area between source face #" << srcFacei
374  << " and target face #" << tgtFacei << (debugWrite ? "\n" : ": ");
375  }
376 
377  // Patch data
378  const pointField& srcPoints = this->srcPatch_.localPoints();
379  const pointField& tgtPoints = this->tgtPatch_.localPoints();
380  const vectorField& srcPointNormals = this->srcPatch_.pointNormals();
381 
382  // Faces
383  const face& srcFace = this->srcPatch_.localFaces()[srcFacei];
384  const face& tgtFace = this->tgtPatch_.localFaces()[tgtFacei];
385 
386  // Write out the faces
387  if (debugWrite)
388  {
389  writeFaceOBJ(srcFace, srcPoints, "source");
390  writeFaceOBJ(tgtFace, tgtPoints, "target");
391  }
392  if (debugTgtFace)
393  {
394  writeFaceOBJ(tgtFace, tgtPoints, "target" + name(tgtFacei));
395  }
396 
397  // Triangulate the faces
398  const faceAreaIntersect::triangulationMode triMode = this->triMode_;
399  faceList srcFaceTris, tgtFaceTris;
400  faceAreaIntersect::triangulate(srcFace, srcPoints, triMode, srcFaceTris);
401  faceAreaIntersect::triangulate(tgtFace, tgtPoints, triMode, tgtFaceTris);
402 
403  // Area sum
404  scalar areaMag = Zero;
405 
406  // Loop the target triangles
407  forAllConstIter(faceList, tgtFaceTris, tgtIter)
408  {
409  const FixedList<point, 3>
410  tgtTri =
411  {
412  tgtPoints[(*tgtIter)[0]],
413  tgtPoints[(*tgtIter)[this->reverseTarget_ ? 2 : 1]],
414  tgtPoints[(*tgtIter)[this->reverseTarget_ ? 1 : 2]]
415  };
416 
417  // Loop the source triangles
418  forAllConstIter(faceList, srcFaceTris, srcIter)
419  {
421  srcTri =
422  {
423  srcPoints[(*srcIter)[0]],
424  srcPoints[(*srcIter)[1]],
425  srcPoints[(*srcIter)[2]],
426  vector::zero
427  };
429  srcNrm =
430  {
431  srcPointNormals[(*srcIter)[0]],
432  srcPointNormals[(*srcIter)[1]],
433  srcPointNormals[(*srcIter)[2]],
434  vector::zero
435  };
436 
437  // Get the source projection modified for any reverse intersections
438  const label srcN = getSourceProjection(srcTri, srcNrm, tgtTri);
439  if (srcN <= 0)
440  {
441  continue;
442  }
443 
444  // Write the source projection
445  if (debugWrite)
446  {
447  writeProjectionOBJ(srcN, srcTri, srcNrm);
448  }
449 
450  // Create the initial cut triangle list
451  cutTriList<8> cutTris;
452  cutTris.append(tgtTri);
453 
454  // Write the initial triangle
455  if (debugWrite)
456  {
457  writeCutTrisVTK(cutTris, "tris0");
458  }
459 
460  // Do all but one of the cuts
461  for
462  (
463  label i = 0;
464  i < srcN - 1 && (debugWrite || cutTris.size());
465  ++ i
466  )
467  {
468  // Do the cut
469  const plane cutPlane =
470  getCutPlane
471  (
472  srcTri[i],
473  srcTri[i+1],
474  srcNrm[i],
475  srcNrm[i+1],
476  tgtTri
477  );
478 
479  cutTriList<8> cutTrisTmp;
480 
481  for (label j = 0; j < cutTris.size(); ++j)
482  {
483  triCut
484  (
485  cutTris[j],
486  cutPlane,
487  cut::noOp(),
488  cut::appendOp<cutTriList<8>>(cutTrisTmp)
489  );
490  }
491  Swap(cutTris, cutTrisTmp);
492 
493  // Write the triangles resulting from the cut
494  if (debugWrite)
495  {
496  writeCutTrisVTK(cutTris, "tris" + name(i + 1));
497  }
498  }
499 
500  // Do the last cut
501  const plane cutPlane =
502  getCutPlane
503  (
504  srcTri[srcN - 1],
505  srcTri[0],
506  srcNrm[srcN - 1],
507  srcNrm[0],
508  tgtTri
509  );
510  cutTriList<8> cutTrisTmp;
511  for (label i = 0; i < cutTris.size(); ++ i)
512  {
513  // Sum the area of the cut triangle
514  areaMag +=
515  mag
516  (
517  triCut
518  (
519  cutTris[i],
520  cutPlane,
521  cut::noOp(),
522  cut::areaOp()
523  )
524  );
525  // Store the cut triangle if it is needed for output
526  if (debugWrite || debugTgtFace)
527  {
528  triCut
529  (
530  cutTris[i],
531  cutPlane,
532  cut::noOp(),
533  cut::appendOp<cutTriList<8>>(cutTrisTmp)
534  );
535  }
536  }
537  Swap(cutTris, cutTrisTmp);
538 
539  // Write the triangles resulting from the cuts
540  if (debugTgtFace && cutTris.size())
541  {
542  static label writeCutTrisVTKIndex = 0;
543  writeCutTrisVTK
544  (
545  cutTris,
546  "target" + name(tgtFacei) + "_"
547  + "tris" + name(writeCutTrisVTKIndex ++)
548  );
549  }
550  if (debugWrite)
551  {
552  writeCutTrisVTK(cutTris, "tris" + name(srcN));
553 
554  Info << "view triangles then press ENTER to continue ...";
555  getchar();
556  }
557  }
558  }
559 
560  // Print the difference between this inter-area and that obtained by the
561  // standard algorithm built around faceAreaIntersect
562  if (debugPrint)
563  {
564  const scalar standardAreaMag =
566  (
567  srcFacei,
568  tgtFacei
569  );
570 
571  Info<<"standard=" << standardAreaMag << ", swept=" << areaMag
572  << endl;
573  }
574 
575  return areaMag;
576 }
577 
578 
579 template<class SourcePatch, class TargetPatch>
580 Foam::scalar
582 {
583  return small;
584 }
585 
586 
587 template<class SourcePatch, class TargetPatch>
588 Foam::scalar
590 {
591  return degToRad(180);
592 }
593 
594 
595 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
596 
597 template<class SourcePatch, class TargetPatch>
599 (
600 )
601 {}
602 
603 
604 // ************************************************************************* //
Swept face area weighted Arbitrary Mesh Interface (AMI) method.
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
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
cut::opAddResult< AboveOp, BelowOp >::type triCut(const FixedList< point, 3 > &tri, const FixedList< scalar, 3 > &level, const AboveOp &aboveOp, const BelowOp &belowOp)
Cut a triangle along the zero plane defined by the given levels. Templated.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensionedScalar j1(const dimensionedScalar &ds)
Face area weighted Arbitrary Mesh Interface (AMI) method.
scalar u0
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dimensionedScalar cos(const dimensionedScalar &ds)
void Swap(T &a, T &b)
Definition: Swap.H:43
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const Cmpt & x() const
Definition: VectorI.H:75
dimensioned< scalar > magSqr(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual scalar interArea(const label srcFacei, const label tgtFacei) const
Area of intersection between source and target faces.
vector point
Point is a vector.
Definition: point.H:41
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & nu
virtual ~sweptFaceAreaWeightAMI()
Destructor.