levelSet.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) 2017-2022 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 Description:
25  Functions to create an averaged field from a discontinuous field defined by
26  a level-set.
27 
28 SourceFiles:
29  levelSet.C
30  levelSetTemplates.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef levelSet_H
35 #define levelSet_H
36 
37 #include "volFields.H"
38 #include "pointFields.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 //- Calculate the average value of two fields, one on each side of a level set
48 // The level set and the fields are given both on the points and cell-centres.
49 template<class Type>
51 (
52  const fvMesh& mesh,
53  const scalarField& levelC,
54  const scalarField& levelP,
55  const Field<Type>& positiveC,
56  const Field<Type>& positiveP,
57  const Field<Type>& negativeC,
58  const Field<Type>& negativeP
59 );
60 
61 //- As the above overload, but on the faces of a patch
62 template<class Type>
64 (
65  const fvPatch& patch,
66  const scalarField& levelF,
67  const scalarField& levelP,
68  const Field<Type>& positiveF,
69  const Field<Type>& positiveP,
70  const Field<Type>& negativeF,
71  const Field<Type>& negativeP
72 );
73 
74 //- As the above overload, but both in cells and on patch faces
75 template<class Type>
77 (
78  const volScalarField& levelC,
79  const pointScalarField& levelP,
80  const VolField<Type>& positiveC,
81  const PointField<Type>& positiveP,
82  const VolField<Type>& negativeC,
83  const PointField<Type>& negativeP
84 );
85 
86 //- Calculate the volume-fraction that a level set occupies. This gives the
87 // same result as levelSetAverage if the fields passed to the latter are
88 // uniformly 0 and 1. The above flag flips the direction.
90 (
91  const fvMesh& mesh,
92  const scalarField& levelC,
93  const scalarField& levelP,
94  const bool above
95 );
96 
97 //- As the above overload, but on the faces of a patch
99 (
100  const fvPatch& patch,
101  const scalarField& levelF,
102  const scalarField& levelP,
103  const bool above
104 );
105 
106 //- As the above overload, but both in cells and on patch faces
108 (
109  const volScalarField& levelC,
110  const pointScalarField& levelP,
111  const bool above
112 );
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
115 
116 } // End namespace Foam
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 #ifdef NoRepository
121  #include "levelSetTemplates.C"
122 #endif
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 #endif
127 
128 // ************************************************************************* //
Generic GeometricField class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
tmp< scalarField > levelSetFraction(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const bool above)
Calculate the volume-fraction that a level set occupies. This gives the.
Definition: levelSet.C:34
tmp< Field< Type > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const Field< Type > &positiveC, const Field< Type > &positiveP, const Field< Type > &negativeC, const Field< Type > &negativeP)
Calculate the average value of two fields, one on each side of a level set.