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-------------------------------------------------------------------------------
10License
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
33template<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
48template<class CloudType>
50:
52 alphaName_(pf.alphaName_),
53 C_(pf.C_),
54 gradInterForceInterpPtr_(pf.gradInterForceInterpPtr_)
55{}
56
57
58// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
59
60template<class CloudType>
62{}
63
64
65// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66
67template<class CloudType>
69{
70 static word fName("gradAlpha");
71
72 bool fieldExists =
73 this->mesh().template foundObject<volVectorField>(fName);
74
75 if (store)
76 {
77 if (!fieldExists)
78 {
79 const volScalarField& alpha = this->mesh().template
80 lookupObject<volScalarField>(alphaName_);
81
82 volVectorField* gradInterForcePtr =
83 new volVectorField(fName, fvc::grad(alpha*(1-alpha)));
84
85 gradInterForcePtr->store();
86 }
87
88 const volVectorField& gradInterForce = this->mesh().template
89 lookupObject<volVectorField>(fName);
90
91 gradInterForceInterpPtr_.reset
92 (
94 (
95 this->owner().solution().interpolationSchemes(),
96 gradInterForce
97 ).ptr()
98 );
99 }
100 else
101 {
102 gradInterForceInterpPtr_.clear();
103
104 if (fieldExists)
105 {
106 volVectorField& gradInterForce =
107 this->mesh().template lookupObjectRef<volVectorField>(fName);
108
109 gradInterForce.checkOut();
110 }
111 }
112}
113
114
115template<class CloudType>
117(
118 const typename CloudType::parcelType& p,
119 const typename CloudType::parcelType::trackingData& td,
120 const scalar dt,
121 const scalar mass,
122 const scalar Re,
123 const scalar muc
124) const
125{
126 forceSuSp value(Zero);
127
128 const interpolation<vector>& interp = gradInterForceInterpPtr_();
129
130 value.Su() =
131 C_
132 *mass
133 *interp.interpolate(p.coordinates(), p.currentTetIndices());
134
135 return value;
136}
137
138
139// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Vector force apply to particles to avoid the crossing of particles from one phase to the other....
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.
virtual ~InterfaceForce()
Destructor.
virtual void cacheFields(const bool store)
Cache fields.
Abstract base class for particle forces.
Definition: ParticleForce.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:67
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:61
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
bool checkOut()
Remove all file watches and remove object from registry.
Definition: regIOobject.C:224
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
Calculate the gradient of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
volScalarField & alpha
dictionary dict