plicRDF.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::reconstruction::plicRDF
28 
29 Description
30  Reconstructs an interface (centre and normal vector) consisting of planes
31  to match the internal fluid distribution in cells. The interface normals
32  are estimated by least square gradient scheme on the RDF function (height).
33  Uses the normal from the previous times step as intial guess.
34 
35  Reference:
36  \verbatim
37  Henning Scheufler, Johan Roenby,
38  Accurate and efficient surface reconstruction from volume
39  fraction data on general meshes,
40  Journal of Computational Physics, 2019,
41  doi 10.1016/j.jcp.2019.01.009
42 
43  \endverbatim
44 
45  Original code supplied by Henning Scheufler, DLR (2019)
46 
47 SourceFiles
48  plicRDF.C
49 
50 \*---------------------------------------------------------------------------*/
51 
52 #ifndef plicRDF_H
53 #define plicRDF_H
54 
55 #include "typeInfo.H"
56 #include "reconstructionSchemes.H"
57 #include "volFields.H"
58 #include "dimensionedScalar.H"
59 #include "autoPtr.H"
60 #include "surfaceIteratorPLIC.H"
62 #include "zoneDistribute.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 namespace reconstruction
69 {
70 
71 /*---------------------------------------------------------------------------*\
72  Class plicRDF Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class plicRDF
76 :
78 {
79  // Private Data
80 
81  //- residuals storage
82  struct normalRes
83  {
84  label celli = {};
85  scalar normalResidual = {};
86  scalar avgAngle = {};
87  };
88 
89  //- Reference to mesh
90  const fvMesh& mesh_;
91 
92  //- Interpolation object from cell centres to points
93  DynamicField<vector> interfaceNormal_;
94 
95 
96  // Switches and tolerances. Tolerances need to go into toleranceSwitches
97 
98  //- Tolerance for search of isoFace giving specified VOF value
99  scalar isoFaceTol_;
100 
101  //- Tolerance for marking of surface cells:
102  // Those with surfCellTol_ < alpha1 < 1 - surfCellTol_
103  scalar surfCellTol_;
104 
105  //- Tolerance
106  scalar tol_;
107 
108  //- Relative tolerance
109  scalar relTol_;
110 
111  //- Number of times that the interface is reconstructed
112  //- has to be bigger than 2
113  label iteration_;
114 
115  //- Interpolated normal from previous time step
116  bool interpolateNormal_;
117 
118  //- Calculates the RDF function
120 
121  //- surfaceIterator finds the plane centre for specified VOF value
122  surfaceIteratorPLIC sIterPLIC_;
123 
124 
125  // Private Member Functions
126 
127  //- Set initial normals by interpolation from the previous
128  //- timestep or with the Young method
129  void setInitNormals(bool interpolate);
130 
131  //- compute gradient at the surfaces
132  void gradSurf(const volScalarField& phi);
133 
134  //- compute the normal residuals
135  void calcResidual
136  (
137  List<normalRes>& normalResidual
138  );
139 
140  //- interpolation of the normals from the previous time step
141  void interpolateNormal();
142 
143  //- No copy construct
144  plicRDF(const plicRDF&) = delete;
145 
146  //- No copy assignment
147  void operator=(const plicRDF&) = delete;
148 
149 
150 public:
151 
152  //- Runtime type information
153  TypeName("plicRDF");
154 
155 
156  //- Construct from components
157  plicRDF
158  (
160  const surfaceScalarField& phi,
161  const volVectorField& U,
162  const dictionary& dict
163  );
164 
165 
166  //- Destructor
167  virtual ~plicRDF() = default;
168 
169 
170  // Member Functions
171 
172  //- Reconstruct interface
173  virtual void reconstruct(bool forceUpdate = true);
174 
175  //- Map VoF Field in case of refinement
176  virtual void mapAlphaField() const;
177 };
178 
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace reconstruction
183 } // End namespace Foam
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 #endif
188 
189 // ************************************************************************* //
volFields.H
typeInfo.H
reconstructionSchemes.H
zoneDistribute.H
surfaceIteratorPLIC.H
Foam::reconstruction::plicRDF::reconstruct
virtual void reconstruct(bool forceUpdate=true)
Reconstruct interface.
Definition: plicRDF.C:383
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::DynamicField
Dynamically sized Field.
Definition: DynamicField.H:49
reconstructedDistanceFunction.H
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::reconstruction::plicRDF::mapAlphaField
virtual void mapAlphaField() const
Map VoF Field in case of refinement.
Definition: plicRDF.C:549
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
dimensionedScalar.H
U
U
Definition: pEqn.H:72
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::reconstructedDistanceFunction
Calculates a reconstructed distance function.
Definition: reconstructedDistanceFunction.H:58
Foam::reconstruction::plicRDF::TypeName
TypeName("plicRDF")
Runtime type information.
Foam::reconstructionSchemes
Original code supplied by Henning Scheufler, DLR (2019)
Definition: reconstructionSchemes.H:58
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::reconstruction::plicRDF::~plicRDF
virtual ~plicRDF()=default
Destructor.
Foam::reconstruction::plicRDF
Reconstructs an interface (centre and normal vector) consisting of planes to match the internal fluid...
Definition: plicRDF.H:74
Foam::surfaceIteratorPLIC
Finds the cutValue that matches the volume fraction.
Definition: surfaceIteratorPLIC.H:66
autoPtr.H