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_),
71 curlUcInterpPtr_(nullptr)
72{}
73
74
75// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
76
77template<class CloudType>
79{}
80
81
82// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83
84template<class CloudType>
86{
87 static word fName("curlUcDt");
88
89 bool fieldExists = this->mesh().template foundObject<volVectorField>(fName);
90
91 if (store)
92 {
93 if (!fieldExists)
94 {
95 const volVectorField& Uc = this->mesh().template
96 lookupObject<volVectorField>(UName_);
97
98 volVectorField* curlUcPtr =
99 new volVectorField(fName, fvc::curl(Uc));
100
101 curlUcPtr->store();
102 }
103
104 const volVectorField& curlUc = this->mesh().template
105 lookupObject<volVectorField>(fName);
106
107 curlUcInterpPtr_.reset
108 (
110 (
111 this->owner().solution().interpolationSchemes(),
112 curlUc
113 ).ptr()
114 );
115 }
116 else
117 {
118 curlUcInterpPtr_.clear();
119
120 if (fieldExists)
121 {
122 volVectorField& curlUc =
123 this->mesh().template lookupObjectRef<volVectorField>(fName);
124
125 curlUc.checkOut();
126 }
127 }
128}
129
130
131template<class CloudType>
133(
134 const typename CloudType::parcelType& p,
135 const typename CloudType::parcelType::trackingData& td,
136 const scalar dt,
137 const scalar mass,
138 const scalar Re,
139 const scalar muc
140) const
141{
142 forceSuSp value(Zero);
143
144 vector curlUc =
145 curlUcInterp().interpolate(p.coordinates(), p.currentTetIndices());
146
147 scalar Cl = this->Cl(p, td, curlUc, Re, muc);
148
149 value.Su() = mass/p.rho()*td.rhoc()*Cl*((td.Uc() - p.U())^curlUc);
150
151 return value;
152}
153
154
155// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Base class for particle lift force models.
Definition: LiftForce.H:60
virtual void cacheFields(const bool store)
Cache fields.
Definition: LiftForce.C:85
virtual ~LiftForce()
Destructor.
Definition: LiftForce.C:78
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:133
Abstract base class for particle forces.
Definition: ParticleForce.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:67
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:61
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
bool checkOut()
Remove all file watches and remove object from registry.
Definition: regIOobject.C:224
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcCurl.C:47
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict