levelSetTemplates.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "levelSet.H"
30#include "cut.H"
32#include "tetIndices.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const fvMesh& mesh,
40 const scalarField& levelC,
41 const scalarField& levelP,
42 const DimensionedField<Type, volMesh>& positiveC,
43 const DimensionedField<Type, pointMesh>& positiveP,
44 const DimensionedField<Type, volMesh>& negativeC,
45 const DimensionedField<Type, pointMesh>& negativeP
46)
47{
48 tmp<DimensionedField<Type, volMesh>> tResult
49 (
50 new DimensionedField<Type, volMesh>
51 (
52 IOobject
53 (
54 positiveC.name() + ":levelSetAverage",
55 mesh.time().timeName(),
56 mesh
57 ),
58 mesh,
59 dimensioned<Type>(positiveC.dimensions(), Zero)
60 )
61 );
62 DimensionedField<Type, volMesh>& result = tResult.ref();
63
64 forAll(result, cI)
65 {
66 const List<tetIndices> cellTetIs =
67 polyMeshTetDecomposition::cellTetIndices(mesh, cI);
68
69 scalar v = 0;
70 Type r = Zero;
71
72 forAll(cellTetIs, cellTetI)
73 {
74 const triFace triIs = cellTetIs[cellTetI].faceTriIs(mesh);
75
76 const FixedList<point, 4>
77 tet =
78 {
79 mesh.cellCentres()[cI],
80 mesh.points()[triIs[0]],
81 mesh.points()[triIs[1]],
82 mesh.points()[triIs[2]]
83 };
84 const FixedList<scalar, 4>
85 level =
86 {
87 levelC[cI],
88 levelP[triIs[0]],
89 levelP[triIs[1]],
90 levelP[triIs[2]]
91 };
92 const cut::volumeIntegrateOp<Type>
93 positive = FixedList<Type, 4>
94 ({
95 positiveC[cI],
96 positiveP[triIs[0]],
97 positiveP[triIs[1]],
98 positiveP[triIs[2]]
99 });
100 const cut::volumeIntegrateOp<Type>
101 negative = FixedList<Type, 4>
102 ({
103 negativeC[cI],
104 negativeP[triIs[0]],
105 negativeP[triIs[1]],
106 negativeP[triIs[2]]
107 });
108
109 v += cut::volumeOp()(tet);
110
111 r += tetCut(tet, level, positive, negative);
112 }
113
114 result[cI] = r/v;
115 }
116
117 return tResult;
118}
119
120
121template<class Type>
123(
124 const fvPatch& patch,
125 const scalarField& levelF,
126 const scalarField& levelP,
127 const Field<Type>& positiveF,
128 const Field<Type>& positiveP,
129 const Field<Type>& negativeF,
130 const Field<Type>& negativeP
131)
132{
133 typedef typename outerProduct<Type, vector>::type sumType;
134
135 tmp<Field<Type>> tResult(new Field<Type>(patch.size(), Zero));
136 Field<Type>& result = tResult.ref();
137
138 forAll(result, fI)
139 {
140 const face& f = patch.patch().localFaces()[fI];
141
142 vector a(Zero);
143 sumType r = Zero;
144
145 for (label edgei = 0; edgei < f.nEdges(); ++edgei)
146 {
147 const edge e = f.edge(edgei);
148
149 const FixedList<point, 3>
150 tri =
151 {
152 patch.patch().faceCentres()[fI],
153 patch.patch().localPoints()[e[0]],
154 patch.patch().localPoints()[e[1]]
155 };
156 const FixedList<scalar, 3>
157 level =
158 {
159 levelF[fI],
160 levelP[e[0]],
161 levelP[e[1]]
162 };
163 const cut::areaIntegrateOp<Type>
164 positive = FixedList<Type, 3>
165 ({
166 positiveF[fI],
167 positiveP[e[0]],
168 positiveP[e[1]]
169 });
170 const cut::areaIntegrateOp<Type>
171 negative = FixedList<Type, 3>
172 ({
173 negativeF[fI],
174 negativeP[e[0]],
175 negativeP[e[1]]
176 });
177
178 a += cut::areaOp()(tri);
179
180 r += triCut(tri, level, positive, negative);
181 }
182
183 result[fI] = a/magSqr(a) & r;
184 }
185
186 return tResult;
187}
188
189
190// ************************************************************************* //
A class for managing temporary objects.
Definition: tmp.H:65
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
dynamicFvMesh & mesh
const std::string patch
OpenFOAM patch number as a std::string.
cut::opAddResult< AboveOp, BelowOp >::type triCut(const FixedList< point, 3 > &tri, const FixedList< scalar, 3 > &level, const AboveOp &aboveOp, const BelowOp &belowOp)
Cut a triangle along the zero plane defined by the given levels.
cut::opAddResult< AboveOp, BelowOp >::type tetCut(const FixedList< point, 4 > &tet, const FixedList< scalar, 4 > &level, const AboveOp &aboveOp, const BelowOp &belowOp)
As triCut, but for a tetrahedron.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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.
face triFace(3)
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333