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 // ************************************************************************* //
volFields.H
Foam::surfaceIteratorIso::VolumeOfFluid
scalar VolumeOfFluid() const
Volume of Fluid for cellI (subCellVolume_/mesh_.V()[cellI])
Definition: surfaceIteratorIso.H:144
Foam::surfaceIteratorIso::surfaceIteratorIso
surfaceIteratorIso(const fvMesh &mesh, scalarField &pointVal, const scalar tol)
Construct from fvMesh and a scalarField.
Definition: surfaceIteratorIso.C:34
Foam::surfaceIteratorIso::facePoints
const DynamicList< point > & facePoints()
The points of the cutting isosurface in sorted order.
Definition: surfaceIteratorIso.H:162
Foam::DynamicList< point >
Foam::cutCellIso::cellStatus
label cellStatus() const noexcept
Returns cellStatus.
Definition: cutCellIso.H:177
Foam::cutCellIso::facePoints
const DynamicList< point > & facePoints()
Returns the points of the cutting PLICface.
Definition: cutCellIso.C:173
surfaceFields.H
Foam::surfaceFields.
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::surfaceIteratorIso::isASurfaceCell
bool isASurfaceCell(const scalar alpha1) const
Determine if a cell is a surface cell.
Definition: surfaceIteratorIso.H:100
Foam::surfaceIteratorIso::subCellCentre
const point & subCellCentre() const
The centre point of cutted volume.
Definition: surfaceIteratorIso.H:120
Foam::cutCellIso::cutValue
scalar cutValue() const noexcept
Returns cutValue.
Definition: cutCellIso.H:189
cutCellIso.H
Foam::Field< scalar >
Foam::surfaceIteratorIso::cutValue
scalar cutValue() const
The cutValue.
Definition: surfaceIteratorIso.H:150
Foam::cutCellIso::subCellVolume
scalar subCellVolume() const noexcept
Returns subCellVolume.
Definition: cutCellIso.H:156
Foam::surfaceIteratorIso::surfaceArea
const vector & surfaceArea() const
The area vector of cutting isosurface.
Definition: surfaceIteratorIso.H:138
Foam::cutCellIso::VolumeOfFluid
scalar VolumeOfFluid() const noexcept
Returns volume of fluid value.
Definition: cutCellIso.H:183
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::surfaceIteratorIso
Finds the isovalue that matches the volume fraction.
Definition: surfaceIteratorIso.H:64
Foam::cutCellIso
Class for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with an isosurface defined ...
Definition: cutCellIso.H:76
Foam::cutCellIso::subCellCentre
const point & subCellCentre() const noexcept
Returns subCellCentre.
Definition: cutCellIso.H:150
Foam::Vector< scalar >
Foam::surfaceIteratorIso::cellStatus
label cellStatus()
The cellStatus.
Definition: surfaceIteratorIso.H:156
Foam::surfaceIteratorIso::surfaceCentre
const point & surfaceCentre() const
The centre of cutting isosurface.
Definition: surfaceIteratorIso.H:132
Foam::cutCellIso::faceCentre
const point & faceCentre() const noexcept
Returns the centre of the cutting PLICface.
Definition: cutCellIso.H:165
Foam::surfaceIteratorIso::vofCutCell
label vofCutCell(const label celli, const scalar alpha1, const scalar tol, const label maxIter)
Definition: surfaceIteratorIso.C:50
Foam::cutCellIso::faceArea
const vector & faceArea() const noexcept
Returns the area normal vector of the cutting PLICface.
Definition: cutCellIso.H:171
Foam::surfaceIteratorIso::subCellVolume
scalar subCellVolume() const
The volume of cutted volume.
Definition: surfaceIteratorIso.H:126