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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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  (
48  (
49  IOobject
50  (
51  "levelSetFraction",
52  mesh.time().timeName(),
53  mesh
54  ),
55  mesh,
57  )
58  );
59  DimensionedField<scalar, volMesh>& result = tResult.ref();
60 
61  forAll(result, cI)
62  {
63  const List<tetIndices> cellTetIs =
64  polyMeshTetDecomposition::cellTetIndices(mesh, cI);
65 
66  scalar v = 0, r = 0;
67 
68  forAll(cellTetIs, cellTetI)
69  {
70  const triFace triIs = cellTetIs[cellTetI].faceTriIs(mesh);
71 
73  tet =
74  {
75  mesh.cellCentres()[cI],
76  mesh.points()[triIs[0]],
77  mesh.points()[triIs[1]],
78  mesh.points()[triIs[2]]
79  };
81  level =
82  {
83  levelC[cI],
84  levelP[triIs[0]],
85  levelP[triIs[1]],
86  levelP[triIs[2]]
87  };
88 
89  v += cut::volumeOp()(tet);
90 
91  if (above)
92  {
93  r += tetCut(tet, level, cut::volumeOp(), cut::noOp());
94  }
95  else
96  {
97  r += tetCut(tet, level, cut::noOp(), cut::volumeOp());
98  }
99  }
100 
101  result[cI] = r/v;
102  }
103 
104  return tResult;
105 }
106 
107 
109 (
110  const fvPatch& patch,
111  const scalarField& levelF,
112  const scalarField& levelP,
113  const bool above
114 )
115 {
116  tmp<scalarField> tResult(new scalarField(patch.size(), Zero));
117  scalarField& result = tResult.ref();
118 
119  forAll(result, fI)
120  {
121  const face& f = patch.patch().localFaces()[fI];
122 
123  vector a(Zero);
124  vector r(Zero);
125 
126  for (label edgei = 0; edgei < f.nEdges(); ++edgei)
127  {
128  const edge e = f.edge(edgei);
129 
130  const FixedList<point, 3>
131  tri =
132  {
133  patch.patch().faceCentres()[fI],
134  patch.patch().localPoints()[e[0]],
135  patch.patch().localPoints()[e[1]]
136  };
137  const FixedList<scalar, 3>
138  level =
139  {
140  levelF[fI],
141  levelP[e[0]],
142  levelP[e[1]]
143  };
144 
145  a += cut::areaOp()(tri);
146 
147  if (above)
148  {
149  r += triCut(tri, level, cut::areaOp(), cut::noOp());
150  }
151  else
152  {
153  r += triCut(tri, level, cut::noOp(), cut::areaOp());
154  }
155  }
156 
157  result[fI] = a/magSqr(a) & r;
158  }
159 
160  return tResult;
161 }
162 
163 // ************************************************************************* //
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
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.
const dimensionSet dimless
Dimensionless.
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:65
Vector< scalar > vector
Definition: vector.H:57
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
labelList f(nPoints)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127