engineValve.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 "engineValve.H"
27 #include "engineTime.H"
28 #include "polyMesh.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 Foam::scalar Foam::engineValve::adjustCrankAngle(const scalar theta) const
33 {
34  if (theta < liftProfileStart_)
35  {
36  scalar adjustedTheta = theta;
37 
38  while (adjustedTheta < liftProfileStart_)
39  {
40  adjustedTheta += liftProfileEnd_ - liftProfileStart_;
41  }
42 
43  return adjustedTheta;
44  }
45  else if (theta > liftProfileEnd_)
46  {
47  scalar adjustedTheta = theta;
48 
49  while (adjustedTheta > liftProfileEnd_)
50  {
51  adjustedTheta -= liftProfileEnd_ - liftProfileStart_;
52  }
53 
54  return adjustedTheta;
55  }
56  else
57  {
58  return theta;
59  }
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 (
67  const word& name,
68  const polyMesh& mesh,
69  const autoPtr<coordinateSystem>& valveCS,
70  const word& bottomPatchName,
71  const word& poppetPatchName,
72  const word& stemPatchName,
73  const word& curtainInPortPatchName,
74  const word& curtainInCylinderPatchName,
75  const word& detachInCylinderPatchName,
76  const word& detachInPortPatchName,
77  const labelList& detachFaces,
78  const Function1s::Table<scalar>& liftProfile,
79  const scalar minLift,
80  const scalar minTopLayer,
81  const scalar maxTopLayer,
82  const scalar minBottomLayer,
83  const scalar maxBottomLayer,
84  const scalar diameter
85 )
86 :
87  name_(name),
88  mesh_(mesh),
89  engineDB_(refCast<const engineTime>(mesh.time())),
90  csPtr_(valveCS),
91  bottomPatch_(bottomPatchName, mesh.boundaryMesh()),
92  poppetPatch_(poppetPatchName, mesh.boundaryMesh()),
93  stemPatch_(stemPatchName, mesh.boundaryMesh()),
94  curtainInPortPatch_(curtainInPortPatchName, mesh.boundaryMesh()),
95  curtainInCylinderPatch_(curtainInCylinderPatchName, mesh.boundaryMesh()),
96  detachInCylinderPatch_(detachInCylinderPatchName, mesh.boundaryMesh()),
97  detachInPortPatch_(detachInPortPatchName, mesh.boundaryMesh()),
98  detachFaces_(detachFaces),
99  liftProfile_(liftProfile),
100  liftProfileStart_(min(liftProfile_.x())),
101  liftProfileEnd_(max(liftProfile_.x())),
102  minLift_(minLift),
103  minTopLayer_(minTopLayer),
104  maxTopLayer_(maxTopLayer),
105  minBottomLayer_(minBottomLayer),
106  maxBottomLayer_(maxBottomLayer),
107  diameter_(diameter)
108 {}
109 
110 
112 (
113  const word& name,
114  const polyMesh& mesh,
115  const dictionary& dict
116 )
117 :
118  name_(name),
119  mesh_(mesh),
120  engineDB_(refCast<const engineTime>(mesh_.time())),
121  csPtr_
122  (
124  (
125  mesh_,
126  dict.subDict("coordinateSystem")
127  )
128  ),
129  bottomPatch_(dict.lookup("bottomPatch"), mesh.boundaryMesh()),
130  poppetPatch_(dict.lookup("poppetPatch"), mesh.boundaryMesh()),
131  stemPatch_(dict.lookup("stemPatch"), mesh.boundaryMesh()),
132  curtainInPortPatch_
133  (
134  dict.lookup("curtainInPortPatch"),
135  mesh.boundaryMesh()
136  ),
137  curtainInCylinderPatch_
138  (
139  dict.lookup("curtainInCylinderPatch"),
140  mesh.boundaryMesh()
141  ),
142  detachInCylinderPatch_
143  (
144  dict.lookup("detachInCylinderPatch"),
145  mesh.boundaryMesh()
146  ),
147  detachInPortPatch_
148  (
149  dict.lookup("detachInPortPatch"),
150  mesh.boundaryMesh()
151  ),
152  detachFaces_(dict.lookup("detachFaces")),
153  liftProfile_("liftProfile", dict),
154  liftProfileStart_(min(liftProfile_.x())),
155  liftProfileEnd_(max(liftProfile_.x())),
156  minLift_(dict.lookup<scalar>("minLift")),
157  minTopLayer_(dict.lookup<scalar>("minTopLayer")),
158  maxTopLayer_(dict.lookup<scalar>("maxTopLayer")),
159  minBottomLayer_(dict.lookup<scalar>("minBottomLayer")),
160  maxBottomLayer_(dict.lookup<scalar>("maxBottomLayer")),
161  diameter_(dict.lookup<scalar>("diameter"))
162 {}
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
167 Foam::scalar Foam::engineValve::lift(const scalar theta) const
168 {
169  return liftProfile_.value(adjustCrankAngle(theta));
170 }
171 
172 
174 {
175  return lift(engineDB_.theta()) >= minLift_;
176 }
177 
178 
179 Foam::scalar Foam::engineValve::curLift() const
180 {
181  return max
182  (
183  lift(engineDB_.theta()),
184  minLift_
185  );
186 }
187 
188 
189 Foam::scalar Foam::engineValve::curVelocity() const
190 {
191  return
192  -(
193  curLift()
194  - max
195  (
196  lift(engineDB_.theta() - engineDB_.deltaTheta()),
197  minLift_
198  )
199  )/(engineDB_.deltaTValue() + vSmall);
200 }
201 
202 
204 {
205  labelList mpIDs(2);
206  label nMpIDs = 0;
207 
208  if (bottomPatch_.active())
209  {
210  mpIDs[nMpIDs] = bottomPatch_.index();
211  nMpIDs++;
212  }
213 
214  if (poppetPatch_.active())
215  {
216  mpIDs[nMpIDs] = poppetPatch_.index();
217  nMpIDs++;
218  }
219 
220  mpIDs.setSize(nMpIDs);
221 
222  return mpIDs;
223 }
224 
225 
227 {
228  os << nl << name() << nl << token::BEGIN_BLOCK;
229 
230  cs().writeDict(os);
231 
232  os << "bottomPatch " << bottomPatch_.name() << token::END_STATEMENT << nl
233  << "poppetPatch " << poppetPatch_.name() << token::END_STATEMENT << nl
234  << "stemPatch " << stemPatch_.name() << token::END_STATEMENT << nl
235  << "curtainInPortPatch " << curtainInPortPatch_.name()
237  << "curtainInCylinderPatch " << curtainInCylinderPatch_.name()
239  << "detachInCylinderPatch " << detachInCylinderPatch_.name()
241  << "detachInPortPatch " << detachInPortPatch_.name()
243  << "detachFaces " << detachFaces_ << token::END_STATEMENT << nl
244  << "liftProfile " << nl << token::BEGIN_BLOCK
245  << liftProfile_ << token::END_BLOCK << token::END_STATEMENT << nl
246  << "minLift " << minLift_ << token::END_STATEMENT << nl
247  << "minTopLayer " << minTopLayer_ << token::END_STATEMENT << nl
248  << "maxTopLayer " << maxTopLayer_ << token::END_STATEMENT << nl
249  << "minBottomLayer " << minBottomLayer_ << token::END_STATEMENT << nl
250  << "maxBottomLayer " << maxBottomLayer_ << token::END_STATEMENT << nl
251  << "diameter " << diameter_ << token::END_STATEMENT << nl
252  << token::END_BLOCK << endl;
253 }
254 
255 
256 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool isOpen() const
Is the valve open?
Definition: engineValve.C:173
virtual Type value(const scalar x) const
Return Table value as a function of scalar x.
Definition: Table.C:204
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
engineValve(const word &name, const polyMesh &mesh, const autoPtr< coordinateSystem > &valveCS, const word &bottomPatchName, const word &poppetPatchName, const word &stemPatchName, const word &curtainInPortPatchName, const word &curtainInCylinderPatchName, const word &detachInCylinderPatchName, const word &detachInPortPatchName, const labelList &detachFaces, const Function1s::Table< scalar > &liftProfile, const scalar minLift, const scalar minTopLayer, const scalar maxTopLayer, const scalar minBottomLayer, const scalar maxBottomLayer, const scalar diameter)
Construct from components.
Definition: engineValve.C:66
void writeDict(Ostream &) const
Write dictionary.
Definition: engineValve.C:226
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
scalar curVelocity() const
Return valve velocity for current time-step.
Definition: engineValve.C:189
A class for handling words, derived from string.
Definition: word.H:59
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
label index() const
Return index of first matching zone.
Definition: DynamicID.H:107
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
const Time & time() const
Return time.
scalar curLift() const
Return current lift.
Definition: engineValve.C:179
scalar lift(const scalar theta) const
Return valve lift given crank angle in degrees.
Definition: engineValve.C:167
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual scalar deltaTheta() const =0
Return engine time increment.
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual scalar theta() const =0
Return current engine time.
static autoPtr< coordinateSystem > New(const objectRegistry &obr, const dictionary &dict)
Select constructed from dictionary and objectRegistry.
bool active() const
Has the zone been found.
Definition: DynamicID.H:113
const word & name() const
Return name.
Definition: engineValve.H:182
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const coordinateSystem & cs() const
Return coordinate system.
Definition: engineValve.H:188
const wordRe & name() const
Return name.
Definition: DynamicID.H:95
labelList movingPatchIDs() const
Return list of active patch labels for the valve head.
Definition: engineValve.C:203
void writeDict(Ostream &, bool subDict=true) const
Write dictionary.
virtual tmp< scalarField > x() const
Return the reference values.
Definition: Table.C:266
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844