fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / engine / engineTime / engineTime.C
blob1f3bfa9834ac87cb976db6e7845bf93f639fa78a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "engineTime.H"
28 #include "mathematicalConstants.H"
30 // * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * * * //
32 defineTypeNameAndDebug(Foam::engineTime, 0);
35 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
37 void Foam::engineTime::timeAdjustment()
39     deltaT_  = degToTime(deltaT_);
40     endTime_ = degToTime(endTime_);
42     if
43     (
44         writeControl_ == wcRunTime
45      || writeControl_ == wcAdjustableRunTime
46     )
47     {
48         writeInterval_ = degToTime(writeInterval_);
49     }
53 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
55 //- Construct from objectRegistry arguments
56 Foam::engineTime::engineTime
58     const word& name,
59     const fileName& rootPath,
60     const fileName& caseName,
61     const fileName& systemName,
62     const fileName& constantName,
63     const fileName& dictName
66     Time
67     (
68         name,
69         rootPath,
70         caseName,
71         systemName,
72         constantName
73     ),
74     dict_
75     (
76         IOobject
77         (
78             dictName,
79             constant(),
80             *this,
81             IOobject::MUST_READ,
82             IOobject::NO_WRITE,
83             false
84         )
85     ),
86     rpm_(dict_.lookup("rpm")),
87     conRodLength_(dimensionedScalar("conRodLength", dimLength, 0)),
88     bore_(dimensionedScalar("bore", dimLength, 0)),
89     stroke_(dimensionedScalar("stroke", dimLength, 0)),
90     clearance_(dimensionedScalar("clearance", dimLength, 0))
92     // the geometric parameters are not strictly required for Time
93     if (dict_.found("conRodLength"))
94     {
95         dict_.lookup("conRodLength") >> conRodLength_;
96     }
97     if (dict_.found("bore"))
98     {
99         dict_.lookup("bore") >> bore_;
100     }
101     if (dict_.found("stroke"))
102     {
103         dict_.lookup("stroke") >> stroke_;
104     }
105     if (dict_.found("clearance"))
106     {
107         dict_.lookup("clearance") >> clearance_;
108     }
110     timeAdjustment();
112     startTime_ = degToTime(startTime_);
113     value()    = degToTime(value());
114     deltaT0_   = deltaT_;
118 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
120 // Read the controlDict and set all the parameters
121 void Foam::engineTime::readDict()
123     Time::readDict();
124     timeAdjustment();
128 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
130 // Read the controlDict and set all the parameters
131 bool Foam::engineTime::read()
133     if (Time::read())
134     {
135         timeAdjustment();
136         return true;
137     }
138     else
139     {
140         return false;
141     }
145 Foam::scalar Foam::engineTime::degToRad(const scalar deg) const
147     return mathematicalConstant::pi*deg/180.0;
151 Foam::scalar Foam::engineTime::degToTime(const scalar theta) const
153     // 6 * rpm => deg/s
154     return theta/(6.0*rpm_.value());
158 Foam::scalar Foam::engineTime::timeToDeg(const scalar t) const
160     // 6 * rpm => deg/s
161     return t*(6.0*rpm_.value());
165 Foam::scalar Foam::engineTime::theta() const
167     return timeToDeg(value());
171 // Return current crank-angle translated to a single revolution
172 // (value between -180 and 180 with 0 = top dead centre)
173 Foam::scalar Foam::engineTime::thetaRevolution() const
175     scalar t = theta();
177     while (t > 180.0)
178     {
179         t -= 360.0;
180     }
182     while (t < -180.0)
183     {
184         t += 360.0;
185     }
187     return t;
191 Foam::scalar Foam::engineTime::deltaTheta() const
193     return timeToDeg(deltaT().value());
197 Foam::scalar Foam::engineTime::pistonPosition(const scalar theta) const
199     return
200     (
201         conRodLength_.value()
202       + stroke_.value()/2.0
203       + clearance_.value()
204     )
205   - (
206         stroke_.value()*::cos(degToRad(theta))/2.0
207       + ::sqrt
208         (
209             sqr(conRodLength_.value())
210             - sqr(stroke_.value()*::sin(degToRad(theta))/2.0)
211         )
212     );
216 Foam::dimensionedScalar Foam::engineTime::pistonPosition() const
218     return dimensionedScalar
219     (
220         "pistonPosition",
221         dimLength,
222         pistonPosition(theta())
223     );
227 Foam::dimensionedScalar Foam::engineTime::pistonDisplacement() const
229     return dimensionedScalar
230     (
231         "pistonDisplacement",
232         dimLength,
233         pistonPosition(theta() - deltaTheta()) - pistonPosition().value()
234     );
238 Foam::dimensionedScalar Foam::engineTime::pistonSpeed() const
240     return dimensionedScalar
241     (
242         "pistonSpeed",
243         dimVelocity,
244         pistonDisplacement().value()/(deltaT().value() + VSMALL)
245     );
249 Foam::scalar Foam::engineTime::userTimeToTime(const scalar theta) const
251     return degToTime(theta);
255 Foam::scalar Foam::engineTime::timeToUserTime(const scalar t) const
257     return timeToDeg(t);
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 // ************************************************************************* //