Loading...
Searching...
No Matches
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 | 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\*---------------------------------------------------------------------------*/
28
29#include "LiftForce.H"
30#include "fvcCurl.H"
31
32// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33
34template<class CloudType>
36(
37 const typename CloudType::parcelType& p,
38 const typename CloudType::parcelType::trackingData& td,
39 const vector& curlUc,
40 const scalar Re,
41 const scalar muc
42) const
43{
44 // dummy
45 return 0.0;
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
51template<class CloudType>
53(
54 CloudType& owner,
55 const fvMesh& mesh,
56 const dictionary& dict,
57 const word& forceType
58)
59:
60 ParticleForce<CloudType>(owner, mesh, dict, forceType, true),
61 UName_(this->coeffs().template getOrDefault<word>("U", "U")),
62 curlUcInterpPtr_(nullptr)
63{}
64
65
66template<class CloudType>
68:
70 UName_(lf.UName_),
72{}
73
74
75// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
76
77template<class CloudType>
79{}
80
81
82// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83
84template<class CloudType>
85void Foam::LiftForce<CloudType>::cacheFields(const bool store)
86{
87 static word resultName("curlUcDt");
88
89 volVectorField* resultPtr =
90 this->mesh().template getObjectPtr<volVectorField>(resultName);
91
92 if (store)
93 {
94 if (!resultPtr)
95 {
96 const volVectorField& Uc = this->mesh().template
97 lookupObject<volVectorField>(UName_);
98
99 resultPtr = new volVectorField(resultName, fvc::curl(Uc));
100
101 resultPtr->store();
102 }
103
104 const volVectorField& curlUc = *resultPtr;
105
106 curlUcInterpPtr_.reset
107 (
108 interpolation<vector>::New
109 (
110 this->owner().solution().interpolationSchemes(),
111 curlUc
112 ).ptr()
113 );
114 }
115 else
116 {
117 curlUcInterpPtr_.clear();
118
119 if (resultPtr)
120 {
121 resultPtr->checkOut();
122 }
123 }
124}
125
126
127template<class CloudType>
129(
130 const typename CloudType::parcelType& p,
131 const typename CloudType::parcelType::trackingData& td,
132 const scalar dt,
133 const scalar mass,
134 const scalar Re,
135 const scalar muc
136) const
137{
138 forceSuSp value(Zero);
139
140 vector curlUc =
141 curlUcInterp().interpolate(p.coordinates(), p.currentTetIndices());
142
143 scalar Cl = this->Cl(p, td, curlUc, Re, muc);
144
145 value.Su() = mass/p.rho()*td.rhoc()*Cl*((td.Uc() - p.U())^curlUc);
146
147 return value;
148}
149
150
151// ************************************************************************* //
Base class for particle lift force models.
Definition LiftForce.H:56
autoPtr< interpolation< vector > > curlUcInterpPtr_
Curl of carrier phase velocity interpolator.
Definition LiftForce.H:69
const interpolation< vector > & curlUcInterp() const
Return the curl of the carrier phase velocity interpolator.
Definition LiftForceI.H:26
virtual void cacheFields(const bool store)
Cache fields.
Definition LiftForce.C:78
LiftForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType)
Construct from mesh.
Definition LiftForce.C:46
virtual ~LiftForce()
Destructor.
Definition LiftForce.C:71
virtual scalar Cl(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const vector &curlUc, const scalar Re, const scalar muc) const
Calculate the lift coefficient.
const word UName_
Name of velocity field.
Definition LiftForce.H:64
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
Definition LiftForce.C:122
Abstract base class for particle forces.
const fvMesh & mesh() const noexcept
Return the mesh database.
const CloudType & owner() const noexcept
Return const access to the cloud owner.
ParticleForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType, const bool readCoeffs)
Construct from mesh.
const dictionary & coeffs() const noexcept
Return the force coefficients dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Helper container for force Su and Sp terms.
Definition forceSuSp.H:63
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition forceSuSpI.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
bool store()
Register object with its registry and transfer ownership to the registry.
bool checkOut()
Remove object from registry, and remove all file watches.
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvcCurl.C:40
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Vector< scalar > vector
Definition vector.H:57
dictionary dict