isoAlpha.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) 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 \*---------------------------------------------------------------------------*/
27 
28 #include "isoAlpha.H"
30 #include "cutCellPLIC.H"
31 #include "profiling.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace reconstruction
38 {
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 Foam::reconstruction::isoAlpha::isoAlpha
48 (
50  const surfaceScalarField& phi,
51  const volVectorField& U,
52  const dictionary& dict
53 )
54 :
56  (
57  typeName,
58  alpha1,
59  phi,
60  U,
61  dict
62  ),
63  mesh_(alpha1.mesh()),
64  // Interpolation data
65  ap_(mesh_.nPoints()),
66 
67  // Tolerances and solution controls
68  isoFaceTol_(modelDict().getOrDefault<scalar>("isoFaceTol", 1e-8)),
69  surfCellTol_(modelDict().getOrDefault<scalar>("surfCellTol", 1e-8)),
70  sIterIso_(mesh_, ap_, surfCellTol_)
71 {
72  reconstruct();
73 }
74 
75 
76 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
77 
79 {
80  addProfilingInFunction(geometricVoF);
81  const bool uptodate = alreadyReconstructed(forceUpdate);
82 
83  if (uptodate && !forceUpdate)
84  {
85  return;
86  }
87 
88  // Interpolating alpha1 cell centre values to mesh points (vertices)
89  if (mesh_.topoChanging())
90  {
91  // Introduced resizing to cope with changing meshes
92  if (ap_.size() != mesh_.nPoints())
93  {
94  ap_.resize(mesh_.nPoints());
95 
96  }
97  if (interfaceCell_.size() != mesh_.nCells())
98  {
99  interfaceCell_.resize(mesh_.nCells());
100  }
101  }
103 
104  DynamicList<List<point>> facePts;
105 
107 
108  forAll(alpha1_,cellI)
109  {
110  if (sIterIso_.isASurfaceCell(alpha1_[cellI]))
111  {
112  interfaceLabels_.append(cellI);
113 
114  sIterIso_.vofCutCell
115  (
116  cellI,
117  alpha1_[cellI],
118  isoFaceTol_,
119  100
120  );
121 
122  if (sIterIso_.cellStatus() == 0)
123  {
124  normal_[cellI] = sIterIso_.surfaceArea();
125  centre_[cellI] = sIterIso_.surfaceCentre();
126  if (mag(normal_[cellI]) != 0)
127  {
128  interfaceCell_[cellI] = true;
129  }
130  else
131  {
132  interfaceCell_[cellI] = false;
133  }
134  }
135  else
136  {
137  normal_[cellI] = Zero;
138  centre_[cellI] = Zero;
139  interfaceCell_[cellI] = false;
140  }
141  }
142  else
143  {
144  normal_[cellI] = Zero;
145  centre_[cellI] = Zero;
146  interfaceCell_[cellI] = false;
147  }
148  }
149 }
150 
151 
152 // ************************************************************************* //
profiling.H
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:56
Foam::reconstruction::defineTypeNameAndDebug
defineTypeNameAndDebug(isoAlpha, 0)
Foam::List::resize
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::New
static const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::reconstructionSchemes::normal_
volVectorField normal_
Interface area normals.
Definition: reconstructionSchemes.H:82
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Foam::reconstructionSchemes::interfaceCell_
boolList interfaceCell_
Is interface cell?
Definition: reconstructionSchemes.H:88
Foam::surfaceIteratorIso::isASurfaceCell
bool isASurfaceCell(const scalar alpha1) const
Determine if a cell is a surface cell.
Definition: surfaceIteratorIso.H:100
Foam::polyMesh::topoChanging
bool topoChanging() const noexcept
Is mesh topology changing.
Definition: polyMesh.H:536
Foam::primitiveMesh::nPoints
label nPoints() const noexcept
Number of mesh points.
Definition: primitiveMeshI.H:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::reconstructionSchemes::centre_
volVectorField centre_
Interface centres.
Definition: reconstructionSchemes.H:85
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::reconstruction::isoAlpha
Reconstructs an interface (centre and normal vectors) consisting of isosurfaces to match the internal...
Definition: isoAlpha.H:75
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::surfaceIteratorIso::surfaceArea
const vector & surfaceArea() const
The area vector of cutting isosurface.
Definition: surfaceIteratorIso.H:138
Foam::reconstructionSchemes::alpha1_
volScalarField & alpha1_
Reference to the VoF Field.
Definition: reconstructionSchemes.H:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
addProfilingInFunction
#define addProfilingInFunction(name)
Definition: profilingTrigger.H:124
Foam::reconstruction::addToRunTimeSelectionTable
addToRunTimeSelectionTable(reconstructionSchemes, isoAlpha, components)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
U
U
Definition: pEqn.H:72
Foam::reconstructionSchemes::alreadyReconstructed
bool alreadyReconstructed(bool forceUpdate=true) const
Is the interface already up-to-date?
Definition: reconstructionSchemes.C:42
Foam::surfaceIteratorIso::cellStatus
label cellStatus()
The cellStatus.
Definition: surfaceIteratorIso.H:156
cutCellPLIC.H
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::surfaceIteratorIso::surfaceCentre
const point & surfaceCentre() const
The centre of cutting isosurface.
Definition: surfaceIteratorIso.H:132
Foam::reconstructionSchemes
Original code supplied by Henning Scheufler, DLR (2019)
Definition: reconstructionSchemes.H:58
Foam::reconstruction::isoAlpha::reconstruct
virtual void reconstruct(bool forceUpdate=true)
Reconstructs the interface.
Definition: isoAlpha.C:78
Foam::GeometricField< scalar, fvPatchField, volMesh >
isoAlpha.H
Foam::DynamicField::append
DynamicField< T, SizeMin > & append(const T &val)
Append an element at the end of the list.
Definition: DynamicFieldI.H:587
Foam::DynamicField::clear
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicFieldI.H:431
Foam::surfaceIteratorIso::vofCutCell
label vofCutCell(const label celli, const scalar alpha1, const scalar tol, const label maxIter)
Definition: surfaceIteratorIso.C:50
Foam::reconstructionSchemes::interfaceLabels_
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
Definition: reconstructionSchemes.H:91