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 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 namespace reconstruction
68 {
69 
70 /*---------------------------------------------------------------------------*\
71  Class plicRDF Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 class plicRDF
75 :
77 {
78  // Private Data
79 
80  //- residuals storage
81  struct normalRes
82  {
83  label celli = {};
84  scalar normalResidual = {};
85  scalar avgAngle = {};
86  };
87 
88  //- Reference to mesh
89  const fvMesh& mesh_;
90 
91  //- Interpolation object from cell centres to points
92  DynamicField<vector> interfaceNormal_;
93 
94 
95  // Switches and tolerances. Tolerances need to go into toleranceSwitches
96 
97  //- Tolerance for search of isoFace giving specified VOF value
98  scalar isoFaceTol_;
99 
100  //- Tolerance for marking of surface cells:
101  // Those with surfCellTol_ < alpha1 < 1 - surfCellTol_
102  scalar surfCellTol_;
103 
104  //- Tolerance
105  scalar tol_;
106 
107  //- Relative tolerance
108  scalar relTol_;
109 
110  //- Number of times that the interface is reconstructed
111  //- has to be bigger than 2
112  label iteration_;
113 
114  //- Interpolated normal from previous time step
115  bool interpolateNormal_;
116 
117  //- Calculates the RDF function
119 
120  //- surfaceIterator finds the plane centre for specified VOF value
121  surfaceIteratorPLIC sIterPLIC_;
122 
123 
124  // Private Member Functions
125 
126  //- Set initial normals by interpolation from the previous
127  //- timestep or with the Young method
128  void setInitNormals(bool interpolate);
129 
130  //- compute gradient at the surfaces
131  void gradSurf(const volScalarField& phi);
132 
133  //- compute the normal residuals
134  void calcResidual
135  (
136  List<normalRes>& normalResidual
137  );
138 
139  //- interpolation of the normals from the previous time step
140  void interpolateNormal();
141 
142  //- No copy construct
143  plicRDF(const plicRDF&) = delete;
144 
145  //- No copy assignment
146  void operator=(const plicRDF&) = delete;
147 
148 
149 public:
150 
151  //- Runtime type information
152  TypeName("plicRDF");
153 
154 
155  //- Construct from components
156  plicRDF
157  (
159  const surfaceScalarField& phi,
160  const volVectorField& U,
161  const dictionary& dict
162  );
163 
164 
165  //- Destructor
166  virtual ~plicRDF() = default;
167 
168 
169  // Member Functions
170 
171  //- Reconstruct interface
172  virtual void reconstruct(bool forceUpdate = true);
173 
174  //- Map VoF Field in case of refinement
175  virtual void mapAlphaField() const;
176 };
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace reconstruction
182 } // End namespace Foam
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #endif
187 
188 // ************************************************************************* //
Original code supplied by Henning Scheufler, DLR (2019)
dictionary dict
Reconstructs an interface (centre and normal vector) consisting of planes to match the internal fluid...
Definition: plicRDF.H:69
Finds the cutValue that matches the volume fraction.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
virtual void reconstruct(bool forceUpdate=true)
Reconstruct interface.
Definition: plicRDF.C:377
virtual ~plicRDF()=default
Destructor.
Dynamically sized Field.
Definition: DynamicField.H:45
Calculates a reconstructed distance function.
U
Definition: pEqn.H:72
virtual void mapAlphaField() const
Map VoF Field in case of refinement.
Definition: plicRDF.C:543
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Namespace for OpenFOAM.
const volScalarField & alpha1
TypeName("plicRDF")
Runtime type information.