Loading...
Searching...
No Matches
gradAlpha.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 "gradAlpha.H"
29#include "fvc.H"
30#include "leastSquareGrad.H"
31#include "zoneDistribute.H"
33#include "profiling.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace reconstruction
40{
43}
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::reconstruction::gradAlpha::gradSurf(const volScalarField& phi)
50{
51 addProfilingInFunction(geometricVoF);
52 leastSquareGrad<scalar> lsGrad("polyDegree1",mesh_.geometricD());
53
54 zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
55 exchangeFields.setUpCommforZone(interfaceCell_, true);
56
57 Map<vector> mapCC
58 (
59 exchangeFields.getDatafromOtherProc(interfaceCell_, mesh_.C())
60 );
61 Map<scalar> mapPhi
62 (
64 );
65
66 DynamicField<vector> cellCentre(100);
67 DynamicField<scalar> phiValues(100);
68
69 const labelListList& stencil = exchangeFields.getStencil();
70
72 {
73 const label celli = interfaceLabels_[i];
74
75 cellCentre.clear();
76 phiValues.clear();
77
78 for (const label gblIdx : stencil[celli])
79 {
80 cellCentre.append
81 (
82 exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
83 );
84 phiValues.append
85 (
86 exchangeFields.getValue(phi, mapPhi, gblIdx)
87 );
88 }
89
90 cellCentre -= mesh_.C()[celli];
91 interfaceNormal_[i] = lsGrad.grad(cellCentre, phiValues);
92 }
93}
94
95
96// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97
98Foam::reconstruction::gradAlpha::gradAlpha
99(
101 const surfaceScalarField& phi,
102 const volVectorField& U,
103 const dictionary& dict
104)
105:
107 (
108 typeName,
109 alpha1,
110 phi,
111 U,
112 dict
113 ),
114 mesh_(alpha1.mesh()),
115 interfaceNormal_(fvc::grad(alpha1)),
116 isoFaceTol_(modelDict().getOrDefault<scalar>("isoFaceTol", 1e-8)),
117 surfCellTol_(modelDict().getOrDefault<scalar>("surfCellTol", 1e-8)),
118 sIterPLIC_(mesh_,surfCellTol_)
120 reconstruct();
121}
122
123
124// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125
127{
128 addProfilingInFunction(geometricVoF);
129 const bool uptodate = alreadyReconstructed(forceUpdate);
130
131 if (uptodate && !forceUpdate)
132 {
133 return;
134 }
135
136 if (mesh_.topoChanging())
137 {
138 // Introduced resizing to cope with changing meshes
139 interfaceCell_.resize_nocopy(mesh_.nCells());
140 }
141 interfaceCell_ = false;
142
143 interfaceLabels_.clear();
144
145 forAll(alpha1_, celli)
146 {
147 if (sIterPLIC_.isASurfaceCell(alpha1_[celli]))
148 {
149 interfaceCell_[celli] = true; // is set to false earlier
150 interfaceLabels_.append(celli);
151 }
152 }
153 interfaceNormal_.resize(interfaceLabels_.size());
154 centre_ = Zero;
155 normal_ = Zero;
156
157 gradSurf(alpha1_);
158
159 forAll(interfaceLabels_, i)
160 {
161 const label celli = interfaceLabels_[i];
162 if (mag(interfaceNormal_[i]) == 0)
163 {
164 continue;
165 }
166
167 sIterPLIC_.vofCutCell
168 (
169 celli,
170 alpha1_[celli],
171 isoFaceTol_,
172 100,
173 interfaceNormal_[i]
174 );
175
176 if (sIterPLIC_.cellStatus() == 0)
177 {
178 normal_[celli] = sIterPLIC_.surfaceArea();
179 centre_[celli] = sIterPLIC_.surfaceCentre();
180 if (mag(normal_[celli]) == 0)
181 {
182 normal_[celli] = Zero;
183 centre_[celli] = Zero;
184 }
185 }
186 else
187 {
188 normal_[celli] = Zero;
189 centre_[celli] = Zero;
190 }
191 }
192}
193
194
196{
197 // Without this line, we seem to get a race condition
198 mesh_.C();
199
200 cutCellPLIC cutCell(mesh_);
201
202 forAll(normal_, celli)
203 {
204 if (mag(normal_[celli]) != 0)
205 {
206 vector n = normal_[celli]/mag(normal_[celli]);
207 scalar cutValue = (centre_[celli] - mesh_.C()[celli]) & (n);
208 cutCell.calcSubCell
209 (
210 celli,
211 cutValue,
212 n
213 );
214 alpha1_[celli] = cutCell.VolumeOfFluid();
215 }
216 }
217
218 alpha1_.correctBoundaryConditions();
219 alpha1_.oldTime () = alpha1_;
220 alpha1_.oldTime().correctBoundaryConditions();
221}
222
223
224// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha1
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
Class for cutting a cell, cellI, of an fvMesh, mesh_, at its intersection with an surface defined by ...
Definition cutCellPLIC.H:70
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
Definition cutCell.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const volVectorField & C() const
Return cell centres as volVectorField.
Estimates the gradient with a least square scheme in a cell.
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition polyMesh.C:850
Original code supplied by Henning Scheufler, DLR (2019).
boolList interfaceCell_
Is interface cell?
volVectorField normal_
Interface area normals.
bool alreadyReconstructed(bool forceUpdate=true) const
Is the interface already up-to-date?
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
reconstructionSchemes(const reconstructionSchemes &)=delete
No copy construct.
volScalarField & alpha1_
Reference to the VoF Field.
volVectorField centre_
Interface centres.
dictionary & modelDict() noexcept
Access to the model dictionary.
Reconstructs an interface (centre and normal vector) consisting of planes to match the internal fluid...
Definition gradAlpha.H:61
virtual void reconstruct(bool forceUpdate=true)
Reconstruct interface.
Definition gradAlpha.C:119
virtual void mapAlphaField() const
map VoF Field in case of refinement
Definition gradAlpha.C:188
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Class for parallel communication in a narrow band. It either provides a Map with the neighbouring val...
void setUpCommforZone(const boolList &zone, bool updateStencil=true)
Update stencil with boolList the size has to match mesh nCells.
Type getValue(const VolumeField< Type > &phi, const Map< Type > &valuesFromOtherProc, const label gblIdx) const
Gives patchNumber and patchFaceNumber for a given Geometric volume field.
static zoneDistribute & New(const fvMesh &)
Selector.
const labelListList & getStencil() noexcept
Stencil reference.
Map< Type > getDatafromOtherProc(const boolList &zone, const VolumeField< Type > &phi)
Returns stencil and provides a Map with globalNumbering.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
Namespace of functions to calculate explicit derivatives.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
#define addProfilingInFunction(Name)
Define profiling trigger with specified name and description corresponding to the compiler-defined fu...
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299