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-------------------------------------------------------------------------------
10License
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
35namespace Foam
36{
37namespace reconstruction
38{
41}
42}
43
44
45// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46
48(
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{
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 }
102 ap_ = volPointInterpolation::New(mesh_).interpolate(alpha1_);
103
105
106 interfaceLabels_.clear();
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
const volScalarField & alpha1
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Original code supplied by Henning Scheufler, DLR (2019)
Reconstructs an interface (centre and normal vectors) consisting of isosurfaces to match the internal...
Definition: isoAlpha.H:78
virtual void reconstruct(bool forceUpdate=true)
Reconstructs the interface.
Definition: isoAlpha.C:78
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
label nPoints
Namespace for OpenFOAM.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
#define addProfilingInFunction(name)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333