All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PLIC.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::PLIC
26 
27 Description
28  Piecewise-Linear Interface Calculation (PLIC) corrected scheme is a surface
29  interpolation scheme for flux calculation in advection of a bounded
30  variable, e.g. phase fraction and for interface capturing in the volume of
31  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  Example:
42  \verbatim
43  divSchemes
44  {
45  .
46  .
47  div(phi,alpha) Gauss PLIC interfaceCompression vanLeer 1;
48  .
49  .
50  }
51  \endverbatim
52 
53 See also
54  Foam::PLICU
55  Foam::MPLIC
56  Foam::MPLICU
57  Foam::interfaceCompression
58 
59 SourceFiles
60  PLIC.C
61 
62 \*---------------------------------------------------------------------------*/
63 
64 #ifndef PLIC_H
65 #define PLIC_H
66 
67 #include "MPLIC.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 
74 /*---------------------------------------------------------------------------*\
75  Class PLIC Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class PLIC
79 :
80  public MPLIC
81 {
82 
83 protected:
84 
85  // Protected member data
86 
87  const surfaceScalarField& phi_;
88 
89  //- Base scheme to which the compression is applied
91 
92 
93 public:
94 
95  //- Runtime type information
96  TypeName("PLIC");
97 
98 
99  // Constructors
100 
101  //- Construct from faceFlux and Istream
103  (
104  const fvMesh& mesh,
105  const surfaceScalarField& faceFlux,
106  Istream& is
107  )
108  :
109  MPLIC(mesh, faceFlux, is),
110  phi_(faceFlux),
111  tScheme_
112  (
113  surfaceInterpolationScheme<scalar>::New(mesh, faceFlux, is)
114  )
115  {}
116 
117 
118  // Member Functions
119 
120  //- Return the face-interpolate of the given cell field
122  (
123  const volScalarField& vf
124  ) const;
125 
126 
127  // Member Operators
128 
129  //- Disallow default bitwise assignment
130  void operator=(const PLIC&) = delete;
131 };
132 
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 } // End namespace Foam
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 #endif
141 
142 // ************************************************************************* //
void operator=(const PLIC &)=delete
Disallow default bitwise assignment.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< scalar > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
const surfaceScalarField & phi_
Definition: PLIC.H:86
TypeName("PLIC")
Runtime type information.
PLIC(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from faceFlux and Istream.
Definition: PLIC.H:102
Piecewise-Linear Interface Calculation (PLIC) corrected scheme is a surface interpolation scheme for ...
Definition: PLIC.H:77
virtual tmp< surfaceScalarField > interpolate(const volScalarField &vf) const
Return the face-interpolate of the given cell field.
Definition: PLIC.C:43
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Multicut Piecewise-Linear Interface Calculation (MPLIC) corrected scheme is a surface interpolation s...
Definition: MPLIC.H:83
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceInterpolationScheme< scalar > > tScheme_
Base scheme to which the compression is applied.
Definition: PLIC.H:89
Abstract base class for surface interpolation schemes.
Namespace for OpenFOAM.