FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

engineTime.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "engineTime.H"
00027 #include <OpenFOAM/mathematicalConstants.H>
00028 
00029 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00030 
00031 void Foam::engineTime::timeAdjustment()
00032 {
00033     deltaT_  = degToTime(deltaT_);
00034     endTime_ = degToTime(endTime_);
00035 
00036     if
00037     (
00038         writeControl_ == wcRunTime
00039      || writeControl_ == wcAdjustableRunTime
00040     )
00041     {
00042         writeInterval_ = degToTime(writeInterval_);
00043     }
00044 }
00045 
00046 
00047 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00048 
00049 //- Construct from objectRegistry arguments
00050 Foam::engineTime::engineTime
00051 (
00052     const word& name,
00053     const fileName& rootPath,
00054     const fileName& caseName,
00055     const fileName& systemName,
00056     const fileName& constantName,
00057     const fileName& dictName
00058 )
00059 :
00060     Time
00061     (
00062         name,
00063         rootPath,
00064         caseName,
00065         systemName,
00066         constantName
00067     ),
00068     dict_
00069     (
00070         IOobject
00071         (
00072             "engineGeometry",
00073             constant(),
00074             *this,
00075             IOobject::MUST_READ,
00076             IOobject::NO_WRITE,
00077             false
00078         )
00079     ),
00080     rpm_(dict_.lookup("rpm")),
00081     conRodLength_(dimensionedScalar("conRodLength", dimLength, 0)),
00082     bore_(dimensionedScalar("bore", dimLength, 0)),
00083     stroke_(dimensionedScalar("stroke", dimLength, 0)),
00084     clearance_(dimensionedScalar("clearance", dimLength, 0))
00085 {
00086     // geometric parameters are not strictly required for Time
00087     dict_.readIfPresent("conRodLength", conRodLength_);
00088     dict_.readIfPresent("bore", bore_);
00089     dict_.readIfPresent("stroke", stroke_);
00090     dict_.readIfPresent("clearance", clearance_);
00091 
00092     timeAdjustment();
00093 
00094     startTime_ = degToTime(startTime_);
00095     value()    = degToTime(value());
00096     deltaT0_   = deltaT_;
00097 }
00098 
00099 
00100 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00101 
00102 // Read the controlDict and set all the parameters
00103 void Foam::engineTime::readDict()
00104 {
00105     Time::readDict();
00106     timeAdjustment();
00107 }
00108 
00109 
00110 // Read the controlDict and set all the parameters
00111 bool Foam::engineTime::read()
00112 {
00113     if (Time::read())
00114     {
00115         timeAdjustment();
00116         return true;
00117     }
00118     else
00119     {
00120         return false;
00121     }
00122 }
00123 
00124 
00125 Foam::scalar Foam::engineTime::degToRad(const scalar deg) const
00126 {
00127     return mathematicalConstant::pi*deg/180.0;
00128 }
00129 
00130 
00131 Foam::scalar Foam::engineTime::degToTime(const scalar theta) const
00132 {
00133     // 6 * rpm => deg/s
00134     return theta/(6.0*rpm_.value());
00135 }
00136 
00137 
00138 Foam::scalar Foam::engineTime::timeToDeg(const scalar t) const
00139 {
00140     // 6 * rpm => deg/s
00141     return t*(6.0*rpm_.value());
00142 }
00143 
00144 
00145 Foam::scalar Foam::engineTime::theta() const
00146 {
00147     return timeToDeg(value());
00148 }
00149 
00150 
00151 // Return current crank-angle translated to a single revolution
00152 // (value between -180 and 180 with 0 = top dead centre)
00153 Foam::scalar Foam::engineTime::thetaRevolution() const
00154 {
00155     scalar t = theta();
00156 
00157     while (t > 180.0)
00158     {
00159         t -= 360.0;
00160     }
00161 
00162     while (t < -180.0)
00163     {
00164         t += 360.0;
00165     }
00166 
00167     return t;
00168 }
00169 
00170 
00171 Foam::scalar Foam::engineTime::deltaTheta() const
00172 {
00173     return timeToDeg(deltaT().value());
00174 }
00175 
00176 
00177 Foam::scalar Foam::engineTime::pistonPosition(const scalar theta) const
00178 {
00179     return
00180     (
00181         conRodLength_.value()
00182       + stroke_.value()/2.0
00183       + clearance_.value()
00184     )
00185   - (
00186         stroke_.value()*::cos(degToRad(theta))/2.0
00187       + ::sqrt
00188         (
00189             sqr(conRodLength_.value())
00190             - sqr(stroke_.value()*::sin(degToRad(theta))/2.0)
00191         )
00192     );
00193 }
00194 
00195 
00196 Foam::dimensionedScalar Foam::engineTime::pistonPosition() const
00197 {
00198     return dimensionedScalar
00199     (
00200         "pistonPosition",
00201         dimLength,
00202         pistonPosition(theta())
00203     );
00204 }
00205 
00206 
00207 Foam::dimensionedScalar Foam::engineTime::pistonDisplacement() const
00208 {
00209     return dimensionedScalar
00210     (
00211         "pistonDisplacement",
00212         dimLength,
00213         pistonPosition(theta() - deltaTheta()) - pistonPosition().value()
00214     );
00215 }
00216 
00217 
00218 Foam::dimensionedScalar Foam::engineTime::pistonSpeed() const
00219 {
00220     return dimensionedScalar
00221     (
00222         "pistonSpeed",
00223         dimVelocity,
00224         pistonDisplacement().value()/(deltaT().value() + VSMALL)
00225     );
00226 }
00227 
00228 
00229 Foam::scalar Foam::engineTime::userTimeToTime(const scalar theta) const
00230 {
00231     return degToTime(theta);
00232 }
00233 
00234 
00235 Foam::scalar Foam::engineTime::timeToUserTime(const scalar t) const
00236 {
00237     return timeToDeg(t);
00238 }
00239 
00240 
00241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00242 
00243 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines