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 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 {
39  defineTypeNameAndDebug(edgeInterpolation, 0);
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  orthogonal_(false),
63  correctionVectors_(nullptr),
64  skew_(true),
65  skewCorrectionVectors_(nullptr)
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"
129  << abort(FatalError);
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"
154  << abort(FatalError);
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 void Foam::edgeInterpolation::makeLPN() const
178 {
180  << "Constructing geodesic distance between points P and N"
181  << endl;
182 
183 
184  lPN_ = new edgeScalarField
185  (
186  IOobject
187  (
188  "lPN",
189  faMesh_.time().constant(),
190  faMesh_(),
193  false
194  ),
195  mesh(),
196  dimLength
197  );
198  edgeScalarField& lPN = *lPN_;
199 
200  // Set local references to mesh data
201  const edgeVectorField& edgeCentres = mesh().edgeCentres();
202  const areaVectorField& faceCentres = mesh().areaCentres();
203  const labelUList& owner = mesh().owner();
204  const labelUList& neighbour = mesh().neighbour();
205 
206  scalarField& lPNIn = lPN.primitiveFieldRef();
207 
208  forAll(owner, edgeI)
209  {
210  vector curSkewCorrVec(Zero);
211 
212  if (skew())
213  {
214  curSkewCorrVec = skewCorrectionVectors()[edgeI];
215  }
216 
217  scalar lPE =
218  mag
219  (
220  edgeCentres[edgeI]
221  - curSkewCorrVec
222  - faceCentres[owner[edgeI]]
223  );
224 
225  scalar lEN =
226  mag
227  (
228  faceCentres[neighbour[edgeI]]
229  - edgeCentres[edgeI]
230  + curSkewCorrVec
231  );
232 
233  lPNIn[edgeI] = (lPE + lEN);
234  }
235 
236 
237  forAll(lPN.boundaryField(), patchI)
238  {
239  mesh().boundary()[patchI].makeDeltaCoeffs
240  (
241  lPN.boundaryFieldRef()[patchI]
242  );
243 
244  lPN.boundaryFieldRef()[patchI] = 1.0/lPN.boundaryField()[patchI];
245  }
246 
247 
249  << "Finished constructing geodesic distance PN"
250  << endl;
251 }
252 
253 
254 void Foam::edgeInterpolation::makeWeights() const
255 {
257  << "Constructing weighting factors for edge interpolation"
258  << endl;
259 
260 
261  weightingFactors_ = new edgeScalarField
262  (
263  IOobject
264  (
265  "weightingFactors",
266  mesh()().pointsInstance(),
267  mesh()(),
270  false
271  ),
272  mesh(),
273  dimless
274  );
275  edgeScalarField& weightingFactors = *weightingFactors_;
276 
277 
278  // Set local references to mesh data
279  const edgeVectorField& edgeCentres = mesh().edgeCentres();
280  const areaVectorField& faceCentres = mesh().areaCentres();
281  const labelUList& owner = mesh().owner();
282  const labelUList& neighbour = mesh().neighbour();
283 
284  scalarField& weightingFactorsIn = weightingFactors.primitiveFieldRef();
285 
286  forAll(owner, edgeI)
287  {
288  vector curSkewCorrVec(Zero);
289 
290  if (skew())
291  {
292  curSkewCorrVec = skewCorrectionVectors()[edgeI];
293  }
294 
295  scalar lPE =
296  mag
297  (
298  edgeCentres[edgeI]
299  - curSkewCorrVec
300  - faceCentres[owner[edgeI]]
301  );
302 
303  scalar lEN =
304  mag
305  (
306  faceCentres[neighbour[edgeI]]
307  - edgeCentres[edgeI]
308  + curSkewCorrVec
309  );
310 
311  weightingFactorsIn[edgeI] =
312  lEN
313  /(
314  lPE
315 # ifdef BAD_MESH_STABILISATION
316  + VSMALL
317 # endif
318  + lEN
319  );
320  }
321 
322  forAll(mesh().boundary(), patchI)
323  {
324  mesh().boundary()[patchI].makeWeights
325  (
326  weightingFactors.boundaryFieldRef()[patchI]
327  );
328  }
329 
331  << "Finished constructing weighting factors for face interpolation"
332  << endl;
333 }
334 
335 
336 void Foam::edgeInterpolation::makeDeltaCoeffs() const
337 {
339  << "Constructing differencing factors array for edge gradient"
340  << endl;
341 
342  // Force the construction of the weighting factors
343  // needed to make sure deltaCoeffs are calculated for parallel runs.
344  weights();
345 
346  differenceFactors_ = new edgeScalarField
347  (
348  IOobject
349  (
350  "differenceFactors_",
351  mesh()().pointsInstance(),
352  mesh()(),
355  false
356  ),
357  mesh(),
359  );
360  edgeScalarField& DeltaCoeffs = *differenceFactors_;
361  scalarField& dc = DeltaCoeffs.primitiveFieldRef();
362 
363  // Set local references to mesh data
364  const edgeVectorField& edgeCentres = mesh().edgeCentres();
365  const areaVectorField& faceCentres = mesh().areaCentres();
366  const labelUList& owner = mesh().owner();
367  const labelUList& neighbour = mesh().neighbour();
368  const edgeVectorField& lengths = mesh().Le();
369 
370  const edgeList& edges = mesh().edges();
371  const pointField& points = mesh().points();
372 
373 
374  forAll(owner, edgeI)
375  {
376  // Edge normal - area normal
377  vector edgeNormal =
378  normalised(lengths[edgeI] ^ edges[edgeI].vec(points));
379 
380  // Unit delta vector
381  vector unitDelta =
382  faceCentres[neighbour[edgeI]]
383  - faceCentres[owner[edgeI]];
384 
385  unitDelta -= edgeNormal*(edgeNormal & unitDelta);
386  unitDelta.normalise();
387 
388 
389  // Calc PN arc length
390  vector curSkewCorrVec(Zero);
391 
392  if (skew())
393  {
394  curSkewCorrVec = skewCorrectionVectors()[edgeI];
395  }
396 
397  scalar lPE =
398  mag
399  (
400  edgeCentres[edgeI]
401  - curSkewCorrVec
402  - faceCentres[owner[edgeI]]
403  );
404 
405  scalar lEN =
406  mag
407  (
408  faceCentres[neighbour[edgeI]]
409  - edgeCentres[edgeI]
410  + curSkewCorrVec
411  );
412 
413  scalar lPN = lPE + lEN;
414 
415 
416  // Edge normal - area tangent
417  edgeNormal = normalised(lengths[edgeI]);
418 
419  // Delta coefficient as per definition
420 // dc[edgeI] = 1.0/(lPN*(unitDelta & edgeNormal));
421 
422  // Stabilised form for bad meshes. HJ, 23/Jul/2009
423  dc[edgeI] = 1.0/max((lPN*(unitDelta & edgeNormal)), 0.05*lPN);
424  }
425 
426 
427  forAll(DeltaCoeffs.boundaryField(), patchI)
428  {
429  mesh().boundary()[patchI].makeDeltaCoeffs
430  (
431  DeltaCoeffs.boundaryFieldRef()[patchI]
432  );
433  }
434 }
435 
436 
437 void Foam::edgeInterpolation::makeCorrectionVectors() const
438 {
440  << "Constructing non-orthogonal correction vectors"
441  << endl;
442 
443  correctionVectors_ = new edgeVectorField
444  (
445  IOobject
446  (
447  "correctionVectors",
448  mesh()().pointsInstance(),
449  mesh()(),
452  false
453  ),
454  mesh(),
455  dimless
456  );
457  edgeVectorField& CorrVecs = *correctionVectors_;
458 
459  // Set local references to mesh data
460  const areaVectorField& faceCentres = mesh().areaCentres();
461 
462  const labelUList& owner = mesh().owner();
463  const labelUList& neighbour = mesh().neighbour();
464 
465  const edgeVectorField& lengths = mesh().Le();
466 
467  const edgeList& edges = mesh().edges();
468  const pointField& points = mesh().points();
469 
470  scalarField deltaCoeffs(owner.size());
471 
472  vectorField& CorrVecsIn = CorrVecs.primitiveFieldRef();
473 
474  forAll(owner, edgeI)
475  {
476  // Edge normal - area normal
477  vector edgeNormal =
478  normalised(lengths[edgeI] ^ edges[edgeI].vec(points));
479 
480  // Unit delta vector
481  vector unitDelta =
482  faceCentres[neighbour[edgeI]]
483  - faceCentres[owner[edgeI]];
484 
485  unitDelta -= edgeNormal*(edgeNormal & unitDelta);
486  unitDelta.normalise();
487 
488  // Edge normal - area tangent
489  edgeNormal = normalised(lengths[edgeI]);
490 
491  // Delta coeffs
492  deltaCoeffs[edgeI] = 1.0/(unitDelta & edgeNormal);
493 
494  // Edge correction vector
495  CorrVecsIn[edgeI] =
496  edgeNormal
497  - deltaCoeffs[edgeI]*unitDelta;
498  }
499 
500 
501  edgeVectorField::Boundary& CorrVecsbf = CorrVecs.boundaryFieldRef();
502 
503  forAll(CorrVecs.boundaryField(), patchI)
504  {
505  mesh().boundary()[patchI].makeCorrectionVectors(CorrVecsbf[patchI]);
506  }
507 
508  scalar NonOrthogCoeff = 0.0;
509 
510  if (owner.size() > 0)
511  {
512  scalarField sinAlpha(deltaCoeffs*mag(CorrVecs.internalField()));
513 
514  forAll(sinAlpha, edgeI)
515  {
516  sinAlpha[edgeI] = max(-1, min(sinAlpha[edgeI], 1));
517  }
518 
519  NonOrthogCoeff = max(Foam::asin(sinAlpha)*180.0/M_PI);
520  }
521 
522  reduce(NonOrthogCoeff, maxOp<scalar>());
523 
525  << "non-orthogonality coefficient = " << NonOrthogCoeff << " deg."
526  << endl;
527 
528  if (NonOrthogCoeff < 0.1)
529  {
530  orthogonal_ = true;
531  deleteDemandDrivenData(correctionVectors_);
532  }
533  else
534  {
535  orthogonal_ = false;
536  }
537 
539  << "Finished constructing non-orthogonal correction vectors"
540  << endl;
541 }
542 
543 
544 void Foam::edgeInterpolation::makeSkewCorrectionVectors() const
545 {
547  << "Constructing skew correction vectors"
548  << endl;
549 
550  skewCorrectionVectors_ = new edgeVectorField
551  (
552  IOobject
553  (
554  "skewCorrectionVectors",
555  mesh()().pointsInstance(),
556  mesh()(),
559  false
560  ),
561  mesh(),
562  dimless
563  );
564  edgeVectorField& SkewCorrVecs = *skewCorrectionVectors_;
565 
566  // Set local references to mesh data
567  const areaVectorField& C = mesh().areaCentres();
568  const edgeVectorField& Ce = mesh().edgeCentres();
569 
570  const labelUList& owner = mesh().owner();
571  const labelUList& neighbour = mesh().neighbour();
572 
573  const pointField& points = mesh().points();
574  const edgeList& edges = mesh().edges();
575 
576 
577  forAll(neighbour, edgeI)
578  {
579  const vector& P = C[owner[edgeI]];
580  const vector& N = C[neighbour[edgeI]];
581  const vector& S = points[edges[edgeI].start()];
582  vector e = edges[edgeI].vec(points);
583 
584  scalar alpha =
585  -(((N - P)^(S - P))&((N - P)^e ))/(((N - P)^e)&((N - P)^e));
586 
587  vector E(S + alpha*e);
588 
589  SkewCorrVecs[edgeI] = Ce[edgeI] - E;
590  }
591 
592 
593  edgeVectorField::Boundary& bSkewCorrVecs =
594  SkewCorrVecs.boundaryFieldRef();
595 
596  forAll(SkewCorrVecs.boundaryField(), patchI)
597  {
598  faePatchVectorField& patchSkewCorrVecs = bSkewCorrVecs[patchI];
599 
600  if (patchSkewCorrVecs.coupled())
601  {
602  const labelUList& edgeFaces =
603  mesh().boundary()[patchI].edgeFaces();
604 
605  const edgeList::subList patchEdges =
606  mesh().boundary()[patchI].patchSlice(edges);
607 
608  vectorField ngbC(C.boundaryField()[patchI].patchNeighbourField());
609 
610  forAll(patchSkewCorrVecs, edgeI)
611  {
612  const vector& P = C[edgeFaces[edgeI]];
613  const vector& N = ngbC[edgeI];
614  const vector& S = points[patchEdges[edgeI].start()];
615  vector e = patchEdges[edgeI].vec(points);
616 
617  scalar alpha =
618  - (((N - P)^(S - P))&((N - P)^e))
619  /(((N - P)^e)&((N - P)^e));
620 
621  vector E(S + alpha*e);
622 
623  patchSkewCorrVecs[edgeI] =
624  Ce.boundaryField()[patchI][edgeI] - E;
625  }
626  }
627  else
628  {
629  patchSkewCorrVecs = vector::zero;
630  }
631  }
632 
633 
634  scalar skewCoeff = 0.0;
635 
636  // Calculating PN arc length
637  scalarField lPN(owner.size());
638 
639  forAll(owner, edgeI)
640  {
641  lPN[edgeI] =
642  mag
643  (
644  Ce[edgeI]
645  - SkewCorrVecs[edgeI]
646  - C[owner[edgeI]]
647  )
648  + mag
649  (
650  C[neighbour[edgeI]]
651  - Ce[edgeI]
652  + SkewCorrVecs[edgeI]
653  );
654  }
655 
656  if (lPN.size() > 0)
657  {
658  skewCoeff = max(mag(SkewCorrVecs.internalField())/mag(lPN));
659  }
660 
662  << "skew coefficient = " << skewCoeff << endl;
663 
664  if (skewCoeff < 0.1)
665  {
666  skew_ = false;
667  deleteDemandDrivenData(skewCorrectionVectors_);
668  }
669  else
670  {
671  skew_ = true;
672  }
673 
675  << "Finished constructing skew correction vectors"
676  << endl;
677 }
678 
679 
680 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:138
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::edgeVectorField
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
Definition: edgeFieldsFwd.H:57
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::edgeInterpolation::skew
bool skew() const
Return whether mesh is skew or not.
Definition: edgeInterpolation.C:136
Foam::edgeInterpolation::skewCorrectionVectors
const edgeVectorField & skewCorrectionVectors() const
Return reference to skew vectors array.
Definition: edgeInterpolation.C:148
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
faMesh.H
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
Foam::faePatchVectorField
faePatchField< vector > faePatchVectorField
Definition: faePatchFieldsFwd.H:47
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
C
volScalarField & C
Definition: readThermalProperties.H:102
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::areaVectorField
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:58
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::Field< scalar >
Foam::edgeInterpolation::deltaCoeffs
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
Definition: edgeInterpolation.C:101
Foam::edgeInterpolation::~edgeInterpolation
~edgeInterpolation()
Destructor.
Definition: edgeInterpolation.C:71
Foam::edgeInterpolation::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: edgeInterpolation.C:44
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
edgeFields.H
Foam::edgeInterpolation::lPN
const edgeScalarField & lPN() const
Return reference to PN geodesic distance.
Definition: edgeInterpolation.C:79
Foam::edgeScalarField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:52
Foam::edgeInterpolation::orthogonal
bool orthogonal() const
Return whether mesh is orthogonal or not.
Definition: edgeInterpolation.C:112
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
areaFields.H
Foam::List::subList
SubList< T > subList
Declare type of subList.
Definition: List.H:112
Foam::FatalError
error FatalError
Foam::fvMesh::neighbour
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:413
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::skewCorrectionVectors
Skew-correction vectors for the skewness-corrected interpolation scheme.
Definition: skewCorrectionVectors.H:54
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::edgeInterpolation::movePoints
bool movePoints() const
Do what is necessary if the mesh has moved.
Definition: edgeInterpolation.C:161
Foam::edgeInterpolation::weights
const edgeScalarField & weights() const
Return reference to weighting factors array.
Definition: edgeInterpolation.C:90
Foam::fvMesh::owner
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:407
Foam::edgeInterpolation::correctionVectors
const edgeVectorField & correctionVectors() const
Return reference to non-orthogonality correction vectors array.
Definition: edgeInterpolation.C:123
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
faPatchFields.H
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:82
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
N
const Vector< label > N(dict.get< Vector< label >>("N"))
fam
Calculate the matrix for the second temporal derivative.
Foam::edgeInterpolation::edgeInterpolation
edgeInterpolation(const faMesh &)
Construct given an faMesh.
Definition: edgeInterpolation.C:56
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::GeometricField< scalar, faePatchField, edgeMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189