solidParticleI.H
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 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30
32(
33 const solidParticleCloud& spc,
34 const interpolationCellPoint<scalar>& rhoInterp,
35 const interpolationCellPoint<vector>& UInterp,
36 const interpolationCellPoint<scalar>& nuInterp,
37 const vector& g
38)
39:
41 rhoInterp_(rhoInterp),
42 UInterp_(UInterp),
43 nuInterp_(nuInterp),
44 g_(g)
45{}
46
47
49(
50 const polyMesh& mesh,
51 const vector& position,
52 const label celli
53)
54:
55 particle(mesh, position, celli),
56 d_(0),
57 U_(Zero)
58{}
59
60
62(
63 const polyMesh& mesh,
65 const label celli,
66 const label tetFacei,
67 const label tetPti,
68 const scalar d,
69 const vector& U
70)
71:
72 particle(mesh, coordinates, celli, tetFacei, tetPti),
73 d_(d),
74 U_(U)
75{}
76
77
78// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79
82{
83 return rhoInterp_;
84}
85
86
89{
90 return UInterp_;
91}
92
93
96{
97 return nuInterp_;
98}
99
101{
102 return g_;
103}
104
105
106inline Foam::scalar Foam::solidParticle::d() const
107{
108 return d_;
109}
110
111
113{
114 return U_;
115}
116
117
118// ************************************************************************* //
const uniformDimensionedVectorField & g
Class used to pass data into container.
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Base particle class.
Definition: particle.H:79
vector position() const
Return current particle position.
Definition: particleI.H:314
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:137
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A Cloud of solid particles.
Class used to pass tracking data to the trackToFace function.
Definition: solidParticle.H:86
const interpolationCellPoint< scalar > & rhoInterp() const
const interpolationCellPoint< vector > & UInterp() const
const interpolationCellPoint< scalar > & nuInterp() const
Simple solid spherical particle class with one-way coupling with the continuous phase.
Definition: solidParticle.H:68
const vector & U() const
Return velocity.
scalar d() const
Return diameter.
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131