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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "PressureGradientForce.H"
29 #include "fvcDdt.H"
30 #include "fvcGrad.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  CloudType& owner,
38  const fvMesh& mesh,
39  const dictionary& dict,
40  const word& forceType
41 )
42 :
43  ParticleForce<CloudType>(owner, mesh, dict, forceType, true),
44  UName_(this->coeffs().template lookupOrDefault<word>("U", "U")),
45  DUcDtInterpPtr_(nullptr)
46 {}
47 
48 
49 template<class CloudType>
51 (
52  const PressureGradientForce& pgf
53 )
54 :
56  UName_(pgf.UName_),
57  DUcDtInterpPtr_(nullptr)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
62 
63 template<class CloudType>
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class CloudType>
72 {
73  static word fName("DUcDt");
74 
75  bool fieldExists = this->mesh().template foundObject<volVectorField>(fName);
76 
77  if (store)
78  {
79  if (!fieldExists)
80  {
81  const volVectorField& Uc = this->mesh().template
82  lookupObject<volVectorField>(UName_);
83 
84  volVectorField* DUcDtPtr = new volVectorField
85  (
86  fName,
87  fvc::ddt(Uc) + (Uc & fvc::grad(Uc))
88  );
89 
90  DUcDtPtr->store();
91  }
92 
93  const volVectorField& DUcDt = this->mesh().template
94  lookupObject<volVectorField>(fName);
95 
96  DUcDtInterpPtr_.reset
97  (
99  (
100  this->owner().solution().interpolationSchemes(),
101  DUcDt
102  ).ptr()
103  );
104  }
105  else
106  {
107  DUcDtInterpPtr_.clear();
108 
109  if (fieldExists)
110  {
111  volVectorField& DUcDt =
112  this->mesh().template lookupObjectRef<volVectorField>(fName);
113 
114  DUcDt.checkOut();
115  }
116  }
117 }
118 
119 
120 template<class CloudType>
122 (
123  const typename CloudType::parcelType& p,
124  const typename CloudType::parcelType::trackingData& td,
125  const scalar dt,
126  const scalar mass,
127  const scalar Re,
128  const scalar muc
129 ) const
130 {
131  forceSuSp value(Zero);
132 
133  vector DUcDt =
134  DUcDtInterp().interpolate(p.coordinates(), p.currentTetIndices());
135 
136  value.Su() = mass*td.rhoc()/p.rho()*DUcDt;
137 
138  return value;
139 }
140 
141 
142 template<class CloudType>
144 (
145  const typename CloudType::parcelType&,
146  const typename CloudType::parcelType::trackingData& td,
147  const scalar
148 ) const
149 {
150  return 0.0;
151 }
152 
153 
154 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::PressureGradientForce::PressureGradientForce
PressureGradientForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType=typeName)
Construct from mesh.
Definition: PressureGradientForce.C:36
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
PressureGradientForce.H
Foam::PressureGradientForce::~PressureGradientForce
virtual ~PressureGradientForce()
Destructor.
Definition: PressureGradientForce.C:64
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
Foam::PressureGradientForce::cacheFields
virtual void cacheFields(const bool store)
Cache fields.
Definition: PressureGradientForce.C:71
Foam::forceSuSp::Su
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:64
Foam::PressureGradientForce
Calculates particle pressure gradient force.
Definition: PressureGradientForce.H:57
Foam::PressureGradientForce::calcCoupled
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: PressureGradientForce.C:122
Foam::forceSuSp
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:63
Foam::ParticleForce
Abstract base class for particle forces.
Definition: ParticleForce.H:59
Foam::PressureGradientForce::UName_
const word UName_
Name of velocity field.
Definition: PressureGradientForce.H:66
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:95
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
Foam::Vector< scalar >
fvcGrad.H
Calculate the gradient of the given field.
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
fvcDdt.H
Calculate the first temporal derivative.
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::PressureGradientForce::massAdd
virtual scalar massAdd(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar mass) const
Return the added mass.
Definition: PressureGradientForce.C:144