faMeshDemandDrivenData.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2018-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "faMesh.H"
30 #include "faMeshLduAddressing.H"
31 #include "areaFields.H"
32 #include "edgeFields.H"
33 #include "fac.H"
34 #include "cartesianCS.H"
35 #include "scalarMatrices.H"
36 #include "processorFaPatch.H"
37 #include "processorFaPatchFields.H"
38 #include "emptyFaPatchFields.H"
39 #include "wedgeFaPatch.H"
40 #include "triangle.H"
41 
42 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // Define an area-weighted normal for three points (defining a triangle)
48 // (p0, p1, p2) are the base, first, second points respectively
49 //
50 // From the original Tukovic code:
51 //
52 // vector n = (d1^d2)/mag(d1^d2);
53 // scalar sinAlpha = mag(d1^d2)/(mag(d1)*mag(d2));
54 // scalar w = sinAlpha/(mag(d1)*mag(d2));
55 //
56 // ie, normal weighted by area, sine angle and inverse distance squared.
57 // - area : larger weight for larger areas
58 // - sin : lower weight for narrow angles (eg, shards)
59 // - inv distance squared : lower weights for distant points
60 //
61 // The above refactored, with 0.5 for area:
62 //
63 // (d1 ^ d2) / (2 * magSqr(d1) * magSqr(d2))
64 
66 (
67  const vector& a,
68  const vector& b
69 )
70 {
71  // TBD: Use volatile to avoid aggressive branch optimization
72  const scalar s(2*magSqr(a)*magSqr(b));
73 
74  return s < ROOTVSMALL ? Zero : (a ^ b) / s;
75 }
76 
77 
78 // The area normal for the face dual (around base-point)
79 // described by the right/left edge points and the centre point
80 //
81 // The adjustment for 1/2 edge point (the dual point) is done internally
83 (
84  const point& basePoint,
85  const point& rightPoint,
86  const point& leftPoint,
87  const point& centrePoint
88 )
89 {
90  const vector mid(centrePoint - basePoint);
91 
92  return
93  (
95  (
96  0.5*(rightPoint - basePoint), // vector to right 1/2 edge
97  mid
98  )
100  (
101  mid,
102  0.5*(leftPoint - basePoint) // vector to left 1/2 edge
103  )
104  );
105 }
106 
107 
108 // Weighted area normal for the face triangle (around base-point)
109 // to the edge points and the centre point.
110 // Used for example for area-weighted edge normals
111 // ---------------------------------------------------------------------------
112 // own |
113 // * | Don't bother trying to split the triangle into
114 // /.\ | sub-triangles. Planar anyhow and skew weighting
115 // / . \ | is less relevant here.
116 // / . \ |
117 // / . \ | Note: need negative for neighbour side orientation.
118 // e0 *----:----* e1 |
119 // ---------------------------------------------------------------------------
120 
122 (
123  const linePointRef& edgeLine,
124  const point& faceCentre
125 )
126 {
127  return
128  (
130  (
131  (edgeLine.a() - faceCentre), // From centre to right edge
132  (edgeLine.b() - faceCentre) // From centre to left edge
133  )
134  );
135 }
136 
137 
138 // Simple area normal calculation for specified edge vector and face centre
139 // The face centre comes last since it is less accurate
140 static inline vector areaNormalFaceTriangle
141 (
142  const linePointRef& edgeLine,
143  const point& faceCentre
144 )
145 {
146  return triPointRef::areaNormal(edgeLine.a(), edgeLine.b(), faceCentre);
147 }
148 
149 } // End namespace Foam
150 
151 
152 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
153 
154 namespace Foam
155 {
156 
157 // Calculate transform tensor with reference vector (unitAxis1)
158 // and direction vector (axis2).
159 //
160 // This is nearly identical to the meshTools axesRotation
161 // with an E3_E1 transformation with the following exceptions:
162 //
163 // - axis1 (e3 == unitAxis1): is already normalized (unit vector)
164 // - axis2 (e1 == dirn): no difference
165 // - transformation is row-vectors, not column-vectors
166 static inline tensor rotation_e3e1
167 (
168  const vector& unitAxis1,
169  vector dirn
170 )
171 {
172  dirn.removeCollinear(unitAxis1);
173  dirn.normalise();
174 
175  // Set row vectors
176  return tensor
177  (
178  dirn,
179  (unitAxis1^dirn),
180  unitAxis1
181  );
182 }
183 
184 } // End namespace Foam
185 
187 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
188 
189 namespace Foam
190 {
191 
192 // A bitSet (size patch nPoints()) with boundary points marked
194 {
195  // Initially all unmarked
196  bitSet markPoints(p.nPoints());
197  for (const edge& e : p.boundaryEdges())
198  {
199  // Mark boundary points
200  markPoints.set(e.first());
201  markPoints.set(e.second());
202  }
203 
204  return markPoints;
205 }
206 
207 
208 } // End namespace Foam
209 
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
213 void Foam::faMesh::calcLduAddressing() const
214 {
216  << "Calculating addressing" << endl;
217 
218  if (lduPtr_)
219  {
221  << "lduPtr_ already allocated"
222  << abort(FatalError);
223  }
224 
225  lduPtr_ = new faMeshLduAddressing(*this);
226 }
227 
228 
229 void Foam::faMesh::calcPatchStarts() const
230 {
232  << "Calculating patch starts" << endl;
233 
234  if (patchStartsPtr_)
235  {
237  << "patchStarts already allocated"
238  << abort(FatalError);
239  }
240 
241  patchStartsPtr_ = new labelList(boundary().patchStarts());
242 }
243 
244 
245 void Foam::faMesh::calcWhichPatchFaces() const
246 {
247  // Usually need both together
248  if (polyPatchFacesPtr_ || polyPatchIdsPtr_)
249  {
251  << "Already allocated polyPatchFaces/polyPatchIds"
252  << abort(FatalError);
253  }
254 
255  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
256 
257  polyPatchFacesPtr_.reset
258  (
259  new List<labelPair>(pbm.whichPatchFace(faceLabels_))
260  );
261 
262  labelHashSet ids;
263 
264  // Extract patch ids from (patch, facei) tuples
265  for (const labelPair& tup : *polyPatchFacesPtr_)
266  {
267  ids.insert(tup.first());
268  }
269 
270  ids.erase(-1); // Without internal faces (patchi == -1)
271 
272  // parSync
274  (
275  ids,
276  bitOrOp<labelHashSet>(),
278  this->comm()
279  );
280 
281  polyPatchIdsPtr_.reset(new labelList(ids.sortedToc()));
282 }
283 
284 
285 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
286 
287 namespace Foam
288 {
289 
290 //- Calculate the 'Le' vector from faceCentre to edge centre
291 // using the edge normal to correct for curvature
292 //
293 // Normalise and rescaled to the edge length
294 static inline vector calcLeVector
295 (
296  const point& faceCentre,
297  const linePointRef& edgeLine,
298  const vector& edgeNormal // (unit or area normal)
299 )
300 {
301  const vector centreToEdge(edgeLine.centre() - faceCentre);
302 
303  vector leVector(edgeLine.vec() ^ edgeNormal);
304 
305  scalar s(mag(leVector));
306 
307  if (s < ROOTVSMALL)
308  {
309  // The calculated edgeNormal somehow degenerate and thus a
310  // bad cross-product?
311  // Revert to basic centre -> edge
312 
313  leVector = centreToEdge;
314  leVector.removeCollinear(edgeLine.unitVec());
315  s = mag(leVector);
316 
317  if (s < ROOTVSMALL)
318  {
319  // Unlikely that this should happen
320  return Zero;
321  }
322 
323  leVector *= edgeLine.mag()/s;
324  }
325  else
326  {
327  // The additional orientation is probably unnecessary
328  leVector *= edgeLine.mag()/s * sign(leVector & centreToEdge);
329  }
330 
331  return leVector;
332 }
333 
334 } // End namespace Foam
335 
336 
337 namespace Foam
338 {
339 
340 // Impose a minimum (edge) vector length
341 // - disallow any mag(vec) < SMALL
342 static inline void imposeMinVectorLength(vectorField& fld)
343 {
344  // If it is small, no orientation information reasonably possible
345  // so use a vector(1,1,1) * sqrt(1/3)*min-length
346 
347  // sqrt(1/3) = 0.5773502691896257, but slightly rounded down
348  const vector minVector(vector::uniform(0.57735*SMALL));
349  const scalar minLenSqr(SMALL*SMALL);
350 
351  for (vector& v : fld)
352  {
353  if (v.magSqr() < minLenSqr)
354  {
355  v = minVector;
356  }
357  }
358 }
359 
360 } // End namespace Foam
361 
362 
363 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
364 
365 Foam::tmp<Foam::vectorField> Foam::faMesh::calcRawEdgeNormals(int order) const
366 {
367  // Return edge normals with flat boundary addressing
368  auto tedgeNormals = tmp<vectorField>::New(nEdges_);
369  auto& edgeNormals = tedgeNormals.ref();
370 
371  // Need face centres
372  const areaVectorField& fCentres = areaCentres();
373 
374  // Also need local points
375  const pointField& localPoints = points();
376 
377 
378  if (order < 0)
379  {
380  // geometryOrder (-1): no communication
381 
382  // Simple (primitive) edge normal calculations.
383  // These are primarly designed to avoid any communication
384  // but are thus necessarily inconsistent across processor boundaries!
385 
387  << "Using geometryOrder < 0 : "
388  "simplified edge area-normals, without processor connectivity"
389  << endl;
390 
391  // Internal (edge normals) - contributions from owner/neighbour
392  for (label edgei = 0; edgei < nInternalEdges_; ++edgei)
393  {
394  const linePointRef edgeLine(edges_[edgei].line(localPoints));
395 
396  edgeNormals[edgei] =
397  (
399  (
400  edgeLine,
401  fCentres[edgeOwner()[edgei]]
402  )
403  // NB: reversed sign (edge orientation flipped for neighbour)
405  (
406  edgeLine,
407  fCentres[edgeNeighbour()[edgei]]
408  )
409  );
410  }
411 
412  // Boundary (edge normals). Like above, but only has owner
413  for (label edgei = nInternalEdges_; edgei < nEdges_; ++edgei)
414  {
415  const linePointRef edgeLine(edges_[edgei].line(localPoints));
416 
417  edgeNormals[edgei] =
418  (
420  (
421  edgeLine,
422  fCentres[edgeOwner()[edgei]]
423  )
424  );
425  }
426  }
427  else
428  {
429  // geometryOrder (0) - no communication
430  // otherwise with communication
431 
432  if (order < 1)
433  {
435  << "Using geometryOrder == 0 : "
436  "weighted edge normals, without processor connectivity"
437  << endl;
438  }
439 
440  // Internal (edge normals)
441  // - area-weighted contributions from owner/neighbour
442  for (label edgei = 0; edgei < nInternalEdges_; ++edgei)
443  {
444  const linePointRef edgeLine(edges_[edgei].line(localPoints));
445 
446  edgeNormals[edgei] =
447  (
449  (
450  edgeLine,
451  fCentres[edgeOwner()[edgei]]
452  )
453  // NB: reversed sign (edge orientation flipped for neighbour)
455  (
456  edgeLine,
457  fCentres[edgeNeighbour()[edgei]]
458  )
459  );
460  }
461 
462  // Boundary (edge normals). Like above, but only has owner
463  for (label edgei = nInternalEdges_; edgei < nEdges_; ++edgei)
464  {
465  const linePointRef edgeLine(edges_[edgei].line(localPoints));
466 
467  edgeNormals[edgei] =
468  (
470  (
471  edgeLine,
472  fCentres[edgeOwner()[edgei]]
473  )
474  );
475  }
476  }
477 
478 
479  // Boundary edge corrections
480  bitSet nbrBoundaryAdjust;
481 
482  // Processor-processor first (for convenience)
483  if (order >= 1)
484  {
485  DynamicList<UPstream::Request> requests
486  (
487  2*(boundary().size() - boundary().nNonProcessor())
488  );
489 
490  PtrList<vectorField> nbrEdgeNormals(boundary().size());
491 
492  // Prefer our own slicing into the raw list (instead of patchSlice).
493  // Although patchSlice does properly handle the start() offset in
494  // the presence of emptyPaPatch, we are slightly paranoia,
495  // and want to avoid kicking off the patchStarts() data...
496 
497  label patchOffset = nInternalEdges_;
498 
499  forAll(boundary(), patchi)
500  {
501  const faPatch& fap = boundary()[patchi];
502 
503  auto edgeNorms = edgeNormals.slice(patchOffset, fap.nEdges());
504  patchOffset += edgeNorms.size();
505 
506  const auto* fapp = isA<processorFaPatch>(fap);
507 
508  if (fapp)
509  {
510  const label nbrProci = fapp->neighbProcNo();
511 
512  if (UPstream::parRun() && !edgeNorms.empty())
513  {
514  nbrEdgeNormals.emplace(patchi, edgeNorms.size());
515 
516  // Recv accumulated weighted edge normals
518  (
519  requests.emplace_back(),
520  nbrProci,
521  nbrEdgeNormals[patchi]
522  );
523 
524  // Send accumulated weighted edge normals
526  (
527  requests.emplace_back(),
528  nbrProci,
529  edgeNorms
530  );
531  }
532  }
533  else if (isA<wedgeFaPatch>(fap))
534  {
535  // Correct wedge edges
536  const auto& wedgePatch = refCast<const wedgeFaPatch>(fap);
537 
538  const vector wedgeNorm
539  (
540  transform
541  (
542  wedgePatch.edgeT(),
543  wedgePatch.centreNormal()
544  ).normalise()
545  );
546 
547  for (vector& e : edgeNorms)
548  {
549  e.removeCollinear(wedgeNorm);
550  }
551  }
552  // TBD:
557  else if (correctPatchPointNormals(patchi) && !fap.coupled())
558  {
559  // Neighbour correction requested
560  nbrBoundaryAdjust.set(patchi);
561  }
562  }
563 
564  // Wait for outstanding requests
565  UPstream::waitRequests(requests);
566 
567  if (!requests.empty())
568  {
569  // Add in received values
570 
571  patchOffset = nInternalEdges_;
572 
573  forAll(boundary(), patchi)
574  {
575  const faPatch& fap = boundary()[patchi];
576 
577  auto edgeNorms = edgeNormals.slice(patchOffset, fap.nEdges());
578  patchOffset += edgeNorms.size();
579 
580  if (nbrEdgeNormals.set(patchi))
581  {
582  const vectorField& nbrNorms = nbrEdgeNormals[patchi];
583 
584  forAll(edgeNorms, patchEdgei)
585  {
586  edgeNorms[patchEdgei] += nbrNorms[patchEdgei];
587  }
588  }
589  }
590  }
591  }
592 
593 
594  // Apply boundary edge corrections
595  if
596  (
597  order >= 1
598  && hasPatchPointNormalsCorrection()
599  && returnReduceOr(nbrBoundaryAdjust.any())
600  )
601  {
603  << "Apply " << nbrBoundaryAdjust.count()
604  << " boundary neighbour corrections" << nl;
605 
606  // Collect face normals, per boundary edge
607 
608  (void)this->haloFaceNormals();
609 
610  // Using our own slicing (instead of patchSlice) - see above
611  label patchOffset = nInternalEdges_;
612 
613  for (const label patchi : nbrBoundaryAdjust)
614  {
615  const faPatch& fap = boundary()[patchi];
616 
617  auto edgeNorms = edgeNormals.slice(patchOffset, fap.nEdges());
618  patchOffset += edgeNorms.size();
619 
620  if (fap.ngbPolyPatchIndex() < 0)
621  {
623  << "Neighbour polyPatch index is not defined "
624  << "for faPatch " << fap.name()
625  << abort(FatalError);
626  }
627 
628  // Apply the correction
629 
630  // Note from Zeljko Tukovic:
631  //
632  // This posibility is used for free-surface tracking
633  // calculations to enforce 90 deg contact angle between
634  // finite-area mesh and neighbouring polyPatch. It is very
635  // important for curvature calculation to have correct normal
636  // at contact line points.
637 
638  vectorField nbrNorms(this->haloFaceNormals(patchi));
639 
640  const label nPatchEdges = min(edgeNorms.size(), nbrNorms.size());
641 
642  for (label patchEdgei = 0; patchEdgei < nPatchEdges; ++patchEdgei)
643  {
644  edgeNorms[patchEdgei].removeCollinear(nbrNorms[patchEdgei]);
645  }
646  }
647  }
648 
649 
650  // Remove collinear components and normalise
651 
652  forAll(edgeNormals, edgei)
653  {
654  const linePointRef edgeLine(edges_[edgei].line(localPoints));
655 
656  edgeNormals[edgei].removeCollinear(edgeLine.unitVec());
657  edgeNormals[edgei].normalise();
658  }
659 
660  imposeMinVectorLength(edgeNormals);
661 
662  return tedgeNormals;
663 }
664 
665 
666 void Foam::faMesh::calcLe() const
667 {
669  << "Calculating local edges" << endl;
670 
671  if (LePtr_)
672  {
674  << "LePtr_ already allocated"
675  << abort(FatalError);
676  }
677 
678  LePtr_ =
679  new edgeVectorField
680  (
681  IOobject
682  (
683  "Le",
684  mesh().pointsInstance(),
686  faMesh::thisDb(),
690  ),
691  *this,
692  dimLength
693  // -> calculatedType()
694  );
695  edgeVectorField& Le = *LePtr_;
696 
697  // Need face centres
698  const areaVectorField& fCentres = areaCentres();
699 
700  // Also need local points
701  const pointField& localPoints = points();
702 
703 
704  // The edgeAreaNormals _may_ use communication (depends on geometryOrder)
705 
706  {
707  const edgeVectorField& edgeNormals = edgeAreaNormals();
708 
709  // Internal (edge vector)
710  {
711  vectorField& fld = Le.primitiveFieldRef();
712  for (label edgei = 0; edgei < nInternalEdges_; ++edgei)
713  {
714  fld[edgei] = calcLeVector
715  (
716  fCentres[edgeOwner()[edgei]],
717  edges_[edgei].line(localPoints),
718  edgeNormals[edgei]
719  );
720  }
721  }
722 
723  // Boundary (edge vector)
724  forAll(boundary(), patchi)
725  {
726  const faPatch& fap = boundary()[patchi];
727  vectorField& pfld = Le.boundaryFieldRef()[patchi];
728 
729  const vectorField& bndEdgeNormals =
730  edgeNormals.boundaryField()[patchi];
731 
732  label edgei = fap.start();
733 
734  forAll(pfld, patchEdgei)
735  {
736  pfld[patchEdgei] = calcLeVector
737  (
738  fCentres[edgeOwner()[edgei]],
739  edges_[edgei].line(localPoints),
740  bndEdgeNormals[patchEdgei]
741  );
742 
743  ++edgei;
744  }
745  }
746  }
747 
748 
749  // Impose a minimum (edge) vector length
750 
751  // Internal (edge vector)
752  {
753  imposeMinVectorLength(Le.primitiveFieldRef());
754  }
755 
756  // Boundary (edge vector)
757  for (vectorField& pfld : Le.boundaryFieldRef())
758  {
759  imposeMinVectorLength(pfld);
760  }
761 }
762 
763 
764 void Foam::faMesh::calcMagLe() const
765 {
767  << "Calculating local edge magnitudes" << endl;
768 
769  if (magLePtr_)
770  {
772  << "magLe() already allocated"
773  << abort(FatalError);
774  }
775 
776  magLePtr_ =
777  new edgeScalarField
778  (
779  IOobject
780  (
781  "magLe",
782  mesh().pointsInstance(),
784  faMesh::thisDb(),
788  ),
789  *this,
790  dimLength
791  );
792 
793  edgeScalarField& magLe = *magLePtr_;
794 
795  const pointField& localPoints = points();
796 
797  // Internal (edge length)
798  {
799  auto iter = magLe.primitiveFieldRef().begin();
800 
801  for (const edge& e : internalEdges())
802  {
803  *iter = e.mag(localPoints);
804 
805  // Do not allow any mag(val) < SMALL
806  if (mag(*iter) < SMALL)
807  {
808  *iter = SMALL;
809  }
810 
811  ++iter;
812  }
813  }
814 
815  // Boundary (edge length)
816  {
817  auto& bfld = magLe.boundaryFieldRef();
818 
819  forAll(boundary(), patchi)
820  {
821  auto iter = bfld[patchi].begin();
822 
823  for (const edge& e : boundary()[patchi].patchSlice(edges_))
824  {
825  *iter = e.mag(localPoints);
826 
827  // Do not allow any mag(val) < SMALL
828  if (mag(*iter) < SMALL)
829  {
830  *iter = SMALL;
831  }
832 
833  ++iter;
834  }
835  }
836  }
837 }
838 
839 
840 void Foam::faMesh::calcFaceCentres() const
841 {
843  << "Calculating face centres" << endl;
844 
845  if (faceCentresPtr_)
846  {
848  << "areaCentres already allocated"
849  << abort(FatalError);
850  }
851 
852  faceCentresPtr_ =
853  new areaVectorField
854  (
855  IOobject
856  (
857  "centres",
858  mesh().pointsInstance(),
860  faMesh::thisDb(),
864  ),
865  *this,
866  dimLength
867  // -> calculatedType()
868  );
869 
870  areaVectorField& centres = *faceCentresPtr_;
871 
872  // Need local points
873  const pointField& localPoints = points();
874 
875 
876  // Internal (face centres)
877  {
878  if (mesh().hasFaceCentres())
879  {
880  // The volume mesh has faceCentres, can reuse them
881 
882  centres.primitiveFieldRef()
883  = UIndirectList<vector>(mesh().faceCentres(), faceLabels());
884  }
885  else
886  {
887  // Calculate manually
888  auto iter = centres.primitiveFieldRef().begin();
889 
890  for (const face& f : faces())
891  {
892  *iter = f.centre(localPoints);
893  ++iter;
894  }
895  }
896  }
897 
898  // Boundary (edge centres or neighbour face centres)
899  {
900  auto& bfld = centres.boundaryFieldRef();
901 
902  forAll(boundary(), patchi)
903  {
904  auto iter = bfld[patchi].begin();
905 
906  for (const edge& e : boundary()[patchi].patchSlice(edges_))
907  {
908  *iter = e.centre(localPoints);
909  ++iter;
910  }
911  }
912  }
913 
914  // Parallel consistency, exchange on processor patches
915  if (UPstream::parRun())
916  {
917  centres.boundaryFieldRef().evaluateCoupled<processorFaPatch>();
918  }
919 }
920 
921 
922 void Foam::faMesh::calcEdgeCentres() const
923 {
925  << "Calculating edge centres" << endl;
926 
927  if (edgeCentresPtr_)
928  {
930  << "edgeCentres already allocated"
931  << abort(FatalError);
932  }
933 
934  edgeCentresPtr_ =
935  new edgeVectorField
936  (
937  IOobject
938  (
939  "edgeCentres",
940  mesh().pointsInstance(),
942  faMesh::thisDb(),
946  ),
947  *this,
948  dimLength
949  // -> calculatedType()
950  );
951 
952  edgeVectorField& centres = *edgeCentresPtr_;
953 
954  // Need local points
955  const pointField& localPoints = points();
956 
957 
958  // Internal (edge centres)
959  {
960  auto iter = centres.primitiveFieldRef().begin();
961 
962  for (const edge& e : internalEdges())
963  {
964  *iter = e.centre(localPoints);
965  ++iter;
966  }
967  }
968 
969  // Boundary (edge centres)
970  {
971  auto& bfld = centres.boundaryFieldRef();
972 
973  forAll(boundary(), patchi)
974  {
975  auto iter = bfld[patchi].begin();
976 
977  for (const edge& e : boundary()[patchi].patchSlice(edges_))
978  {
979  *iter = e.centre(localPoints);
980  ++iter;
981  }
982  }
983  }
984 }
985 
986 
987 void Foam::faMesh::calcS() const
988 {
990  << "Calculating areas" << endl;
991 
992  if (SPtr_)
993  {
995  << "S() already allocated"
996  << abort(FatalError);
997  }
998 
999  SPtr_ = new DimensionedField<scalar, areaMesh>
1000  (
1001  IOobject
1002  (
1003  "S",
1004  time().timeName(),
1005  faMesh::thisDb(),
1009  ),
1010  *this,
1011  dimArea
1012  );
1013  auto& areas = *SPtr_;
1014 
1015 
1016  // No access to fvMesh::magSf(), only polyMesh::faceAreas()
1017  if (mesh().hasFaceAreas())
1018  {
1019  // The volume mesh has faceAreas, can reuse them
1020  UIndirectList<vector> meshFaceAreas(mesh().faceAreas(), faceLabels());
1021 
1022  auto& fld = areas.field();
1023 
1024  forAll(fld, facei)
1025  {
1026  fld[facei] = Foam::mag(meshFaceAreas[facei]);
1027 
1028  // Do not allow any mag(val) < SMALL
1029  if (mag(fld[facei]) < SMALL)
1030  {
1031  fld[facei] = SMALL;
1032  }
1033  }
1034  }
1035  else
1036  {
1037  // Calculate manually
1038 
1039  const pointField& localPoints = points();
1040 
1041  auto iter = areas.field().begin();
1042 
1043  for (const face& f : faces())
1044  {
1045  *iter = f.mag(localPoints);
1046 
1047  // Do not allow any mag(val) < SMALL
1048  if (mag(*iter) < SMALL)
1049  {
1050  *iter = SMALL;
1051  }
1052 
1053  ++iter;
1054  }
1055  }
1056 }
1057 
1058 
1059 void Foam::faMesh::calcFaceAreaNormals() const
1060 {
1062  << "Calculating face area normals" << endl;
1063 
1064  if (faceAreaNormalsPtr_)
1065  {
1067  << "faceAreaNormals already allocated"
1068  << abort(FatalError);
1069  }
1070 
1071  faceAreaNormalsPtr_ =
1072  new areaVectorField
1073  (
1074  IOobject
1075  (
1076  "faceAreaNormals",
1077  mesh().pointsInstance(),
1079  faMesh::thisDb(),
1083  ),
1084  *this,
1085  dimless
1086  );
1087 
1088  areaVectorField& faceNormals = *faceAreaNormalsPtr_;
1089 
1090  const pointField& localPoints = points();
1091 
1092  // Internal
1093  {
1094  auto& fld = faceNormals.primitiveFieldRef();
1095 
1096  if (mesh().hasFaceAreas())
1097  {
1098  // The volume mesh has faceAreas, can reuse them
1099  fld = UIndirectList<vector>(mesh().faceAreas(), faceLabels());
1100  }
1101  else
1102  {
1103  // Calculate manually
1104 
1105  auto iter = fld.begin();
1106 
1107  for (const face& f : faces())
1108  {
1109  *iter = f.areaNormal(localPoints);
1110  ++iter;
1111  }
1112  }
1113 
1114  // Make unit normals
1115  fld.normalise();
1116 
1118  }
1119 
1120 
1121  // Boundary - copy from edges
1122  {
1123  const auto& edgeNormalsBoundary = edgeAreaNormals().boundaryField();
1124 
1125  forAll(boundary(), patchi)
1126  {
1127  faceNormals.boundaryFieldRef()[patchi]
1128  = edgeNormalsBoundary[patchi];
1129  }
1130  }
1131 
1132  // Parallel consistency, exchange on processor patches
1133  if (UPstream::parRun())
1134  {
1135  faceNormals.boundaryFieldRef().evaluateCoupled<processorFaPatch>();
1136  }
1137 }
1138 
1139 
1140 void Foam::faMesh::calcEdgeAreaNormals() const
1141 {
1143  << "Calculating edge area normals" << endl;
1144 
1145  if (edgeAreaNormalsPtr_)
1146  {
1148  << "edgeAreaNormalsPtr_ already allocated"
1149  << abort(FatalError);
1150  }
1151 
1152  edgeAreaNormalsPtr_ =
1153  new edgeVectorField
1154  (
1155  IOobject
1156  (
1157  "edgeAreaNormals",
1158  mesh().pointsInstance(),
1160  faMesh::thisDb(),
1164  ),
1165  *this,
1166  dimless
1167  // -> calculatedType()
1168  );
1169  edgeVectorField& edgeAreaNormals = *edgeAreaNormalsPtr_;
1170 
1171 
1172  if (faMesh::geometryOrder() < 2)
1173  {
1174  // The edge normals with flat boundary addressing
1175  // (which _may_ use communication)
1176  vectorField edgeNormals
1177  (
1178  calcRawEdgeNormals(faMesh::geometryOrder())
1179  );
1180 
1181  // Copy internal field
1182  edgeAreaNormals.primitiveFieldRef()
1183  = vectorField::subList(edgeNormals, nInternalEdges_);
1184 
1185  // Using our own slicing (instead of patchSlice) - see above
1186 
1187  auto& bfld = edgeAreaNormals.boundaryFieldRef();
1188 
1189  label patchOffset = nInternalEdges_;
1190 
1191  forAll(boundary(), patchi)
1192  {
1193  const faPatch& fap = boundary()[patchi];
1194 
1195  bfld[patchi] = vectorField::subList
1196  (
1197  edgeNormals,
1198  bfld[patchi].size(),
1199  patchOffset
1200  );
1201 
1202  patchOffset += fap.nEdges();
1203  }
1204 
1205  return;
1206  }
1207 
1208 
1209  // ------------------------------------------------------------------------
1210 
1211  // This is the original approach using an average of the
1212  // point normals. May be removed in the future (2022-09)
1213 
1214  // Starting from point area normals
1215  const vectorField& pointNormals = pointAreaNormals();
1216 
1217  // Also need local points
1218  const pointField& localPoints = points();
1219 
1220  // Internal edges
1221  {
1222  auto& fld = edgeAreaNormals.primitiveFieldRef();
1223 
1224  forAll(fld, edgei)
1225  {
1226  const edge& e = edges_[edgei];
1227  const linePointRef edgeLine(e.line(localPoints));
1228 
1229  // Average of both ends
1230  fld[edgei] = (pointNormals[e.first()] + pointNormals[e.second()]);
1231 
1232  fld[edgei].removeCollinear(edgeLine.unitVec());
1233  fld[edgei].normalise();
1234  }
1235 
1237  }
1238 
1239  // Boundary
1240  {
1241  auto& bfld = edgeAreaNormals.boundaryFieldRef();
1242 
1243  forAll(boundary(), patchi)
1244  {
1245  const faPatch& fap = boundary()[patchi];
1246 
1247  auto& pfld = bfld[patchi];
1248 
1249  const edgeList::subList bndEdges = fap.patchSlice(edges_);
1250 
1251  forAll(bndEdges, patchEdgei)
1252  {
1253  const edge& e = bndEdges[patchEdgei];
1254  const linePointRef edgeLine(e.line(localPoints));
1255 
1256  // Average of both ends
1257  pfld[patchEdgei] =
1258  (pointNormals[e.first()] + pointNormals[e.second()]);
1259 
1260  pfld[patchEdgei].removeCollinear(edgeLine.unitVec());
1261  pfld[patchEdgei].normalise();
1262  }
1263 
1264  imposeMinVectorLength(pfld);
1265  }
1266  }
1267 }
1268 
1269 
1270 void Foam::faMesh::calcFaceCurvatures() const
1271 {
1273  << "Calculating face curvatures" << endl;
1274 
1275  if (faceCurvaturesPtr_)
1276  {
1278  << "faceCurvaturesPtr_ already allocated"
1279  << abort(FatalError);
1280  }
1281 
1282  faceCurvaturesPtr_ =
1283  new areaScalarField
1284  (
1285  IOobject
1286  (
1287  "faceCurvatures",
1288  mesh().pointsInstance(),
1290  faMesh::thisDb(),
1294  ),
1295  *this,
1297  );
1298 
1299  areaScalarField& faceCurvatures = *faceCurvaturesPtr_;
1300 
1301 
1302 // faceCurvatures =
1303 // fac::edgeIntegrate(Le()*edgeLengthCorrection())
1304 // &faceAreaNormals();
1305 
1306  areaVectorField kN(fac::edgeIntegrate(Le()*edgeLengthCorrection()));
1307 
1308  faceCurvatures = sign(kN&faceAreaNormals())*mag(kN);
1309 }
1310 
1311 
1312 void Foam::faMesh::calcEdgeTransformTensors() const
1313 {
1315  << "Calculating edge transformation tensors" << endl;
1316 
1317  if (edgeTransformTensorsPtr_)
1318  {
1320  << "edgeTransformTensorsPtr_ already allocated"
1321  << abort(FatalError);
1322  }
1323 
1324  edgeTransformTensorsPtr_ = new FieldField<Field, tensor>(nEdges_);
1325  auto& edgeTransformTensors = *edgeTransformTensorsPtr_;
1326 
1327  // Initialize all transformation tensors
1328  for (label edgei = 0; edgei < nEdges_; ++edgei)
1329  {
1330  edgeTransformTensors.set(edgei, new Field<tensor>(3, tensor::I));
1331  }
1332 
1333  const areaVectorField& Nf = faceAreaNormals();
1334  const areaVectorField& Cf = areaCentres();
1335 
1336  const edgeVectorField& Ne = edgeAreaNormals();
1337  const edgeVectorField& Ce = edgeCentres();
1338 
1339  // Internal edges transformation tensors
1340  for (label edgei = 0; edgei < nInternalEdges_; ++edgei)
1341  {
1342  const label ownFacei = owner()[edgei];
1343  const label neiFacei = neighbour()[edgei];
1344  auto& tensors = edgeTransformTensors[edgei];
1345 
1346  vector edgeCtr = Ce.internalField()[edgei];
1347 
1348  if (skew())
1349  {
1350  edgeCtr -= skewCorrectionVectors().internalField()[edgei];
1351  }
1352 
1353  // Edge transformation tensor
1354  tensors[0] =
1356  (
1357  Ne.internalField()[edgei],
1358  (edgeCtr - Cf.internalField()[ownFacei])
1359  );
1360 
1361  // Owner transformation tensor
1362  tensors[1] =
1364  (
1365  Nf.internalField()[ownFacei],
1366  (edgeCtr - Cf.internalField()[ownFacei])
1367  );
1368 
1369  // Neighbour transformation tensor
1370  tensors[2] =
1372  (
1373  Nf.internalField()[neiFacei],
1374  (Cf.internalField()[neiFacei] - edgeCtr)
1375  );
1376  }
1377 
1378  // Boundary edges transformation tensors
1379  forAll(boundary(), patchi)
1380  {
1381  const labelUList& edgeFaces = boundary()[patchi].edgeFaces();
1382 
1383  if (boundary()[patchi].coupled())
1384  {
1385  vectorField nbrCf(Cf.boundaryField()[patchi].patchNeighbourField());
1386  vectorField nbrNf(Nf.boundaryField()[patchi].patchNeighbourField());
1387 
1388  forAll(edgeFaces, bndEdgei)
1389  {
1390  const label ownFacei = edgeFaces[bndEdgei];
1391  const label meshEdgei(boundary()[patchi].start() + bndEdgei);
1392 
1393  auto& tensors = edgeTransformTensors[meshEdgei];
1394 
1395  vector edgeCtr = Ce.boundaryField()[patchi][bndEdgei];
1396 
1397  if (skew())
1398  {
1399  edgeCtr -= skewCorrectionVectors()
1400  .boundaryField()[patchi][bndEdgei];
1401  }
1402 
1403  // Edge transformation tensor
1404  tensors[0] =
1406  (
1407  Ne.boundaryField()[patchi][bndEdgei],
1408  (edgeCtr - Cf.internalField()[ownFacei])
1409  );
1410 
1411  // Owner transformation tensor
1412  tensors[1] =
1414  (
1415  Nf.internalField()[ownFacei],
1416  (edgeCtr - Cf.internalField()[ownFacei])
1417  );
1418 
1419  // Neighbour transformation tensor
1420  tensors[2] =
1422  (
1423  nbrNf[bndEdgei],
1424  (nbrCf[bndEdgei] - edgeCtr)
1425  );
1426  }
1427  }
1428  else
1429  {
1430  forAll(edgeFaces, bndEdgei)
1431  {
1432  const label ownFacei = edgeFaces[bndEdgei];
1433  const label meshEdgei(boundary()[patchi].start() + bndEdgei);
1434 
1435  auto& tensors = edgeTransformTensors[meshEdgei];
1436 
1437  vector edgeCtr = Ce.boundaryField()[patchi][bndEdgei];
1438 
1439  if (skew())
1440  {
1441  edgeCtr -= skewCorrectionVectors()
1442  .boundaryField()[patchi][bndEdgei];
1443  }
1444 
1445  // Edge transformation tensor
1446  tensors[0] =
1448  (
1449  Ne.boundaryField()[patchi][bndEdgei],
1450  (edgeCtr - Cf.internalField()[ownFacei])
1451  );
1452 
1453  // Owner transformation tensor
1454  tensors[1] =
1456  (
1457  Nf.internalField()[ownFacei],
1458  (edgeCtr - Cf.internalField()[ownFacei])
1459  );
1460 
1461  // Neighbour transformation tensor
1462  tensors[2] = tensor::I;
1463  }
1464  }
1465  }
1466 }
1467 
1468 
1470 {
1472  << "Calculating internal points" << endl;
1473 
1474  bitSet markPoints(markupBoundaryPoints(this->patch()));
1475  markPoints.flip();
1476 
1477  return markPoints.sortedToc();
1478 }
1479 
1480 
1482 {
1484  << "Calculating boundary points" << endl;
1485 
1486  bitSet markPoints(markupBoundaryPoints(this->patch()));
1487 
1488  return markPoints.sortedToc();
1489 }
1490 
1491 
1492 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1493 // Point normal calculations
1494 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1495 
1496 // Revised method (general)
1497 // ------------------------
1498 //
1499 // - For each patch face obtain a centre point (mathematical avg)
1500 // and use that to define the face dual as a pair of triangles:
1501 // - tri1: base-point / mid-side of right edge / face centre
1502 // - tri2: base-point / face centre / mid-side of left edge
1503 //
1504 // - Walk all face points, inserting directly into the corresponding
1505 // locations. No distinction between internal or boundary points (yet).
1506 //
1507 // Revised method (boundary correction)
1508 // ------------------------------------
1509 //
1510 // - correct wedge directly, use processor patch information to exchange
1511 // the current summed values. [No change from original].
1512 //
1513 // - explicit correction of other boundaries.
1514 // Use the new boundary halo information for the face normals.
1515 // Calculate the equivalent face-point normals locally and apply
1516 // correction as before.
1517 
1518 void Foam::faMesh::calcPointAreaNormals(vectorField& result) const
1519 {
1521  << "Calculating pointAreaNormals : face dual method" << endl;
1522 
1523  result.resize_nocopy(nPoints());
1524  result = Zero;
1525 
1526  const pointField& points = patch().localPoints();
1527  const faceList& faces = patch().localFaces();
1528 
1529 
1530  // Loop over all faces
1531 
1532  for (const face& f : faces)
1533  {
1534  const label nVerts(f.size());
1535 
1536  point centrePoint(Zero);
1537  for (const label fp : f)
1538  {
1539  centrePoint += points[fp];
1540  }
1541  centrePoint /= nVerts;
1542 
1543  for (label i = 0; i < nVerts; ++i)
1544  {
1545  const label pt0 = f.thisLabel(i); // base
1546 
1547  result[pt0] +=
1549  (
1550  points[pt0], // base
1551  points[f.nextLabel(i)], // right
1552  points[f.prevLabel(i)], // left
1553  centrePoint
1554  );
1555  }
1556  }
1557 
1558 
1559  // Boundary edge corrections
1560  bitSet nbrBoundaryAdjust;
1561 
1562  forAll(boundary(), patchi)
1563  {
1564  const faPatch& fap = boundary()[patchi];
1565 
1566  if (isA<wedgeFaPatch>(fap))
1567  {
1568  // Correct wedge points
1569 
1570  const auto& wedgePatch = refCast<const wedgeFaPatch>(fap);
1571 
1572  const labelList& patchPoints = wedgePatch.pointLabels();
1573 
1574  const vector N
1575  (
1576  transform
1577  (
1578  wedgePatch.edgeT(),
1579  wedgePatch.centreNormal()
1580  ).normalise()
1581  );
1582 
1583  for (const label pti : patchPoints)
1584  {
1585  result[pti].removeCollinear(N);
1586  }
1587 
1588  // Axis point correction
1589  if (wedgePatch.axisPoint() > -1)
1590  {
1591  result[wedgePatch.axisPoint()] =
1592  wedgePatch.axis()
1593  *(
1594  wedgePatch.axis()
1595  &result[wedgePatch.axisPoint()]
1596  );
1597  }
1598  }
1599  else if (Pstream::parRun() && isA<processorFaPatch>(fap))
1600  {
1601  // Correct processor patch points
1602  const auto& procPatch = refCast<const processorFaPatch>(fap);
1603 
1604  const labelList& patchPoints = procPatch.pointLabels();
1605  const labelList& nbrPatchPoints = procPatch.neighbPoints();
1606 
1607  const labelList& nonGlobalPatchPoints =
1608  procPatch.nonGlobalPatchPoints();
1609 
1610  // Send my values
1611 
1612  vectorField patchPointNormals
1613  (
1614  UIndirectList<vector>(result, patchPoints)
1615  );
1616 
1617  {
1619  (
1621  procPatch.neighbProcNo(),
1622  patchPointNormals
1623  );
1624  }
1625 
1626  // Receive neighbour values
1627  patchPointNormals.resize_nocopy(nbrPatchPoints.size());
1628 
1629  {
1631  (
1633  procPatch.neighbProcNo(),
1634  patchPointNormals
1635  );
1636  }
1637 
1638  for (const label pti : nonGlobalPatchPoints)
1639  {
1640  result[patchPoints[pti]] +=
1641  patchPointNormals[nbrPatchPoints[pti]];
1642  }
1643  }
1644  // TBD:
1649  else if (correctPatchPointNormals(patchi) && !fap.coupled())
1650  {
1651  // Neighbour correction requested
1652  nbrBoundaryAdjust.set(patchi);
1653  }
1654  }
1655 
1656 
1657  // Correct global points
1658  if (globalData().nGlobalPoints())
1659  {
1660  const labelList& spLabels = globalData().sharedPointLabels();
1661  const labelList& addr = globalData().sharedPointAddr();
1662 
1663  vectorField spNormals
1664  (
1665  UIndirectList<vector>(result, spLabels)
1666  );
1667 
1668  vectorField gpNormals
1669  (
1670  globalData().nGlobalPoints(),
1671  Zero
1672  );
1673 
1674  forAll(addr, i)
1675  {
1676  gpNormals[addr[i]] += spNormals[i];
1677  }
1678 
1679  Pstream::combineReduce(gpNormals, plusEqOp<vectorField>());
1680 
1681  // Extract local data
1682  forAll(addr, i)
1683  {
1684  spNormals[i] = gpNormals[addr[i]];
1685  }
1686 
1687  forAll(spNormals, pointI)
1688  {
1689  result[spLabels[pointI]] = spNormals[pointI];
1690  }
1691  }
1692 
1693 
1694  if
1695  (
1696  hasPatchPointNormalsCorrection()
1697  && returnReduceOr(nbrBoundaryAdjust.any())
1698  )
1699  {
1701  << "Apply " << nbrBoundaryAdjust.count()
1702  << " boundary neighbour corrections" << nl;
1703 
1704  // Apply boundary points correction
1705 
1706  // Collect face normals as point normals
1707 
1708  const auto& haloNormals = this->haloFaceNormals();
1709 
1710  Map<vector> fpNormals(4*nBoundaryEdges());
1711 
1712  for (const label patchi : nbrBoundaryAdjust)
1713  {
1714  const faPatch& fap = boundary()[patchi];
1715 
1716  if (fap.ngbPolyPatchIndex() < 0)
1717  {
1719  << "Neighbour polyPatch index is not defined "
1720  << "for faPatch " << fap.name()
1721  << abort(FatalError);
1722  }
1723 
1724  // NB: haloFaceNormals uses primitivePatch edge indexing
1725  for (const label edgei : fap.edgeLabels())
1726  {
1727  const edge& e = patch().edges()[edgei];
1728 
1729  // Halo face unitNormal at boundary edge (starts as 0)
1730  const vector& fnorm = haloNormals[edgei - nInternalEdges_];
1731 
1732  fpNormals(e.first()) += fnorm;
1733  fpNormals(e.second()) += fnorm;
1734  }
1735  }
1736 
1737  // Apply the correction
1738 
1739  // Note from Zeljko Tukovic:
1740  //
1741  // This posibility is used for free-surface tracking
1742  // calculations to enforce 90 deg contact angle between
1743  // finite-area mesh and neighbouring polyPatch. It is very
1744  // important for curvature calculation to have correct normal
1745  // at contact line points.
1746 
1747  forAllConstIters(fpNormals, iter)
1748  {
1749  const label pointi = iter.key();
1750  result[pointi].removeCollinear(normalised(iter.val()));
1751  }
1752  }
1753 
1754  result.normalise();
1755 }
1756 
1757 
1758 void Foam::faMesh::calcPointAreaNormalsByQuadricsFit(vectorField& result) const
1759 {
1760  const labelList intPoints(internalPoints());
1761  const labelList bndPoints(boundaryPoints());
1762 
1763  const pointField& points = patch().localPoints();
1764  const faceList& faces = patch().localFaces();
1765  const labelListList& pointFaces = patch().pointFaces();
1766 
1767  forAll(intPoints, pointI)
1768  {
1769  label curPoint = intPoints[pointI];
1770 
1771  const labelHashSet faceSet(pointFaces[curPoint]);
1772  const labelList curFaces(faceSet.toc());
1773 
1774  labelHashSet pointSet;
1775 
1776  for (const label facei : curFaces)
1777  {
1778  const labelList& facePoints = faces[facei];
1779  pointSet.insert(facePoints);
1780  }
1781  pointSet.erase(curPoint);
1782  labelList curPoints(pointSet.toc());
1783 
1784  if (pointSet.size() < 5)
1785  {
1786  DebugInfo
1787  << "WARNING: Extending point set for fitting." << endl;
1788 
1789  labelHashSet faceSet(pointFaces[curPoint]);
1790  labelList curFaces(faceSet.toc());
1791  for (const label facei : curFaces)
1792  {
1793  const labelList& curFaceFaces = patch().faceFaces()[facei];
1794  faceSet.insert(curFaceFaces);
1795  }
1796  curFaces = faceSet.toc();
1797 
1798  pointSet.clear();
1799 
1800  for (const label facei : curFaces)
1801  {
1802  const labelList& facePoints = faces[facei];
1803  pointSet.insert(facePoints);
1804  }
1805  pointSet.erase(curPoint);
1806  curPoints = pointSet.toc();
1807  }
1808 
1809  pointField allPoints(curPoints.size());
1810  scalarField W(curPoints.size(), 1.0);
1811  for (label i=0; i<curPoints.size(); ++i)
1812  {
1813  allPoints[i] = points[curPoints[i]];
1814  W[i] = 1.0/magSqr(allPoints[i] - points[curPoint]);
1815  }
1816 
1817  // Transform points
1818  coordSystem::cartesian cs
1819  (
1820  points[curPoint], // origin
1821  result[curPoint], // axis [e3] (normalized by constructor)
1822  allPoints[0] - points[curPoint] // direction [e1]
1823  );
1824 
1825  for (point& p : allPoints)
1826  {
1827  p = cs.localPosition(p);
1828  }
1829 
1831  (
1832  allPoints.size(),
1833  5,
1834  0.0
1835  );
1836 
1837  for (label i = 0; i < allPoints.size(); ++i)
1838  {
1839  M[i][0] = sqr(allPoints[i].x());
1840  M[i][1] = sqr(allPoints[i].y());
1841  M[i][2] = allPoints[i].x()*allPoints[i].y();
1842  M[i][3] = allPoints[i].x();
1843  M[i][4] = allPoints[i].y();
1844  }
1845 
1846  scalarSquareMatrix MtM(5, Zero);
1847 
1848  for (label i = 0; i < MtM.n(); ++i)
1849  {
1850  for (label j = 0; j < MtM.m(); ++j)
1851  {
1852  for (label k = 0; k < M.n(); ++k)
1853  {
1854  MtM[i][j] += M[k][i]*M[k][j]*W[k];
1855  }
1856  }
1857  }
1858 
1859  scalarField MtR(5, Zero);
1860 
1861  for (label i=0; i<MtR.size(); ++i)
1862  {
1863  for (label j=0; j<M.n(); ++j)
1864  {
1865  MtR[i] += M[j][i]*allPoints[j].z()*W[j];
1866  }
1867  }
1868 
1869  LUsolve(MtM, MtR);
1870 
1871  vector curNormal = vector(MtR[3], MtR[4], -1);
1872 
1873  curNormal = cs.globalVector(curNormal);
1874 
1875  curNormal *= sign(curNormal&result[curPoint]);
1876 
1877  result[curPoint] = curNormal;
1878  }
1879 
1880 
1881  forAll(boundary(), patchI)
1882  {
1883  const faPatch& fap = boundary()[patchI];
1884 
1885  if (Pstream::parRun() && isA<processorFaPatch>(fap))
1886  {
1887  const processorFaPatch& procPatch =
1888  refCast<const processorFaPatch>(boundary()[patchI]);
1889 
1890  const labelList& patchPointLabels = procPatch.pointLabels();
1891 
1892  labelList toNgbProcLsPointStarts(patchPointLabels.size(), Zero);
1893  vectorField toNgbProcLsPoints
1894  (
1895  10*patchPointLabels.size(),
1896  Zero
1897  );
1898  label nPoints = 0;
1899 
1900  for (label pointI=0; pointI<patchPointLabels.size(); ++pointI)
1901  {
1902  label curPoint = patchPointLabels[pointI];
1903 
1904  toNgbProcLsPointStarts[pointI] = nPoints;
1905 
1906  const labelHashSet faceSet(pointFaces[curPoint]);
1907  const labelList curFaces(faceSet.toc());
1908 
1909  labelHashSet pointSet;
1910 
1911  for (const label facei : curFaces)
1912  {
1913  const labelList& facePoints = faces[facei];
1914  pointSet.insert(facePoints);
1915  }
1916  pointSet.erase(curPoint);
1917  labelList curPoints = pointSet.toc();
1918 
1919  for (label i=0; i<curPoints.size(); ++i)
1920  {
1921  toNgbProcLsPoints[nPoints++] = points[curPoints[i]];
1922  }
1923  }
1924 
1925  toNgbProcLsPoints.setSize(nPoints);
1926 
1927  {
1928  OPstream toNeighbProc
1929  (
1931  procPatch.neighbProcNo(),
1932  toNgbProcLsPoints.size_bytes()
1933  + toNgbProcLsPointStarts.size_bytes()
1934  + 10*sizeof(label)
1935  );
1936 
1937  toNeighbProc
1938  << toNgbProcLsPoints
1939  << toNgbProcLsPointStarts;
1940  }
1941  }
1942  }
1943 
1944  for (const faPatch& fap : boundary())
1945  {
1946  if (Pstream::parRun() && isA<processorFaPatch>(fap))
1947  {
1948  const processorFaPatch& procPatch =
1949  refCast<const processorFaPatch>(fap);
1950 
1951  const labelList& patchPointLabels = procPatch.pointLabels();
1952 
1953  labelList fromNgbProcLsPointStarts(patchPointLabels.size(), Zero);
1954  vectorField fromNgbProcLsPoints;
1955 
1956  {
1957  IPstream fromNeighbProc
1958  (
1960  procPatch.neighbProcNo(),
1961  10*patchPointLabels.size()*sizeof(vector)
1962  + fromNgbProcLsPointStarts.size_bytes()
1963  + 10*sizeof(label)
1964  );
1965 
1966  fromNeighbProc
1967  >> fromNgbProcLsPoints
1968  >> fromNgbProcLsPointStarts;
1969  }
1970 
1971  const labelList& nonGlobalPatchPoints =
1972  procPatch.nonGlobalPatchPoints();
1973 
1974  forAll(nonGlobalPatchPoints, pointI)
1975  {
1976  label curPoint =
1977  patchPointLabels[nonGlobalPatchPoints[pointI]];
1978  label curNgbPoint =
1979  procPatch.neighbPoints()[nonGlobalPatchPoints[pointI]];
1980 
1981  labelHashSet faceSet;
1982  faceSet.insert(pointFaces[curPoint]);
1983 
1984  labelList curFaces = faceSet.toc();
1985 
1986  labelHashSet pointSet;
1987 
1988  for (const label facei : curFaces)
1989  {
1990  const labelList& facePoints = faces[facei];
1991  pointSet.insert(facePoints);
1992  }
1993  pointSet.erase(curPoint);
1994  labelList curPoints = pointSet.toc();
1995 
1996  label nAllPoints = curPoints.size();
1997 
1998  if (curNgbPoint == fromNgbProcLsPointStarts.size() - 1)
1999  {
2000  nAllPoints +=
2001  fromNgbProcLsPoints.size()
2002  - fromNgbProcLsPointStarts[curNgbPoint];
2003  }
2004  else
2005  {
2006  nAllPoints +=
2007  fromNgbProcLsPointStarts[curNgbPoint + 1]
2008  - fromNgbProcLsPointStarts[curNgbPoint];
2009  }
2010 
2011  vectorField allPointsExt(nAllPoints);
2012  label counter = 0;
2013  for (label i=0; i<curPoints.size(); ++i)
2014  {
2015  allPointsExt[counter++] = points[curPoints[i]];
2016  }
2017 
2018  if (curNgbPoint == fromNgbProcLsPointStarts.size() - 1)
2019  {
2020  for
2021  (
2022  label i=fromNgbProcLsPointStarts[curNgbPoint];
2023  i<fromNgbProcLsPoints.size();
2024  ++i
2025  )
2026  {
2027  allPointsExt[counter++] = fromNgbProcLsPoints[i];
2028  }
2029  }
2030  else
2031  {
2032  for
2033  (
2034  label i=fromNgbProcLsPointStarts[curNgbPoint];
2035  i<fromNgbProcLsPointStarts[curNgbPoint+1];
2036  ++i
2037  )
2038  {
2039  allPointsExt[counter++] = fromNgbProcLsPoints[i];
2040  }
2041  }
2042 
2043  // Remove duplicate points
2044  vectorField allPoints(nAllPoints, Zero);
2045  const scalar tol = 0.001*treeBoundBox(allPointsExt).mag();
2046 
2047  nAllPoints = 0;
2048  for (const point& pt : allPointsExt)
2049  {
2050  bool duplicate = false;
2051  for (label i = 0; i < nAllPoints; ++i)
2052  {
2053  if (pt.dist(allPoints[i]) < tol)
2054  {
2055  duplicate = true;
2056  break;
2057  }
2058  }
2059 
2060  if (!duplicate)
2061  {
2062  allPoints[nAllPoints++] = pt;
2063  }
2064  }
2065 
2066  allPoints.setSize(nAllPoints);
2067 
2068  if (nAllPoints < 5)
2069  {
2071  << "There are no enough points for quadrics "
2072  << "fitting for a point at processor patch"
2073  << abort(FatalError);
2074  }
2075 
2076  // Transform points
2077  const vector& origin = points[curPoint];
2078  const vector axis = normalised(result[curPoint]);
2079  vector dir(allPoints[0] - points[curPoint]);
2080  dir.removeCollinear(axis);
2081  dir.normalise();
2082 
2083  const coordinateSystem cs("cs", origin, axis, dir);
2084 
2085  scalarField W(allPoints.size(), 1.0);
2086 
2087  forAll(allPoints, pI)
2088  {
2089  W[pI] = 1.0/magSqr(allPoints[pI] - points[curPoint]);
2090 
2091  allPoints[pI] = cs.localPosition(allPoints[pI]);
2092  }
2093 
2095  (
2096  allPoints.size(),
2097  5,
2098  0.0
2099  );
2100 
2101  for (label i=0; i<allPoints.size(); ++i)
2102  {
2103  M[i][0] = sqr(allPoints[i].x());
2104  M[i][1] = sqr(allPoints[i].y());
2105  M[i][2] = allPoints[i].x()*allPoints[i].y();
2106  M[i][3] = allPoints[i].x();
2107  M[i][4] = allPoints[i].y();
2108  }
2109 
2110  scalarSquareMatrix MtM(5, Zero);
2111 
2112  for (label i = 0; i < MtM.n(); ++i)
2113  {
2114  for (label j = 0; j < MtM.m(); ++j)
2115  {
2116  for (label k = 0; k < M.n(); ++k)
2117  {
2118  MtM[i][j] += M[k][i]*M[k][j]*W[k];
2119  }
2120  }
2121  }
2122 
2123  scalarField MtR(5, Zero);
2124 
2125  for (label i = 0; i < MtR.size(); ++i)
2126  {
2127  for (label j = 0; j < M.n(); ++j)
2128  {
2129  MtR[i] += M[j][i]*allPoints[j].z()*W[j];
2130  }
2131  }
2132 
2133  LUsolve(MtM, MtR);
2134 
2135  vector curNormal = vector(MtR[3], MtR[4], -1);
2136 
2137  curNormal = cs.globalVector(curNormal);
2138 
2139  curNormal *= sign(curNormal&result[curPoint]);
2140 
2141  result[curPoint] = curNormal;
2142  }
2143  }
2144  }
2145 
2146  // Correct global points
2147  if (globalData().nGlobalPoints() > 0)
2148  {
2149  const labelList& spLabels = globalData().sharedPointLabels();
2150 
2151  const labelList& addr = globalData().sharedPointAddr();
2152 
2153  for (label k=0; k<globalData().nGlobalPoints(); ++k)
2154  {
2155  List<List<vector>> procLsPoints(Pstream::nProcs());
2156 
2157  const label curSharedPointIndex = addr.find(k);
2158 
2159  scalar tol = 0.0;
2160 
2161  if (curSharedPointIndex != -1)
2162  {
2163  label curPoint = spLabels[curSharedPointIndex];
2164 
2165  const labelHashSet faceSet(pointFaces[curPoint]);
2166  const labelList curFaces(faceSet.toc());
2167 
2168  labelHashSet pointSet;
2169  for (const label facei : curFaces)
2170  {
2171  const labelList& facePoints = faces[facei];
2172  pointSet.insert(facePoints);
2173  }
2174  pointSet.erase(curPoint);
2175  labelList curPoints = pointSet.toc();
2176 
2177  vectorField locPoints(points, curPoints);
2178 
2179  procLsPoints[Pstream::myProcNo()] = locPoints;
2180 
2181  tol = 0.001*treeBoundBox(locPoints).mag();
2182  }
2183 
2184  Pstream::allGatherList(procLsPoints);
2185 
2186  if (curSharedPointIndex != -1)
2187  {
2188  label curPoint = spLabels[curSharedPointIndex];
2189 
2190  label nAllPoints = 0;
2191  forAll(procLsPoints, procI)
2192  {
2193  nAllPoints += procLsPoints[procI].size();
2194  }
2195 
2196  vectorField allPoints(nAllPoints, Zero);
2197 
2198  nAllPoints = 0;
2199  forAll(procLsPoints, procI)
2200  {
2201  forAll(procLsPoints[procI], pointI)
2202  {
2203  bool duplicate = false;
2204  for (label i=0; i<nAllPoints; ++i)
2205  {
2206  if
2207  (
2208  mag
2209  (
2210  allPoints[i]
2211  - procLsPoints[procI][pointI]
2212  )
2213  < tol
2214  )
2215  {
2216  duplicate = true;
2217  break;
2218  }
2219  }
2220 
2221  if (!duplicate)
2222  {
2223  allPoints[nAllPoints++] =
2224  procLsPoints[procI][pointI];
2225  }
2226  }
2227  }
2228 
2229  allPoints.setSize(nAllPoints);
2230 
2231  if (nAllPoints < 5)
2232  {
2234  << "There are no enough points for quadratic "
2235  << "fitting for a global processor point "
2236  << abort(FatalError);
2237  }
2238 
2239  // Transform points
2240  const vector& origin = points[curPoint];
2241  const vector axis = normalised(result[curPoint]);
2242  vector dir(allPoints[0] - points[curPoint]);
2243  dir.removeCollinear(axis);
2244  dir.normalise();
2245 
2246  coordinateSystem cs("cs", origin, axis, dir);
2247 
2248  scalarField W(allPoints.size(), 1.0);
2249 
2250  forAll(allPoints, pointI)
2251  {
2252  W[pointI]=
2253  1.0/magSqr(allPoints[pointI] - points[curPoint]);
2254 
2255  allPoints[pointI] = cs.localPosition(allPoints[pointI]);
2256  }
2257 
2259  (
2260  allPoints.size(),
2261  5,
2262  0.0
2263  );
2264 
2265  for (label i=0; i<allPoints.size(); ++i)
2266  {
2267  M[i][0] = sqr(allPoints[i].x());
2268  M[i][1] = sqr(allPoints[i].y());
2269  M[i][2] = allPoints[i].x()*allPoints[i].y();
2270  M[i][3] = allPoints[i].x();
2271  M[i][4] = allPoints[i].y();
2272  }
2273 
2274  scalarSquareMatrix MtM(5, Zero);
2275  for (label i = 0; i < MtM.n(); ++i)
2276  {
2277  for (label j = 0; j < MtM.m(); ++j)
2278  {
2279  for (label k = 0; k < M.n(); ++k)
2280  {
2281  MtM[i][j] += M[k][i]*M[k][j]*W[k];
2282  }
2283  }
2284  }
2285 
2286  scalarField MtR(5, Zero);
2287  for (label i = 0; i < MtR.size(); ++i)
2288  {
2289  for (label j = 0; j < M.n(); ++j)
2290  {
2291  MtR[i] += M[j][i]*allPoints[j].z()*W[j];
2292  }
2293  }
2294 
2295  LUsolve(MtM, MtR);
2296 
2297  vector curNormal = vector(MtR[3], MtR[4], -1);
2298 
2299  curNormal = cs.globalVector(curNormal);
2300 
2301  curNormal *= sign(curNormal&result[curPoint]);
2302 
2303  result[curPoint] = curNormal;
2304  }
2305  }
2306  }
2307 
2308  for (vector& n : result)
2309  {
2310  n.normalise();
2311  }
2312 }
2313 
2314 
2316 {
2318  << "Calculating edge length correction" << endl;
2319 
2320  auto tcorrection = tmp<edgeScalarField>::New
2321  (
2322  IOobject
2323  (
2324  "edgeLengthCorrection",
2325  mesh().pointsInstance(),
2327  faMesh::thisDb(),
2331  ),
2332  *this,
2333  dimless
2334  );
2335  auto& correction = tcorrection.ref();
2336 
2337  const vectorField& pointNormals = pointAreaNormals();
2338 
2339  const auto angleCorrection =
2340  [](const vector& a, const vector& b) -> scalar
2341  {
2342  return Foam::cos(0.5*Foam::asin(Foam::mag(a ^ b)));
2343  };
2344 
2345 
2346  // Internal
2347  {
2348  auto& fld = correction.primitiveFieldRef();
2349 
2350  forAll(fld, edgei)
2351  {
2352  fld[edgei] = angleCorrection
2353  (
2354  pointNormals[edges_[edgei].start()],
2355  pointNormals[edges_[edgei].end()]
2356  );
2357  }
2358  }
2359 
2360  // Boundary
2361  {
2362  auto& bfld = correction.boundaryFieldRef();
2363 
2364  forAll(boundary(), patchi)
2365  {
2366  const faPatch& fap = boundary()[patchi];
2367  scalarField& pfld = bfld[patchi];
2368 
2369  label edgei = fap.start();
2370 
2371  forAll(pfld, patchEdgei)
2372  {
2373  pfld[patchEdgei] = angleCorrection
2374  (
2375  pointNormals[edges_[edgei].start()],
2376  pointNormals[edges_[edgei].end()]
2377  );
2378 
2379  ++edgei;
2380  }
2381  }
2382  }
2383 
2384  return tcorrection;
2385 }
2386 
2387 
2389 {
2390  auto tunitVectors = tmp<edgeVectorField>::New
2391  (
2392  IOobject
2393  (
2394  "unit(Le)",
2395  mesh().pointsInstance(),
2397  faMesh::thisDb(),
2401  ),
2402  *this,
2403  dimless,
2404  (this->Le() / this->magLe())
2405  );
2406 
2407  // The above is not quite correct when a min-length limiter
2408  // has been imposed on both Le() and magLe() fields
2409 
2410  // tunitVectors.ref().oriented() = this->Le().oriented();
2411  return tunitVectors;
2412 }
2413 
2414 
2415 // ************************************************************************* //
faceListList boundary
dimensionedScalar sign(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
labelList internalPoints() const
Return internal point labels.
Point unitVec() const
Return the unit vector (start-to-end)
Definition: lineI.H:122
const polyBoundaryMesh & pbm
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
"blocking" : (MPI_Bsend, MPI_Recv)
A line primitive.
Definition: line.H:52
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:504
tmp< edgeScalarField > edgeLengthCorrection() const
Return edge length correction.
static label read(const UPstream::commsTypes commsType, const int fromProcNo, char *buf, const std::streamsize bufSize, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm, UPstream::Request *req=nullptr)
Read buffer contents from given processor.
Definition: UIPstreamRead.C:35
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar mag() const
The magnitude (length) of the line.
Definition: lineI.H:96
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedTensor skew(const dimensionedTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
labelList sortedToc() const
The indices of the on bits as a sorted labelList.
Definition: bitSetI.H:447
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
labelPair whichPatchFace(const label meshFacei) const
Lookup mesh face index and return (patchi, patchFacei) tuple or (-1, meshFacei) for internal faces...
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
SubList< vector > subList
Declare type of subList.
Definition: List.H:144
static vector areaInvDistSqrWeightedNormalFaceTriangle(const linePointRef &edgeLine, const point &faceCentre)
label k
Boltzmann constant.
Ignore writing from objectRegistry::writeObject()
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const dimensionSet dimless
Dimensionless.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
labelList faceLabels(nFaceLabels)
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Namespace of functions to calculate explicit derivatives.
static vector areaNormalFaceTriangle(const linePointRef &edgeLine, const point &faceCentre)
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:815
dimensionedScalar asin(const dimensionedScalar &ds)
static void waitRequests()
Wait for all requests to finish.
Definition: UPstream.H:1538
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
static void imposeMinVectorLength(vectorField &fld)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static int geometryOrder() noexcept
Return the current geometry treatment.
Definition: faMesh.H:797
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
word timeName
Definition: getTimeIndex.H:3
scalar y
A list of faces which address into the list of points.
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
line< point, const point & > linePointRef
A line using referred points.
Definition: line.H:66
const pointField & points
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
Definition: SubList.H:246
static const Identity< scalar > I
Definition: Identity.H:100
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
Definition: edgeFieldsFwd.H:58
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
PointRef a() const noexcept
The first point.
Definition: line.H:222
Vector< scalar > vector
Definition: vector.H:57
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:46
PointRef b() const noexcept
The second point.
Definition: line.H:227
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
static void combineReduce(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
#define DebugInfo
Report an information message using Foam::Info.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
static Foam::bitSet markupBoundaryPoints(const uindirectPrimitivePatch &p)
labelList f(nPoints)
bool hasFaceAreas() const noexcept
tmp< GeometricField< Type, faPatchField, areaMesh > > edgeIntegrate(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
const Vector< label > N(dict.get< Vector< label >>("N"))
Point centre() const
Return centre (centroid)
Definition: lineI.H:83
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Point vec() const
Return start-to-end vector.
Definition: lineI.H:109
const vectorField & faceCentres() const
labelList boundaryPoints() const
Return boundary point labels.
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:701
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
static vector calcLeVector(const point &faceCentre, const linePointRef &edgeLine, const vector &edgeNormal)
Calculate the &#39;Le&#39; vector from faceCentre to edge centre.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
const vectorField & faceAreas() const
Nothing to be read.
const std::string patch
OpenFOAM patch number as a std::string.
static bool write(const UPstream::commsTypes commsType, const int toProcNo, const char *buf, const std::streamsize bufSize, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm, UPstream::Request *req=nullptr, const UPstream::sendModes sendMode=UPstream::sendModes::normal)
Write buffer contents to given processor.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
label n
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool coupled
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Vector< Cmpt > & removeCollinear(const Vector< Cmpt > &unitVec)
Inplace removal of components that are collinear to the given unit vector.
Definition: VectorI.H:136
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Tensor of scalars, i.e. Tensor<scalar>.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< edgeVectorField > unitLe() const
Return normalised edge length vectors.
SquareMatrix< scalar > scalarSquareMatrix
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:81
static vector areaInvDistSqrWeightedNormalDualEdge(const point &basePoint, const point &rightPoint, const point &leftPoint, const point &centrePoint)
#define M(I)
Do not request registration (bool: false)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
static tensor rotation_e3e1(const vector &unitAxis1, vector dirn)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting returning the LU form of the matrix and the sol...
RectangularMatrix< scalar > scalarRectangularMatrix
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:151
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:80
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
static vector areaInvDistSqrWeightedNormal(const vector &a, const vector &b)