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-2024 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 "fvmD2dt2.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().alpha();
68 }
69 
70 
71 template<class CloudType>
73 (
74  const Implicit<CloudType>& cm
75 )
76 :
78  alpha_(cm.alpha_),
79  phiCorrect_
80  (
81  cm.phiCorrect_.valid()
82  ? cm.phiCorrect_().clone()
83  : tmp<surfaceScalarField>(nullptr)
84  ),
85  uCorrect_
86  (
87  cm.uCorrect_.valid()
88  ? cm.uCorrect_().clone()
89  : tmp<volVectorField>(nullptr)
90  ),
91  applyLimiting_(cm.applyLimiting_),
92  applyGravity_(cm.applyGravity_),
93  alphaMin_(cm.alphaMin_),
94  rhoMin_(cm.rhoMin_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
99 
100 template<class CloudType>
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
107 template<class CloudType>
109 {
111 
112  if (store)
113  {
114  const fvMesh& mesh = this->owner().mesh();
115  const dimensionedScalar deltaT = this->owner().db().time().deltaT();
116  const word& cloudName = this->owner().name();
117 
118  const dimensionedVector& g = this->owner().g();
119  const volScalarField& rhoc = this->owner().rho();
120 
121  const AveragingMethod<scalar>& rhoAverage =
123  (
124  cloudName + ":rhoAverage"
125  );
126  const AveragingMethod<vector>& uAverage =
128  (
129  cloudName + ":uAverage"
130  );
131  const AveragingMethod<scalar>& uSqrAverage =
133  (
134  cloudName + ":uSqrAverage"
135  );
136 
137  mesh.schemes().setFluxRequired(alpha_.name());
138 
139 
140  // Property fields
141  // ~~~~~~~~~~~~~~~
142 
143  // volume fraction field
144  alpha_ = max(this->owner().alpha(), alphaMin_);
145  alpha_.correctBoundaryConditions();
146 
147  // average density
149  (
150  IOobject
151  (
152  cloudName + ":rho",
153  this->owner().db().time().name(),
154  mesh,
157  ),
158  mesh,
161  );
162  rho.primitiveFieldRef() = max(rhoAverage.primitiveField(), rhoMin_);
163  rho.correctBoundaryConditions();
164 
165 
166  // Stress field
167  // ~~~~~~~~~~~~
168 
169  // stress derivative wrt volume fraction
170  volScalarField tauPrime
171  (
172  IOobject
173  (
174  cloudName + ":tauPrime",
175  this->owner().db().time().name(),
176  mesh,
179  ),
180  mesh,
183  );
184 
185  tauPrime.primitiveFieldRef() =
186  this->particleStressModel_->dTaudTheta
187  (
188  alpha_.primitiveField(),
189  rho.primitiveField(),
190  uSqrAverage.primitiveField()
191  )();
192 
193  tauPrime.correctBoundaryConditions();
194 
195 
196  // Gravity flux
197  // ~~~~~~~~~~~~
198 
199  tmp<surfaceScalarField> phiGByA;
200 
201  if (applyGravity_)
202  (
203  phiGByA = surfaceScalarField::New
204  (
205  "phiGByA",
206  (g & mesh.Sf())*fvc::interpolate(1.0 - rhoc/rho)
207  )
208  );
209 
210 
211  // Implicit solution for the volume fraction
212  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
213 
215  tauPrimeByRhoAf
216  (
217  "tauPrimeByRhoAf",
218  fvc::interpolate(tauPrime/rho)
219  );
220 
221  fvScalarMatrix alphaEqn
222  (
223  fvm::d2dt2(alpha_)
224  - fvm::laplacian(tauPrimeByRhoAf, alpha_)
225  );
226 
227  if (applyGravity_)
228  {
229  alphaEqn += fvm::div(phiGByA(), alpha_);
230  }
231 
232  alphaEqn.solve();
233 
234 
235  // Generate correction fields
236  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
237 
238  // correction volumetric flux
239  phiCorrect_ = surfaceScalarField::New
240  (
241  cloudName + ":phiCorrect",
242  deltaT*alphaEqn.flux()/fvc::interpolate(alpha_)
243  );
244 
245  // limit the correction flux
246  if (applyLimiting_)
247  {
249  (
250  IOobject
251  (
252  cloudName + ":U",
253  this->owner().db().time().name(),
254  mesh,
257  ),
258  mesh,
261  );
262  U.primitiveFieldRef() = uAverage.primitiveField();
263  U.correctBoundaryConditions();
264 
266  (
267  cloudName + ":phi",
268  linearInterpolate(U) & mesh.Sf()
269  );
270 
271  if (applyGravity_)
272  {
273  phiCorrect_.ref() -= deltaT*phiGByA();
274  }
275 
276  forAll(phiCorrect_(), facei)
277  {
278  // Current and correction fluxes
279  const scalar phiCurr = phi[facei];
280  scalar& phiCorr = phiCorrect_.ref()[facei];
281 
282  // Don't limit if the correction is in the opposite direction to
283  // the flux. We need all the help we can get in this state.
284  if (phiCurr*phiCorr < 0)
285  {}
286 
287  // If the correction and the flux are in the same direction then
288  // don't apply any more correction than is already present in
289  // the flux.
290  else if (phiCorr > 0)
291  {
292  phiCorr = max(phiCorr - phiCurr, 0);
293  }
294  else
295  {
296  phiCorr = min(phiCorr - phiCurr, 0);
297  }
298  }
299 
300  if (applyGravity_)
301  {
302  phiCorrect_.ref() += deltaT*phiGByA();
303  }
304  }
305 
306  // correction velocity
307  uCorrect_ = volVectorField::New
308  (
309  cloudName + ":uCorrect",
310  fvc::reconstruct(phiCorrect_())
311  );
312  uCorrect_->correctBoundaryConditions();
313 
314  // Info<< endl;
315  // Info<< " alpha: " << alpha_.primitiveField() << endl;
316  // Info<< "phiCorrect: " << phiCorrect_->primitiveField() << endl;
317  // Info<< " uCorrect: " << uCorrect_->primitiveField() << endl;
318  // Info<< endl;
319  }
320  else
321  {
322  phiCorrect_.clear();
323  uCorrect_.clear();
324  }
325 }
326 
327 
328 template<class CloudType>
330 (
331  typename CloudType::parcelType& p,
332  const scalar deltaT
333 ) const
334 {
335  const fvMesh& mesh = this->owner().mesh();
336 
337  // containing tetrahedron and parcel coordinates within
338  const label celli = p.cell();
339  const label facei = p.tetFace();
340 
341  // cell velocity
342  const vector U = uCorrect_()[celli];
343 
344  // face geometry
345  vector nHat = mesh.faces()[facei].area(mesh.points());
346  const scalar nMag = mag(nHat);
347  nHat /= nMag;
348 
349  // get face flux
350  scalar phi;
351  const label patchi = mesh.boundaryMesh().whichPatch(facei);
352  if (patchi == -1)
353  {
354  phi = phiCorrect_()[facei];
355  }
356  else
357  {
358  phi =
359  phiCorrect_().boundaryField()[patchi]
360  [
361  mesh.boundaryMesh()[patchi].whichFace(facei)
362  ];
363  }
364 
365  // interpolant equal to 1 at the cell centre and 0 at the face
366  const scalar t = p.coordinates()[0];
367 
368  // the normal component of the velocity correction is interpolated linearly
369  // the tangential component is equal to that at the cell centre
370  return U + (1.0 - t)*nHat*(phi/nMag - (U & nHat));
371 }
372 
373 
374 // ************************************************************************* //
#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:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
void correctBoundaryConditions()
Correct boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
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:101
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Implicit.C:108
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition: Implicit.C:330
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
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:964
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
const surfaceVectorField & Sf() const
Return cell face area vectors.
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
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:1326
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
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 second-order 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
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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< 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 > > d2dt2(const VolField< Type > &vf)
Definition: fvmD2dt2.C:46
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p
const word cloudName(propsDict.lookup("cloudName"))