All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MPLICU.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::MPLICU
26 
27 Description
28  Velocity-weighted Multicut Piecewise-Linear Interface Calculation (MPLICU)
29  corrected scheme is a surface interpolation scheme for flux calculation in
30  advection of a bounded variable, e.g. phase fraction and for interface
31  capturing in the volume of fluid (VoF) method.
32 
33  The interface is represented by multiple 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 cuts.
38 
39  Three progressively more complex algorithms are used to ensure the cell
40  volume fraction is accurately reproduced:
41  -# single cut: cuts all the cell faces regardless the order
42  -# multi cut: topological face-edge-face walk which can split cell into
43  multiple sub-volumes
44  -# tetrahedron cut: decomposes cell into tetrahedrons which are cut
45 
46  Additionally the face point velocity values are used to calculate the face
47  flux which is likely to be more accurate in the presence of high shear.
48 
49  Example:
50  \verbatim
51  divSchemes
52  {
53  .
54  .
55  div(phi,alpha1) Gauss MPLICU;
56  .
57  .
58  }
59  \endverbatim
60 
61 See also
62  Foam::MPLIC
63  Foam::PLIC
64  Foam::PLICU
65  Foam::interfaceCompression
66 
67 SourceFiles
68  MPLICU.C
69 
70 \*---------------------------------------------------------------------------*/
71 
72 #ifndef MPLICU_H
73 #define MPLICU_H
74 
75 #include "MPLIC.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 namespace Foam
80 {
81 
82 /*---------------------------------------------------------------------------*\
83  Class MPLICU Declaration
84 \*---------------------------------------------------------------------------*/
85 
86 class MPLICU
87 :
88  public MPLIC
89 {
90 
91 public:
92 
93  //- Runtime type information
94  TypeName("MPLICU");
95 
96 
97  // Constructors
98 
99  //- Construct from faceFlux and Istream
101  (
102  const fvMesh& mesh,
103  const surfaceScalarField& faceFlux,
104  Istream& is
105  )
106  :
107  MPLIC(mesh, faceFlux, is)
108  {}
109 
110 
111  // Member Functions
112 
113  //- Return the interpolation weighting factors
115  (
117  ) const
118  {
120 
121  return tmp<surfaceScalarField>(nullptr);
122  }
123 
124  //- Return the face-interpolate of the given cell field
126  (
128  ) const;
129 
130 
131  // Member Operators
132 
133  //- Disallow default bitwise assignment
134  void operator=(const MPLICU&) = delete;
135 };
136 
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace Foam
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
144 #endif
145 
146 // ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void operator=(const MPLICU &)=delete
Disallow default bitwise assignment.
const fvMesh & mesh() const
Return mesh reference.
TypeName("MPLICU")
Runtime type information.
Velocity-weighted Multicut Piecewise-Linear Interface Calculation (MPLICU) corrected scheme is a surf...
Definition: MPLICU.H:85
virtual tmp< surfaceScalarField > interpolate(const GeometricField< scalar, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Definition: MPLICU.C:44
MPLICU(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from faceFlux and Istream.
Definition: MPLICU.H:100
virtual tmp< surfaceScalarField > weights(const GeometricField< scalar, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: MPLICU.H:114
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
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
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
Namespace for OpenFOAM.