CollidingParcelI.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) 2011-2016 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 
28 template<class ParcelType>
31 :
32  ParcelType::constantProperties(),
33  youngsModulus_(this->dict_, 0.0),
34  poissonsRatio_(this->dict_, 0.0)
35 {}
36 
37 
38 template<class ParcelType>
40 (
41  const constantProperties& cp
42 )
43 :
44  ParcelType::constantProperties(cp),
45  youngsModulus_(cp.youngsModulus_),
46  poissonsRatio_(cp.poissonsRatio_)
47 {}
48 
49 
50 template<class ParcelType>
52 (
53  const dictionary& parentDict
54 )
55 :
56  ParcelType::constantProperties(parentDict),
57  youngsModulus_(this->dict_, "youngsModulus"),
58  poissonsRatio_(this->dict_, "poissonsRatio")
59 {}
60 
61 
62 template<class ParcelType>
64 (
65  const polyMesh& owner,
66  const vector& position,
67  const label celli,
68  const label tetFacei,
69  const label tetPtI
70 )
71 :
72  ParcelType(owner, position, celli, tetFacei, tetPtI),
73  f_(Zero),
75  torque_(Zero),
77 {}
78 
79 
80 template<class ParcelType>
82 (
83  const polyMesh& owner,
84  const vector& position,
85  const label celli,
86  const label tetFacei,
87  const label tetPtI,
88  const label typeId,
89  const scalar nParticle0,
90  const scalar d0,
91  const scalar dTarget0,
92  const vector& U0,
93  const vector& f0,
94  const vector& angularMomentum0,
95  const vector& torque0,
96  const typename ParcelType::constantProperties& constProps
97 )
98 :
99  ParcelType
100  (
101  owner,
102  position,
103  celli,
104  tetFacei,
105  tetPtI,
106  typeId,
107  nParticle0,
108  d0,
109  dTarget0,
110  torque0,
111  constProps
112  ),
113  f_(f0),
114  angularMomentum_(angularMomentum0),
115  torque_(torque0),
117 {}
118 
119 
120 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
121 
122 template<class ParcelType>
123 inline Foam::scalar
125 {
126  return youngsModulus_.value();
127 }
128 
129 
130 template<class ParcelType>
131 inline Foam::scalar
133 {
134  return poissonsRatio_.value();
135 }
136 
137 
138 // * * * * * * * * * * CollidingParcel Member Functions * * * * * * * * * * //
139 
140 template<class ParcelType>
142 {
143  return f_;
144 }
145 
146 
147 template<class ParcelType>
148 inline const Foam::vector&
150 {
151  return angularMomentum_;
152 }
153 
154 
155 template<class ParcelType>
157 {
158  return torque_;
159 }
160 
161 
162 template<class ParcelType>
163 inline const Foam::collisionRecordList&
165 {
166  return collisionRecords_;
167 }
168 
169 
170 template<class ParcelType>
172 {
173  return f_;
174 }
175 
176 
177 template<class ParcelType>
179 {
180  return angularMomentum_;
181 }
182 
183 
184 template<class ParcelType>
186 {
187  return torque_;
188 }
189 
190 
191 template<class ParcelType>
194 {
195  return collisionRecords_;
196 }
197 
198 
199 template<class ParcelType>
201 {
202  return angularMomentum_/this->momentOfInertia();
203 }
204 
205 
206 // ************************************************************************* //
const vector & torque() const
Return const access to torque.
vector f_
Force on particle due to collisions [N].
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
scalar poissonsRatio() const
Return const access to Poisson&#39;s ratio.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Class to hold thermo particle constant properties.
const Type & value() const
Return the value.
const vector & f() const
Return const access to force.
collisionRecordList collisionRecords_
Particle collision records.
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
static const zero Zero
Definition: zero.H:91
CollidingParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
vector omega() const
Particle angular velocity.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
vector torque_
Torque on particle due to collisions in global.
const vector & angularMomentum() const
Return const access to angular momentum.
scalar youngsModulus() const
Return const access to Young&#39;s Modulus.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Calculates the inertia tensor and principal axes and moments of a polyhedra/cells/triSurfaces. Inertia can either be of the solid body or of a thin shell.