reconstructedDistanceFunction.H
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) 2019-2020 DLR
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::reconstructedDistanceFunction
28 
29 Description
30  Calculates a reconstructed distance function
31 
32  Original code supplied by Henning Scheufler, DLR (2019)
33 
34 SourceFiles
35  reconstructedDistanceFunction.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef reconstructedDistanceFunction_H
40 #define reconstructedDistanceFunction_H
41 
42 #include "fvMesh.H"
43 #include "volFields.H"
44 #include "surfaceFields.H"
45 #include "globalIndex.H"
46 #include "Map.H"
47 #include "zoneDistribute.H"
48 #include "dimensionedScalar.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class reconstructedDistanceFunction Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 :
61  public volScalarField
62 {
63  // Private Data
64 
65  //- Reference to mesh
66  const fvMesh& mesh_;
67 
68  //- Stores the coupled boundary points which have to be synced
69  labelList coupledBoundaryPoints_;
70 
71  //- Distance of the interface band to the interface
72  volScalarField cellDistLevel_;
73 
74  //- Is the cell in the interface band?
75  boolList nextToInterface_;
76 
77  //- Return patch of all coupled faces.
78  autoPtr<indirectPrimitivePatch> coupledFacesPatch() const;
79 
80 
81 public:
82 
83  //- Construct from fvMesh
84  explicit reconstructedDistanceFunction(const fvMesh& mesh);
85 
86 
87  // Member Functions
88 
90  (
91  const boolList& interfaceCells,
92  const label neiRingLevel
93  );
94 
96  (
98  const volVectorField& centre,
99  const volVectorField& normal,
100  zoneDistribute& distribute,
101  bool updateStencil=true
102  );
103 
104  void updateContactAngle
105  (
106  const volScalarField& alpha,
107  const volVectorField& U,
109  );
110 
111  const volScalarField& cellDistLevel() const noexcept
112  {
113  return cellDistLevel_;
114  }
115 
116  const boolList& nextToInterface() const noexcept
117  {
118  return nextToInterface_;
119  }
120 };
121 
122 
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124 
125 } // End namespace Foam
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 #endif
130 
131 // ************************************************************************* //
Foam::surfaceFields.
const volScalarField & cellDistLevel() const noexcept
reconstructedDistanceFunction(const fvMesh &mesh)
Construct from fvMesh.
const boolList & nextToInterface() const noexcept
void updateContactAngle(const volScalarField &alpha, const volVectorField &U, surfaceVectorField::Boundary &nHatb)
void markCellsNearSurf(const boolList &interfaceCells, const label neiRingLevel)
const volScalarField & constructRDF(const boolList &nextToInterface, const volVectorField &centre, const volVectorField &normal, zoneDistribute &distribute, bool updateStencil=true)
Calculates a reconstructed distance function.
const Mesh & mesh() const noexcept
Return mesh.
const direction noexcept
Definition: Scalar.H:258
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Class for parallel communication in a narrow band. It either provides a Map with the neighbouring val...
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.