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 | 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ParcelType>
33 :
34  ParcelType::constantProperties(),
35  youngsModulus_(this->dict_, 0.0),
36  poissonsRatio_(this->dict_, 0.0)
37 {}
38 
39 
40 template<class ParcelType>
42 (
43  const constantProperties& cp
44 )
45 :
46  ParcelType::constantProperties(cp),
47  youngsModulus_(cp.youngsModulus_),
48  poissonsRatio_(cp.poissonsRatio_)
49 {}
50 
51 
52 template<class ParcelType>
54 (
55  const dictionary& parentDict
56 )
57 :
58  ParcelType::constantProperties(parentDict),
59  youngsModulus_(this->dict_, "youngsModulus"),
60  poissonsRatio_(this->dict_, "poissonsRatio")
61 {}
62 
63 
64 template<class ParcelType>
66 (
67  const polyMesh& owner,
68  const barycentric& coordinates,
69  const label celli,
70  const label tetFacei,
71  const label tetPti
72 )
73 :
74  ParcelType(owner, coordinates, celli, tetFacei, tetPti),
75  f_(Zero),
77  torque_(Zero),
79 {}
80 
81 
82 template<class ParcelType>
84 (
85  const polyMesh& owner,
86  const vector& position,
87  const label celli
88 )
89 :
90  ParcelType(owner, position, celli),
91  f_(Zero),
93  torque_(Zero),
95 {}
96 
97 
98 template<class ParcelType>
100 (
101  const polyMesh& owner,
102  const barycentric& coordinates,
103  const label celli,
104  const label tetFacei,
105  const label tetPti,
106  const label typeId,
107  const scalar nParticle0,
108  const scalar d0,
109  const scalar dTarget0,
110  const vector& U0,
111  const vector& f0,
112  const vector& angularMomentum0,
113  const vector& torque0,
114  const typename ParcelType::constantProperties& constProps
115 )
116 :
117  ParcelType
118  (
119  owner,
120  coordinates,
121  celli,
122  tetFacei,
123  tetPti,
124  typeId,
125  nParticle0,
126  d0,
127  dTarget0,
128  U0,
129  constProps
130  ),
131  f_(f0),
132  angularMomentum_(angularMomentum0),
133  torque_(torque0),
135 {}
136 
137 
138 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
139 
140 template<class ParcelType>
141 inline Foam::scalar
143 {
144  return youngsModulus_.value();
145 }
146 
147 
148 template<class ParcelType>
149 inline Foam::scalar
151 {
152  return poissonsRatio_.value();
153 }
154 
155 
156 // * * * * * * * * * * CollidingParcel Member Functions * * * * * * * * * * //
157 
158 template<class ParcelType>
160 {
161  return f_;
162 }
163 
164 
165 template<class ParcelType>
166 inline const Foam::vector&
168 {
169  return angularMomentum_;
170 }
171 
172 
173 template<class ParcelType>
175 {
176  return torque_;
177 }
178 
179 
180 template<class ParcelType>
181 inline const Foam::collisionRecordList&
183 {
184  return collisionRecords_;
185 }
186 
187 
188 template<class ParcelType>
190 {
191  return f_;
192 }
193 
194 
195 template<class ParcelType>
197 {
198  return angularMomentum_;
199 }
200 
201 
202 template<class ParcelType>
204 {
205  return torque_;
206 }
207 
208 
209 template<class ParcelType>
212 {
213  return collisionRecords_;
214 }
215 
216 
217 template<class ParcelType>
219 {
220  return angularMomentum_/this->momentOfInertia();
221 }
222 
223 
224 // ************************************************************************* //
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::CollidingParcel::constantProperties::constantProperties
constantProperties()
Null constructor.
Definition: CollidingParcelI.H:32
Foam::CollidingParcel::collisionRecords_
collisionRecordList collisionRecords_
Particle collision records.
Definition: CollidingParcel.H:145
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::CollidingParcel::constantProperties
Class to hold thermo particle constant properties.
Definition: CollidingParcel.H:88
Foam::CollidingParcel::torque
const vector & torque() const
Return const access to torque.
Definition: CollidingParcelI.H:174
Foam::CollisionRecordList< vector, vector >
Foam::CollidingParcel::constantProperties::youngsModulus
scalar youngsModulus() const
Return const access to Young's Modulus.
Definition: CollidingParcelI.H:142
Foam::CollidingParcel::omega
vector omega() const
Particle angular velocity.
Definition: CollidingParcelI.H:218
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::CollidingParcel::angularMomentum
const vector & angularMomentum() const
Return const access to angular momentum.
Definition: CollidingParcelI.H:167
Foam::Barycentric< scalar >
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::CollidingParcel::constantProperties::poissonsRatio
scalar poissonsRatio() const
Return const access to Poisson's ratio.
Definition: CollidingParcelI.H:150
Foam::CollidingParcel::collisionRecords
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
Definition: CollidingParcelI.H:182
Foam::momentOfInertia
Calculates the inertia tensor and principal axes and moments of a polyhedra/cells/triSurfaces....
Definition: momentOfInertia.H:65
Foam::cp
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: MSwindows.C:802
Foam::Vector< scalar >
Foam::CollidingParcel::torque_
vector torque_
Torque on particle due to collisions in global.
Definition: CollidingParcel.H:142
Foam::CollidingParcel::f_
vector f_
Force on particle due to collisions [N].
Definition: CollidingParcel.H:135
Foam::CollidingParcel::angularMomentum_
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
Definition: CollidingParcel.H:138
Foam::CollidingParcel::f
const vector & f() const
Return const access to force.
Definition: CollidingParcelI.H:159
Foam::CollidingParcel::CollidingParcel
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: CollidingParcelI.H:66