InterfaceForce.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) 2016 OpenCFD Ltd.
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 "InterfaceForce.H"
29 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  CloudType& owner,
37  const fvMesh& mesh,
38  const dictionary& dict
39 )
40 :
41  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
42  alphaName_(this->coeffs().lookup("alpha")),
43  C_(this->coeffs().getScalar("C")),
44  gradInterForceInterpPtr_(nullptr)
45 {}
46 
47 
48 template<class CloudType>
50 :
52  alphaName_(pf.alphaName_),
53  C_(pf.C_),
54  gradInterForceInterpPtr_(pf.gradInterForceInterpPtr_)
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
59 
60 template<class CloudType>
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 template<class CloudType>
69 {
70  static word resultName("gradAlpha");
71 
72  volVectorField* resultPtr =
73  this->mesh().template getObjectPtr<volVectorField>(resultName);
74 
75  if (store)
76  {
77  if (!resultPtr)
78  {
79  const volScalarField& alpha = this->mesh().template
80  lookupObject<volScalarField>(alphaName_);
81 
82  resultPtr = new volVectorField
83  (
84  resultName,
85  fvc::grad(alpha*(1-alpha))
86  );
87 
88  resultPtr->store();
89  }
90 
91  const volVectorField& gradInterForce = *resultPtr;
92 
93  gradInterForceInterpPtr_.reset
94  (
96  (
97  this->owner().solution().interpolationSchemes(),
98  gradInterForce
99  ).ptr()
100  );
101  }
102  else
103  {
104  gradInterForceInterpPtr_.clear();
105 
106  if (resultPtr)
107  {
108  resultPtr->checkOut();
109  }
110  }
111 }
112 
113 
114 template<class CloudType>
116 (
117  const typename CloudType::parcelType& p,
118  const typename CloudType::parcelType::trackingData& td,
119  const scalar dt,
120  const scalar mass,
121  const scalar Re,
122  const scalar muc
123 ) const
124 {
125  forceSuSp value(Zero);
126 
127  const interpolation<vector>& interp = gradInterForceInterpPtr_();
128 
129  value.Su() =
130  C_
131  *mass
132  *interp.interpolate(p.coordinates(), p.currentTetIndices());
133 
134  return value;
135 }
136 
137 
138 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Vector force apply to particles to avoid the crossing of particles from one phase to the other...
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Abstract base class for particle forces.
Definition: ParticleForce.H:54
Lookup type of boundary radiation properties.
Definition: lookup.H:57
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:60
virtual ~InterfaceForce()
Destructor.
dynamicFvMesh & mesh
Calculate the gradient of the given field.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
A class for handling words, derived from Foam::string.
Definition: word.H:63
InterfaceForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual forceSuSp calcNonCoupled(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.
volScalarField & p
virtual void cacheFields(const bool store)
Cache fields.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127