Loading...
Searching...
No Matches
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-2023 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 auto tresult = DimensionedField<Type, volMesh>::New
49 (
50 IOobject::scopedName(positiveC.name(), "levelSetAverage"),
51 mesh,
52 Foam::zero{}, // value
53 positiveC.dimensions()
54 );
55 auto& result = tresult.ref();
56
57 forAll(result, cI)
58 {
59 const List<tetIndices> cellTetIs =
60 polyMeshTetDecomposition::cellTetIndices(mesh, cI);
61
62 scalar v = 0;
63 Type r = Zero;
64
65 forAll(cellTetIs, cellTetI)
66 {
67 const triFace triIs = cellTetIs[cellTetI].faceTriIs(mesh);
68
69 const FixedList<point, 4>
70 tet =
71 {
72 mesh.cellCentres()[cI],
73 mesh.points()[triIs[0]],
74 mesh.points()[triIs[1]],
75 mesh.points()[triIs[2]]
76 };
77 const FixedList<scalar, 4>
78 level =
79 {
80 levelC[cI],
81 levelP[triIs[0]],
82 levelP[triIs[1]],
83 levelP[triIs[2]]
84 };
85 const cut::volumeIntegrateOp<Type>
86 positive = FixedList<Type, 4>
87 ({
88 positiveC[cI],
89 positiveP[triIs[0]],
90 positiveP[triIs[1]],
91 positiveP[triIs[2]]
92 });
93 const cut::volumeIntegrateOp<Type>
94 negative = FixedList<Type, 4>
95 ({
96 negativeC[cI],
97 negativeP[triIs[0]],
98 negativeP[triIs[1]],
99 negativeP[triIs[2]]
100 });
101
102 v += cut::volumeOp()(tet);
103
104 r += tetCut(tet, level, positive, negative);
105 }
106
107 result[cI] = r/v;
108 }
109
110 return tresult;
111}
112
113
114template<class Type>
116(
117 const fvPatch& patch,
118 const scalarField& levelF,
119 const scalarField& levelP,
120 const Field<Type>& positiveF,
121 const Field<Type>& positiveP,
122 const Field<Type>& negativeF,
123 const Field<Type>& negativeP
124)
125{
126 typedef typename outerProduct<Type, vector>::type sumType;
127
128 auto tResult = tmp<Field<Type>>::New(patch.size(), Zero);
129 auto& result = tResult.ref();
130
131 forAll(result, fI)
132 {
133 const face& f = patch.patch().localFaces()[fI];
134
135 vector a(Zero);
136 sumType r = Zero;
137
138 for (label edgei = 0; edgei < f.nEdges(); ++edgei)
139 {
140 const edge e = f.edge(edgei);
141
142 const FixedList<point, 3>
143 tri =
144 {
145 patch.patch().faceCentres()[fI],
146 patch.patch().localPoints()[e[0]],
147 patch.patch().localPoints()[e[1]]
148 };
149 const FixedList<scalar, 3>
150 level =
151 {
152 levelF[fI],
153 levelP[e[0]],
154 levelP[e[1]]
155 };
156 const cut::areaIntegrateOp<Type>
157 positive = FixedList<Type, 3>
158 ({
159 positiveF[fI],
160 positiveP[e[0]],
161 positiveP[e[1]]
162 });
163 const cut::areaIntegrateOp<Type>
164 negative = FixedList<Type, 3>
165 ({
166 negativeF[fI],
167 negativeP[e[0]],
168 negativeP[e[1]]
169 });
170
171 a += cut::areaOp()(tri);
172
173 r += triCut(tri, level, positive, negative);
174 }
175
176 result[fI] = a/magSqr(a) & r;
177 }
178
179 return tResult;
180}
181
182
183// ************************************************************************* //
A class for managing temporary objects.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
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
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299