PLICU.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::PLICU
26 
27 Description
28  Velocity-weighted Piecewise-Linear Interface Calculation (PLICU) corrected
29  scheme is a surface interpolation scheme for flux calculation in advection
30  of a bounded variable, e.g. phase fraction and for interface capturing in
31  the volume of fluid (VoF) method.
32 
33  The interface is represented by single cuts which split each cell to match
34  the volume fraction of the phase in the cell. The cut planes are oriented
35  according to the point field of the local phase fraction. The phase
36  fraction at each cell face - the interpolated value - is then calculated
37  from the face area on either side of the cut. For cases where the
38  single-cut does not accurately represent the cell volume fraction the
39  specified default scheme is used, e.g. interfaceCompression.
40 
41  Additionally the face point velocity values are used to calculate the face
42  flux which is likely to be more accurate in the presence of high shear.
43 
44  Example:
45  \verbatim
46  divSchemes
47  {
48  .
49  .
50  div(phi,alpha) Gauss PLICU interfaceCompression vanLeer 1;
51  .
52  .
53  }
54  \endverbatim
55 
56 See also
57  Foam::PLIC
58  Foam::MPLIC
59  Foam::MPLICU
60  Foam::interfaceCompression
61 
62 SourceFiles
63  PLICU.C
64 
65 \*---------------------------------------------------------------------------*/
66 
67 #ifndef PLICU_H
68 #define PLICU_H
69 
70 #include "PLIC.H"
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 namespace Foam
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class PLIC Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class PLICU
82 :
83  public PLIC
84 {
85 
86 public:
87 
88  //- Runtime type information
89  TypeName("PLICU");
90 
91 
92  // Constructors
93 
94  //- Construct from faceFlux and Istream
96  (
97  const fvMesh& mesh,
98  const surfaceScalarField& faceFlux,
99  Istream& is
100  )
101  :
102  PLIC(mesh, faceFlux, is)
103  {}
104 
105 
106  // Member Functions
107 
108  //- Return the face-interpolate of the given cell field
110  (
111  const volScalarField& vf
112  ) const;
113 
114 
115  // Member Operators
116 
117  //- Disallow default bitwise assignment
118  void operator=(const PLICU&) = delete;
119 };
120 
121 
122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 
124 } // End namespace Foam
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 #endif
129 
130 // ************************************************************************* //
Velocity-weighted Piecewise-Linear Interface Calculation (PLICU) corrected scheme is a surface interp...
Definition: PLICU.H:80
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
PLICU(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from faceFlux and Istream.
Definition: PLICU.H:95
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< surfaceScalarField > interpolate(const volScalarField &vf) const
Return the face-interpolate of the given cell field.
Definition: PLICU.C:43
TypeName("PLICU")
Runtime type information.
Piecewise-Linear Interface Calculation (PLIC) corrected scheme is a surface interpolation scheme for ...
Definition: PLIC.H:77
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
void operator=(const PLICU &)=delete
Disallow default bitwise assignment.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.