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  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "crankConRod.H"
30 #include "unitConversion.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(crankConRod, 0);
38  addToRunTimeSelectionTable(engineTime, crankConRod, dictionary);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::crankConRod::timeAdjustment()
45 {
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  deltaT_ = degToTime(deltaT_);
97  deltaTSave_ = deltaT_;
98  deltaT0_ = deltaT_;
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 {
106  Time::readDict();
107  timeAdjustment();
108 }
109 
110 
112 {
113  if (Time::read())
114  {
115  timeAdjustment();
116  return true;
117  }
118 
119  return false;
120 }
121 
122 
123 Foam::scalar Foam::crankConRod::degToTime(const scalar theta) const
124 {
125  // 6 * rpm => deg/s
126  return theta/(6.0*rpm_.value());
127 }
128 
129 
130 Foam::scalar Foam::crankConRod::timeToDeg(const scalar t) const
131 {
132  // 6 * rpm => deg/s
133  return t*(6.0*rpm_.value());
134 }
135 
136 
137 Foam::scalar Foam::crankConRod::theta() const
138 {
139  return timeToDeg(value());
140 }
141 
142 
144 {
145  return " CAD";
146 }
147 
148 
150 {
151  scalar t = theta();
152 
153  while (t > 180.0)
154  {
155  t -= 360.0;
156  }
157 
158  while (t < -180.0)
159  {
160  t += 360.0;
161  }
162 
163  return t;
164 }
165 
166 
167 Foam::scalar Foam::crankConRod::deltaTheta() const
168 {
169  return timeToDeg(deltaTValue());
170 }
171 
172 
173 Foam::scalar Foam::crankConRod::pistonPosition(const scalar theta) const
174 {
175  return
176  (
177  conRodLength_.value()
178  + stroke_.value()/2.0
179  + clearance_.value()
180  )
181  - (
182  stroke_.value()*::cos(degToRad(theta))/2.0
183  + ::sqrt
184  (
185  sqr(conRodLength_.value())
186  - sqr(stroke_.value()*::sin(degToRad(theta))/2.0)
187  )
188  );
189 }
190 
191 
192 
193 Foam::scalar Foam::crankConRod::userTimeToTime(const scalar theta) const
194 {
195  return degToTime(theta);
196 }
197 
198 
199 Foam::scalar Foam::crankConRod::timeToUserTime(const scalar t) const
200 {
201  return timeToDeg(t);
202 }
203 
204 
205 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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:167
Foam::crankConRod::unit
virtual word unit() const
Return time unit.
Definition: crankConRod.C:143
dictName
const word dictName("faMeshDefinition")
unitConversion.H
Unit conversion functions.
Foam::crankConRod::read
virtual bool read()
Read the controlDict and set all the parameters.
Definition: crankConRod.C:111
Foam::Time::read
virtual bool read()
Read control dictionary, update controls and time.
Definition: TimeIO.C:442
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:193
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &dict)
Definition: dimensionedType.C:483
Foam::crankConRod::thetaRevolution
scalar thetaRevolution() const
Return current crank-angle translated to a single revolution.
Definition: crankConRod.C:149
Foam::crankConRod::theta
virtual scalar theta() const
Return current crank-angle.
Definition: crankConRod.C:137
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:199
Foam::crankConRod::timeToDeg
scalar timeToDeg(const scalar t) const
Convert seconds to degrees (for given engine speed in RPM)
Definition: crankConRod.C:130
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::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::Time::endTime_
scalar endTime_
Definition: Time.H:153
Foam::Time::writeControl_
writeControls writeControl_
Definition: Time.H:157
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::crankConRod::degToTime
scalar degToTime(const scalar theta) const
Convert degrees to seconds (for given engine speed in RPM)
Definition: crankConRod.C:123
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:104
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265