Loading...
Searching...
No Matches
cutCellPLIC.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) 2016-2017 DHI
9 Copyright (C) 2018-2019 Johan Roenby
10 Copyright (C) 2019-2020 DLR
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28\*---------------------------------------------------------------------------*/
29
30#include "cutCellPLIC.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
35:
37 mesh_(mesh),
38 cellI_(-1),
39 normal_(Zero),
40 cutValue_(0),
41 cutFace_(mesh_),
42 cutFaceCentres_(10),
43 cutFaceAreas_(10),
44 plicFaceEdges_(10),
45 facePoints_(10),
46 faceCentre_(Zero),
47 faceArea_(Zero),
48 subCellCentre_(Zero),
49 subCellVolume_(-10),
50 VOF_(-10),
51 cellStatus_(-1)
54}
55
56
57// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58
60(
61 const label celli,
62 const scalar cutValue,
63 const vector& normal
64)
65{
67 normal_ = normal;
68 cellI_ = celli;
69 cutValue_ = cutValue;
70 const cell& c = mesh_.cells()[celli];
71
72 vector base = mesh_.C()[cellI_] + normal_ * cutValue_;
73 bool fullyBelow = true;
74 bool fullyAbove = true;
75
76 label nFaceBelowInterface = 0;
77
78 // loop over all cell faces
79 for (const label faceI : c)
80 {
81 const label faceStatus = cutFace_.calcSubFace(faceI, normal_, base);
82
83 if (faceStatus == 0) // face is cut
84 {
85 cutFaceCentres_.append(cutFace_.subFaceCentre());
86 cutFaceAreas_.append(cutFace_.subFaceArea());
87 plicFaceEdges_.append(cutFace_.surfacePoints());
88 fullyBelow = false;
89 fullyAbove = false;
90 }
91 else if (faceStatus == -1) // face fully below
92 {
93 cutFaceCentres_.append(cutFace_.subFaceCentre());
94 cutFaceAreas_.append(cutFace_.subFaceArea());
95 fullyAbove = false;
96 nFaceBelowInterface++;
97 }
98 else
99 {
100 fullyBelow = false;
101 }
102 }
103
104 if (!fullyBelow && !fullyAbove) // cell cut at least at one face
105 {
106 cellStatus_ = 0;
107
108 // calc faceArea and faceCentre
110 (
111 plicFaceEdges_,
112 average(cutFaceCentres_),
113 faceArea_,
114 faceCentre_
115 );
116
117 // In the rare but occuring cases where a cell is only touched at a
118 // point or a line the isoFaceArea_ will have zero length and here the
119 // cell should be treated as either completely empty or full.
120 if (mag(faceArea_) < ROOTVSMALL)
121 {
122 if (nFaceBelowInterface == 0)
123 {
124 // Cell fully above isosurface
125 cellStatus_ = 1;
126 subCellCentre_ = Zero;
127 subCellVolume_ = 0;
128 VOF_ = 0;
129 return cellStatus_;
130 }
131 else
132 {
133 // Cell fully below isosurface
134 cellStatus_ = -1;
135 subCellCentre_ = mesh_.C()[cellI_];
136 subCellVolume_ = mesh_.V()[cellI_];
137 VOF_ = 1;
138 return cellStatus_;
139 }
140 }
141
142 cutFaceCentres_.append(faceCentre_);
143 cutFaceAreas_.append(faceArea_);
144
145 // calc volume and sub cell centre
147 (
148 cutFaceCentres_,
149 cutFaceAreas_,
150 subCellCentre_,
151 subCellVolume_
152 );
153
154 VOF_ = subCellVolume_ / mesh_.V()[cellI_];
155 }
156 else if (fullyAbove) // cell fully above isosurface
157 {
158 cellStatus_ = 1;
159 subCellCentre_ = Zero;
160 subCellVolume_ = 0;
161 VOF_ = 0;
162 }
163 else if (fullyBelow) // cell fully below isosurface
164 {
165 cellStatus_ = -1;
166 subCellCentre_ = mesh_.C()[cellI_];
167 subCellVolume_ = mesh_.V()[cellI_];
168 VOF_ = 1;
169 }
170
171 return cellStatus_;
172}
173
174
176{
177 if (facePoints_.empty())
178 {
179 // get face points in sorted order
180 calcIsoFacePointsFromEdges
181 (
182 faceArea_,
183 faceCentre_,
184 plicFaceEdges_,
185 facePoints_
186 );
187 }
188
189 return facePoints_;
190}
191
192
194{
195 cellI_ = -1;
196 cutValue_ = 0;
197 cutFaceCentres_.clear();
198 cutFaceAreas_.clear();
199 plicFaceEdges_.clear();
200 facePoints_.clear();
201 faceCentre_ = Zero;
202 faceArea_ = Zero;
203 subCellCentre_ = Zero;
204 subCellVolume_ = -10;
205 VOF_ = -10;
206 cellStatus_ = -1;
207}
208
209
210// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
label calcSubCell(const label celli, const scalar cutValue, const vector &normal)
Sets internal values and returns face status.
Definition cutCellPLIC.C:53
scalar cutValue() const noexcept
Returns cutValue.
const DynamicList< point > & facePoints()
Returns the points of the cutting PLICface.
void clearStorage()
Resets internal values.
cutCellPLIC(const fvMesh &mesh)
Construct from fvMesh.
Definition cutCellPLIC.C:27
static void calcCellData(const DynamicList< point > &cutFaceCentres, const DynamicList< vector > &cutFaceAreas, vector &subCellCentre, scalar &subCellVolume)
Calculates volume and centre of the cutted cell.
Definition cutCell.C:49
cutCell(const fvMesh &mesh)
Construct from fvMesh.
Definition cutCell.C:34
static void calcIsoFacePointsFromEdges(const vector &faceArea, const vector &faceCentre, const DynamicList< DynamicList< point > > &faceEdges, DynamicList< point > &facePoints)
Calculates the point of the cutting face.
Definition cutCell.C:158
static void calcGeomDataCutFace(const DynamicList< DynamicList< point > > &faceEdges, const vector &subCellCentre, vector &faceArea, vector &faceCentre)
Calculates area and centre of the cutting face.
Definition cutCell.C:88
label calcSubFace(const label faceI, const vector &normal, const vector &base)
Calculate cut points along edges of faceI.
Definition cutFacePLIC.C:46
const vector & subFaceArea() const noexcept
Returns area vector of cutted face.
const DynamicList< point > & surfacePoints() const noexcept
Returns point of the face in sorted of cutted face.
const point & subFaceCentre() const noexcept
Returns centre of cutted face.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const volVectorField & C() const
Return cell centres as volVectorField.
const cellList & cells() const
dynamicFvMesh & mesh
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)