PressureGradientForce.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "PressureGradientForce.H"
27 #include "fvcDdt.H"
28 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 (
35  CloudType& owner,
36  const fvMesh& mesh,
37  const dictionary& dict,
38  const word& forceType
39 )
40 :
41  ParticleForce<CloudType>(owner, mesh, dict, forceType, true),
42  UName_(this->coeffs().template lookupOrDefault<word>("U", "U")),
43  DUcDtInterpPtr_(nullptr)
44 {}
45 
46 
47 template<class CloudType>
49 (
50  const PressureGradientForce& pgf
51 )
52 :
54  UName_(pgf.UName_),
55  DUcDtInterpPtr_(nullptr)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
60 
61 template<class CloudType>
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
68 template<class CloudType>
70 {
71  static word fName("DUcDt");
72 
73  bool fieldExists = this->mesh().template foundObject<volVectorField>(fName);
74 
75  if (store)
76  {
77  if (!fieldExists)
78  {
79  const volVectorField& Uc = this->mesh().template
80  lookupObject<volVectorField>(UName_);
81 
82  volVectorField* DUcDtPtr = new volVectorField
83  (
84  fName,
85  fvc::ddt(Uc) + (Uc & fvc::grad(Uc))
86  );
87 
88  DUcDtPtr->store();
89  }
90 
91  const volVectorField& DUcDt = this->mesh().template
92  lookupObject<volVectorField>(fName);
93 
94  DUcDtInterpPtr_.reset
95  (
97  (
98  this->owner().solution().interpolationSchemes(),
99  DUcDt
100  ).ptr()
101  );
102  }
103  else
104  {
105  DUcDtInterpPtr_.clear();
106 
107  if (fieldExists)
108  {
109  const volVectorField& DUcDt = this->mesh().template
110  lookupObject<volVectorField>(fName);
111 
112  const_cast<volVectorField&>(DUcDt).checkOut();
113  }
114  }
115 }
116 
117 
118 template<class CloudType>
120 (
121  const typename CloudType::parcelType& p,
122  const typename CloudType::parcelType::trackingData& td,
123  const scalar dt,
124  const scalar mass,
125  const scalar Re,
126  const scalar muc
127 ) const
128 {
129  forceSuSp value(Zero, 0.0);
130 
131  const vector DUcDt =
132  DUcDtInterp().interpolate
133  (
134  p.coordinates(),
135  p.currentTetIndices(td.mesh)
136  );
137 
138  value.Su() = mass*td.rhoc()/p.rho()*DUcDt;
139 
140  return value;
141 }
142 
143 
144 template<class CloudType>
146 (
147  const typename CloudType::parcelType&,
148  const typename CloudType::parcelType::trackingData& td,
149  const scalar
150 ) const
151 {
152  return 0.0;
153 }
154 
155 
156 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Generic GeometricField class.
Abstract base class for particle forces.
Definition: ParticleForce.H:54
Calculates particle pressure gradient force.
virtual void cacheFields(const bool store)
Cache fields.
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.
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 any number of values (e....
Definition: dictionary.H:162
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:64
const vector & Su() const
Return const access to the explicit contribution [kg m/s^2].
Definition: forceSuSpI.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Abstract base class for interpolation.
Definition: interpolation.H:55
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
A class for handling words, derived from string.
Definition: word.H:62
Calculate the first temporal derivative.
Calculate the gradient of the given field.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
static const zero Zero
Definition: zero.H:97
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict
volScalarField & p