initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / intermediate / particleForces / particleForces.C
blob5320566c84382e91f744263062c1e2c7b493dd3f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "particleForces.H"
28 #include "fvMesh.H"
29 #include "volFields.H"
30 #include "fvcGrad.H"
32 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
34 Foam::particleForces::particleForces
36     const fvMesh& mesh,
37     const dictionary& dict,
38     const vector& g
41     mesh_(mesh),
42     dict_(dict.subDict("particleForces")),
43     g_(g),
44     gradUPtr_(NULL),
45     gravity_(dict_.lookup("gravity")),
46     virtualMass_(dict_.lookup("virtualMass")),
47     Cvm_(0.0),
48     pressureGradient_(dict_.lookup("pressureGradient")),
49     UName_(dict_.lookupOrDefault<word>("U", "U"))
51     if (virtualMass_)
52     {
53         dict_.lookup("Cvm") >> Cvm_;
54     }
58 Foam::particleForces::particleForces(const particleForces& f)
60     mesh_(f.mesh_),
61     dict_(f.dict_),
62     g_(f.g_),
63     gradUPtr_(f.gradUPtr_),
64     gravity_(f.gravity_),
65     virtualMass_(f.virtualMass_),
66     Cvm_(f.Cvm_),
67     pressureGradient_(f.pressureGradient_),
68     UName_(f.UName_)
72 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
74 Foam::particleForces::~particleForces()
76     cacheFields(false);
80 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
82 const Foam::dictionary& Foam::particleForces::dict() const
84     return dict_;
88 const Foam::vector& Foam::particleForces::g() const
90     return g_;
94 Foam::Switch Foam::particleForces::gravity() const
96     return gravity_;
100 Foam::Switch Foam::particleForces::virtualMass() const
102     return virtualMass_;
106 Foam::Switch Foam::particleForces::pressureGradient() const
108     return pressureGradient_;
112 const Foam::word& Foam::particleForces::UName() const
114     return UName_;
118 void Foam::particleForces::cacheFields(const bool store)
120     if (store && pressureGradient_)
121     {
122         const volVectorField U = mesh_.lookupObject<volVectorField>(UName_);
123         gradUPtr_ = fvc::grad(U).ptr();
124     }
125     else
126     {
127         if (gradUPtr_)
128         {
129             delete gradUPtr_;
130         }
131     }
135 Foam::vector Foam::particleForces::calcCoupled
137     const label cellI,
138     const scalar dt,
139     const scalar rhoc,
140     const scalar rho,
141     const vector& Uc,
142     const vector& U
143 ) const
145     vector Ftot = vector::zero;
147     // Virtual mass force
148     if (virtualMass_)
149     {
150         notImplemented
151         (
152             "Foam::particleForces::calcCoupled(...) - virtual mass force"
153         );
154 //        Ftot += Cvm_*rhoc/rho*d(Uc - U)/dt;
155     }
157     // Pressure gradient force
158     if (pressureGradient_)
159     {
160         const volTensorField& gradU = *gradUPtr_;
161         Ftot += rhoc/rho*(U & gradU[cellI]);
162     }
164     return Ftot;
168 Foam::vector Foam::particleForces::calcNonCoupled
170     const label cellI,
171     const scalar dt,
172     const scalar rhoc,
173     const scalar rho,
174     const vector& Uc,
175     const vector& U
176 ) const
178     vector Ftot = vector::zero;
180     // Gravity force
181     if (gravity_)
182     {
183         Ftot += g_*(1.0 - rhoc/rho);
184     }
186     return Ftot;
190 // ************************************************************************* //