CollidingParcelI.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-2020 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 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
97 
98 template<class ParcelType>
99 inline Foam::scalar
101 {
102  return youngsModulus_.value();
103 }
104 
105 
106 template<class ParcelType>
107 inline Foam::scalar
109 {
110  return poissonsRatio_.value();
111 }
112 
113 
114 // * * * * * * * * * * CollidingParcel Member Functions * * * * * * * * * * //
115 
116 template<class ParcelType>
118 {
119  return f_;
120 }
121 
122 
123 template<class ParcelType>
124 inline const Foam::vector&
126 {
127  return angularMomentum_;
128 }
129 
130 
131 template<class ParcelType>
133 {
134  return torque_;
135 }
136 
137 
138 template<class ParcelType>
139 inline const Foam::collisionRecordList&
141 {
142  return collisionRecords_;
143 }
144 
145 
146 template<class ParcelType>
148 {
149  return f_;
150 }
151 
152 
153 template<class ParcelType>
155 {
156  return angularMomentum_;
157 }
158 
159 
160 template<class ParcelType>
162 {
163  return torque_;
164 }
165 
166 
167 template<class ParcelType>
170 {
171  return collisionRecords_;
172 }
173 
174 
175 template<class ParcelType>
177 {
178  return angularMomentum_/this->momentOfInertia();
179 }
180 
181 
182 // ************************************************************************* //
vector f_
Force on particle due to collisions [N].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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:97
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:76
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.