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 -------------------------------------------------------------------------------
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 "gradAlpha.H"
29 #include "fvc.H"
30 #include "leastSquareGrad.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace reconstruction
38 {
39  defineTypeNameAndDebug(gradAlpha, 0);
40  addToRunTimeSelectionTable(reconstructionSchemes, gradAlpha, components);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::reconstruction::gradAlpha::gradSurf(const volScalarField& phi)
48 {
49  leastSquareGrad<scalar> lsGrad("polyDegree1",mesh_.geometricD());
50 
51  exchangeFields_.setUpCommforZone(interfaceCell_,true);
52 
53  Map<vector> mapCC
54  (
55  exchangeFields_.getDatafromOtherProc(interfaceCell_, mesh_.C())
56  );
57  Map<scalar> mapPhi
58  (
59  exchangeFields_.getDatafromOtherProc(interfaceCell_, phi)
60  );
61 
62  DynamicField<vector> cellCentre(100);
63  DynamicField<scalar> phiValues(100);
64 
65  const labelListList& stencil = exchangeFields_.getStencil();
66 
68  {
69  const label celli = interfaceLabels_[i];
70 
71  cellCentre.clear();
72  phiValues.clear();
73 
74  for (const label gblIdx : stencil[celli])
75  {
76  cellCentre.append
77  (
78  exchangeFields_.getValue(mesh_.C(), mapCC, gblIdx)
79  );
80  phiValues.append
81  (
82  exchangeFields_.getValue(phi, mapPhi, gblIdx)
83  );
84  }
85 
86  cellCentre -= mesh_.C()[celli];
87  interfaceNormal_[i] = lsGrad.grad(cellCentre, phiValues);
88  }
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
94 Foam::reconstruction::gradAlpha::gradAlpha
95 (
97  const surfaceScalarField& phi,
98  const volVectorField& U,
99  const dictionary& dict
100 )
101 :
103  (
104  typeName,
105  alpha1,
106  phi,
107  U,
108  dict
109  ),
110  mesh_(alpha1.mesh()),
111  interfaceNormal_(fvc::grad(alpha1)),
112  isoFaceTol_(modelDict().getOrDefault<scalar>("isoFaceTol", 1e-8)),
113  surfCellTol_(modelDict().getOrDefault<scalar>("surfCellTol", 1e-8)),
114  exchangeFields_(zoneDistribute::New(mesh_)),
115  sIterPLIC_(mesh_,surfCellTol_)
116 {
117  reconstruct();
118 }
119 
120 
121 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
122 
124 {
125  const bool uptodate = alreadyReconstructed(forceUpdate);
126 
127  if (uptodate && !forceUpdate)
128  {
129  return;
130  }
131 
132  if (mesh_.topoChanging())
133  {
134  // Introduced resizing to cope with changing meshes
135  if (interfaceCell_.size() != mesh_.nCells())
136  {
137  interfaceCell_.resize(mesh_.nCells());
138  }
139  }
140  interfaceCell_ = false;
141 
142  interfaceLabels_.clear();
143 
144  forAll(alpha1_, celli)
145  {
146  if (sIterPLIC_.isASurfaceCell(alpha1_[celli]))
147  {
148  interfaceCell_[celli] = true; // is set to false earlier
149  interfaceLabels_.append(celli);
150  }
151  }
152  interfaceNormal_.resize(interfaceLabels_.size());
153  centre_ = dimensionedVector("centre", dimLength, Zero);
154  normal_ = dimensionedVector("normal", dimArea, Zero);
155 
156  gradSurf(alpha1_);
157 
158  forAll(interfaceLabels_, i)
159  {
160  const label celli = interfaceLabels_[i];
161  if (mag(interfaceNormal_[i]) == 0)
162  {
163  continue;
164  }
165 
166  sIterPLIC_.vofCutCell
167  (
168  celli,
169  alpha1_[celli],
170  isoFaceTol_,
171  100,
172  interfaceNormal_[i]
173  );
174 
175  if (sIterPLIC_.cellStatus() == 0)
176  {
177  normal_[celli] = sIterPLIC_.surfaceArea();
178  centre_[celli] = sIterPLIC_.surfaceCentre();
179  if (mag(normal_[celli]) == 0)
180  {
181  normal_[celli] = Zero;
182  centre_[celli] = Zero;
183  }
184  }
185  else
186  {
187  normal_[celli] = Zero;
188  centre_[celli] = Zero;
189  }
190  }
191 }
192 
193 
195 {
196  // Without this line, we seem to get a race condition
197  mesh_.C();
198 
199  cutCellPLIC cutCell(mesh_);
200 
201  forAll(normal_, celli)
202  {
203  if (mag(normal_[celli]) != 0)
204  {
205  vector n = normal_[celli]/mag(normal_[celli]);
206  scalar cutValue = (centre_[celli] - mesh_.C()[celli]) & (n);
207  cutCell.calcSubCell
208  (
209  celli,
210  cutValue,
211  n
212  );
213  alpha1_[celli] = cutCell.VolumeOfFluid();
214  }
215  }
216 
217  alpha1_.correctBoundaryConditions();
218  alpha1_.oldTime () = alpha1_;
219  alpha1_.oldTime().correctBoundaryConditions();
220 }
221 
222 
223 // ************************************************************************* //
Foam::zoneDistribute::getValue
Type getValue(const GeometricField< Type, fvPatchField, volMesh > &phi, const Map< Type > &valuesFromOtherProc, const label gblIdx) const
Definition: zoneDistributeI.H:85
Foam::DynamicField::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicFieldI.H:350
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::zoneDistribute::getDatafromOtherProc
Map< Type > getDatafromOtherProc(const boolList &zone, const GeometricField< Type, fvPatchField, volMesh > &phi)
Returns stencil and provides a Map with globalNumbering.
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
leastSquareGrad.H
Foam::cutCellPLIC
Class for cutting a cell, cellI, of an fvMesh, mesh_, at its intersection with an surface defined by ...
Definition: cutCellPLIC.H:71
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::polyMesh::geometricD
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:858
Foam::reconstructionSchemes::interfaceCell_
boolList interfaceCell_
Is interface cell?
Definition: reconstructionSchemes.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::reconstruction::gradAlpha::mapAlphaField
virtual void mapAlphaField() const
map VoF Field in case of refinement
Definition: gradAlpha.C:194
Foam::reconstruction::gradAlpha::reconstruct
virtual void reconstruct(bool forceUpdate=true)
Reconstruct interface.
Definition: gradAlpha.C:123
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::zoneDistribute::getStencil
const labelListList & getStencil()
Stencil reference.
Definition: zoneDistribute.H:129
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
Foam::cutCell
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
Definition: cutCell.H:59
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:121
gradAlpha.H
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::reconstruction::addToRunTimeSelectionTable
addToRunTimeSelectionTable(reconstructionSchemes, isoAlpha, components)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::zoneDistribute::New
static zoneDistribute & New(const fvMesh &)
Definition: zoneDistribute.C:104
Foam::zoneDistribute::setUpCommforZone
void setUpCommforZone(const boolList &zone, bool updateStencil=true)
Update stencil with boolList the size has to match mesh nCells.
Definition: zoneDistribute.C:131
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
U
U
Definition: pEqn.H:72
Foam::Vector< scalar >
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
fvc.H
Foam::reconstructionSchemes
Original code supplied by Henning Scheufler, DLR (2019)
Definition: reconstructionSchemes.H:58
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::reconstructionSchemes::interfaceLabels_
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
Definition: reconstructionSchemes.H:91