DSMCParcelI.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 -------------------------------------------------------------------------------
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 "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class ParcelType>
34 :
35  mass_(0),
36  d_(0)
37 {}
38 
39 
40 template<class ParcelType>
42 (
43  const dictionary& dict
44 )
45 :
46  mass_(dict.get<scalar>("mass")),
47  d_(dict.get<scalar>("diameter")),
48  internalDegreesOfFreedom_(dict.get<int>("internalDegreesOfFreedom")),
49  omega_(dict.get<scalar>("omega"))
50 {}
51 
52 
53 template<class ParcelType>
55 (
56  const polyMesh& mesh,
57  const barycentric& coordinates,
58  const label celli,
59  const label tetFacei,
60  const label tetPti,
61  const vector& U,
62  const scalar Ei,
63  const label typeId
64 )
65 :
66  ParcelType(mesh, coordinates, celli, tetFacei, tetPti),
67  U_(U),
68  Ei_(Ei),
69  typeId_(typeId)
70 {}
71 
72 
73 template<class ParcelType>
75 (
76  const polyMesh& mesh,
77  const vector& position,
78  const label celli,
79  const vector& U,
80  const scalar Ei,
81  const label typeId
82 )
83 :
84  ParcelType(mesh, position, celli),
85  U_(U),
86  Ei_(Ei),
87  typeId_(typeId)
88 {}
89 
90 
91 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
92 
93 template<class ParcelType>
94 inline Foam::scalar
96 {
97  return mass_;
98 }
99 
100 
101 template<class ParcelType>
102 inline Foam::scalar Foam::DSMCParcel<ParcelType>::constantProperties::d() const
103 {
104  return d_;
105 }
106 
107 
108 template<class ParcelType>
109 inline Foam::scalar
111 {
112  return constant::mathematical::pi*d_*d_;
113 }
114 
115 
116 template<class ParcelType>
117 inline Foam::direction
119 const
120 {
121  return internalDegreesOfFreedom_;
122 }
123 
124 
125 template<class ParcelType>
126 inline Foam::scalar
128 {
129  return omega_;
130 }
131 
132 
133 // * * * * * * * * * * DSMCParcel Member Functions * * * * * * * * * * //
134 
135 template<class ParcelType>
136 inline Foam::label Foam::DSMCParcel<ParcelType>::typeId() const
137 {
138  return typeId_;
139 }
140 
141 
142 template<class ParcelType>
144 {
145  return U_;
146 }
147 
148 
149 template<class ParcelType>
150 inline Foam::scalar Foam::DSMCParcel<ParcelType>::Ei() const
151 {
152  return Ei_;
153 }
154 
155 
156 template<class ParcelType>
158 {
159  return U_;
160 }
161 
162 
163 template<class ParcelType>
164 inline Foam::scalar& Foam::DSMCParcel<ParcelType>::Ei()
165 {
166  return Ei_;
167 }
168 
169 
170 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
uint8_t direction
Definition: direction.H:46
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label typeId_
Parcel type id.
Definition: DSMCParcel.H:180
scalar sigmaT() const
Return the reference total collision cross section.
Definition: DSMCParcelI.H:103
constantProperties()
Null constructor, allows List of constantProperties to be.
Definition: DSMCParcelI.H:26
scalar omega() const
Return the viscosity index.
Definition: DSMCParcelI.H:120
dynamicFvMesh & mesh
direction internalDegreesOfFreedom() const
Return the internalDegreesOfFreedom.
Definition: DSMCParcelI.H:111
scalar mass() const
Return const access to the particle mass [kg].
Definition: DSMCParcelI.H:88
constexpr scalar pi(M_PI)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
scalar d() const
Return const access to the hard sphere diameter [m].
Definition: DSMCParcelI.H:95
U
Definition: pEqn.H:72
PtrList< coordinateSystem > coordinates(solidRegions.size())
scalar Ei_
Internal energy of the Parcel, covering all non-translational.
Definition: DSMCParcel.H:175
scalar Ei() const
Return const access to internal energy.
Definition: DSMCParcelI.H:143
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
vector U_
Velocity of Parcel [m/s].
Definition: DSMCParcel.H:168
label typeId() const
Return type id.
Definition: DSMCParcelI.H:129
DSMCParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const vector &U, const scalar Ei, const label typeId)
Construct from components.
Definition: DSMCParcelI.H:48
const vector & U() const
Return const access to velocity.
Definition: DSMCParcelI.H:136