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 | Copyright (C) 2013-2014 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_(NULL),
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 =
84  mesh.lookupObject<AveragingMethod<scalar> >
85  (
86  cloudName + ":volumeAverage"
87  );
88  const AveragingMethod<scalar>& rhoAverage =
89  mesh.lookupObject<AveragingMethod<scalar> >
90  (
91  cloudName + ":rhoAverage"
92  );
93  const AveragingMethod<vector>& uAverage =
94  mesh.lookupObject<AveragingMethod<vector> >
95  (
96  cloudName + ":uAverage"
97  );
98  const AveragingMethod<scalar>& uSqrAverage =
99  mesh.lookupObject<AveragingMethod<scalar> >
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_ = NULL;
133  uAverage_ = NULL;
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 fvMesh& mesh = this->owner().mesh();
147  const tetIndices tetIs(p.cell(), p.tetFace(), p.tetPt(), mesh);
148 
149  // interpolated quantities
150  const scalar alpha =
151  this->volumeAverage_->interpolate(p.position(), tetIs);
152  const vector alphaGrad =
153  this->volumeAverage_->interpolateGrad(p.position(), tetIs);
154  const vector uMean =
155  this->uAverage_->interpolate(p.position(), tetIs);
156 
157  // stress gradient
158  const vector tauGrad =
159  stressAverage_->interpolateGrad(p.position(), tetIs);
160 
161  // parcel relative velocity
162  const vector uRelative = p.U() - uMean;
163 
164  // correction velocity
165  vector dU = vector::zero;
166 
168  //const scalar Re = p.Re(p.U(), p.d(), p.rhoc(), p.muc());
169  //const typename CloudType::forceType& forces = this->owner().forces();
170  //const forceSuSp F =
171  // forces.calcCoupled(p, deltaT, p.mass(), Re, p.muc())
172  // + forces.calcNonCoupled(p, deltaT, p.mass(), Re, p.muc());
173 
174  // correction velocity
175  if ((uRelative & alphaGrad) > 0)
176  {
177  dU = - deltaT*tauGrad/(p.rho()*alpha/* + deltaT*F.Sp()*/);
178  }
179 
180  // apply the velocity limiters
181  return
182  correctionLimiting_->limitedVelocity
183  (
184  p.U(),
185  dU,
186  uMean
187  );
188 }
189 
190 
191 // ************************************************************************* //
Foam::AveragingMethod< scalar >
Foam::PackingModels::Explicit::cacheFields
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition: Explicit.C:74
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::PackingModels::Explicit
Explicit model for applying an inter-particle stress to the particles.
Definition: Explicit.H:65
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::PackingModels::Explicit::stressAverage_
autoPtr< AveragingMethod< scalar > > stressAverage_
Stress average field.
Definition: Explicit.H:80
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Foam::PackingModels::Explicit::correctionLimiting_
autoPtr< CorrectionLimitingMethod > correctionLimiting_
Correction limiter.
Definition: Explicit.H:83
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Foam::PackingModels::Explicit::~Explicit
virtual ~Explicit()
Destructor.
Definition: Explicit.C:67
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::PackingModel< CloudType >
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Explicit.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::PackingModels::Explicit::velocityCorrection
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition: Explicit.C:141
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:73
Foam::PackingModels::Explicit::Explicit
Explicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition: Explicit.C:32
Foam::Vector< scalar >
cloudName
const word cloudName(propsDict.lookup("cloudName"))
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217