Loading...
Searching...
No Matches
levelSet.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
38(
39 const fvMesh& mesh,
40 const scalarField& levelC,
41 const scalarField& levelP,
42 const bool above
43)
44{
46 (
47 "levelSetFraction",
49 mesh,
51 );
52 auto& result = tResult.ref();
53
54 forAll(result, cI)
55 {
56 const List<tetIndices> cellTetIs =
58
59 scalar v = 0, r = 0;
60
61 forAll(cellTetIs, cellTetI)
62 {
63 const triFace triIs = cellTetIs[cellTetI].faceTriIs(mesh);
64
66 tet =
67 {
68 mesh.cellCentres()[cI],
69 mesh.points()[triIs[0]],
70 mesh.points()[triIs[1]],
71 mesh.points()[triIs[2]]
72 };
74 level =
75 {
76 levelC[cI],
77 levelP[triIs[0]],
78 levelP[triIs[1]],
79 levelP[triIs[2]]
80 };
81
82 v += cut::volumeOp()(tet);
83
84 if (above)
85 {
86 r += tetCut(tet, level, cut::volumeOp(), cut::noOp());
87 }
88 else
89 {
90 r += tetCut(tet, level, cut::noOp(), cut::volumeOp());
91 }
92 }
93
94 result[cI] = r/v;
95 }
96
97 return tResult;
98}
99
100
102(
103 const fvPatch& patch,
104 const scalarField& levelF,
105 const scalarField& levelP,
106 const bool above
107)
108{
109 auto tResult = tmp<scalarField>::New(patch.size(), Zero);
110 auto& result = tResult.ref();
111
112 forAll(result, fI)
113 {
114 const face& f = patch.patch().localFaces()[fI];
115
116 vector a(Zero);
117 vector r(Zero);
118
119 for (label edgei = 0; edgei < f.nEdges(); ++edgei)
120 {
121 const edge e = f.edge(edgei);
122
123 const FixedList<point, 3>
124 tri =
125 {
126 patch.patch().faceCentres()[fI],
127 patch.patch().localPoints()[e[0]],
128 patch.patch().localPoints()[e[1]]
129 };
130 const FixedList<scalar, 3>
131 level =
132 {
133 levelF[fI],
134 levelP[e[0]],
135 levelP[e[1]]
136 };
137
138 a += cut::areaOp()(tri);
139
140 if (above)
141 {
142 r += triCut(tri, level, cut::areaOp(), cut::noOp());
143 }
144 else
145 {
146 r += triCut(tri, level, cut::noOp(), cut::areaOp());
147 }
148 }
149
150 result[fI] = a/magSqr(a) & r;
151 }
152
153 return tResult;
154}
155
156// ************************************************************************* //
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field....
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
@ NO_REGISTER
Do not request registration (bool: false).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
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.
const dimensionSet dimless
Dimensionless.
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:31
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299