surfaceIteratorIso.H
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) 2019-2020 DLR
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 Class
27  Foam::surfaceIteratorIso
28 
29 Description
30  Finds the isovalue that matches the volume fraction
31 
32  Reference:
33  \verbatim
34  Roenby, J., Bredmose, H. and Jasak, H. (2016).
35  A computational method for sharp interface advection
36  Royal Society Open Science, 3
37  doi 10.1098/rsos.160405
38  \endverbatim
39 
40 Author
41  Johan Roenby, DHI, all rights reserved.
42 
43 SourceFiles
44  surfaceIteratorIso.C
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #ifndef surfaceIteratorIso_H
49 #define surfaceIteratorIso_H
50 
51 #include "fvMesh.H"
52 #include "volFields.H"
53 #include "surfaceFields.H"
54 #include "cutCellIso.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class surfaceIteratorIso Declaration
63 \*---------------------------------------------------------------------------*/
64 
66 {
67  // Private Data
68 
69  //- Mesh whose cells and faces to cut at their intersection with an
70  //- isosurface.
71  const fvMesh& mesh_;
72 
73  //- Reference to the point values
74  scalarField& ap_;
75 
76  //- Cuts the cell and returns volume centre and normal
77  cutCellIso cutCell_;
78 
79  //- Tolerance for marking of surface cells:
80  // Those with surfCellTol_ < alpha1 < 1 - surfCellTol_
81  scalar surfCellTol_;
82 
83 
84 public:
85 
86  // Constructors
87 
88  //- Construct from fvMesh and a scalarField
89  // The scalarField size should equal the number of mesh points
91  (
92  const fvMesh& mesh,
93  scalarField& pointVal,
94  const scalar tol
95  );
96 
97 
98  // Member Functions
99 
100  //- Determine if a cell is a surface cell
101  bool isASurfaceCell(const scalar alpha1) const
102  {
103  return
104  (
105  surfCellTol_ < alpha1
106  && alpha1 < 1 - surfCellTol_
107  );
108  }
109 
110  //- finds matching isoValue for the given value fraction
111  //- returns the cellStatus
112  label vofCutCell
113  (
114  const label celli,
115  const scalar alpha1,
116  const scalar tol,
117  const label maxIter
118  );
119 
120  //- The centre point of cutted volume
121  const point& subCellCentre() const
122  {
123  return cutCell_.subCellCentre();
124  }
125 
126  //- The volume of cutted volume
127  scalar subCellVolume() const
128  {
129  return cutCell_.subCellVolume();
130  }
131 
132  //- The centre of cutting isosurface
133  const point& surfaceCentre() const
134  {
135  return cutCell_.faceCentre();
136  }
137 
138  //- The area vector of cutting isosurface
139  const vector& surfaceArea() const
140  {
141  return cutCell_.faceArea();
142  }
143 
144  //- Volume of Fluid for cellI (subCellVolume_/mesh_.V()[cellI])
145  scalar VolumeOfFluid() const
146  {
147  return cutCell_.VolumeOfFluid();
148  }
149 
150  //- The cutValue
151  scalar cutValue() const
152  {
153  return cutCell_.cutValue();
154  }
155 
156  //- The cellStatus
157  label cellStatus()
158  {
159  return cutCell_.cellStatus();
160  }
161 
162  //- The points of the cutting isosurface in sorted order
164  {
165  return cutCell_.facePoints();
166  }
167 };
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 } // End namespace Foam
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 #endif
177 
178 // ************************************************************************* //
Foam::surfaceFields.
const point & subCellCentre() const
The centre point of cutted volume.
const vector & faceArea() const noexcept
Returns the area normal vector of the cutting PLICface.
Definition: cutCellIso.H:216
const point & subCellCentre() const noexcept
Returns subCellCentre.
Definition: cutCellIso.H:187
Class for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with an isosurface defined ...
Definition: cutCellIso.H:72
const point & faceCentre() const noexcept
Returns the centre of the cutting PLICface.
Definition: cutCellIso.H:208
scalar cutValue() const noexcept
Returns cutValue.
Definition: cutCellIso.H:240
const DynamicList< point > & facePoints()
The points of the cutting isosurface in sorted order.
label cellStatus() const noexcept
Returns cellStatus.
Definition: cutCellIso.H:224
scalar VolumeOfFluid() const
Volume of Fluid for cellI (subCellVolume_/mesh_.V()[cellI])
const vector & surfaceArea() const
The area vector of cutting isosurface.
surfaceIteratorIso(const fvMesh &mesh, scalarField &pointVal, const scalar tol)
Construct from fvMesh and a scalarField.
scalar subCellVolume() const noexcept
Returns subCellVolume.
Definition: cutCellIso.H:195
dynamicFvMesh & mesh
label vofCutCell(const label celli, const scalar alpha1, const scalar tol, const label maxIter)
finds matching isoValue for the given value fraction returns the cellStatus
scalar VolumeOfFluid() const noexcept
Returns volume of fluid value.
Definition: cutCellIso.H:232
const point & surfaceCentre() const
The centre of cutting isosurface.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Finds the isovalue that matches the volume fraction.
scalar cutValue() const
The cutValue.
vector point
Point is a vector.
Definition: point.H:37
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const DynamicList< point > & facePoints()
Returns the points of the cutting PLICface.
Definition: cutCellIso.C:166
bool isASurfaceCell(const scalar alpha1) const
Determine if a cell is a surface cell.
label cellStatus()
The cellStatus.
scalar subCellVolume() const
The volume of cutted volume.
Namespace for OpenFOAM.
const volScalarField & alpha1