1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "engineTime.H"
27 #include "unitConversion.H"
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 void Foam::engineTime::timeAdjustment()
33 deltaT_ = degToTime(deltaT_);
34 endTime_ = degToTime(endTime_);
38 writeControl_ == wcRunTime
39 || writeControl_ == wcAdjustableRunTime
42 writeInterval_ = degToTime(writeInterval_);
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 //- Construct from objectRegistry arguments
50 Foam::engineTime::engineTime
53 const fileName& rootPath,
54 const fileName& caseName,
55 const fileName& systemName,
56 const fileName& constantName,
57 const fileName& dictName
75 IOobject::MUST_READ_IF_MODIFIED,
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))
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_);
94 startTime_ = degToTime(startTime_);
95 value() = degToTime(value());
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 // Read the controlDict and set all the parameters
103 void Foam::engineTime::readDict()
110 // Read the controlDict and set all the parameters
111 bool Foam::engineTime::read()
125 Foam::scalar Foam::engineTime::degToTime(const scalar theta) const
128 return theta/(6.0*rpm_.value());
132 Foam::scalar Foam::engineTime::timeToDeg(const scalar t) const
135 return t*(6.0*rpm_.value());
139 Foam::scalar Foam::engineTime::theta() const
141 return timeToDeg(value());
145 // Return current crank-angle translated to a single revolution
146 // (value between -180 and 180 with 0 = top dead centre)
147 Foam::scalar Foam::engineTime::thetaRevolution() const
165 Foam::scalar Foam::engineTime::deltaTheta() const
167 return timeToDeg(deltaTValue());
171 Foam::scalar Foam::engineTime::pistonPosition(const scalar theta) const
175 conRodLength_.value()
176 + stroke_.value()/2.0
180 stroke_.value()*::cos(degToRad(theta))/2.0
183 sqr(conRodLength_.value())
184 - sqr(stroke_.value()*::sin(degToRad(theta))/2.0)
190 Foam::dimensionedScalar Foam::engineTime::pistonPosition() const
192 return dimensionedScalar
196 pistonPosition(theta())
201 Foam::dimensionedScalar Foam::engineTime::pistonDisplacement() const
203 return dimensionedScalar
205 "pistonDisplacement",
207 pistonPosition(theta() - deltaTheta()) - pistonPosition().value()
212 Foam::dimensionedScalar Foam::engineTime::pistonSpeed() const
214 return dimensionedScalar
218 pistonDisplacement().value()/(deltaTValue() + VSMALL)
223 Foam::scalar Foam::engineTime::userTimeToTime(const scalar theta) const
225 return degToTime(theta);
229 Foam::scalar Foam::engineTime::timeToUserTime(const scalar t) const
235 // ************************************************************************* //