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-2020 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_->clone()),
57  correctionLimiting_
58  (
59  cm.correctionLimiting_->clone()
60  )
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
66 template<class CloudType>
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 template<class CloudType>
75 {
77 
78  if (store)
79  {
80  const fvMesh& mesh = this->owner().mesh();
81  const word& cloudName = this->owner().name();
82 
83  const AveragingMethod<scalar>& volumeAverage =
85  (
86  cloudName + ":volumeAverage"
87  );
88  const AveragingMethod<scalar>& rhoAverage =
90  (
91  cloudName + ":rhoAverage"
92  );
93  const AveragingMethod<vector>& uAverage =
95  (
96  cloudName + ":uAverage"
97  );
98  const AveragingMethod<scalar>& uSqrAverage =
100  (
101  cloudName + ":uSqrAverage"
102  );
103 
104  volumeAverage_ = &volumeAverage;
105  uAverage_ = &uAverage;
106 
107  stressAverage_.reset
108  (
110  (
111  IOobject
112  (
113  cloudName + ":stressAverage",
114  this->owner().db().time().timeName(),
115  mesh
116  ),
117  this->owner().solution().dict(),
118  mesh
119  ).ptr()
120  );
121 
122  stressAverage_() =
123  this->particleStressModel_->tau
124  (
125  *volumeAverage_,
126  rhoAverage,
127  uSqrAverage
128  )();
129  }
130  else
131  {
132  volumeAverage_ = nullptr;
133  uAverage_ = nullptr;
134  stressAverage_.clear();
135  }
136 }
137 
138 
139 template<class CloudType>
141 (
142  typename CloudType::parcelType& p,
143  const scalar deltaT
144 ) const
145 {
146  const tetIndices tetIs(p.currentTetIndices());
147 
148  // interpolated quantities
149  const scalar alpha =
150  this->volumeAverage_->interpolate(p.coordinates(), tetIs);
151  const vector alphaGrad =
152  this->volumeAverage_->interpolateGrad(p.coordinates(), tetIs);
153  const vector uMean =
154  this->uAverage_->interpolate(p.coordinates(), tetIs);
155 
156  // stress gradient
157  const vector tauGrad =
158  stressAverage_->interpolateGrad(p.coordinates(), tetIs);
159 
160  // parcel relative velocity
161  const vector uRelative = p.U() - uMean;
162 
163  // correction velocity
164  vector dU = Zero;
165 
167  // const scalar Re = p.Re(td);
168  // const typename CloudType::forceType& forces = this->owner().forces();
169  // const forceSuSp F =
170  // forces.calcCoupled(p, td, deltaT, p.mass(), Re, td.muc())
171  // + forces.calcNonCoupled(p, td, deltaT, p.mass(), Re, td.muc());
172 
173  // correction velocity
174  if ((uRelative & alphaGrad) > 0)
175  {
176  dU = - deltaT*tauGrad/(p.rho()*alpha/* + deltaT*F.Sp()*/);
177  }
178 
179  // apply the velocity limiters
180  return
181  correctionLimiting_->limitedVelocity
182  (
183  p.U(),
184  dU,
185  uMean
186  );
187 }
188 
189 
190 // ************************************************************************* //
dictionary dict
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition: Explicit.C:141
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Base class for packing models.
Definition: MPPICCloud.H:53
fvMesh & mesh
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Explicit model for applying an inter-particle stress to the particles.
Definition: Explicit.H:65
A class for handling words, derived from string.
Definition: word.H:59
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
Explicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition: Explicit.C:32
virtual ~Explicit()
Destructor.
Definition: Explicit.C:67
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:174
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Explicit.C:74
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
const word cloudName(propsDict.lookup("cloudName"))
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75