Loading...
Searching...
No Matches
Implicit.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) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "Implicit.H"
31#include "fvmDdt.H"
32#include "fvmDiv.H"
33#include "fvmLaplacian.H"
37
38// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39
40template<class CloudType>
42(
43 const dictionary& dict,
45)
46:
48 alpha_
49 (
51 (
52 IOobject::scopedName(this->owner().name(), "alpha"),
53 this->owner().db().time().timeName(),
54 this->owner().mesh(),
55 IOobject::NO_READ,
56 IOobject::NO_WRITE
57 ),
58 this->owner().mesh(),
60 fvPatchFieldBase::zeroGradientType()
61 ),
62 phiCorrect_(nullptr),
63 uCorrect_(nullptr),
64 applyLimiting_(this->coeffDict().lookup("applyLimiting")),
65 applyGravity_(this->coeffDict().lookup("applyGravity")),
66 alphaMin_(this->coeffDict().getScalar("alphaMin")),
67 rhoMin_(this->coeffDict().getScalar("rhoMin"))
69 alpha_ = this->owner().theta();
70 alpha_.oldTime();
71}
72
73
74template<class CloudType>
76(
77 const Implicit<CloudType>& cm
78)
79:
81 alpha_(cm.alpha_),
82 phiCorrect_(cm.phiCorrect_()),
83 uCorrect_(cm.uCorrect_()),
84 applyLimiting_(cm.applyLimiting_),
85 applyGravity_(cm.applyGravity_),
86 alphaMin_(cm.alphaMin_),
87 rhoMin_(cm.rhoMin_)
88{
89 alpha_.oldTime();
90}
91
92
93// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94
95template<class CloudType>
97{}
98
99
100// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101
102template<class CloudType>
104{
106
107 if (store)
108 {
109 const fvMesh& mesh = this->owner().mesh();
110 const dimensionedScalar deltaT = this->owner().db().time().deltaT();
111 const word& cloudName = this->owner().name();
112
113 const dimensionedVector& g = this->owner().g();
114 const volScalarField& rhoc = this->owner().rho();
115
116 const auto& rhoAverage =
117 mesh.lookupObject<AveragingMethod<scalar>>
118 (
119 IOobject::scopedName(cloudName, "rhoAverage")
120 );
121 const auto& uAverage =
122 mesh.lookupObject<AveragingMethod<vector>>
123 (
125 );
126 const auto& uSqrAverage =
127 mesh.lookupObject<AveragingMethod<scalar>>
128 (
129 IOobject::scopedName(cloudName, "uSqrAverage")
130 );
131
132 mesh.setFluxRequired(alpha_.name());
133
134 // Property fields
135 // ~~~~~~~~~~~~~~~
136
137 // volume fraction field
138 alpha_ = max(this->owner().theta(), alphaMin_);
139 alpha_.correctBoundaryConditions();
140
141 // average density
143 (
145 mesh,
148 );
149 auto& rho = trho.ref();
150
151 rho.primitiveFieldRef() = max(rhoAverage.primitiveField(), rhoMin_);
152 rho.correctBoundaryConditions();
153
154 // Stress field
155 // ~~~~~~~~~~~~
156
157 // stress derivative wrt volume fraction
158 auto ttauPrime = volScalarField::New
159 (
160 IOobject::scopedName(cloudName, "tauPrime"),
161 mesh,
164 );
165 auto& tauPrime = ttauPrime.ref();
166
167 tauPrime.primitiveFieldRef() =
168 this->particleStressModel_->dTaudTheta
169 (
170 alpha_.primitiveField(),
171 rho.primitiveField(),
172 uSqrAverage.primitiveField()
173 )();
174
175 tauPrime.correctBoundaryConditions();
176
177
178 // Gravity flux
179 // ~~~~~~~~~~~~
180
182
183 if (applyGravity_)
184 {
186 (
187 "phiGByA",
189 deltaT*(g & mesh.Sf())*fvc::interpolate(1.0 - rhoc/rho)
190 );
191 }
192
193
194 // Implicit solution for the volume fraction
195 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
196
197 surfaceScalarField tauPrimeByRhoAf
198 (
199 "tauPrimeByRhoAf",
200 fvc::interpolate(deltaT*tauPrime/rho)
201 );
202
203 fvScalarMatrix alphaEqn
204 (
205 fvm::ddt(alpha_)
206 - fvc::ddt(alpha_)
207 - fvm::laplacian(tauPrimeByRhoAf, alpha_)
208 );
209
210 if (applyGravity_)
211 {
212 alphaEqn += fvm::div(phiGByA(), alpha_);
213 }
214
215 alphaEqn.solve();
216
217
218 // Generate correction fields
219 // ~~~~~~~~~~~~~~~~~
220
221 // correction volumetric flux
222 phiCorrect_ = surfaceScalarField::New
223 (
224 IOobject::scopedName(cloudName, "phiCorrect"),
226 alphaEqn.flux()/fvc::interpolate(alpha_)
227 );
228
229 // limit the correction flux
230 if (applyLimiting_)
231 {
232 auto tU = volVectorField::New
233 (
236 mesh,
239 );
240 auto& U = tU.ref();
241
242 U.primitiveFieldRef() = uAverage.primitiveField();
243 U.correctBoundaryConditions();
244
246 (
248 linearInterpolate(U) & mesh.Sf()
249 );
250
251 if (applyGravity_)
252 {
253 phiCorrect_.ref() -= phiGByA();
254 }
255
256 forAll(phiCorrect_(), facei)
257 {
258 // Current and correction fluxes
259 const scalar phiCurr = phi[facei];
260 scalar& phiCorr = phiCorrect_.ref()[facei];
261
262 // Don't limit if the correction is in the opposite direction to
263 // the flux. We need all the help we can get in this state.
264 if (phiCurr*phiCorr < 0)
265 {}
266
267 // If the correction and the flux are in the same direction then
268 // don't apply any more correction than is already present in
269 // the flux.
270 else if (phiCorr > 0)
271 {
272 phiCorr = max(phiCorr - phiCurr, 0);
273 }
274 else
275 {
276 phiCorr = min(phiCorr - phiCurr, 0);
277 }
278 }
279
280 if (applyGravity_)
281 {
282 phiCorrect_.ref() += phiGByA();
283 }
284 }
285
286 // correction velocity
287 uCorrect_ = volVectorField::New
288 (
289 IOobject::scopedName(cloudName, "uCorrect"),
291 fvc::reconstruct(phiCorrect_())
292 );
293 uCorrect_->correctBoundaryConditions();
294
295 //Info << endl;
296 //Info << " alpha: " << alpha_.primitiveField() << endl;
297 //Info << "phiCorrect: " << phiCorrect_->primitiveField() << endl;
298 //Info << " uCorrect: " << uCorrect_->primitiveField() << endl;
299 //Info << endl;
300 }
301 else
302 {
303 alpha_.oldTime();
304 phiCorrect_.clear();
305 uCorrect_.clear();
306 }
307}
308
309
310template<class CloudType>
312(
313 typename CloudType::parcelType& p,
314 const scalar deltaT
315) const
316{
317 const fvMesh& mesh = this->owner().mesh();
318
319 // containing tetrahedron and parcel coordinates within
320 const label celli = p.cell();
321 const label facei = p.tetFace();
322
323 // cell velocity
324 const vector U = uCorrect_()[celli];
325
326 // face geometry
327 vector nHat = mesh.faces()[facei].areaNormal(mesh.points());
328 const scalar nMag = mag(nHat);
329 nHat /= nMag;
330
331 // get face flux
332 scalar phi;
333 const label patchi = mesh.boundaryMesh().whichPatch(facei);
334 if (patchi == -1)
335 {
336 phi = phiCorrect_()[facei];
337 }
338 else
339 {
340 phi =
341 phiCorrect_().boundaryField()[patchi]
342 [
343 mesh.boundaryMesh()[patchi].whichFace(facei)
344 ];
345 }
346
347 // interpolant equal to 1 at the cell centre and 0 at the face
348 const scalar t = p.coordinates()[0];
349
350 // the normal component of the velocity correction is interpolated linearly
351 // the tangential component is equal to that at the cell centre
352 return U + (1.0 - t)*nHat*(phi/nMag - (U & nHat));
353}
354
355
356// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
const uniformDimensionedVectorField & g
Base class for lagrangian averaging methods.
virtual tmp< Field< Type > > primitiveField() const =0
Return an internal field of the average.
const CloudType & owner() const
Return const access to the owner cloud.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Base class for packing models.
autoPtr< ParticleStressModel > particleStressModel_
Protected data.
PackingModel(CloudType &owner)
Construct null from owner.
Implicit model for applying an inter-particle stress to the particles.
Definition Implicit.H:60
Implicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition Implicit.C:35
virtual ~Implicit()
Destructor.
Definition Implicit.C:89
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition Implicit.C:96
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition Implicit.C:305
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition fvMatrix.C:1415
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Template invariant parts for fvPatchField.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Lookup type of boundary radiation properties.
Definition lookup.H:60
virtual void cacheFields(const bool store)
Cache dependent sub-model fields.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
Reconstruct volField from a face flux field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
word timeName
Definition getTimeIndex.H:3
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
const dimensionSet dimPressure
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
const dimensionSet dimless
Dimensionless.
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition linear.H:120
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
const dimensionSet dimDensity
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< volScalarField > trho
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299