FreeFOAM The Cross-Platform CFD Toolkit
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 | Copyright (C) 2009-2010 OpenCFD Ltd.
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 
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_(readScalar(dict.lookup("mass"))),
45  d_(readScalar(dict.lookup("diameter"))),
46  internalDegreesOfFreedom_
47  (
48  readScalar(dict.lookup("internalDegreesOfFreedom"))
49  ),
50  omega_(readScalar(dict.lookup("omega")))
51 {}
52 
53 
54 template <class ParcelType>
56 (
58 )
59 :
61  cloud_(cloud)
62 {}
63 
64 
65 template <class ParcelType>
67 (
68  DsmcCloud<ParcelType>& owner,
69  const vector& position,
70  const vector& U,
71  const scalar Ei,
72  const label celli,
73  const label typeId
74 )
75 :
76  Particle<ParcelType>(owner, position, celli),
77  U_(U),
78  Ei_(Ei),
79  typeId_(typeId)
80 {}
81 
82 
83 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
84 
85 template <class ParcelType>
86 inline Foam::scalar
88 {
89  return mass_;
90 }
91 
92 
93 template <class ParcelType>
94 inline Foam::scalar
96 {
97  return d_;
98 }
99 
100 
101 template <class ParcelType>
102 inline Foam::scalar
104 {
105  return mathematicalConstant::pi*d_*d_;
106 }
107 
108 
109 template <class ParcelType>
110 inline Foam::scalar
112 const
113 {
114  return internalDegreesOfFreedom_;
115 }
116 
117 
118 template <class ParcelType>
119 inline Foam::scalar
121 {
122  return omega_;
123 }
124 
125 
126 // * * * * * * * * * * * trackData Member Functions * * * * * * * * * * * * //
127 
128 template <class ParcelType>
131 {
132  return cloud_;
133 }
134 
135 
136 // * * * * * * * * * * DsmcParcel Member Functions * * * * * * * * * * //
137 
138 template <class ParcelType>
139 inline Foam::label Foam::DsmcParcel<ParcelType>::typeId() const
140 {
141  return typeId_;
142 }
143 
144 
145 template <class ParcelType>
147 {
148  return U_;
149 }
150 
151 
152 template <class ParcelType>
153 inline Foam::scalar Foam::DsmcParcel<ParcelType>::Ei() const
154 {
155  return Ei_;
156 }
157 
158 
159 template <class ParcelType>
161 {
162  return U_;
163 }
164 
165 
166 template <class ParcelType>
167 inline Foam::scalar& Foam::DsmcParcel<ParcelType>::Ei()
168 {
169  return Ei_;
170 }
171 
172 
173 // ************************ vim: set sw=4 sts=4 et: ************************ //