Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39}
40
41
42// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43
44void Foam::crankConRod::timeAdjustment()
45{
47
48 if
49 (
52 )
53 {
55 }
56}
57
58
59// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60
61Foam::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:
71 engineTime
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
94 value() = degToTime(value());
95
99}
100
101
102// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103
105{
107 timeAdjustment();
108}
109
110
112{
113 if (Time::read())
114 {
115 timeAdjustment();
116 return true;
117 }
118
119 return false;
120}
121
122
123Foam::scalar Foam::crankConRod::degToTime(const scalar theta) const
124{
125 // 6 * rpm => deg/s
126 return theta/(6.0*rpm_.value());
127}
128
129
130Foam::scalar Foam::crankConRod::timeToDeg(const scalar t) const
131{
132 // 6 * rpm => deg/s
133 return t*(6.0*rpm_.value());
134}
135
137Foam::scalar Foam::crankConRod::theta() const
138{
139 return timeToDeg(value());
140}
141
144{
145 return " CAD";
146}
147
148
149Foam::scalar Foam::crankConRod::thetaRevolution() const
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
167Foam::scalar Foam::crankConRod::deltaTheta() const
168{
169 return timeToDeg(deltaTValue());
170}
171
172
173Foam::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
193Foam::scalar Foam::crankConRod::userTimeToTime(const scalar theta) const
194{
195 return degToTime(theta);
196}
197
198
199Foam::scalar Foam::crankConRod::timeToUserTime(const scalar t) const
200{
201 return timeToDeg(t);
202}
203
204
205// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
scalar deltaTSave_
Definition TimeState.H:69
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
scalar deltaT0_
Definition TimeState.H:68
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition TimeIO.C:84
@ wcAdjustableRunTime
"adjustable" / "adjustableRunTime"
Definition Time.H:86
@ wcRunTime
"runTime"
Definition Time.H:85
scalar writeInterval_
Definition Time.H:163
const fileName & rootPath() const noexcept
The root path.
Definition TimePathsI.H:66
writeControls writeControl_
Definition Time.H:161
scalar endTime_
Definition Time.H:157
scalar startTime_
Definition Time.H:155
virtual bool read()
Read control dictionary, update controls and time.
Definition TimeIO.C:434
const fileName & caseName() const noexcept
The case name.
Definition TimePathsI.H:78
Manage time in terms of engine RPM and crank-angle.
Definition crankConRod.H:67
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition crankConRod.C:97
virtual scalar timeToUserTime(const scalar t) const
Convert the real-time (s) into user-time (CA deg).
virtual scalar theta() const
Return current crank-angle.
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (CA deg) to real-time (s).
scalar thetaRevolution() const
Return current crank-angle translated to a single revolution.
scalar degToTime(const scalar theta) const
Convert degrees to seconds (for given engine speed in RPM).
scalar timeToDeg(const scalar t) const
Convert seconds to degrees (for given engine speed in RPM).
virtual scalar deltaTheta() const
Return crank-angle increment.
virtual word unit() const
Return time unit.
virtual bool read()
Read the controlDict and set all the parameters.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
const scalar & value() const noexcept
An abstract class for the time description of the piston motion.
Definition engineTime.H:53
engineTime(const word &name, const fileName &rootPath, const fileName &caseName, const fileName &systemName="system", const fileName &constantName="constant", const fileName &dictName="engineGeometry")
Construct from objectRegistry arguments.
Definition engineTime.C:36
const IOdictionary dict_
Definition engineTime.H:57
dimensionedScalar pistonPosition() const
Return current piston position.
Definition engineTime.C:87
A class for handling file names.
Definition fileName.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const word dictName("faMeshDefinition")
auto & name
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
dimensionedScalar sin(const dimensionedScalar &ds)
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensionedScalar cos(const dimensionedScalar &ds)
Unit conversion functions.