solidParticleI.H
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-2018 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 (
30  const solidParticleCloud& spc,
31  const interpolationCellPoint<scalar>& rhoInterp,
32  const interpolationCellPoint<vector>& UInterp,
33  const interpolationCellPoint<scalar>& nuInterp,
34  const vector& g
35 )
36 :
38  rhoInterp_(rhoInterp),
39  UInterp_(UInterp),
40  nuInterp_(nuInterp),
41  g_(g)
42 {}
43 
44 
46 (
47  const polyMesh& mesh,
48  const barycentric& coordinates,
49  const label celli,
50  const label tetFacei,
51  const label tetPti,
52  const scalar d,
53  const vector& U
54 )
55 :
56  particle(mesh, coordinates, celli, tetFacei, tetPti),
57  d_(d),
58  U_(U)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
66 {
67  return rhoInterp_;
68 }
69 
70 
73 {
74  return UInterp_;
75 }
76 
77 
80 {
81  return nuInterp_;
82 }
83 
85 {
86  return g_;
87 }
88 
89 
90 inline Foam::scalar Foam::solidParticle::d() const
91 {
92  return d_;
93 }
94 
95 
96 inline const Foam::vector& Foam::solidParticle::U() const
97 {
98  return U_;
99 }
100 
101 
102 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
trackingData(const solidParticleCloud &spc, const interpolationCellPoint< scalar > &rhoInterp, const interpolationCellPoint< vector > &UInterp, const interpolationCellPoint< scalar > &nuInterp, const vector &g)
scalar d() const
Return diameter.
Base particle class.
Definition: particle.H:84
const vector & U() const
Return velocity.
A Cloud of solid particles.
const interpolationCellPoint< scalar > & rhoInterp() const
const interpolationCellPoint< scalar > & nuInterp() const
const interpolationCellPoint< vector > & UInterp() const
solidParticle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const scalar d, const vector &U)
Construct from components.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74