edgeInterpolation.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) 2020-2022 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 "areaFields.H"
31 #include "edgeFields.H"
32 #include "demandDrivenData.H"
33 #include "faPatchFields.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
43 
45 {
47  deleteDemandDrivenData(weightingFactors_);
48  deleteDemandDrivenData(differenceFactors_);
49  deleteDemandDrivenData(correctionVectors_);
50  deleteDemandDrivenData(skewCorrectionVectors_);
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 :
58  faMesh_(fam),
59  lPN_(nullptr),
60  weightingFactors_(nullptr),
61  differenceFactors_(nullptr),
62  correctionVectors_(nullptr),
63  skewCorrectionVectors_(nullptr),
64  orthogonal_(false),
65  skew_(true)
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70 
72 {
73  clearOut();
74 }
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  if (!lPN_)
82  {
83  makeLPN();
84  }
85 
86  return (*lPN_);
87 }
88 
89 
91 {
92  if (!weightingFactors_)
93  {
94  makeWeights();
95  }
96 
97  return (*weightingFactors_);
98 }
99 
100 
102 {
103  if (!differenceFactors_)
104  {
105  makeDeltaCoeffs();
106  }
107 
108  return (*differenceFactors_);
109 }
110 
111 
113 {
114  if (orthogonal_ == false && !correctionVectors_)
115  {
116  makeCorrectionVectors();
117  }
118 
119  return orthogonal_;
120 }
121 
122 
124 {
125  if (orthogonal())
126  {
128  << "cannot return correctionVectors; mesh is orthogonal"
130  }
131 
132  return (*correctionVectors_);
133 }
134 
135 
137 {
138  if (skew_ == true && !skewCorrectionVectors_)
139  {
140  makeSkewCorrectionVectors();
141  }
142 
143  return skew_;
144 }
145 
146 
149 {
150  if (!skew())
151  {
153  << "cannot return skewCorrectionVectors; mesh is now skewed"
155  }
156 
157  return (*skewCorrectionVectors_);
158 }
159 
160 
162 {
164  deleteDemandDrivenData(weightingFactors_);
165  deleteDemandDrivenData(differenceFactors_);
166 
167  orthogonal_ = false;
168  deleteDemandDrivenData(correctionVectors_);
169 
170  skew_ = true;
171  deleteDemandDrivenData(skewCorrectionVectors_);
172 
173  return true;
174 }
175 
176 
177 const Foam::vector& Foam::edgeInterpolation::skewCorr(const label edgeI) const
178 {
179  #ifdef FA_SKEW_CORRECTION
180 
181  return
182  (
183  skewCorrectionVectors_
184  ? (*skewCorrectionVectors_)[edgeI]
185  : pTraits<vector>::zero
186  );
187 
188  #else
189 
190  return (*skewCorrectionVectors_)[edgeI];
191 
192  #endif
193 }
194 
195 
196 void Foam::edgeInterpolation::makeLPN() const
197 {
199  << "Constructing geodesic distance between points P and N"
200  << endl;
201 
202 
203  lPN_ = new edgeScalarField
204  (
205  IOobject
206  (
207  "lPN",
208  mesh().time().constant(),
209  mesh().thisDb(),
213  ),
214  mesh(),
215  dimLength
216  );
217  edgeScalarField& lPN = *lPN_;
218 
219  // Set local references to mesh data
220  const edgeVectorField& edgeCentres = mesh().edgeCentres();
221  const areaVectorField& faceCentres = mesh().areaCentres();
222  const labelUList& owner = mesh().owner();
223  const labelUList& neighbour = mesh().neighbour();
224 
225  scalarField& lPNIn = lPN.primitiveFieldRef();
226 
227  // Calculate skewness correction vectors if necessary
228  (void) skew();
229 
230  forAll(owner, edgeI)
231  {
232  const vector& skewCorrEdge = skewCorr(edgeI);
233 
234  scalar lPE =
235  mag
236  (
237  edgeCentres[edgeI]
238  - skewCorrEdge
239  - faceCentres[owner[edgeI]]
240  );
241 
242  scalar lEN =
243  mag
244  (
245  faceCentres[neighbour[edgeI]]
246  - edgeCentres[edgeI]
247  + skewCorrEdge
248  );
249 
250  lPNIn[edgeI] = (lPE + lEN);
251 
252  // Do not allow any mag(val) < SMALL
253  if (mag(lPNIn[edgeI]) < SMALL)
254  {
255  lPNIn[edgeI] = SMALL;
256  }
257  }
258 
259 
260  forAll(lPN.boundaryField(), patchI)
261  {
262  mesh().boundary()[patchI].makeDeltaCoeffs
263  (
264  lPN.boundaryFieldRef()[patchI]
265  );
266 
267  lPN.boundaryFieldRef()[patchI] = 1.0/lPN.boundaryField()[patchI];
268  }
269 
270 
272  << "Finished constructing geodesic distance PN"
273  << endl;
274 }
275 
276 
277 void Foam::edgeInterpolation::makeWeights() const
278 {
280  << "Constructing weighting factors for edge interpolation"
281  << endl;
282 
283 
284  weightingFactors_ = new edgeScalarField
285  (
286  IOobject
287  (
288  "weightingFactors",
289  mesh().pointsInstance(),
290  mesh().thisDb(),
294  ),
295  mesh(),
297  );
298  edgeScalarField& weightingFactors = *weightingFactors_;
299 
300 
301  // Set local references to mesh data
302  const edgeVectorField& edgeCentres = mesh().edgeCentres();
303  const areaVectorField& faceCentres = mesh().areaCentres();
304  const labelUList& owner = mesh().owner();
305  const labelUList& neighbour = mesh().neighbour();
306 
307  scalarField& weightingFactorsIn = weightingFactors.primitiveFieldRef();
308 
309  // Calculate skewness correction vectors if necessary
310  (void) skew();
311 
312  forAll(owner, edgeI)
313  {
314  const vector& skewCorrEdge = skewCorr(edgeI);
315 
316  scalar lPE =
317  mag
318  (
319  edgeCentres[edgeI]
320  - skewCorrEdge
321  - faceCentres[owner[edgeI]]
322  );
323 
324  scalar lEN =
325  mag
326  (
327  faceCentres[neighbour[edgeI]]
328  - edgeCentres[edgeI]
329  + skewCorrEdge
330  );
331 
332  // weight = (0,1]
333  const scalar lPN = lPE + lEN;
334  if (mag(lPN) > SMALL)
335  {
336  weightingFactorsIn[edgeI] = lEN/lPN;
337  }
338  }
339 
340  forAll(mesh().boundary(), patchI)
341  {
342  mesh().boundary()[patchI].makeWeights
343  (
344  weightingFactors.boundaryFieldRef()[patchI]
345  );
346  }
347 
349  << "Finished constructing weighting factors for face interpolation"
350  << endl;
351 }
352 
353 
354 void Foam::edgeInterpolation::makeDeltaCoeffs() const
355 {
357  << "Constructing differencing factors array for edge gradient"
358  << endl;
359 
360  // Force the construction of the weighting factors
361  // needed to make sure deltaCoeffs are calculated for parallel runs.
362  weights();
363 
364  differenceFactors_ = new edgeScalarField
365  (
366  IOobject
367  (
368  "differenceFactors_",
369  mesh().pointsInstance(),
370  mesh().thisDb(),
374  ),
375  mesh(),
377  );
378  edgeScalarField& DeltaCoeffs = *differenceFactors_;
379  scalarField& dc = DeltaCoeffs.primitiveFieldRef();
380 
381  // Set local references to mesh data
382  const edgeVectorField& edgeCentres = mesh().edgeCentres();
383  const areaVectorField& faceCentres = mesh().areaCentres();
384  const labelUList& owner = mesh().owner();
385  const labelUList& neighbour = mesh().neighbour();
386  const edgeVectorField& lengths = mesh().Le();
387 
388  const edgeList& edges = mesh().edges();
389  const pointField& points = mesh().points();
390 
391  // Calculate skewness correction vectors if necessary
392  (void) skew();
393 
394  forAll(owner, edgeI)
395  {
396  // Edge normal - area normal
397  vector edgeNormal =
398  normalised(lengths[edgeI] ^ edges[edgeI].vec(points));
399 
400  // Unit delta vector
401  vector unitDelta =
402  faceCentres[neighbour[edgeI]]
403  - faceCentres[owner[edgeI]];
404 
405  unitDelta.removeCollinear(edgeNormal);
406  unitDelta.normalise();
407 
408 
409  const vector& skewCorrEdge = skewCorr(edgeI);
410 
411  scalar lPE =
412  mag
413  (
414  edgeCentres[edgeI]
415  - skewCorrEdge
416  - faceCentres[owner[edgeI]]
417  );
418 
419  scalar lEN =
420  mag
421  (
422  faceCentres[neighbour[edgeI]]
423  - edgeCentres[edgeI]
424  + skewCorrEdge
425  );
426 
427  scalar lPN = lPE + lEN;
428 
429 
430  // Edge normal - area tangent
431  edgeNormal = normalised(lengths[edgeI]);
432 
433  // Do not allow any mag(val) < SMALL
434  const scalar alpha = lPN*(unitDelta & edgeNormal);
435  if (mag(alpha) > SMALL)
436  {
437  dc[edgeI] = scalar(1)/max(alpha, 0.05*lPN);
438  }
439  }
440 
441 
442  forAll(DeltaCoeffs.boundaryField(), patchI)
443  {
444  mesh().boundary()[patchI].makeDeltaCoeffs
445  (
446  DeltaCoeffs.boundaryFieldRef()[patchI]
447  );
448  }
449 }
450 
451 
452 void Foam::edgeInterpolation::makeCorrectionVectors() const
453 {
455  << "Constructing non-orthogonal correction vectors"
456  << endl;
457 
458  correctionVectors_ = new edgeVectorField
459  (
460  IOobject
461  (
462  "correctionVectors",
463  mesh().pointsInstance(),
464  mesh().thisDb(),
468  ),
469  mesh(),
470  dimless
471  );
472  edgeVectorField& CorrVecs = *correctionVectors_;
473 
474  // Set local references to mesh data
475  const areaVectorField& faceCentres = mesh().areaCentres();
476 
477  const labelUList& owner = mesh().owner();
478  const labelUList& neighbour = mesh().neighbour();
479 
480  const edgeVectorField& lengths = mesh().Le();
481 
482  const edgeList& edges = mesh().edges();
483  const pointField& points = mesh().points();
484 
485  scalarField deltaCoeffs(owner.size(), SMALL);
486 
487  vectorField& CorrVecsIn = CorrVecs.primitiveFieldRef();
488 
489  forAll(owner, edgeI)
490  {
491  // Edge normal - area normal
492  vector edgeNormal =
493  normalised(lengths[edgeI] ^ edges[edgeI].vec(points));
494 
495  // Unit delta vector
496  vector unitDelta =
497  faceCentres[neighbour[edgeI]]
498  - faceCentres[owner[edgeI]];
499 
500  unitDelta.removeCollinear(edgeNormal);
501  unitDelta.normalise();
502 
503  // Edge normal - area tangent
504  edgeNormal = normalised(lengths[edgeI]);
505 
506  // Do not allow any mag(val) < SMALL
507  const scalar alpha = unitDelta & edgeNormal;
508  if (mag(alpha) > SMALL)
509  {
510  deltaCoeffs[edgeI] = scalar(1)/alpha;
511  }
512 
513  // Edge correction vector
514  CorrVecsIn[edgeI] =
515  edgeNormal
516  - deltaCoeffs[edgeI]*unitDelta;
517  }
518 
519 
520  edgeVectorField::Boundary& CorrVecsbf = CorrVecs.boundaryFieldRef();
521 
522  forAll(CorrVecs.boundaryField(), patchI)
523  {
524  mesh().boundary()[patchI].makeCorrectionVectors(CorrVecsbf[patchI]);
525  }
526 
527  scalar NonOrthogCoeff = 0.0;
528 
529  if (owner.size() > 0)
530  {
531  scalarField sinAlpha(deltaCoeffs*mag(CorrVecs.internalField()));
532  sinAlpha.clamp_range(-1, 1);
533 
534  NonOrthogCoeff = max(Foam::asin(sinAlpha)*180.0/M_PI);
535  }
536 
537  reduce(NonOrthogCoeff, maxOp<scalar>());
538 
540  << "non-orthogonality coefficient = " << NonOrthogCoeff << " deg."
541  << endl;
542 
543  if (NonOrthogCoeff < 0.1)
544  {
545  orthogonal_ = true;
546  deleteDemandDrivenData(correctionVectors_);
547  }
548  else
549  {
550  orthogonal_ = false;
551  }
552 
554  << "Finished constructing non-orthogonal correction vectors"
555  << endl;
556 }
557 
558 
559 void Foam::edgeInterpolation::makeSkewCorrectionVectors() const
560 {
562  << "Constructing skew correction vectors"
563  << endl;
564 
565  skewCorrectionVectors_ = new edgeVectorField
566  (
567  IOobject
568  (
569  "skewCorrectionVectors",
570  mesh().pointsInstance(),
571  mesh().thisDb(),
575  ),
576  mesh(),
578  );
579  edgeVectorField& SkewCorrVecs = *skewCorrectionVectors_;
580 
581  // Set local references to mesh data
582  const areaVectorField& C = mesh().areaCentres();
583  const edgeVectorField& Ce = mesh().edgeCentres();
584 
585  const labelUList& owner = mesh().owner();
586  const labelUList& neighbour = mesh().neighbour();
587 
588  const pointField& points = mesh().points();
589  const edgeList& edges = mesh().edges();
590 
591 
592  forAll(neighbour, edgeI)
593  {
594  const vector& P = C[owner[edgeI]];
595  const vector& N = C[neighbour[edgeI]];
596  const vector& S = points[edges[edgeI].start()];
597 
598  // (T:Eq. 5.4)
599  const vector d(N - P);
600  const vector e(edges[edgeI].vec(points));
601  const vector de(d^e);
602  const scalar alpha = magSqr(de);
603 
604  if (alpha < SMALL)
605  {
606  // Too small - skew correction remains zero
607  continue;
608  }
609  const scalar beta = -((d^(S - P)) & de)/alpha;
610 
611  // (T:Eq. 5.3)
612  const vector E(S + beta*e);
613 
614  SkewCorrVecs[edgeI] = Ce[edgeI] - E;
615  }
616 
617 
618  edgeVectorField::Boundary& bSkewCorrVecs =
619  SkewCorrVecs.boundaryFieldRef();
620 
621  forAll(SkewCorrVecs.boundaryField(), patchI)
622  {
623  faePatchVectorField& patchSkewCorrVecs = bSkewCorrVecs[patchI];
624 
625  if (patchSkewCorrVecs.coupled())
626  {
627  const labelUList& edgeFaces =
628  mesh().boundary()[patchI].edgeFaces();
629 
630  const edgeList::subList patchEdges =
631  mesh().boundary()[patchI].patchSlice(edges);
632 
633  vectorField ngbC(C.boundaryField()[patchI].patchNeighbourField());
634 
635  forAll(patchSkewCorrVecs, edgeI)
636  {
637  const vector& P = C[edgeFaces[edgeI]];
638  const vector& N = ngbC[edgeI];
639  const vector& S = points[patchEdges[edgeI].start()];
640 
641  // (T:Eq. 5.4)
642  const vector d(N - P);
643  const vector e(patchEdges[edgeI].vec(points));
644  const vector de(d^e);
645  const scalar alpha = magSqr(de);
646 
647  if (alpha < SMALL)
648  {
649  // Too small - skew correction remains zero
650  continue;
651  }
652  const scalar beta = -((d^(S - P)) & de)/alpha;
653 
654  const vector E(S + beta*e);
655 
656  patchSkewCorrVecs[edgeI] =
657  Ce.boundaryField()[patchI][edgeI] - E;
658  }
659  }
660  }
661 
662  #ifdef FA_SKEW_CORRECTION
663 
664  constexpr scalar maxSkewRatio = 0.1;
665  scalar skewCoeff = 0;
666 
667  forAll(own, edgeI)
668  {
669  const scalar magSkew = mag(skewCorrVecs[edgeI]);
670 
671  const scalar lPN =
672  mag
673  (
674  Ce[edgeI]
675  - skewCorrVecs[edgeI]
676  - C[owner[edgeI]]
677  )
678  + mag
679  (
680  C[neighbour[edgeI]]
681  - Ce[edgeI]
682  + skewCorrVecs[edgeI]
683  );
684 
685  const scalar ratio = magSkew/lPN;
686 
687  if (skewCoeff < ratio)
688  {
689  skewCoeff = ratio;
690 
691  if (skewCoeff > maxSkewRatio)
692  {
693  break;
694  }
695  }
696  }
697 
699  << "skew coefficient = " << skewCoeff << endl;
700 
701  if (skewCoeff < maxSkewRatio)
702  {
703  deleteDemandDrivenData(skewCorrectionVectors_);
704  }
705 
706  #endif
707 
708  skew_ = bool(skewCorrectionVectors_);
709 
710 
712  << "Finished constructing skew correction vectors"
713  << endl;
714 }
715 
716 
717 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
faceListList boundary
const edgeVectorField & skewCorrectionVectors() const
Return reference to skew vectors array.
const edgeScalarField & weights() const
Return reference to weighting factors array.
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
const edgeVectorField & correctionVectors() const
Return reference to non-orthogonality correction vectors array.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
dimensionedTensor skew(const dimensionedTensor &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Face to edge interpolation scheme. Included in faMesh.
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition: VectorI.H:114
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of boundary fields.
SubList< edge > subList
Declare type of subList.
Definition: List.H:144
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.
bool movePoints() const
Do what is necessary if the mesh has moved.
const edgeScalarField & lPN() const
Return reference to PN geodesic distance.
dimensionedScalar asin(const dimensionedScalar &ds)
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Declarations for faPatchField types.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:580
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
const pointField & points
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
void clearOut()
Clear all geometry and addressing.
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
bool orthogonal() const
Return whether mesh is orthogonal or not.
Vector< scalar > vector
Definition: vector.H:57
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:46
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
~edgeInterpolation()
Destructor.
bool skew() const
Return whether mesh is skew or not.
const Vector< label > N(dict.get< Vector< label >>("N"))
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:572
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Template functions to aid in the implementation of demand driven data.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
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.
faePatchField< vector > faePatchVectorField
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Field< vector > vectorField
Specialisation of Field<T> for vector.
edgeInterpolation(const faMesh &)
Construct given an faMesh.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Vector< Cmpt > & removeCollinear(const Vector< Cmpt > &unitVec)
Inplace removal of components that are collinear to the given unit vector.
Definition: VectorI.H:136
void deleteDemandDrivenData(DataPtr &dataPtr)
Calculate the matrix for the second temporal derivative.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:81
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Do not request registration (bool: false)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127