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 -------------------------------------------------------------------------------
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 
35 template<class Type>
37 (
38  const fvMesh& mesh,
39  const scalarField& levelC,
40  const scalarField& levelP,
41  const DimensionedField<Type, volMesh>& positiveC,
42  const DimensionedField<Type, pointMesh>& positiveP,
43  const DimensionedField<Type, volMesh>& negativeC,
44  const DimensionedField<Type, pointMesh>& negativeP
45 )
46 {
47  tmp<DimensionedField<Type, volMesh>> tResult
48  (
49  new DimensionedField<Type, volMesh>
50  (
51  IOobject
52  (
53  positiveC.name() + ":levelSetAverage",
54  mesh.time().timeName(),
55  mesh
56  ),
57  mesh,
58  dimensioned<Type>(positiveC.dimensions(), Zero)
59  )
60  );
61  DimensionedField<Type, volMesh>& result = tResult.ref();
62 
63  forAll(result, cI)
64  {
65  const List<tetIndices> cellTetIs =
66  polyMeshTetDecomposition::cellTetIndices(mesh, cI);
67 
68  scalar v = 0;
69  Type r = Zero;
70 
71  forAll(cellTetIs, cellTetI)
72  {
73  const triFace triIs = cellTetIs[cellTetI].faceTriIs(mesh);
74 
75  const FixedList<point, 4>
76  tet =
77  {
78  mesh.cellCentres()[cI],
79  mesh.points()[triIs[0]],
80  mesh.points()[triIs[1]],
81  mesh.points()[triIs[2]]
82  };
83  const FixedList<scalar, 4>
84  level =
85  {
86  levelC[cI],
87  levelP[triIs[0]],
88  levelP[triIs[1]],
89  levelP[triIs[2]]
90  };
91  const cut::volumeIntegrateOp<Type>
92  positive = FixedList<Type, 4>
93  ({
94  positiveC[cI],
95  positiveP[triIs[0]],
96  positiveP[triIs[1]],
97  positiveP[triIs[2]]
98  });
99  const cut::volumeIntegrateOp<Type>
100  negative = FixedList<Type, 4>
101  ({
102  negativeC[cI],
103  negativeP[triIs[0]],
104  negativeP[triIs[1]],
105  negativeP[triIs[2]]
106  });
107 
108  v += cut::volumeOp()(tet);
109 
110  r += tetCut(tet, level, positive, negative);
111  }
112 
113  result[cI] = r/v;
114  }
115 
116  return tResult;
117 }
118 
119 
120 template<class Type>
122 (
123  const fvPatch& patch,
124  const scalarField& levelF,
125  const scalarField& levelP,
126  const Field<Type>& positiveF,
127  const Field<Type>& positiveP,
128  const Field<Type>& negativeF,
129  const Field<Type>& negativeP
130 )
131 {
132  typedef typename outerProduct<Type, vector>::type sumType;
133 
134  tmp<Field<Type>> tResult(new Field<Type>(patch.size(), Zero));
135  Field<Type>& result = tResult.ref();
136 
137  forAll(result, fI)
138  {
139  const face& f = patch.patch().localFaces()[fI];
140 
141  vector a(Zero);
142  sumType r = Zero;
143 
144  for (label eI = 0; eI < f.size(); ++eI)
145  {
146  const edge e = f.faceEdge(eI);
147 
148  const FixedList<point, 3>
149  tri =
150  {
151  patch.patch().faceCentres()[fI],
152  patch.patch().localPoints()[e[0]],
153  patch.patch().localPoints()[e[1]]
154  };
155  const FixedList<scalar, 3>
156  level =
157  {
158  levelF[fI],
159  levelP[e[0]],
160  levelP[e[1]]
161  };
162  const cut::areaIntegrateOp<Type>
163  positive = FixedList<Type, 3>
164  ({
165  positiveF[fI],
166  positiveP[e[0]],
167  positiveP[e[1]]
168  });
169  const cut::areaIntegrateOp<Type>
170  negative = FixedList<Type, 3>
171  ({
172  negativeF[fI],
173  negativeP[e[0]],
174  negativeP[e[1]]
175  });
176 
177  a += cut::areaOp()(tri);
178 
179  r += triCut(tri, level, positive, negative);
180  }
181 
182  result[fI] = a/magSqr(a) & r;
183  }
184 
185  return tResult;
186 }
187 
188 
189 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
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
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)
levelSet.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
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.
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
triFace
face triFace(3)
tetIndices.H
Foam::levelSetAverage
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.
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.