piston.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) 2024 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 "multiValveEngine.H"
27 #include "pointDist.H"
28 
29 /* * * * * * * * * * * * * Static Private Member Data * * * * * * * * * * * */
30 
32 (
33  "pistonBowl"
34 );
35 
36 
37 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
38 
39 Foam::scalar Foam::fvMeshMovers::multiValveEngine::pistonObject::bore() const
40 {
41  const polyBoundaryMesh& pbm = meshMover_.mesh().boundaryMesh();
42 
43  // Find the maximum and minimum coordinates of the piston patch-set
44  vector pistonMax(vector::min);
45  vector pistonMin(vector::max);
46 
48  {
49  const label patchi = iter.key();
50  if (pbm[patchi].localPoints().size())
51  {
52  pistonMax = max(pistonMax, max(pbm[patchi].localPoints()));
53  pistonMin = min(pistonMin, min(pbm[patchi].localPoints()));
54  }
55  }
56 
57  reduce(pistonMax, maxOp<point>());
58  reduce(pistonMin, minOp<point>());
59 
60  // Assuming the piston moves in the positive axis direction
61  // remove the axis_ component to find the lateral extent of the piston
62  return mag
63  (
64  (pistonMax - (axis& pistonMax)*pistonMax)
65  - (pistonMin - (axis& pistonMin)*pistonMin)
66  )/sqrt(2.0);
67 }
68 
69 
70 void Foam::fvMeshMovers::multiValveEngine::pistonObject::correctClearance
71 (
72  pointDist& pDist
73 )
74 {
75  clearance_ = great;
76 
77  forAllConstIter(labelHashSet, staticPatchSet_, iter)
78  {
79  const polyPatch& pp = meshMover_.mesh().boundaryMesh()[iter.key()];
80  const labelList& meshPoints = pp.meshPoints();
81 
82  forAll(meshPoints, pointi)
83  {
84  clearance_ = min(clearance_, pDist[meshPoints[pointi]]);
85  }
86  }
87 
88  reduce(clearance_, minOp<scalar>());
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
95 (
96  const word& name,
97  const multiValveEngine& engine,
98  const dictionary& dict
99 )
100 :
101  movingObject(name, engine, dict),
102  bore_(bore()),
103  clearance_(0)
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 (
111  const scalar theta
112 ) const
113 {
114  return motion_->value(theta);
115 }
116 
117 
118 Foam::scalar
120 {
121  return position(meshMover_.userTime());
122 }
123 
124 
125 Foam::scalar
127 {
128  return
129  position(meshMover_.userTime() - meshMover_.userDeltaT())
130  - position();
131 }
132 
133 
135 {
136  return displacement()/(meshMover_.mesh().time().deltaTValue() + vSmall);
137 }
138 
139 
140 Foam::scalar
142 {
143  if (mag(position() - position0_) > travelInterval_)
144  {
145  return clearance_;
146  }
147  else
148  {
149  // Note, valve movement is not considered as valveSet may include
150  // other valves than ones related to ports. Furthermore, this value
151  // is only an estimate and updated rather frequently anyway.
152  return clearance_ - displacement();
153  }
154 }
155 
156 
158 (
159  pointField& newPoints
160 )
161 {
162  const scalar position = this->position();
163 
164  // Update a cached scale_ field if needed
165  if
166  (
167  executionCount_ == 0
168  || mag(position - position0_) > travelInterval_
169  )
170  {
171  Info << " Updating scale field" << endl;
172 
173  const pointMesh& pMesh = pointMesh::New(meshMover_.mesh());
174  const pointField& points(meshMover_.mesh().points());
175 
176  pointDist pDistMoving
177  (
178  pMesh,
179  patchSet,
180  movingPointZones(),
181  points,
182  maxMotionDistance_
183  );
184 
185  pointDist pDistStatic
186  (
187  pMesh,
188  staticPatchSet_,
189  staticPointZones(),
190  points,
191  maxMotionDistance_
192  );
193 
194  // Update the clearance from the distance to piston field
195  correctClearance(pDistMoving);
196 
197  calcScale
198  (
199  pMesh,
200  pDistMoving,
201  pDistStatic,
202  movingFrozenLayerThickness_,
203  maxMotionDistance_,
204  staticFrozenLayerThickness_
205  );
206 
207  position0_ = position;
208  }
209 
210  const vector translationVector(displacement()*axis);
211  transformPoints(newPoints, translationVector);
212 }
213 
214 
216 (
217  const polyMeshMap& map
218 )
219 {
221 }
222 
223 
224 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static const Form max
Definition: VectorSpace.H:120
static const Form min
Definition: VectorSpace.H:121
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dictionary & dict() const
Return the dynamicMeshDict/mover dict.
Definition: fvMeshMover.H:144
fvMesh & mesh()
Return the fvMesh.
Definition: fvMeshMover.H:132
const multiValveEngine & meshMover_
Reference to engine mesh mover.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: movingObject.C:337
const labelHashSet & patchSet
Object patchSet.
pistonObject(const word &name, const multiValveEngine &engine, const dictionary &dict)
Construct from dictionary.
Definition: piston.C:95
static word pistonBowlName
Name of the piston bowl pointZone.
scalar position() const
Return the current piston position.
Definition: piston.C:119
scalar clearance() const
Return clearance estimate.
Definition: piston.C:141
scalar displacement() const
Return piston displacement for current time-step.
Definition: piston.C:126
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: piston.C:216
void updatePoints(pointField &)
update points due to piston motion
Definition: piston.C:158
scalar speed() const
Return piston position for current time-step.
Definition: piston.C:134
A mesh mover using explicit node translation based on scaled distance functions per moving object....
Calculates the distance to the specified sets of patch and pointZone points or for all points.
Definition: pointDist.H:54
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
A class for handling words, derived from string.
Definition: word.H:62
label patchi
const pointField & points
List< label > labelList
A List of labels.
Definition: labelList.H:56
void transformPoints(vectorField &, const spatialTransform &, const vectorField &)
Transform given vectorField of coordinates with the given spatialTransform.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213