Loading...
Searching...
No Matches
ParticleTrap.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2020 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\*---------------------------------------------------------------------------*/
29#include "ParticleTrap.H"
30#include "fvcGrad.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class CloudType>
36(
37 const dictionary& dict,
39 const word& modelName
40)
41:
43 alphaName_
44 (
45 this->coeffDict().template getOrDefault<word>("alpha", "alpha")
46 ),
47 alphaPtr_(nullptr),
48 gradAlphaPtr_(nullptr),
49 threshold_(this->coeffDict().getScalar("threshold"))
50{}
51
52
53template<class CloudType>
55(
57)
58:
60 alphaName_(pt.alphaName_),
61 alphaPtr_(pt.alphaPtr_),
62 gradAlphaPtr_(nullptr),
63 threshold_(pt.threshold_)
64{}
65
66
67// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68
69template<class CloudType>
71(
72 const typename parcelType::trackingData& td
73)
74{
75 if (alphaPtr_ == nullptr)
76 {
77 const fvMesh& mesh = this->owner().mesh();
78 const volScalarField& alpha =
79 mesh.lookupObject<volScalarField>(alphaName_);
80
81 alphaPtr_ = &alpha;
82 }
83
84 if (gradAlphaPtr_)
85 {
86 *gradAlphaPtr_ == fvc::grad(*alphaPtr_);
87 }
88 else
89 {
90 gradAlphaPtr_.reset(new volVectorField(fvc::grad(*alphaPtr_)));
91 }
92}
93
94
95template<class CloudType>
97(
98 const typename parcelType::trackingData& td
100{
101 gradAlphaPtr_.clear();
102}
103
104
105template<class CloudType>
107(
108 parcelType& p,
109 const scalar,
110 const point&,
111 const typename parcelType::trackingData& td
112)
113{
114 if (alphaPtr_->primitiveField()[p.cell()] < threshold_)
115 {
116 const vector& gradAlpha = gradAlphaPtr_()[p.cell()];
117 vector nHat = gradAlpha/mag(gradAlpha);
118 scalar nHatU = nHat & p.U();
119
120 if (nHatU < 0)
121 {
122 p.U() -= 2*nHat*nHatU;
123 }
124 }
125
126 return true;
127}
128
129
130// ************************************************************************* //
Templated cloud function object base class.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
particle::trackingData trackingData
Definition DSMCParcel.H:155
Traps particles within a given phase fraction for multi-phase cases.
ParticleTrap(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void postEvolve(const typename parcelType::trackingData &td)
Post-evolve hook.
virtual void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve hook.
virtual bool postMove(parcelType &p, const scalar dt, const point &position0, const typename parcelType::trackingData &td)
Post-move hook.
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
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
Calculate the gradient of the given field.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
Definition point.H:37
Vector< scalar > vector
Definition vector.H:57
volScalarField & alpha