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-2022 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& mesh,
66  const barycentric& coordinates,
67  const label celli,
68  const label tetFacei,
69  const label tetPti,
70  const label facei
71 )
72 :
73  ParcelType(mesh, coordinates, celli, tetFacei, tetPti, facei),
74  f_(Zero),
76  torque_(Zero),
78 {}
79 
80 
81 template<class ParcelType>
83 (
84  const polyMesh& mesh,
85  const vector& position,
86  const label celli
87 )
88 :
89  ParcelType(mesh, position, celli),
90  f_(Zero),
91  angularMomentum_(Zero),
92  torque_(Zero),
93  collisionRecords_()
94 {}
95 
96 
97 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
98 
99 template<class ParcelType>
100 inline Foam::scalar
102 {
103  return youngsModulus_.value();
104 }
105 
106 
107 template<class ParcelType>
108 inline Foam::scalar
110 {
111  return poissonsRatio_.value();
112 }
113 
114 
115 // * * * * * * * * * * CollidingParcel Member Functions * * * * * * * * * * //
116 
117 template<class ParcelType>
119 {
120  return f_;
121 }
122 
123 
124 template<class ParcelType>
125 inline const Foam::vector&
127 {
128  return angularMomentum_;
129 }
130 
131 
132 template<class ParcelType>
134 {
135  return torque_;
136 }
137 
138 
139 template<class ParcelType>
140 inline const Foam::collisionRecordList&
142 {
143  return collisionRecords_;
144 }
145 
146 
147 template<class ParcelType>
149 {
150  return f_;
151 }
152 
153 
154 template<class ParcelType>
156 {
157  return angularMomentum_;
158 }
159 
160 
161 template<class ParcelType>
163 {
164  return torque_;
165 }
166 
167 
168 template<class ParcelType>
171 {
172  return collisionRecords_;
173 }
174 
175 
176 template<class ParcelType>
178 {
179  return angularMomentum_/this->momentOfInertia();
180 }
181 
182 
183 // ************************************************************************* //
Class to hold thermo particle constant properties.
scalar youngsModulus() const
Return const access to Young's Modulus.
scalar poissonsRatio() const
Return const access to Poisson's ratio.
const vector & angularMomentum() const
Return const access to angular momentum.
vector f_
Force on particle due to collisions [N].
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
vector omega() const
Particle angular velocity.
const vector & f() const
Return const access to force.
collisionRecordList collisionRecords_
Particle collision records.
const vector & torque() const
Return const access to torque.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
vector torque_
Torque on particle due to collisions in global.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Calculates the inertia tensor and principal axes and moments of a polyhedra/cells/triSurfaces....
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static const zero Zero
Definition: zero.H:97
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
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753