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 | Copyright (C) 2012-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 "ParticleTrap.H"
27 #include "fvcGrad.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class CloudType>
33 (
34  const dictionary& dict,
35  CloudType& owner,
36  const word& modelName
37 )
38 :
39  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
40  alphaName_
41  (
42  this->coeffDict().template lookupOrDefault<word>("alphaName", "alpha")
43  ),
44  alphaPtr_(NULL),
45  gradAlphaPtr_(NULL),
46  threshold_(readScalar(this->coeffDict().lookup("threshold")))
47 {}
48 
49 
50 template<class CloudType>
52 (
53  const ParticleTrap<CloudType>& pt
54 )
55 :
57  alphaName_(pt.alphaName_),
58  alphaPtr_(pt.alphaPtr_),
59  gradAlphaPtr_(NULL),
60  threshold_(pt.threshold_)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
66 template<class CloudType>
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 template<class CloudType>
75 {
76  if (alphaPtr_ == NULL)
77  {
78  const fvMesh& mesh = this->owner().mesh();
79  const volScalarField& alpha =
80  mesh.lookupObject<volScalarField>(alphaName_);
81 
82  alphaPtr_ = &alpha;
83  }
84 
85  if (gradAlphaPtr_.valid())
86  {
87  gradAlphaPtr_() == fvc::grad(*alphaPtr_);
88  }
89  else
90  {
91  gradAlphaPtr_.reset(new volVectorField(fvc::grad(*alphaPtr_)));
92  }
93 }
94 
95 
96 template<class CloudType>
98 {
99  gradAlphaPtr_.clear();
100 }
101 
102 
103 template<class CloudType>
105 (
106  parcelType& p,
107  const label cellI,
108  const scalar,
109  const point&,
110  bool&
111 )
112 {
113  if (alphaPtr_->internalField()[cellI] < threshold_)
114  {
115  const vector& gradAlpha = gradAlphaPtr_()[cellI];
116  vector nHat = gradAlpha/mag(gradAlpha);
117  scalar nHatU = nHat & p.U();
118 
119  if (nHatU < 0)
120  {
121  p.U() -= 2*nHat*nHatU;
122  }
123  }
124 }
125 
126 
127 // ************************************************************************* //
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::ParticleTrap
Traps particles within a given phase fraction for multi-phase cases.
Definition: ParticleTrap.H:60
Foam::ParticleTrap::postEvolve
virtual void postEvolve()
Post-evolve hook.
Definition: ParticleTrap.C:97
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::ParticleTrap::~ParticleTrap
virtual ~ParticleTrap()
Destructor.
Definition: ParticleTrap.C:67
Foam::ParticleTrap::parcelType
CloudType::parcelType parcelType
Convenience typedef for parcel type.
Definition: ParticleTrap.H:69
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::ParticleTrap::alphaPtr_
const volScalarField * alphaPtr_
Pointer to the volume fraction field.
Definition: ParticleTrap.H:76
Foam::ParticleTrap::postMove
virtual void postMove(typename CloudType::parcelType &p, const label cellI, const scalar dt, const point &position0, bool &keepParticle)
Post-move hook.
Definition: ParticleTrap.C:105
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
ParticleTrap.H
Foam::Vector< scalar >
Foam::ParticleTrap::alphaName_
const word alphaName_
Name of vol fraction field.
Definition: ParticleTrap.H:73
Foam::CloudFunctionObject
Templated cloud function object base class.
Definition: CloudFunctionObject.H:56
fvcGrad.H
Calculate the gradient of the given field.
Foam::ParticleTrap::threshold_
scalar threshold_
Threshold beyond which model is active.
Definition: ParticleTrap.H:82
Foam::ParticleTrap::preEvolve
virtual void preEvolve()
Pre-evolve hook.
Definition: ParticleTrap.C:74
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::ParticleTrap::ParticleTrap
ParticleTrap(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ParticleTrap.C:33
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress