Loading...
Searching...
No Matches
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\*---------------------------------------------------------------------------*/
28#include "InterfaceForce.H"
29#include "fvcGrad.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CloudType>
35(
37 const fvMesh& mesh,
38 const dictionary& dict
39)
40:
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 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,
86 );
87
88 resultPtr->store();
89 }
90
91 const volVectorField& gradInterForce = *resultPtr;
92
93 gradInterForceInterpPtr_.reset
94 (
95 interpolation<vector>::New
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
114template<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// ************************************************************************* //
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.
InterfaceForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
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.
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
Abstract base class for volume field interpolation.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
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:60
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 class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
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
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
volScalarField & alpha
dictionary dict