crankConRod.C
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) 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 #include "crankConRod.H"
29 #include "unitConversion.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(crankConRod, 0);
37  addToRunTimeSelectionTable(engineTime, crankConRod, dictionary);
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::crankConRod::timeAdjustment()
44 {
47 
48  if
49  (
52  )
53  {
55  }
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 Foam::crankConRod::crankConRod
62 (
63  const word& name,
64  const fileName& rootPath,
65  const fileName& caseName,
66  const fileName& systemName,
67  const fileName& constantName,
68  const fileName& dictName
69 )
70 :
72  (
73  name,
74  rootPath,
75  caseName,
76  systemName,
77  constantName
78  ),
79  rpm_("rpm", dimless/dimTime, dict_),
80  conRodLength_("conRodLength", dimLength, Zero),
81  bore_("bore", dimLength, Zero),
82  stroke_("stroke", dimLength, Zero),
83  clearance_("clearance", dimLength, Zero)
84 {
85  // geometric parameters are not strictly required for Time
86  dict_.readIfPresent("conRodLength", conRodLength_);
87  dict_.readIfPresent("bore", bore_);
88  dict_.readIfPresent("stroke", stroke_);
89  dict_.readIfPresent("clearance", clearance_);
90 
91  timeAdjustment();
92 
93  startTime_ = degToTime(startTime_);
94  value() = degToTime(value());
95 
96  deltaTSave_ = deltaT_;
97  deltaT0_ = deltaT_;
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  Time::readDict();
106  timeAdjustment();
107 }
108 
109 
111 {
112  if (Time::read())
113  {
114  timeAdjustment();
115  return true;
116  }
117 
118  return false;
119 }
120 
121 
122 Foam::scalar Foam::crankConRod::degToTime(const scalar theta) const
123 {
124  // 6 * rpm => deg/s
125  return theta/(6.0*rpm_.value());
126 }
127 
128 
129 Foam::scalar Foam::crankConRod::timeToDeg(const scalar t) const
130 {
131  // 6 * rpm => deg/s
132  return t*(6.0*rpm_.value());
133 }
134 
135 
136 Foam::scalar Foam::crankConRod::theta() const
137 {
138  return timeToDeg(value());
139 }
140 
141 
143 {
144  return " CAD";
145 }
146 
147 
149 {
150  scalar t = theta();
151 
152  while (t > 180.0)
153  {
154  t -= 360.0;
155  }
156 
157  while (t < -180.0)
158  {
159  t += 360.0;
160  }
161 
162  return t;
163 }
164 
165 
166 Foam::scalar Foam::crankConRod::deltaTheta() const
167 {
168  return timeToDeg(deltaTValue());
169 }
170 
171 
172 Foam::scalar Foam::crankConRod::pistonPosition(const scalar theta) const
173 {
174  return
175  (
176  conRodLength_.value()
177  + stroke_.value()/2.0
178  + clearance_.value()
179  )
180  - (
181  stroke_.value()*::cos(degToRad(theta))/2.0
182  + ::sqrt
183  (
184  sqr(conRodLength_.value())
185  - sqr(stroke_.value()*::sin(degToRad(theta))/2.0)
186  )
187  );
188 }
189 
190 
191 
192 Foam::scalar Foam::crankConRod::userTimeToTime(const scalar theta) const
193 {
194  return degToTime(theta);
195 }
196 
197 
198 Foam::scalar Foam::crankConRod::timeToUserTime(const scalar t) const
199 {
200  return timeToDeg(t);
201 }
202 
203 
204 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::crankConRod::deltaTheta
virtual scalar deltaTheta() const
Return crank-angle increment.
Definition: crankConRod.C:166
Foam::crankConRod::unit
virtual word unit() const
Return time unit.
Definition: crankConRod.C:142
dictName
const word dictName("blockMeshDict")
unitConversion.H
Unit conversion functions.
Foam::crankConRod::read
virtual bool read()
Read the controlDict and set all the parameters.
Definition: crankConRod.C:110
Foam::Time::read
virtual bool read()
Read control dictionary, update controls and time.
Definition: TimeIO.C:445
crankConRod.H
Foam::crankConRod::userTimeToTime
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (CA deg) to real-time (s).
Definition: crankConRod.C:192
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &dict)
Definition: dimensionedType.C:483
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::crankConRod::thetaRevolution
scalar thetaRevolution() const
Return current crank-angle translated to a single revolution.
Definition: crankConRod.C:148
Foam::crankConRod::theta
virtual scalar theta() const
Return current crank-angle.
Definition: crankConRod.C:136
Foam::engineTime
An abstract class for the time description of the piston motion.
Definition: engineTime.H:54
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::Time::wcAdjustableRunTime
"adjustable" / "adjustableRunTime"
Definition: Time.H:89
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::Time::readDict
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition: TimeIO.C:90
Foam::crankConRod::timeToUserTime
virtual scalar timeToUserTime(const scalar t) const
Convert the real-time (s) into user-time (CA deg)
Definition: crankConRod.C:198
Foam::crankConRod::timeToDeg
scalar timeToDeg(const scalar t) const
Convert seconds to degrees (for given engine speed in RPM)
Definition: crankConRod.C:129
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::Time::wcRunTime
"runTime"
Definition: Time.H:88
Foam::Time::writeInterval_
scalar writeInterval_
Definition: Time.H:159
Foam::Time::endTime_
scalar endTime_
Definition: Time.H:153
Foam::Time::writeControl_
writeControls writeControl_
Definition: Time.H:157
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::TimeState::deltaT_
scalar deltaT_
Definition: TimeState.H:60
Foam::crankConRod::degToTime
scalar degToTime(const scalar theta) const
Convert degrees to seconds (for given engine speed in RPM)
Definition: crankConRod.C:122
Foam::engineTime::pistonPosition
dimensionedScalar pistonPosition() const
Return current piston position.
Definition: engineTime.C:94
Foam::crankConRod::readDict
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition: crankConRod.C:103
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265