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