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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "levelSet.H"
29 #include "cut.H"
31 #include "tetIndices.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
37 (
38  const fvMesh& mesh,
39  const scalarField& levelC,
40  const scalarField& levelP,
41  const bool above
42 )
43 {
45  (
47  (
48  IOobject
49  (
50  "levelSetFraction",
51  mesh.time().timeName(),
52  mesh
53  ),
54  mesh,
56  )
57  );
58  DimensionedField<scalar, volMesh>& result = tResult.ref();
59 
60  forAll(result, cI)
61  {
62  const List<tetIndices> cellTetIs =
63  polyMeshTetDecomposition::cellTetIndices(mesh, cI);
64 
65  scalar v = 0, r = 0;
66 
67  forAll(cellTetIs, cellTetI)
68  {
69  const triFace triIs = cellTetIs[cellTetI].faceTriIs(mesh);
70 
72  tet =
73  {
74  mesh.cellCentres()[cI],
75  mesh.points()[triIs[0]],
76  mesh.points()[triIs[1]],
77  mesh.points()[triIs[2]]
78  };
80  level =
81  {
82  levelC[cI],
83  levelP[triIs[0]],
84  levelP[triIs[1]],
85  levelP[triIs[2]]
86  };
87 
88  v += cut::volumeOp()(tet);
89 
90  if (above)
91  {
92  r += tetCut(tet, level, cut::volumeOp(), cut::noOp());
93  }
94  else
95  {
96  r += tetCut(tet, level, cut::noOp(), cut::volumeOp());
97  }
98  }
99 
100  result[cI] = r/v;
101  }
102 
103  return tResult;
104 }
105 
106 
108 (
109  const fvPatch& patch,
110  const scalarField& levelF,
111  const scalarField& levelP,
112  const bool above
113 )
114 {
115  tmp<scalarField> tResult(new scalarField(patch.size(), Zero));
116  scalarField& result = tResult.ref();
117 
118  forAll(result, fI)
119  {
120  const face& f = patch.patch().localFaces()[fI];
121 
122  vector a(Zero);
123  vector r(Zero);
124 
125  for (label eI = 0; eI < f.size(); ++eI)
126  {
127  const edge e = f.faceEdge(eI);
128 
129  const FixedList<point, 3>
130  tri =
131  {
132  patch.patch().faceCentres()[fI],
133  patch.patch().localPoints()[e[0]],
134  patch.patch().localPoints()[e[1]]
135  };
137  level =
138  {
139  levelF[fI],
140  levelP[e[0]],
141  levelP[e[1]]
142  };
143 
144  a += cut::areaOp()(tri);
145 
146  if (above)
147  {
148  r += triCut(tri, level, cut::areaOp(), cut::noOp());
149  }
150  else
151  {
152  r += triCut(tri, level, cut::noOp(), cut::areaOp());
153  }
154  }
155 
156  result[fI] = a/magSqr(a) & r;
157  }
158 
159  return tResult;
160 }
161 
162 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::cut::volumeOp
Definition: cut.H:183
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::cut::areaOp
Definition: cut.H:138
Foam::levelSetFraction
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:37
polyMeshTetDecomposition.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
levelSet.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::cut::noOp
Definition: cut.H:92
Foam::tetCut
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.
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::FixedList< point, 4 >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
tetIndices.H
Foam::triCut
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.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54