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-2017 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 barycentric& coordinates,
67  const label celli,
68  const label tetFacei,
69  const label tetPti
70 )
71 :
72  ParcelType(owner, coordinates, 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 )
87 :
88  ParcelType(owner, position, celli),
89  f_(Zero),
91  torque_(Zero),
93 {}
94 
95 
96 template<class ParcelType>
98 (
99  const polyMesh& owner,
100  const barycentric& coordinates,
101  const label celli,
102  const label tetFacei,
103  const label tetPti,
104  const label typeId,
105  const scalar nParticle0,
106  const scalar d0,
107  const scalar dTarget0,
108  const vector& U0,
109  const vector& f0,
110  const vector& angularMomentum0,
111  const vector& torque0,
112  const typename ParcelType::constantProperties& constProps
113 )
114 :
115  ParcelType
116  (
117  owner,
118  coordinates,
119  celli,
120  tetFacei,
121  tetPti,
122  typeId,
123  nParticle0,
124  d0,
125  dTarget0,
126  torque0,
127  constProps
128  ),
129  f_(f0),
130  angularMomentum_(angularMomentum0),
131  torque_(torque0),
133 {}
134 
135 
136 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
137 
138 template<class ParcelType>
139 inline Foam::scalar
141 {
142  return youngsModulus_.value();
143 }
144 
145 
146 template<class ParcelType>
147 inline Foam::scalar
149 {
150  return poissonsRatio_.value();
151 }
152 
153 
154 // * * * * * * * * * * CollidingParcel Member Functions * * * * * * * * * * //
155 
156 template<class ParcelType>
158 {
159  return f_;
160 }
161 
162 
163 template<class ParcelType>
164 inline const Foam::vector&
166 {
167  return angularMomentum_;
168 }
169 
170 
171 template<class ParcelType>
173 {
174  return torque_;
175 }
176 
177 
178 template<class ParcelType>
179 inline const Foam::collisionRecordList&
181 {
182  return collisionRecords_;
183 }
184 
185 
186 template<class ParcelType>
188 {
189  return f_;
190 }
191 
192 
193 template<class ParcelType>
195 {
196  return angularMomentum_;
197 }
198 
199 
200 template<class ParcelType>
202 {
203  return torque_;
204 }
205 
206 
207 template<class ParcelType>
210 {
211  return collisionRecords_;
212 }
213 
214 
215 template<class ParcelType>
217 {
218  return angularMomentum_/this->momentOfInertia();
219 }
220 
221 
222 // ************************************************************************* //
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
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 vector & f() const
Return const access to force.
scalar poissonsRatio() const
Return const access to Poisson&#39;s ratio.
scalar youngsModulus() const
Return const access to Young&#39;s Modulus.
collisionRecordList collisionRecords_
Particle collision records.
const Type & value() const
Return the value.
static const zero Zero
Definition: zero.H:91
const vector & torque() const
Return const access to torque.
const vector & angularMomentum() const
Return const access to angular momentum.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
vector torque_
Torque on particle due to collisions in global.
vector omega() const
Particle angular velocity.
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
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.
const collisionRecordList & collisionRecords() const
Return const access to the collision records.