Explicit.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-2022 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 "Explicit.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class CloudType>
32 (
33  const dictionary& dict,
34  CloudType& owner
35 )
36 :
37  PackingModel<CloudType>(dict, owner, typeName),
38  stressAverage_(nullptr),
39  correctionLimiting_
40  (
42  (
43  this->coeffDict().subDict(CorrectionLimitingMethod::typeName)
44  )
45  )
46 {}
47 
48 
49 template<class CloudType>
51 (
52  const Explicit<CloudType>& cm
53 )
54 :
56  stressAverage_(cm.stressAverage_, false),
57  correctionLimiting_(cm.correctionLimiting_, false)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
63 template<class CloudType>
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class CloudType>
72 {
74 
75  if (store)
76  {
77  const fvMesh& mesh = this->owner().mesh();
78  const word& cloudName = this->owner().name();
79 
80  const AveragingMethod<scalar>& volumeAverage =
82  (
83  cloudName + ":volumeAverage"
84  );
85  const AveragingMethod<scalar>& rhoAverage =
87  (
88  cloudName + ":rhoAverage"
89  );
90  const AveragingMethod<vector>& uAverage =
92  (
93  cloudName + ":uAverage"
94  );
95  const AveragingMethod<scalar>& uSqrAverage =
97  (
98  cloudName + ":uSqrAverage"
99  );
100 
101  volumeAverage_ = &volumeAverage;
102  uAverage_ = &uAverage;
103 
104  stressAverage_.reset
105  (
107  (
108  IOobject
109  (
110  cloudName + ":stressAverage",
111  this->owner().db().time().name(),
112  mesh
113  ),
114  this->owner().solution().dict(),
115  mesh
116  ).ptr()
117  );
118 
119  stressAverage_() =
120  this->particleStressModel_->tau
121  (
122  *volumeAverage_,
123  rhoAverage,
124  uSqrAverage
125  )();
126  }
127  else
128  {
129  volumeAverage_ = nullptr;
130  uAverage_ = nullptr;
131  stressAverage_.clear();
132  }
133 }
134 
135 
136 template<class CloudType>
138 (
139  typename CloudType::parcelType& p,
140  const scalar deltaT
141 ) const
142 {
143  const tetIndices tetIs(p.currentTetIndices(this->owner().mesh()));
144 
145  // interpolated quantities
146  const scalar alpha =
147  this->volumeAverage_->interpolate(p.coordinates(), tetIs);
148  const vector alphaGrad =
149  this->volumeAverage_->interpolateGrad(p.coordinates(), tetIs);
150  const vector uMean =
151  this->uAverage_->interpolate(p.coordinates(), tetIs);
152 
153  // stress gradient
154  const vector tauGrad =
155  stressAverage_->interpolateGrad(p.coordinates(), tetIs);
156 
157  // parcel relative velocity
158  const vector uRelative = p.U() - uMean;
159 
160  // correction velocity
161  vector dU = Zero;
162 
164  // const scalar Re = p.Re(td);
165  // const typename CloudType::forceType& forces = this->owner().forces();
166  // const forceSuSp F =
167  // forces.calcCoupled(p, td, deltaT, p.mass(), Re, td.muc())
168  // + forces.calcNonCoupled(p, td, deltaT, p.mass(), Re, td.muc());
169 
170  // correction velocity
171  if ((uRelative & alphaGrad) > 0)
172  {
173  dU = - deltaT*tauGrad/(p.rho()*alpha/* + deltaT*F.Sp()*/);
174  }
175 
176  // apply the velocity limiters
177  return
178  correctionLimiting_->limitedVelocity
179  (
180  p.U(),
181  dU,
182  uMean
183  );
184 }
185 
186 
187 // ************************************************************************* //
Base class for correction limiting methods.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
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
Explicit model for applying an inter-particle stress to the particles.
Definition: Explicit.H:68
virtual ~Explicit()
Destructor.
Definition: Explicit.C:64
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Explicit.C:71
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition: Explicit.C:138
Explicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition: Explicit.C:32
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:198
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
virtual void cacheFields(const bool store)
Cache dependent sub-model fields.
Definition: subModelBase.C:154
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:82
A class for handling words, derived from string.
Definition: word.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static const zero Zero
Definition: zero.H:97
dictionary dict
volScalarField & p
const word cloudName(propsDict.lookup("cloudName"))