engineTime.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 | Copyright (C) 2011-2016 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 #include "engineTime.H"
27 #include "unitConversion.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 void Foam::engineTime::timeAdjustment()
32 {
35 
36  if
37  (
40  )
41  {
43  }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 Foam::engineTime::engineTime
50 (
51  const word& name,
52  const fileName& rootPath,
53  const fileName& caseName,
54  const fileName& systemName,
55  const fileName& constantName,
56  const fileName& dictName
57 )
58 :
59  Time
60  (
61  name,
62  rootPath,
63  caseName,
64  systemName,
65  constantName
66  ),
67  dict_
68  (
69  IOobject
70  (
71  "engineGeometry",
72  constant(),
73  *this,
76  false
77  )
78  ),
79  rpm_(dict_.lookup("rpm")),
80  conRodLength_(dimensionedScalar("conRodLength", dimLength, 0)),
81  bore_(dimensionedScalar("bore", dimLength, 0)),
82  stroke_(dimensionedScalar("stroke", dimLength, 0)),
83  clearance_(dimensionedScalar("clearance", dimLength, 0))
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());
96  deltaT0_ = deltaT_;
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 // Read the controlDict and set all the parameters
104 {
105  Time::readDict();
106  timeAdjustment();
107 }
108 
109 
110 // Read the controlDict and set all the parameters
112 {
113  if (Time::read())
114  {
115  timeAdjustment();
116  return true;
117  }
118  else
119  {
120  return false;
121  }
122 }
123 
124 
125 Foam::scalar Foam::engineTime::degToTime(const scalar theta) const
126 {
127  // 6 * rpm => deg/s
128  return theta/(6.0*rpm_.value());
129 }
130 
131 
132 Foam::scalar Foam::engineTime::timeToDeg(const scalar t) const
133 {
134  // 6 * rpm => deg/s
135  return t*(6.0*rpm_.value());
136 }
137 
138 
139 Foam::scalar Foam::engineTime::theta() const
140 {
141  return timeToDeg(value());
142 }
143 
144 
145 // Return current crank-angle translated to a single revolution
146 // (value between -180 and 180 with 0 = top dead centre)
148 {
149  scalar t = theta();
150 
151  while (t > 180.0)
152  {
153  t -= 360.0;
154  }
155 
156  while (t < -180.0)
157  {
158  t += 360.0;
159  }
160 
161  return t;
162 }
163 
164 
165 Foam::scalar Foam::engineTime::deltaTheta() const
166 {
167  return timeToDeg(deltaTValue());
168 }
169 
170 
171 Foam::scalar Foam::engineTime::pistonPosition(const scalar theta) const
172 {
173  return
174  (
175  conRodLength_.value()
176  + stroke_.value()/2.0
177  + clearance_.value()
178  )
179  - (
180  stroke_.value()*::cos(degToRad(theta))/2.0
181  + ::sqrt
182  (
183  sqr(conRodLength_.value())
184  - sqr(stroke_.value()*::sin(degToRad(theta))/2.0)
185  )
186  );
187 }
188 
189 
191 {
192  return dimensionedScalar
193  (
194  "pistonPosition",
195  dimLength,
197  );
198 }
199 
200 
202 {
203  return dimensionedScalar
204  (
205  "pistonDisplacement",
206  dimLength,
208  );
209 }
210 
211 
213 {
214  return dimensionedScalar
215  (
216  "pistonSpeed",
217  dimVelocity,
218  pistonDisplacement().value()/(deltaTValue() + VSMALL)
219  );
220 }
221 
222 
223 Foam::scalar Foam::engineTime::userTimeToTime(const scalar theta) const
224 {
225  return degToTime(theta);
226 }
227 
228 
229 Foam::scalar Foam::engineTime::timeToUserTime(const scalar t) const
230 {
231  return timeToDeg(t);
232 }
233 
234 
235 // ************************************************************************* //
scalar deltaT_
Definition: TimeState.H:56
A class for handling file names.
Definition: fileName.H:69
dimensionedScalar pistonDisplacement() const
Return piston displacement for current time step.
Definition: engineTime.C:201
scalar degToTime(const scalar theta) const
Convert degrees to seconds (for given engine speed in RPM)
Definition: engineTime.C:125
scalar writeInterval_
Definition: Time.H:132
scalar deltaTSave_
Definition: TimeState.H:57
virtual bool read()
Read the controlDict and set all the parameters.
Definition: engineTime.C:111
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition: engineTime.C:103
dimensionedScalar sqrt(const dimensionedScalar &ds)
scalar timeToDeg(const scalar t) const
Convert seconds to degrees (for given engine speed in RPM)
Definition: engineTime.C:132
const scalar & value() const
Return const reference to value.
scalar startTime_
Definition: Time.H:123
writeControls writeControl_
Definition: Time.H:130
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensionedScalar pistonPosition() const
Return current piston position.
Definition: engineTime.C:190
scalar deltaT0_
Definition: TimeState.H:58
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar pistonSpeed() const
Return piston speed for current time step.
Definition: engineTime.C:212
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual scalar timeToUserTime(const scalar t) const
Convert the real-time (s) into user-time (CA deg)
Definition: engineTime.C:229
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
scalar endTime_
Definition: Time.H:124
Time(const word &name, const argList &args, const word &systemName="system", const word &constantName="constant")
Construct given name of dictionary to read and argument list.
Definition: Time.C:416
dimensionedScalar sin(const dimensionedScalar &ds)
scalar thetaRevolution() const
Return current crank-angle translated to a single revolution.
Definition: engineTime.C:147
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition: TimeIO.C:33
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
scalar theta() const
Return current crank-angle.
Definition: engineTime.C:139
virtual bool read()
Read control dictionary, update controls and time.
Definition: TimeIO.C:363
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (CA deg) to real-time (s).
Definition: engineTime.C:223
scalar deltaTheta() const
Return crank-angle increment.
Definition: engineTime.C:165
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
const dimensionSet dimVelocity