cutFaceAdvect.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) 2016-2017 DHI
9  Copyright (C) 2018-2019 Johan Roenby
10  Copyright (C) 2019-2020 DLR
11 -------------------------------------------------------------------------------
12 License
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/>.
27 
28 Class
29  Foam::cutFaceAdvect
30 
31 Description
32  Calculates the face fluxes
33 
34  Reference:
35  \verbatim
36  Roenby, J., Bredmose, H. and Jasak, H. (2016).
37  A computational method for sharp interface advection
38  Royal Society Open Science, 3
39  doi 10.1098/rsos.160405
40 
41  \endverbatim
42 
43  Original code supplied by Johan Roenby, DHI (2016)
44 
45 SourceFiles
46  cutFaceAdvect.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef cutFaceAdvect_H
51 #define cutFaceAdvect_H
52 
53 #include "cutFace.H"
54 #include "fvMesh.H"
55 #include "surfaceFields.H"
56 #include "volFields.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 
63 /*---------------------------------------------------------------------------*\
64  Class cutFaceAdvect Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class cutFaceAdvect
68 :
69  public cutFace
70 {
71  // Private Data
72 
77 
78  //- Mesh whose cells and faces to cut at their intersection
79  //- with an isoface
80  const fvMesh& mesh_;
81 
82  //- Reference to VoF field
83  const volScalarField& alpha1_;
84 
85  //- Storage for centre of subface
86  point subFaceCentre_;
87 
88  //- Storage for area vector of subface
89  vector subFaceArea_;
90 
91  //- Storage for subFacePoints
92  DynamicList<point> subFacePoints_;
93 
94  //- Storage for subFacePoints
95  DynamicList<point> surfacePoints_;
96 
97  //- Storage for pointStatus_ cuts the cell at 0
98  DynamicList<scalar> pointStatus_;
99 
100  //- Storage of the edge weight
101  DynamicList<scalar> weight_;
102 
103  //- Storage for arrival Times
104  DynamicList<scalar> pTimes_;
105 
106  //- A face status label taking one of the values:
107  //
108  // -1: face is fully below the isosurface
109  // 0: face is cut, i.e. has values larger and smaller than isoValue_
110  // +1: face is fully above the isosurface
111  label faceStatus_;
112 
113 
114  // Private Member Functions
115 
116  //- Write faces to vtk files
117  void isoFacesToFile
118  (
119  const DynamicList<List<point>>& isoFacePts,
120  const word& fileName,
121  const word& fileDir
122  ) const;
123 
124 
125  public:
126 
127  // Constructors
128 
129  //- Construct from fvMesh and a scalarField
130  // Length of scalarField should equal number of mesh points
132 
133 
134  // Member Functions
135 
136  //- Calculates cut centre and cut area (plicReconstruction)
137  // \return face status
138  label calcSubFace
139  (
140  const label faceI,
141  const vector& normal,
142  const vector& base
143  );
144 
145  //- Calculates cut centre and cut area (iso reconstruction)
146  // \return face status
147  label calcSubFace
148  (
149  const face& f,
150  const pointField& points,
151  const scalarField& val,
152  const scalar cutValue
153  );
154 
155  //- Calculates cut centre and cut area (iso reconstruction)
156  // \return face status
157  label calcSubFace
158  (
159  const label faceI,
160  const label sign,
161  const scalar cutValue
162  );
163 
164  //- Calculate time integrated flux for a face
166  (
167  const label faceI,
168  const vector& x0,
169  const vector& n0,
170  const scalar Un0,
171  const scalar dt,
172  const scalar phi,
173  const scalar magSf
174  );
175 
176  //- Calculate time integrated flux for a face
178  (
179  const label faceI,
180  const scalarField& times,
181  const scalar Un0,
182  const scalar dt,
183  const scalar phi,
184  const scalar magSf
185  );
186 
187  //- Calculate time integrated area for a face
188  scalar timeIntegratedArea
189  (
190  const label faceI,
191  const scalar dt,
192  const scalar magSf,
193  const scalar Un0
194  );
195 
196  //- Calculate time integrated area for a face
197  scalar timeIntegratedArea
198  (
199  const pointField& fPts,
200  const scalarField& pTimes,
201  const scalar dt,
202  const scalar magSf,
203  const scalar Un0
204  );
205 
206  //- For face with vertices p calculate its area and integrated area
207  // between isocutting lines with isovalues f0 and f1 given vertex
208  // values f
209  void quadAreaCoeffs
210  (
211  const DynamicPointList& pf0,
212  const DynamicPointList& pf1,
213  scalar& quadArea,
214  scalar& intQuadArea
215  ) const;
216 
217  //- Get cutPoints of face
218  void cutPoints
219  (
220  const label faceI,
221  const scalar f0,
223  );
224 
225  //- Get cutPoints of face
226  void cutPoints
227  (
228  const pointField& pts,
229  const scalarField& f,
230  const scalar f0,
232  );
233 
234  //- Returns centre of cutted face
235  const point& subFaceCentre() const noexcept
236  {
237  return subFaceCentre_;
238  }
239 
240  //- Returns area vector of cutted face
241  const vector& subFaceArea() const noexcept
242  {
243  return subFaceArea_;
244  }
245 
246  //- Returns the cut edge of the cutted face
247  const DynamicList<point>& subFacePoints() const noexcept
248  {
249  return subFacePoints_;
250  }
251 
252  //- Returns point of the face in sorted of cutted face
253  const DynamicList<point>& surfacePoints() const noexcept
254  {
255  return surfacePoints_;
256  }
257 
258 
259  //- Resets internal variables
260  void clearStorage();
261 };
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 } // End namespace Foam
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 #endif
270 
271 // ************************************************************************* //
volFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::cutFaceAdvect::timeIntegratedFaceFlux
scalar timeIntegratedFaceFlux(const label faceI, const vector &x0, const vector &n0, const scalar Un0, const scalar dt, const scalar phi, const scalar magSf)
Calculate time integrated flux for a face.
Definition: cutFaceAdvect.C:188
Foam::cutFaceAdvect::clearStorage
void clearStorage()
Resets internal variables.
Definition: cutFaceAdvect.C:986
Foam::DynamicList< label >
Foam::cutFaceAdvect::quadAreaCoeffs
void quadAreaCoeffs(const DynamicPointList &pf0, const DynamicPointList &pf1, scalar &quadArea, scalar &intQuadArea) const
For face with vertices p calculate its area and integrated area.
Definition: cutFaceAdvect.C:780
Foam::cutFaceAdvect::surfacePoints
const DynamicList< point > & surfacePoints() const noexcept
Returns point of the face in sorted of cutted face.
Definition: cutFaceAdvect.H:252
Foam::cutFaceAdvect::subFacePoints
const DynamicList< point > & subFacePoints() const noexcept
Returns the cut edge of the cutted face.
Definition: cutFaceAdvect.H:246
Foam::cutFaceAdvect::subFaceArea
const vector & subFaceArea() const noexcept
Returns area vector of cutted face.
Definition: cutFaceAdvect.H:240
surfaceFields.H
Foam::surfaceFields.
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::cutFaceAdvect::cutFaceAdvect
cutFaceAdvect(const fvMesh &mesh, const volScalarField &alpha1)
Construct from fvMesh and a scalarField.
Definition: cutFaceAdvect.C:36
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::cutFaceAdvect::cutPoints
void cutPoints(const label faceI, const scalar f0, DynamicList< point > &cutPoints)
Get cutPoints of face.
Definition: cutFaceAdvect.C:882
Foam::cutFaceAdvect::calcSubFace
label calcSubFace(const label faceI, const vector &normal, const vector &base)
Calculates cut centre and cut area (plicReconstruction)
Definition: cutFaceAdvect.C:60
Foam::Field< vector >
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
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::cutFaceAdvect::timeIntegratedArea
scalar timeIntegratedArea(const label faceI, const scalar dt, const scalar magSf, const scalar Un0)
Calculate time integrated area for a face.
Definition: cutFaceAdvect.C:640
Foam::cutFaceAdvect::subFaceCentre
const point & subFaceCentre() const noexcept
Returns centre of cutted face.
Definition: cutFaceAdvect.H:234
Foam::cutFace
Base class for cutting a face, faceI, of an fvMesh, mesh_, at its intersections.
Definition: cutFace.H:59
Foam::cutFaceAdvect
Calculates the face fluxes.
Definition: cutFaceAdvect.H:66
cutFace.H
f
labelList f(nPoints)
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::GeometricField< scalar, fvPatchField, volMesh >