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-2021 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 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 Foam::scalar Foam::engineValve::adjustCrankAngle(const scalar theta) const
31 {
32  if (theta < liftProfileStart_)
33  {
34  scalar adjustedTheta = theta;
35 
36  while (adjustedTheta < liftProfileStart_)
37  {
38  adjustedTheta += liftProfileEnd_ - liftProfileStart_;
39  }
40 
41  return adjustedTheta;
42  }
43  else if (theta > liftProfileEnd_)
44  {
45  scalar adjustedTheta = theta;
46 
47  while (adjustedTheta > liftProfileEnd_)
48  {
49  adjustedTheta -= liftProfileEnd_ - liftProfileStart_;
50  }
51 
52  return adjustedTheta;
53  }
54  else
55  {
56  return theta;
57  }
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
64 (
65  const word& name,
66  const fvMeshMover& meshMover,
67  const autoPtr<coordinateSystem>& valveCS,
68  const word& bottomPatchName,
69  const word& poppetPatchName,
70  const word& stemPatchName,
71  const word& curtainInPortPatchName,
72  const word& curtainInCylinderPatchName,
73  const word& detachInCylinderPatchName,
74  const word& detachInPortPatchName,
75  const labelList& detachFaces,
76  const Function1s::Table<scalar>& liftProfile,
77  const scalar minLift,
78  const scalar minTopLayer,
79  const scalar maxTopLayer,
80  const scalar minBottomLayer,
81  const scalar maxBottomLayer,
82  const scalar diameter
83 )
84 :
85  name_(name),
86  meshMover_(refCast<const fvMeshMovers::engine>(meshMover)),
87  csPtr_(valveCS),
88  bottomPatch_(bottomPatchName, meshMover_.mesh().boundaryMesh()),
89  poppetPatch_(poppetPatchName, meshMover_.mesh().boundaryMesh()),
90  stemPatch_(stemPatchName, meshMover_.mesh().boundaryMesh()),
91  curtainInPortPatch_
92  (
93  curtainInPortPatchName, meshMover_.mesh().boundaryMesh()
94  ),
95  curtainInCylinderPatch_
96  (
97  curtainInCylinderPatchName, meshMover_.mesh().boundaryMesh()
98  ),
99  detachInCylinderPatch_
100  (
101  detachInCylinderPatchName, meshMover_.mesh().boundaryMesh()
102  ),
103  detachInPortPatch_(detachInPortPatchName, meshMover_.mesh().boundaryMesh()),
104  detachFaces_(detachFaces),
105  liftProfile_(liftProfile),
106  liftProfileStart_(min(liftProfile_.x())),
107  liftProfileEnd_(max(liftProfile_.x())),
108  minLift_(minLift),
109  minTopLayer_(minTopLayer),
110  maxTopLayer_(maxTopLayer),
111  minBottomLayer_(minBottomLayer),
112  maxBottomLayer_(maxBottomLayer),
113  diameter_(diameter)
114 {}
115 
116 
118 (
119  const word& name,
120  const fvMeshMover& meshMover,
121  const dictionary& dict
122 )
123 :
124  name_(name),
125  meshMover_(refCast<const fvMeshMovers::engine>(meshMover)),
126  csPtr_
127  (
129  (
130  meshMover_.mesh(),
131  dict.subDict("coordinateSystem")
132  )
133  ),
134  bottomPatch_(dict.lookup("bottomPatch"), meshMover_.mesh().boundaryMesh()),
135  poppetPatch_(dict.lookup("poppetPatch"), meshMover_.mesh().boundaryMesh()),
136  stemPatch_(dict.lookup("stemPatch"), meshMover_.mesh().boundaryMesh()),
137  curtainInPortPatch_
138  (
139  dict.lookup("curtainInPortPatch"),
140  meshMover_.mesh().boundaryMesh()
141  ),
142  curtainInCylinderPatch_
143  (
144  dict.lookup("curtainInCylinderPatch"),
145  meshMover_.mesh().boundaryMesh()
146  ),
147  detachInCylinderPatch_
148  (
149  dict.lookup("detachInCylinderPatch"),
150  meshMover_.mesh().boundaryMesh()
151  ),
152  detachInPortPatch_
153  (
154  dict.lookup("detachInPortPatch"),
155  meshMover_.mesh().boundaryMesh()
156  ),
157  detachFaces_(dict.lookup("detachFaces")),
158  liftProfile_("liftProfile", dict),
159  liftProfileStart_(min(liftProfile_.x())),
160  liftProfileEnd_(max(liftProfile_.x())),
161  minLift_(dict.lookup<scalar>("minLift")),
162  minTopLayer_(dict.lookup<scalar>("minTopLayer")),
163  maxTopLayer_(dict.lookup<scalar>("maxTopLayer")),
164  minBottomLayer_(dict.lookup<scalar>("minBottomLayer")),
165  maxBottomLayer_(dict.lookup<scalar>("maxBottomLayer")),
166  diameter_(dict.lookup<scalar>("diameter"))
167 {}
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
172 Foam::scalar Foam::engineValve::lift(const scalar theta) const
173 {
174  return liftProfile_.value(adjustCrankAngle(theta));
175 }
176 
177 
179 {
180  return lift(meshMover_.theta()) >= minLift_;
181 }
182 
183 
184 Foam::scalar Foam::engineValve::curLift() const
185 {
186  return max
187  (
188  lift(meshMover_.theta()),
189  minLift_
190  );
191 }
192 
193 
194 Foam::scalar Foam::engineValve::curVelocity() const
195 {
196  return
197  -(
198  curLift()
199  - max
200  (
201  lift(meshMover_.theta() - meshMover_.deltaTheta()),
202  minLift_
203  )
204  )/(meshMover_.mesh().time().deltaTValue() + vSmall);
205 }
206 
207 
209 {
210  labelList mpIDs(2);
211  label nMpIDs = 0;
212 
213  if (bottomPatch_.active())
214  {
215  mpIDs[nMpIDs] = bottomPatch_.index();
216  nMpIDs++;
217  }
218 
219  if (poppetPatch_.active())
220  {
221  mpIDs[nMpIDs] = poppetPatch_.index();
222  nMpIDs++;
223  }
224 
225  mpIDs.setSize(nMpIDs);
226 
227  return mpIDs;
228 }
229 
230 
232 {
233  os << nl << name() << nl << token::BEGIN_BLOCK;
234 
235  cs().writeDict(os);
236 
237  os << "bottomPatch " << bottomPatch_.name() << token::END_STATEMENT << nl
238  << "poppetPatch " << poppetPatch_.name() << token::END_STATEMENT << nl
239  << "stemPatch " << stemPatch_.name() << token::END_STATEMENT << nl
240  << "curtainInPortPatch " << curtainInPortPatch_.name()
242  << "curtainInCylinderPatch " << curtainInCylinderPatch_.name()
244  << "detachInCylinderPatch " << detachInCylinderPatch_.name()
246  << "detachInPortPatch " << detachInPortPatch_.name()
248  << "detachFaces " << detachFaces_ << token::END_STATEMENT << nl
249  << "liftProfile " << nl << token::BEGIN_BLOCK
250  << liftProfile_ << token::END_BLOCK << token::END_STATEMENT << nl
251  << "minLift " << minLift_ << token::END_STATEMENT << nl
252  << "minTopLayer " << minTopLayer_ << token::END_STATEMENT << nl
253  << "maxTopLayer " << maxTopLayer_ << token::END_STATEMENT << nl
254  << "minBottomLayer " << minBottomLayer_ << token::END_STATEMENT << nl
255  << "maxBottomLayer " << maxBottomLayer_ << token::END_STATEMENT << nl
256  << "diameter " << diameter_ << token::END_STATEMENT << nl
257  << token::END_BLOCK << endl;
258 }
259 
260 
261 // ************************************************************************* //
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool isOpen() const
Is the valve open?
Definition: engineValve.C:178
virtual Type value(const scalar x) const
Return Table value as a function of scalar x.
Definition: Table.C:206
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Abstract base class for fvMesh movers.
Definition: fvMeshMover.H:52
void writeDict(Ostream &) const
Write dictionary.
Definition: engineValve.C:231
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
scalar curVelocity() const
Return valve velocity for current time-step.
Definition: engineValve.C:194
A class for handling words, derived from string.
Definition: word.H:59
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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
scalar curLift() const
Return current lift.
Definition: engineValve.C:184
scalar lift(const scalar theta) const
Return valve lift given crank angle in degrees.
Definition: engineValve.C:172
void setSize(const label)
Reset size of List.
Definition: List.C:281
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:175
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const coordinateSystem & cs() const
Return coordinate system.
Definition: engineValve.H:181
engineValve(const word &name, const fvMeshMover &meshMover, 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:64
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:208
void writeDict(Ostream &, bool subDict=true) const
Write dictionary.
virtual tmp< scalarField > x() const
Return the reference values.
Definition: Table.C:268
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864