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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "Explicit.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 (
35  const dictionary& dict,
36  CloudType& owner
37 )
38 :
39  PackingModel<CloudType>(dict, owner, typeName),
40  stressAverage_(nullptr),
41  correctionLimiting_
42  (
44  (
45  this->coeffDict().subDict(CorrectionLimitingMethod::typeName)
46  )
47  )
48 {}
49 
50 
51 template<class CloudType>
53 (
54  const Explicit<CloudType>& cm
55 )
56 :
58  stressAverage_(cm.stressAverage_->clone()),
59  correctionLimiting_
60  (
61  cm.correctionLimiting_->clone()
62  )
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
68 template<class CloudType>
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class CloudType>
77 {
79 
80  if (store)
81  {
82  const fvMesh& mesh = this->owner().mesh();
83  const word& cloudName = this->owner().name();
84 
85  const AveragingMethod<scalar>& volumeAverage =
86  mesh.lookupObject<AveragingMethod<scalar>>
87  (
88  cloudName + ":volumeAverage"
89  );
90  const AveragingMethod<scalar>& rhoAverage =
91  mesh.lookupObject<AveragingMethod<scalar>>
92  (
93  cloudName + ":rhoAverage"
94  );
95  const AveragingMethod<vector>& uAverage =
96  mesh.lookupObject<AveragingMethod<vector>>
97  (
98  cloudName + ":uAverage"
99  );
100  const AveragingMethod<scalar>& uSqrAverage =
101  mesh.lookupObject<AveragingMethod<scalar>>
102  (
103  cloudName + ":uSqrAverage"
104  );
105 
106  volumeAverage_ = &volumeAverage;
107  uAverage_ = &uAverage;
108 
109  stressAverage_.reset
110  (
112  (
113  IOobject
114  (
115  cloudName + ":stressAverage",
116  this->owner().db().time().timeName(),
117  mesh
118  ),
119  this->owner().solution().dict(),
120  mesh
121  ).ptr()
122  );
123 
124  stressAverage_() =
125  this->particleStressModel_->tau
126  (
127  *volumeAverage_,
128  rhoAverage,
129  uSqrAverage
130  )();
131  }
132  else
133  {
134  volumeAverage_ = nullptr;
135  uAverage_ = nullptr;
136  stressAverage_.clear();
137  }
138 }
139 
140 
141 template<class CloudType>
143 (
144  typename CloudType::parcelType& p,
145  const scalar deltaT
146 ) const
147 {
148  const tetIndices tetIs(p.currentTetIndices());
149 
150  // interpolated quantities
151  const scalar alpha =
152  this->volumeAverage_->interpolate(p.coordinates(), tetIs);
153  const vector alphaGrad =
154  this->volumeAverage_->interpolateGrad(p.coordinates(), tetIs);
155  const vector uMean =
156  this->uAverage_->interpolate(p.coordinates(), tetIs);
157 
158  // stress gradient
159  const vector tauGrad =
160  stressAverage_->interpolateGrad(p.coordinates(), tetIs);
161 
162  // parcel relative velocity
163  const vector uRelative = p.U() - uMean;
164 
165  // correction velocity
166  vector dU = Zero;
167 
169  //const scalar Re = p.Re(td);
170  //const typename CloudType::forceType& forces = this->owner().forces();
171  //const forceSuSp F =
172  // forces.calcCoupled(p, td, deltaT, p.mass(), Re, td.muc())
173  // + forces.calcNonCoupled(p, td, deltaT, p.mass(), Re, td.muc());
174 
175  // correction velocity
176  if ((uRelative & alphaGrad) > 0)
177  {
178  dU = - deltaT*tauGrad/(p.rho()*alpha/* + deltaT*F.Sp()*/);
179  }
180 
181  // apply the velocity limiters
182  return
183  correctionLimiting_->limitedVelocity
184  (
185  p.U(),
186  dU,
187  uMean
188  );
189 }
190 
191 
192 // ************************************************************************* //
Foam::AveragingMethod< scalar >
Foam::PackingModels::Explicit::cacheFields
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Explicit.C:76
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:50
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
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:69
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::PackingModel
Base class for packing models.
Definition: MPPICCloud.H:58
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Explicit.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
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:143
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:34
Foam::Vector< scalar >
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220