levelSet.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 | Copyright (C) 2017 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 "DimensionedField.H"
38 #include "fvMesh.H"
39 #include "pointMesh.H"
40 #include "volMesh.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 //- Calculate the average value of two fields, one on each side of a level set
50 // The level set and the fields are given both on the points and cell-centres.
51 template<class Type>
52 tmp<DimensionedField<Type, volMesh>> levelSetAverage
53 (
54  const fvMesh& mesh,
55  const scalarField& levelC,
56  const scalarField& levelP,
57  const DimensionedField<Type, volMesh>& positiveC,
58  const DimensionedField<Type, pointMesh>& positiveP,
59  const DimensionedField<Type, volMesh>& negativeC,
60  const DimensionedField<Type, pointMesh>& negativeP
61 );
62 
63 //- As the above overload, but on the faces of a patch
64 template<class Type>
65 tmp<Field<Type>> levelSetAverage
66 (
67  const fvPatch& patch,
68  const scalarField& levelF,
69  const scalarField& levelP,
70  const Field<Type>& positiveF,
71  const Field<Type>& positiveP,
72  const Field<Type>& negativeF,
73  const Field<Type>& negativeP
74 );
75 
76 //- Calculate the volume-fraction that a level set occupies. This gives the the
77 // same result as levelSetAverage if the fields passed to the latter are
78 // uniformly 0 and 1. The above flag flips the direction.
79 tmp<DimensionedField<scalar, volMesh>> levelSetFraction
80 (
81  const fvMesh& mesh,
82  const scalarField& levelC,
83  const scalarField& levelP,
84  const bool above
85 );
86 
87 //- As the above oveload, but on the faces of a patch
88 tmp<scalarField> levelSetFraction
89 (
90  const fvPatch& patch,
91  const scalarField& levelF,
92  const scalarField& levelP,
93  const bool above
94 );
95 
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 
98 } // End namespace Foam
99 
100 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102 #ifdef NoRepository
103  #include "levelSetTemplates.C"
104 #endif
105 
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108 #endif
109 
110 // ************************************************************************* //
tmp< DimensionedField< Type, volMesh > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const DimensionedField< Type, volMesh > &positiveC, const DimensionedField< Type, pointMesh > &positiveP, const DimensionedField< Type, volMesh > &negativeC, const DimensionedField< Type, pointMesh > &negativeP)
Calculate the average value of two fields, one on each side of a level set.
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< DimensionedField< scalar, volMesh > > 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 the.
Definition: levelSet.C:35
Namespace for OpenFOAM.