Loading...
Searching...
No Matches
PressureGradientForce.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) 2011-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
30#include "fvcDdt.H"
31#include "fvcGrad.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class CloudType>
37(
39 const fvMesh& mesh,
40 const dictionary& dict,
41 const word& forceType
42)
43:
44 ParticleForce<CloudType>(owner, mesh, dict, forceType, true),
45 UName_(this->coeffs().template getOrDefault<word>("U", "U")),
46 DUcDtInterpPtr_(nullptr)
47{}
48
49
50template<class CloudType>
52(
53 const PressureGradientForce& pgf
54)
55:
57 UName_(pgf.UName_),
59{}
60
61
62// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
63
64template<class CloudType>
66{}
67
68
69// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70
71template<class CloudType>
73{
74 static word resultName("DUcDt");
75
76 volVectorField* resultPtr =
77 this->mesh().template getObjectPtr<volVectorField>(resultName);
78
79 if (store)
80 {
81 if (!resultPtr)
82 {
83 const volVectorField& Uc = this->mesh().template
84 lookupObject<volVectorField>(UName_);
85
86 resultPtr = new volVectorField
87 (
88 resultName,
89 fvc::ddt(Uc) + (Uc & fvc::grad(Uc))
90 );
91
92 resultPtr->store();
93 }
94 const volVectorField& DUcDt = *resultPtr;
95
96 DUcDtInterpPtr_.reset
97 (
98 interpolation<vector>::New
99 (
100 this->owner().solution().interpolationSchemes(),
101 DUcDt
102 ).ptr()
103 );
104 }
105 else
106 {
107 DUcDtInterpPtr_.clear();
108
109 if (resultPtr)
110 {
111 resultPtr->checkOut();
112 }
113 }
114}
115
116
117template<class CloudType>
119(
120 const typename CloudType::parcelType& p,
121 const typename CloudType::parcelType::trackingData& td,
122 const scalar dt,
123 const scalar mass,
124 const scalar Re,
125 const scalar muc
126) const
127{
128 forceSuSp value(Zero);
129
130 vector DUcDt =
131 DUcDtInterp().interpolate(p.coordinates(), p.currentTetIndices());
132
133 value.Su() = mass*td.rhoc()/p.rho()*DUcDt;
134
135 return value;
136}
137
138
139template<class CloudType>
141(
142 const typename CloudType::parcelType&,
143 const typename CloudType::parcelType::trackingData& td,
144 const scalar
145) const
146{
147 return 0.0;
148}
149
150
151// ************************************************************************* //
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.
Calculates particle pressure gradient force.
autoPtr< interpolation< vector > > DUcDtInterpPtr_
Rate of change of carrier phase velocity interpolator.
virtual void cacheFields(const bool store)
Cache fields.
const interpolation< vector > & DUcDtInterp() const
Return the rate of change of carrier phase velocity interpolator.
virtual scalar massAdd(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar mass) const
Return the added mass.
PressureGradientForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType=typeName)
Construct from mesh.
const word UName_
Name of velocity field.
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.
virtual ~PressureGradientForce()
Destructor.
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 first temporal derivative.
Calculate the gradient of the given field.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.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