Implicit.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Implicit.H"
28 #include "fvmDdt.H"
29 #include "fvmDiv.H"
30 #include "fvmLaplacian.H"
31 #include "fvcReconstruct.H"
32 #include "volPointInterpolation.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class CloudType>
40 (
41  const dictionary& dict,
42  CloudType& owner
43 )
44 :
45  PackingModel<CloudType>(dict, owner, typeName),
46  alpha_
47  (
48  IOobject
49  (
50  this->owner().name() + ":alpha",
51  this->owner().db().time().name(),
52  this->owner().mesh(),
53  IOobject::NO_READ,
54  IOobject::NO_WRITE
55  ),
56  this->owner().mesh(),
58  zeroGradientFvPatchScalarField::typeName
59  ),
60  phiCorrect_(nullptr),
61  uCorrect_(nullptr),
62  applyLimiting_(this->coeffDict().lookup("applyLimiting")),
63  applyGravity_(this->coeffDict().lookup("applyGravity")),
64  alphaMin_(this->coeffDict().template lookup<scalar>("alphaMin")),
65  rhoMin_(this->coeffDict().template lookup<scalar>("rhoMin"))
66 {
67  alpha_ = this->owner().theta();
68  alpha_.oldTime();
69 }
70 
71 
72 template<class CloudType>
74 (
75  const Implicit<CloudType>& cm
76 )
77 :
79  alpha_(cm.alpha_),
80  phiCorrect_
81  (
82  cm.phiCorrect_.valid()
83  ? cm.phiCorrect_().clone()
84  : tmp<surfaceScalarField>(nullptr)
85  ),
86  uCorrect_
87  (
88  cm.uCorrect_.valid()
89  ? cm.uCorrect_().clone()
90  : tmp<volVectorField>(nullptr)
91  ),
92  applyLimiting_(cm.applyLimiting_),
93  applyGravity_(cm.applyGravity_),
94  alphaMin_(cm.alphaMin_),
95  rhoMin_(cm.rhoMin_)
96 {
97  alpha_.oldTime();
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102 
103 template<class CloudType>
105 {}
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
110 template<class CloudType>
112 {
114 
115  if (store)
116  {
117  const fvMesh& mesh = this->owner().mesh();
118  const dimensionedScalar deltaT = this->owner().db().time().deltaT();
119  const word& cloudName = this->owner().name();
120 
121  const dimensionedVector& g = this->owner().g();
122  const volScalarField& rhoc = this->owner().rho();
123 
124  const AveragingMethod<scalar>& rhoAverage =
126  (
127  cloudName + ":rhoAverage"
128  );
129  const AveragingMethod<vector>& uAverage =
131  (
132  cloudName + ":uAverage"
133  );
134  const AveragingMethod<scalar>& uSqrAverage =
136  (
137  cloudName + ":uSqrAverage"
138  );
139 
140  mesh.schemes().setFluxRequired(alpha_.name());
141 
142  // Property fields
143  // ~~~~~~~~~~~~~~~
144 
145  // volume fraction field
146  alpha_ = max(this->owner().theta(), alphaMin_);
147  alpha_.correctBoundaryConditions();
148 
149  // average density
151  (
152  IOobject
153  (
154  cloudName + ":rho",
155  this->owner().db().time().name(),
156  mesh,
159  ),
160  mesh,
163  );
164  rho.primitiveFieldRef() = max(rhoAverage.primitiveField(), rhoMin_);
165  rho.correctBoundaryConditions();
166 
167  // Stress field
168  // ~~~~~~~~~~~~
169 
170  // stress derivative wrt volume fraction
171  volScalarField tauPrime
172  (
173  IOobject
174  (
175  cloudName + ":tauPrime",
176  this->owner().db().time().name(),
177  mesh,
180  ),
181  mesh,
184  );
185 
186  tauPrime.primitiveFieldRef() =
187  this->particleStressModel_->dTaudTheta
188  (
189  alpha_.primitiveField(),
190  rho.primitiveField(),
191  uSqrAverage.primitiveField()
192  )();
193 
194  tauPrime.correctBoundaryConditions();
195 
196 
197  // Gravity flux
198  // ~~~~~~~~~~~~
199 
200  tmp<surfaceScalarField> phiGByA;
201 
202  if (applyGravity_)
203  (
204  phiGByA = surfaceScalarField::New
205  (
206  "phiGByA",
207  deltaT*(g & mesh.Sf())*fvc::interpolate(1.0 - rhoc/rho)
208  )
209  );
210 
211 
212  // Implicit solution for the volume fraction
213  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
214 
216  tauPrimeByRhoAf
217  (
218  "tauPrimeByRhoAf",
219  fvc::interpolate(deltaT*tauPrime/rho)
220  );
221 
222  fvScalarMatrix alphaEqn
223  (
224  fvm::ddt(alpha_)
225  - fvc::ddt(alpha_)
226  - fvm::laplacian(tauPrimeByRhoAf, alpha_)
227  );
228 
229  if (applyGravity_)
230  {
231  alphaEqn += fvm::div(phiGByA(), alpha_);
232  }
233 
234  alphaEqn.solve();
235 
236 
237  // Generate correction fields
238  // ~~~~~~~~~~~~~~~~~
239 
240  // correction volumetric flux
241  phiCorrect_ = surfaceScalarField::New
242  (
243  cloudName + ":phiCorrect",
244  alphaEqn.flux()/fvc::interpolate(alpha_)
245  );
246 
247  // limit the correction flux
248  if (applyLimiting_)
249  {
251  (
252  IOobject
253  (
254  cloudName + ":U",
255  this->owner().db().time().name(),
256  mesh,
259  ),
260  mesh,
263  );
264  U.primitiveFieldRef() = uAverage.primitiveField();
265  U.correctBoundaryConditions();
266 
268  (
269  cloudName + ":phi",
270  linearInterpolate(U) & mesh.Sf()
271  );
272 
273  if (applyGravity_)
274  {
275  phiCorrect_.ref() -= phiGByA();
276  }
277 
278  forAll(phiCorrect_(), facei)
279  {
280  // Current and correction fluxes
281  const scalar phiCurr = phi[facei];
282  scalar& phiCorr = phiCorrect_.ref()[facei];
283 
284  // Don't limit if the correction is in the opposite direction to
285  // the flux. We need all the help we can get in this state.
286  if (phiCurr*phiCorr < 0)
287  {}
288 
289  // If the correction and the flux are in the same direction then
290  // don't apply any more correction than is already present in
291  // the flux.
292  else if (phiCorr > 0)
293  {
294  phiCorr = max(phiCorr - phiCurr, 0);
295  }
296  else
297  {
298  phiCorr = min(phiCorr - phiCurr, 0);
299  }
300  }
301 
302  if (applyGravity_)
303  {
304  phiCorrect_.ref() += phiGByA();
305  }
306  }
307 
308  // correction velocity
309  uCorrect_ = volVectorField::New
310  (
311  cloudName + ":uCorrect",
312  fvc::reconstruct(phiCorrect_())
313  );
314  uCorrect_->correctBoundaryConditions();
315 
316  // Info<< endl;
317  // Info<< " alpha: " << alpha_.primitiveField() << endl;
318  // Info<< "phiCorrect: " << phiCorrect_->primitiveField() << endl;
319  // Info<< " uCorrect: " << uCorrect_->primitiveField() << endl;
320  // Info<< endl;
321  }
322  else
323  {
324  alpha_.oldTime();
325  phiCorrect_.clear();
326  uCorrect_.clear();
327  }
328 }
329 
330 
331 template<class CloudType>
333 (
334  typename CloudType::parcelType& p,
335  const scalar deltaT
336 ) const
337 {
338  const fvMesh& mesh = this->owner().mesh();
339 
340  // containing tetrahedron and parcel coordinates within
341  const label celli = p.cell();
342  const label facei = p.tetFace();
343 
344  // cell velocity
345  const vector U = uCorrect_()[celli];
346 
347  // face geometry
348  vector nHat = mesh.faces()[facei].area(mesh.points());
349  const scalar nMag = mag(nHat);
350  nHat /= nMag;
351 
352  // get face flux
353  scalar phi;
354  const label patchi = mesh.boundaryMesh().whichPatch(facei);
355  if (patchi == -1)
356  {
357  phi = phiCorrect_()[facei];
358  }
359  else
360  {
361  phi =
362  phiCorrect_().boundaryField()[patchi]
363  [
364  mesh.boundaryMesh()[patchi].whichFace(facei)
365  ];
366  }
367 
368  // interpolant equal to 1 at the cell centre and 0 at the face
369  const scalar t = p.coordinates()[0];
370 
371  // the normal component of the velocity correction is interpolated linearly
372  // the tangential component is equal to that at the cell centre
373  return U + (1.0 - t)*nHat*(phi/nMag - (U & nHat));
374 }
375 
376 
377 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
Generic GeometricField class.
Internal & ref()
Return a reference to the dimensioned internal field.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
void correctBoundaryConditions()
Correct boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Base class for packing models.
Definition: PackingModel.H:65
Implicit model for applying an inter-particle stress to the particles.
Definition: Implicit.H:59
Implicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition: Implicit.C:40
virtual ~Implicit()
Destructor.
Definition: Implicit.C:104
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Implicit.C:111
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition: Implicit.C:333
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:963
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
const surfaceVectorField & Sf() const
Return cell face area vectors.
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1374
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1361
virtual void cacheFields(const bool store)
Cache dependent sub-model fields.
Definition: subModelBase.C:154
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
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.
label patchi
U
Definition: pEqn.H:72
bool valid(const PtrList< ModelType > &l)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimPressure
const dimensionSet dimless
tmp< SurfaceField< Type > > linearInterpolate(const VolField< Type > &vf)
Definition: linear.H:108
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimDensity
dimensioned< scalar > mag(const dimensioned< Type > &)
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimVelocity
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dictionary dict
volScalarField & p
const word cloudName(propsDict.lookup("cloudName"))