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 | 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 \*---------------------------------------------------------------------------*/
25 
26 #include "levelSet.H"
27 #include "cut.H"
29 #include "tetIndices.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
35 (
36  const fvMesh& mesh,
37  const scalarField& levelC,
38  const scalarField& levelP,
39  const bool above
40 )
41 {
43  (
45  (
46  IOobject
47  (
48  "levelSetFraction",
49  mesh.time().timeName(),
50  mesh
51  ),
52  mesh,
53  dimensionedScalar("0", dimless, 0)
54  )
55  );
56  DimensionedField<scalar, volMesh>& result = tResult.ref();
57 
58  forAll(result, cI)
59  {
60  const List<tetIndices> cellTetIs =
61  polyMeshTetDecomposition::cellTetIndices(mesh, cI);
62 
63  scalar v = 0, r = 0;
64 
65  forAll(cellTetIs, cellTetI)
66  {
67  const triFace triIs = cellTetIs[cellTetI].faceTriIs(mesh);
68 
70  tet =
71  {
72  mesh.cellCentres()[cI],
73  mesh.points()[triIs[0]],
74  mesh.points()[triIs[1]],
75  mesh.points()[triIs[2]]
76  };
78  level =
79  {
80  levelC[cI],
81  levelP[triIs[0]],
82  levelP[triIs[1]],
83  levelP[triIs[2]]
84  };
85 
86  v += cut::volumeOp()(tet);
87 
88  if (above)
89  {
90  r += tetCut(tet, level, cut::volumeOp(), cut::noOp());
91  }
92  else
93  {
94  r += tetCut(tet, level, cut::noOp(), cut::volumeOp());
95  }
96  }
97 
98  result[cI] = r/v;
99  }
100 
101  return tResult;
102 }
103 
104 
106 (
107  const fvPatch& patch,
108  const scalarField& levelF,
109  const scalarField& levelP,
110  const bool above
111 )
112 {
113  tmp<scalarField> tResult(new scalarField(patch.size(), 0));
114  scalarField& result = tResult.ref();
115 
116  forAll(result, fI)
117  {
118  const face& f = patch.patch().localFaces()[fI];
119 
120  vector a = vector::zero, r = vector::zero;
121 
122  for(label eI = 0; eI < f.size(); ++ eI)
123  {
124  const edge e = f.faceEdge(eI);
125 
126  const FixedList<point, 3>
127  tri =
128  {
129  patch.patch().faceCentres()[fI],
130  patch.patch().localPoints()[e[0]],
131  patch.patch().localPoints()[e[1]]
132  };
134  level =
135  {
136  levelF[fI],
137  levelP[e[0]],
138  levelP[e[1]]
139  };
140 
141  a += cut::areaOp()(tri);
142 
143  if (above)
144  {
145  r += triCut(tri, level, cut::areaOp(), cut::noOp());
146  }
147  else
148  {
149  r += triCut(tri, level, cut::noOp(), cut::areaOp());
150  }
151  }
152 
153  result[fI] = a/magSqr(a) & r;
154  }
155 
156  return tResult;
157 }
158 
159 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
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. Templated.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1011
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
dynamicFvMesh & mesh
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
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
const vectorField & cellCentres() const
virtual label size() const
Return size.
Definition: fvPatch.H:161
dimensioned< scalar > magSqr(const dimensioned< Type > &)
volScalarField scalarField(fieldObject, mesh)
const Field< PointType > & localPoints() const
Return pointField of points in patch.
labelList f(nPoints)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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...
A class for managing temporary objects.
Definition: PtrList.H:53
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:110
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:284