Loading...
Searching...
No Matches
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-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
33template<class CloudType>
35(
36 const dictionary& dict,
38)
39:
41 stressAverage_(nullptr),
42 correctionLimiting_
43 (
45 (
47 )
48 )
49{}
50
51
52template<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
69template<class CloudType>
71{}
72
73
74// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75
76template<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 auto& volumeAverage =
87 mesh.lookupObject<AveragingMethod<scalar>>
88 (
89 IOobject::scopedName(cloudName, "volumeAverage")
90 );
91 const auto& rhoAverage =
92 mesh.lookupObject<AveragingMethod<scalar>>
93 (
94 IOobject::scopedName(cloudName, "rhoAverage")
95 );
96 const auto& uAverage =
97 mesh.lookupObject<AveragingMethod<vector>>
98 (
100 );
101 const auto& uSqrAverage =
102 mesh.lookupObject<AveragingMethod<scalar>>
103 (
104 IOobject::scopedName(cloudName, "uSqrAverage")
105 );
106
107 volumeAverage_ = &volumeAverage;
108 uAverage_ = &uAverage;
109
110 stressAverage_.reset
111 (
113 (
115 (
116 IOobject::scopedName(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
142template<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// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
Base class for lagrangian averaging methods.
virtual TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const =0
Interpolate gradient.
static autoPtr< AveragingMethod< Type > > New(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Selector.
virtual Type interpolate(const barycentric &coordinates, const tetIndices &tetIs) const =0
Interpolate.
const CloudType & owner() const
Return const access to the owner cloud.
Base class for correction limiting methods.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Base class for packing models.
autoPtr< ParticleStressModel > particleStressModel_
Protected data.
static autoPtr< PackingModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector.
PackingModel(CloudType &owner)
Construct null from owner.
Explicit model for applying an inter-particle stress to the particles.
Definition Explicit.H:69
virtual ~Explicit()
Destructor.
Definition Explicit.C:63
virtual autoPtr< PackingModel< CloudType > > clone() const
Construct and return a clone.
Definition Explicit.H:119
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
Definition Explicit.C:70
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
Definition Explicit.C:137
Explicit(const dictionary &dict, CloudType &owner)
Construct from components.
Definition Explicit.C:28
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition PtrListI.H:98
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
virtual void cacheFields(const bool store)
Cache dependent sub-model fields.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
DSMCCloud< dsmcParcel > CloudType
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
volScalarField & alpha
dictionary dict