Explicit.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) 2020 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 "Explicit.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  const dictionary& dict,
37  CloudType& owner
38 )
39 :
40  PackingModel<CloudType>(dict, owner, typeName),
41  stressAverage_(nullptr),
42  correctionLimiting_
43  (
45  (
46  this->coeffDict().subDict(CorrectionLimitingMethod::typeName)
47  )
48  )
49 {}
50 
51 
52 template<class CloudType>
54 (
55  const Explicit<CloudType>& cm
56 )
57 :
59  stressAverage_(cm.stressAverage_->clone()),
60  correctionLimiting_
61  (
62  cm.correctionLimiting_->clone()
63  )
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
69 template<class CloudType>
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
76 template<class CloudType>
78 {
80 
81  if (store)
82  {
83  const fvMesh& mesh = this->owner().mesh();
84  const word& cloudName = this->owner().name();
85 
86  const AveragingMethod<scalar>& volumeAverage =
87  mesh.lookupObject<AveragingMethod<scalar>>
88  (
89  cloudName + ":volumeAverage"
90  );
91  const AveragingMethod<scalar>& rhoAverage =
92  mesh.lookupObject<AveragingMethod<scalar>>
93  (
94  cloudName + ":rhoAverage"
95  );
96  const AveragingMethod<vector>& uAverage =
97  mesh.lookupObject<AveragingMethod<vector>>
98  (
99  cloudName + ":uAverage"
100  );
101  const AveragingMethod<scalar>& uSqrAverage =
102  mesh.lookupObject<AveragingMethod<scalar>>
103  (
104  cloudName + ":uSqrAverage"
105  );
106 
107  volumeAverage_ = &volumeAverage;
108  uAverage_ = &uAverage;
109 
110  stressAverage_.reset
111  (
113  (
114  IOobject
115  (
116  cloudName + ":stressAverage",
117  this->owner().db().time().timeName(),
118  mesh
119  ),
120  this->owner().solution().dict(),
121  mesh
122  ).ptr()
123  );
124 
125  stressAverage_() =
126  this->particleStressModel_->tau
127  (
128  *volumeAverage_,
129  rhoAverage,
130  uSqrAverage
131  )();
132  }
133  else
134  {
135  volumeAverage_ = nullptr;
136  uAverage_ = nullptr;
137  stressAverage_.clear();
138  }
139 }
140 
141 
142 template<class CloudType>
144 (
145  typename CloudType::parcelType& p,
146  const scalar deltaT
147 ) const
148 {
149  const tetIndices tetIs(p.currentTetIndices());
150 
151  // interpolated quantities
152  const scalar alpha =
153  this->volumeAverage_->interpolate(p.coordinates(), tetIs);
154 
155  const vector alphaGrad =
156  this->volumeAverage_->interpolateGrad(p.coordinates(), tetIs);
157 
158  const vector uMean =
159  this->uAverage_->interpolate(p.coordinates(), tetIs);
160 
161  // stress gradient
162  const vector tauGrad =
163  stressAverage_->interpolateGrad(p.coordinates(), tetIs);
164 
165  // parcel relative velocity
166  const vector uRelative = p.U() - uMean;
167 
168  // correction velocity
169  vector dU = Zero;
170 
172  //const scalar Re = p.Re(td);
173  //const typename CloudType::forceType& forces = this->owner().forces();
174  //const forceSuSp F =
175  // forces.calcCoupled(p, td, deltaT, p.mass(), Re, td.muc())
176  // + forces.calcNonCoupled(p, td, deltaT, p.mass(), Re, td.muc());
177 
178  // correction velocity
179  if ((uRelative & alphaGrad) > 0)
180  {
181  dU = - deltaT*tauGrad/(p.rho()*(alpha + SMALL)/* + deltaT*F.Sp()*/);
182  }
183 
184  // apply the velocity limiters
185  return
186  correctionLimiting_->limitedVelocity
187  (
188  p.U(),
189  dU,
190  uMean
191  );
192 }
193 
194 
195 // ************************************************************************* //
Foam::AveragingMethod< scalar >
Foam::PackingModels::Explicit::cacheFields
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Explicit.C:77
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:55
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::PackingModels::Explicit
Explicit model for applying an inter-particle stress to the particles.
Definition: Explicit.H:70
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::PackingModels::Explicit::~Explicit
virtual ~Explicit()
Destructor.
Definition: Explicit.C:70
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::PackingModel
Base class for packing models.
Definition: KinematicCloud.H:98
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Explicit.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::PackingModels::Explicit::velocityCorrection
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition: Explicit.C:144
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::PackingModels::Explicit::Explicit
Explicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition: Explicit.C:35
Foam::Vector< scalar >
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220