LiftForce.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-2016 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 "LiftForce.H"
27 #include "fvcCurl.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class CloudType>
33 (
34  const typename CloudType::parcelType& p,
35  const vector& curlUc,
36  const scalar Re,
37  const scalar muc
38 ) const
39 {
40  // dummy
41  return 0.0;
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 template<class CloudType>
49 (
50  CloudType& owner,
51  const fvMesh& mesh,
52  const dictionary& dict,
53  const word& forceType
54 )
55 :
56  ParticleForce<CloudType>(owner, mesh, dict, forceType, true),
57  UName_(this->coeffs().template lookupOrDefault<word>("U", "U")),
58  curlUcInterpPtr_(NULL)
59 {}
60 
61 
62 template<class CloudType>
64 :
66  UName_(lf.UName_),
67  curlUcInterpPtr_(NULL)
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
72 
73 template<class CloudType>
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class CloudType>
82 {
83  static word fName("curlUcDt");
84 
85  bool fieldExists = this->mesh().template foundObject<volVectorField>(fName);
86 
87  if (store)
88  {
89  if (!fieldExists)
90  {
91  const volVectorField& Uc = this->mesh().template
92  lookupObject<volVectorField>(UName_);
93 
94  volVectorField* curlUcPtr =
95  new volVectorField(fName, fvc::curl(Uc));
96 
97  curlUcPtr->store();
98  }
99 
100  const volVectorField& curlUc = this->mesh().template
101  lookupObject<volVectorField>(fName);
102 
103  curlUcInterpPtr_.reset
104  (
106  (
107  this->owner().solution().interpolationSchemes(),
108  curlUc
109  ).ptr()
110  );
111  }
112  else
113  {
114  curlUcInterpPtr_.clear();
115 
116  if (fieldExists)
117  {
118  const volVectorField& curlUc = this->mesh().template
119  lookupObject<volVectorField>(fName);
120 
121  const_cast<volVectorField&>(curlUc).checkOut();
122  }
123  }
124 }
125 
126 
127 template<class CloudType>
129 (
130  const typename CloudType::parcelType& p,
131  const scalar dt,
132  const scalar mass,
133  const scalar Re,
134  const scalar muc
135 ) const
136 {
137  forceSuSp value(Zero, 0.0);
138 
139  vector curlUc =
140  curlUcInterp().interpolate(p.position(), p.currentTetIndices());
141 
142  scalar Cl = this->Cl(p, curlUc, Re, muc);
143 
144  value.Su() = mass/p.rho()*p.rhoc()*Cl*((p.Uc() - p.U())^curlUc);
145 
146  return value;
147 }
148 
149 
150 // ************************************************************************* //
virtual void cacheFields(const bool store)
Cache fields.
Definition: LiftForce.C:81
dictionary dict
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
Definition: LiftForce.C:129
LiftForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType)
Construct from mesh.
Definition: LiftForce.C:49
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Abstract base class for particle forces.
Definition: ParticleForce.H:53
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:56
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const fvMesh & mesh() const
Return the mesh database.
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcCurl.C:45
autoPtr< interpolation< vector > > curlUcInterpPtr_
Curk of carrier phase velocity interpolator.
Definition: LiftForce.H:64
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:61
dynamicFvMesh & mesh
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
A class for handling words, derived from string.
Definition: word.H:59
virtual scalar Cl(const typename CloudType::parcelType &p, const vector &curlUc, const scalar Re, const scalar muc) const
Calculate the lift coefficient.
static const zero Zero
Definition: zero.H:91
virtual ~LiftForce()
Destructor.
Definition: LiftForce.C:74
const interpolation< vector > & curlUcInterp() const
Return the curl of the carrier phase velocity interpolator.
Definition: LiftForceI.H:30
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Abstract base class for interpolation.
Base class for particle lift force models.
Definition: LiftForce.H:52
const CloudType & owner() const
Return const access to the cloud owner.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
const word UName_
Name of velocity field.
Definition: LiftForce.H:61