valve.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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const word& name,
34  const multiValveEngine& engine,
35  const dictionary& dict
36 )
37 :
38  movingObject(name, engine, dict),
39  minLift_(dict.lookup<scalar>("minLift", dimLength))
40 {}
41 
42 
43 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44 
46 (
47  const scalar theta
48 ) const
49 {
50  return motion_->value(theta);
51 }
52 
53 
55 {
56  return lift(meshMover_.userTime()) >= minLift_;
57 }
58 
59 
61 {
62  return max
63  (
64  lift(meshMover_.userTime()),
65  minLift_
66  );
67 }
68 
69 
71 {
72  return
73  -(
74  lift()
75  - max
76  (
77  lift(meshMover_.userTime() - meshMover_.userDeltaT()),
78  minLift_
79  )
80  )/(meshMover_.mesh().time().deltaTValue() + vSmall);
81 }
82 
83 
84 Foam::scalar
86 {
87  return
88  lift(meshMover_.userTime() - meshMover_.userDeltaT())
89  - lift(meshMover_.userTime());
90 }
91 
92 
94 (
95  pointField& newPoints
96 )
97 {
98  // Update points only if valve is moving
99  if (isOpen())
100  {
101  const scalar position = this->lift();
102 
103  // Update a cached scale_ field if needed
104  if
105  (
106  executionCount_ == 0
107  || mag(position - position0_) > travelInterval_
108  )
109  {
110  const pointMesh& pMesh = pointMesh::New(meshMover_.mesh());
111  const pointField& points(meshMover_.mesh().points());
112 
113  const pointDist pDistMoving
114  (
115  pMesh,
116  patchSet,
117  movingPointZones(),
118  points,
119  maxMotionDistance_
120  );
121 
122  const pointDist pDistStatic
123  (
124  pMesh,
125  staticPatchSet_,
126  staticPointZones(),
127  points,
128  maxMotionDistance_
129  );
130 
131  calcScale
132  (
133  pMesh,
134  pDistMoving,
135  pDistStatic,
136  movingFrozenLayerThickness_,
137  maxMotionDistance_,
138  staticFrozenLayerThickness_
139  );
140 
141  position0_ = position;
142  }
143 
144  const vector translationVector(displacement()*axis);
145  transformPoints(newPoints, translationVector);
146  }
147 }
148 
149 
150 // ************************************************************************* //
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
bool isOpen() const
Is the valve open?
Definition: valve.C:54
valveObject(const word &name, const multiValveEngine &engine, const dictionary &dict)
Construct from dictionary.
Definition: valve.C:32
scalar displacement() const
Return valve displacement for current time-step.
Definition: valve.C:85
void updatePoints(pointField &)
update points due to valve motion
Definition: valve.C:94
scalar speed() const
Return current valve speed.
Definition: valve.C:70
scalar lift() const
Return current valve position.
Definition: valve.C:60
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
A class for handling words, derived from string.
Definition: word.H:62
const pointField & points
void transformPoints(vectorField &, const spatialTransform &, const vectorField &)
Transform given vectorField of coordinates with the given spatialTransform.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict