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-2012 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 
32 {
35 
36  if
37  (
40  )
41  {
43  }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 //- Construct from objectRegistry arguments
51 (
52  const word& name,
53  const fileName& rootPath,
54  const fileName& caseName,
55  const fileName& systemName,
56  const fileName& constantName,
57  const fileName& dictName
58 )
59 :
60  Time
61  (
62  name,
63  rootPath,
64  caseName,
65  systemName,
66  constantName
67  ),
68  dict_
69  (
70  IOobject
71  (
72  "engineGeometry",
73  constant(),
74  *this,
77  false
78  )
79  ),
80  rpm_(dict_.lookup("rpm")),
81  conRodLength_(dimensionedScalar("conRodLength", dimLength, 0)),
82  bore_(dimensionedScalar("bore", dimLength, 0)),
83  stroke_(dimensionedScalar("stroke", dimLength, 0)),
84  clearance_(dimensionedScalar("clearance", dimLength, 0))
85 {
86  // geometric parameters are not strictly required for Time
87  dict_.readIfPresent("conRodLength", conRodLength_);
88  dict_.readIfPresent("bore", bore_);
89  dict_.readIfPresent("stroke", stroke_);
90  dict_.readIfPresent("clearance", clearance_);
91 
92  timeAdjustment();
93 
94  startTime_ = degToTime(startTime_);
95  value() = degToTime(value());
96  deltaTSave_ = deltaT_;
97  deltaT0_ = deltaT_;
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 // Read the controlDict and set all the parameters
105 {
106  Time::readDict();
107  timeAdjustment();
108 }
109 
110 
111 // Read the controlDict and set all the parameters
113 {
114  if (Time::read())
115  {
116  timeAdjustment();
117  return true;
118  }
119  else
120  {
121  return false;
122  }
123 }
124 
125 
126 Foam::scalar Foam::engineTime::degToTime(const scalar theta) const
127 {
128  // 6 * rpm => deg/s
129  return theta/(6.0*rpm_.value());
130 }
131 
132 
133 Foam::scalar Foam::engineTime::timeToDeg(const scalar t) const
134 {
135  // 6 * rpm => deg/s
136  return t*(6.0*rpm_.value());
137 }
138 
139 
140 Foam::scalar Foam::engineTime::theta() const
141 {
142  return timeToDeg(value());
143 }
144 
145 
146 // Return current crank-angle translated to a single revolution
147 // (value between -180 and 180 with 0 = top dead centre)
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::engineTime::deltaTheta() const
167 {
168  return timeToDeg(deltaTValue());
169 }
170 
171 
172 Foam::scalar Foam::engineTime::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 
192 {
193  return dimensionedScalar
194  (
195  "pistonPosition",
196  dimLength,
197  pistonPosition(theta())
198  );
199 }
200 
201 
203 {
204  return dimensionedScalar
205  (
206  "pistonDisplacement",
207  dimLength,
208  pistonPosition(theta() - deltaTheta()) - pistonPosition().value()
209  );
210 }
211 
212 
214 {
215  return dimensionedScalar
216  (
217  "pistonSpeed",
218  dimVelocity,
219  pistonDisplacement().value()/(deltaTValue() + VSMALL)
220  );
221 }
222 
223 
224 Foam::scalar Foam::engineTime::userTimeToTime(const scalar theta) const
225 {
226  return degToTime(theta);
227 }
228 
229 
230 Foam::scalar Foam::engineTime::timeToUserTime(const scalar t) const
231 {
232  return timeToDeg(t);
233 }
234 
235 
236 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::Time::wcRunTime
@ wcRunTime
Definition: Time.H:94
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::engineTime::deltaTheta
scalar deltaTheta() const
Return crank-angle increment.
Definition: engineTime.C:166
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::engineTime::pistonSpeed
dimensionedScalar pistonSpeed() const
Return piston speed for current time step.
Definition: engineTime.C:213
unitConversion.H
Unit conversion functions.
Foam::engineTime::timeToDeg
scalar timeToDeg(const scalar t) const
Convert seconds to degrees (for given engine speed in RPM)
Definition: engineTime.C:133
Foam::Time::wcAdjustableRunTime
@ wcAdjustableRunTime
Definition: Time.H:95
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::Time::read
virtual bool read()
Read control dictionary, update controls and time.
Definition: TimeIO.C:436
Foam::engineTime::timeAdjustment
void timeAdjustment()
Adjust read time values.
Definition: engineTime.C:31
constant
Constant dispersed-phase particle diameter model.
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
dictName
const word dictName("particleTrackDict")
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::engineTime::thetaRevolution
scalar thetaRevolution() const
Return current crank-angle translated to a single revolution.
Definition: engineTime.C:148
Foam::engineTime::timeToUserTime
virtual scalar timeToUserTime(const scalar t) const
Convert the real-time (s) into user-time (CA deg)
Definition: engineTime.C:230
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
Definition: dimensionedType.C:309
Foam::engineTime::readDict
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition: engineTime.C:104
Foam::Time::readDict
virtual void readDict()
Read the control dictionary and set the write controls etc.
Definition: TimeIO.C:33
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
engineTime.H
Foam::Time::writeInterval_
scalar writeInterval_
Definition: Time.H:132
Foam::engineTime::userTimeToTime
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (CA deg) to real-time (s).
Definition: engineTime.C:224
Foam::engineTime::theta
scalar theta() const
Return current crank-angle.
Definition: engineTime.C:140
Foam::Time::endTime_
scalar endTime_
Definition: Time.H:124
Foam::engineTime::read
virtual bool read()
Read the controlDict and set all the parameters.
Definition: engineTime.C:112
Foam::Time::writeControl_
writeControls writeControl_
Definition: Time.H:130
Foam::TimeState::deltaT_
scalar deltaT_
Definition: TimeState.H:57
Foam::engineTime::pistonPosition
dimensionedScalar pistonPosition() const
Return current piston position.
Definition: engineTime.C:191
Foam::engineTime::engineTime
engineTime(const engineTime &)
Disallow default bitwise copy construct.
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::engineTime::degToTime
scalar degToTime(const scalar theta) const
Convert degrees to seconds (for given engine speed in RPM)
Definition: engineTime.C:126
Foam::degToRad
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Definition: unitConversion.H:45
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256
Foam::engineTime::pistonDisplacement
dimensionedScalar pistonDisplacement() const
Return piston displacement for current time step.
Definition: engineTime.C:202