reconstructionSchemes.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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::reconstructionSchemes
28
29Description
30 Original code supplied by Henning Scheufler, DLR (2019)
31
32SourceFiles
33 reconstructionSchemes.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef reconstructionSchemes_H
38#define reconstructionSchemes_H
39
40#include "typeInfo.H"
42#include "dimensionedScalar.H"
43#include "volFields.H"
44#include "surfaceFields.H"
45#include "IOdictionary.H"
46#include "DynamicField.H"
47#include "MeshedSurface.H"
48#include "MeshedSurfacesFwd.H"
49
50// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52namespace Foam
53{
54
55/*---------------------------------------------------------------------------*\
56 Class reconstructionSchemes Declaration
57\*---------------------------------------------------------------------------*/
60:
61 public IOdictionary
62{
63 // Private Data
64
65 //- Coefficients dictionary
66 dictionary reconstructionSchemesCoeffs_;
67
68
69protected:
70
71 // Protected Data
72
73 //- Reference to the VoF Field
75
76 //- Reference to the face fluxes
78
79 //- Reference to the velocity field
80 const volVectorField& U_;
81
82 //- Interface area normals
84
85 //- Interface centres
87
88 //- Is interface cell?
90
91 //- List of cell labels that have an interface
93
94 //- Store timeindex/iteration to avoid multiple reconstruction
96
97
98 // Protected Member Functions
99
100 //- Is the interface already up-to-date?
101 bool alreadyReconstructed(bool forceUpdate = true) const;
102
103 //- No copy construct
105
106 //- No copy assignment
107 void operator=(const reconstructionSchemes&) = delete;
108
109
110public:
111
112 // Public Classes
114 class interface
115 :
116 public meshedSurface
117 {
118 //- For every face the original cell in mesh
119 labelList meshCells_;
120
121 public:
122
124 (
125 const pointField& pointLst,
126 const faceList& faceLst,
127 const labelList& meshCells
128 )
129 :
130 meshedSurface(pointLst, faceLst, surfZoneList()),
131 meshCells_(meshCells)
132 {}
133
135 (
136 pointField&& pointLst,
137 faceList&& faceLst,
139 )
140 :
142 (
143 std::move(pointLst),
144 std::move(faceLst),
146 ),
147 meshCells_(std::move(meshCells))
148 {}
149
150 //- For every face, the original cell in mesh
151 const labelList& meshCells() const
152 {
153 return meshCells_;
154 }
155 };
156
157
158public:
159
160 //- Runtime type information
161 TypeName("reconstructionSchemes");
162
163
164 // Declare run-time constructor selection table
167 (
168 autoPtr,
170 components,
171 (
173 const surfaceScalarField& phi,
174 const volVectorField& U,
175 const dictionary& dict
176 ),
177 (alpha1, phi, U, dict)
178 );
179
180
181 // Selectors
182
183 //- Return a reference to the selected phaseChange model
185 (
187 const surfaceScalarField& phi,
188 const volVectorField& U,
189 const dictionary& dict
190 );
191
192
193 // Constructors
194
195 //- Construct from components
197 (
198 const word& type,
200 const surfaceScalarField& phi,
201 const volVectorField& U,
202 const dictionary& dict
203 );
204
205
206 //- Destructor
207 virtual ~reconstructionSchemes() = default;
208
209
210 // Member Functions
211
212 //- Access to the model dictionary
214 {
215 return reconstructionSchemesCoeffs_;
216 }
217
218 //- Const access to the model dictionary
219 const dictionary& modelDict() const noexcept
220 {
221 return reconstructionSchemesCoeffs_;
222 }
223
224 //- Reconstruct the interface
225 virtual void reconstruct(bool forceUpdate = true) = 0;
226
227 //- const-Reference to interface normals
228 const volVectorField& normal() const noexcept
229 {
230 return normal_;
231 }
232
233 //- const-Reference to interface centres
234 const volVectorField& centre() const noexcept
235 {
236 return centre_;
237 }
238
239 //- Reference to interface normals
241 {
242 return normal_;
243 }
244
245 //- Reference to interface centres
247 {
248 return centre_;
249 }
250
251 //- Does the cell contain interface
252 const boolList& interfaceCell() const noexcept
253 {
254 return interfaceCell_;
255 }
256
257 //- List of cells with an interface
259 {
260 return interfaceLabels_;
261 }
262
263 //- Map VoF Field in case of refinement
264 virtual void mapAlphaField() const
265 {}
266
267 //- Generated interface surface points/faces.
268 // \note the points are disconnected!
269 interface surface();
270};
271
272
273// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274
275} // End namespace Foam
276
277// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278
279#endif
280
281// ************************************************************************* //
surfaceScalarField & phi
const volScalarField & alpha1
Dynamically sized Field.
Definition: DynamicField.H:65
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const labelList & meshCells() const
For every face, the original cell in mesh.
Original code supplied by Henning Scheufler, DLR (2019)
virtual void mapAlphaField() const
Map VoF Field in case of refinement.
boolList interfaceCell_
Is interface cell?
void operator=(const reconstructionSchemes &)=delete
No copy assignment.
const volVectorField & normal() const noexcept
const-Reference to interface normals
volVectorField normal_
Interface area normals.
static autoPtr< reconstructionSchemes > New(volScalarField &alpha1, const surfaceScalarField &phi, const volVectorField &U, const dictionary &dict)
Return a reference to the selected phaseChange model.
Pair< label > timeIndexAndIter_
Store timeindex/iteration to avoid multiple reconstruction.
const surfaceScalarField & phi_
Reference to the face fluxes.
const dictionary & modelDict() const noexcept
Const access to the model dictionary.
TypeName("reconstructionSchemes")
Runtime type information.
const volVectorField & U_
Reference to the velocity field.
declareRunTimeSelectionTable(autoPtr, reconstructionSchemes, components,(volScalarField &alpha1, const surfaceScalarField &phi, const volVectorField &U, const dictionary &dict),(alpha1, phi, U, dict))
bool alreadyReconstructed(bool forceUpdate=true) const
Is the interface already up-to-date?
const volVectorField & centre() const noexcept
const-Reference to interface centres
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
volVectorField & centre() noexcept
Reference to interface centres.
interface surface()
Generated interface surface points/faces.
const DynamicField< label > & interfaceLabels() const noexcept
List of cells with an interface.
const boolList & interfaceCell() const noexcept
Does the cell contain interface.
reconstructionSchemes(const reconstructionSchemes &)=delete
No copy construct.
volScalarField & alpha1_
Reference to the VoF Field.
volVectorField centre_
Interface centres.
volVectorField & normal() noexcept
Reference to interface normals.
virtual ~reconstructionSchemes()=default
Destructor.
virtual void reconstruct(bool forceUpdate=true)=0
Reconstruct the interface.
dictionary & modelDict() noexcept
Access to the model dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
const direction noexcept
Definition: Scalar.H:223
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dictionary dict
Foam::surfaceFields.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73