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 -------------------------------------------------------------------------------
11 License
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"
34 #include "fvcReconstruct.H"
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class CloudType>
42 (
43  const dictionary& dict,
44  CloudType& owner
45 )
46 :
47  PackingModel<CloudType>(dict, owner, typeName),
48  alpha_
49  (
50  IOobject
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(),
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"))
68 {
69  alpha_ = this->owner().theta();
70  alpha_.oldTime();
71 }
72 
73 
74 template<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 
95 template<class CloudType>
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 template<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  (
124  IOobject::scopedName(cloudName, "uAverage")
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 
181  tmp<surfaceScalarField> phiGByA;
182 
183  if (applyGravity_)
184  {
185  phiGByA = surfaceScalarField::New
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 
310 template<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 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
Implicit model for applying an inter-particle stress to the particles.
Definition: Implicit.H:57
dictionary dict
virtual tmp< Field< Type > > primitiveField() const =0
Return an internal field of the average.
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Implicit.C:96
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const word zeroGradientType
A zeroGradient patch field type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:120
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Reconstruct volField from a face flux field.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< volScalarField > trho
Base class for packing models.
virtual ~Implicit()
Destructor.
Definition: Implicit.C:89
Calculate the matrix for the laplacian of the field.
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
const dimensionSet dimless
Dimensionless.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
const CloudType & owner() const
Return const access to the owner cloud.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
word timeName
Definition: getTimeIndex.H:3
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
Implicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition: Implicit.C:35
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
const dimensionSet dimPressure
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
const uniformDimensionedVectorField & g
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
virtual void cacheFields(const bool store)
Cache dependent sub-model fields.
Definition: subModelBase.C:155
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.
const dimensionSet dimDensity
Calculate the matrix for the divergence of the given field and flux.
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:37
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition: Implicit.C:305
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Do not request registration (bool: false)
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:61
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity