Loading...
Searching...
No Matches
reconstructionSchemes.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
29#include "cutCellPLIC.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
37}
38
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
43{
44 const Time& runTime = alpha1_.mesh().time();
45
46 label& curTimeIndex = timeIndexAndIter_.first();
47 label& curIter = timeIndexAndIter_.second();
48
49 // Reset timeIndex and curIter
50 if (curTimeIndex < runTime.timeIndex())
51 {
52 curTimeIndex = runTime.timeIndex();
53 curIter = 0;
54 return false;
55 }
56
57 if (forceUpdate)
58 {
59 curIter = 0;
60 return false;
61 }
62
63 // Always reconstruct when subcycling
64 if (runTime.subCycling() != 0)
65 {
66 return false;
67 }
68
69 ++curIter;
70 if (curIter > 1)
71 {
72 return true;
73 }
75 return false;
76}
77
78
79// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80
82(
83 const word& type,
86 const volVectorField& U,
87 const dictionary& dict
88)
89:
91 (
93 (
94 "reconstructionScheme",
95 alpha1.time().constant(),
96 alpha1.db(),
97 IOobject::NO_READ,
98 IOobject::NO_WRITE
99 )
100 ),
101 reconstructionSchemesCoeffs_(dict),
102 alpha1_(alpha1),
103 phi_(phi),
104 U_(U),
105 normal_
106 (
108 (
109 IOobject::groupName("interfaceNormal", alpha1.group()),
110 alpha1_.mesh().time().timeName(),
111 alpha1_.mesh(),
112 IOobject::NO_READ,
113 dict.getOrDefault("writeFields",false)
114 ? IOobject::AUTO_WRITE
115 : IOobject::NO_WRITE
116 ),
117 alpha1_.mesh(),
119 ),
120 centre_
121 (
123 (
124 IOobject::groupName("interfaceCentre", alpha1.group()),
125 alpha1_.mesh().time().timeName(),
126 alpha1_.mesh(),
127 IOobject::NO_READ,
128 dict.getOrDefault("writeFields",false)
129 ? IOobject::AUTO_WRITE
130 : IOobject::NO_WRITE
131 ),
132 alpha1_.mesh(),
134 ),
135 interfaceCell_(alpha1_.mesh().nCells(), false),
138{}
139
140
141// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142
144{
145 reconstruct(false);
146 const fvMesh& mesh = centre_.mesh();
147
148 cutCellPLIC cellCut(mesh);
149
150 DynamicList<point> dynPts;
151 DynamicList<face> dynFaces(mesh.nCells()*0.1);
152 bitSet interfaceCellAddressing(mesh.nCells());
153
154 forAll(interfaceCell_, celli)
155 {
156 if (interfaceCell_[celli])
157 {
158 if (mag(normal_[celli]) != 0)
159 {
160 interfaceCellAddressing.set(celli);
161 vector n = -normal_[celli]/mag(normal_[celli]);
162
163 scalar cutVal = (centre_[celli] - mesh.C()[celli]) & n;
164 cellCut.calcSubCell(celli, cutVal, n);
165
166 // cellCut.facePoints() are ordered and not connected
167 // to the other face
168 // append face with the starting label: dynPts.size()
169 dynFaces.append
170 (
171 face(identity(cellCut.facePoints().size(), dynPts.size()))
172 );
173 dynPts.append(cellCut.facePoints());
174 }
175 }
176 }
177
178
179 labelList meshCells(interfaceCellAddressing.sortedToc());
180
181 // Transfer to mesh storage
182 pointField pts(std::move(dynPts));
183 faceList faces(std::move(dynFaces));
184
185 return interface(std::move(pts), std::move(faces), std::move(meshCells));
186}
187
188
189// ************************************************************************* //
label n
const volScalarField & alpha1
const Mesh & mesh() const noexcept
Return const reference to mesh.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
static word group(const word &name)
Return group (extension part of name).
Definition IOobject.C:300
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
labelList sortedToc() const
The indices of the on bits as a sorted labelList.
Definition bitSetI.H:441
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
Class for cutting a cell, cellI, of an fvMesh, mesh_, at its intersection with an surface defined by ...
Definition cutCellPLIC.H:70
label calcSubCell(const label celli, const scalar cutValue, const vector &normal)
Sets internal values and returns face status.
Definition cutCellPLIC.C:53
const DynamicList< point > & facePoints()
Returns the points of the cutting PLICface.
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...
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Original code supplied by Henning Scheufler, DLR (2019).
boolList interfaceCell_
Is interface cell?
volVectorField normal_
Interface area normals.
Pair< label > timeIndexAndIter_
Store timeindex/iteration to avoid multiple reconstruction.
const surfaceScalarField & phi_
Reference to the face fluxes.
const volVectorField & U_
Reference to the velocity field.
bool alreadyReconstructed(bool forceUpdate=true) const
Is the interface already up-to-date?
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
interface surface()
Generated interface surface points/faces.
reconstructionSchemes(const reconstructionSchemes &)=delete
No copy construct.
volScalarField & alpha1_
Reference to the VoF Field.
volVectorField centre_
Interface centres.
virtual void reconstruct(bool forceUpdate=true)=0
Reconstruct the interface.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
engineTime & runTime
word timeName
Definition getTimeIndex.H:3
Different types of constants.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
List< face > faceList
List of faces.
Definition faceListFwd.H:41
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
const pointField & pts
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299