multiValveEngine.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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace fvMeshMovers
34 {
37 }
38 }
39 
41 (
42  "cylinderHead"
43 );
44 
45 
46 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
47 
49 Foam::fvMeshMovers::multiValveEngine::findLinerPatchSet() const
50 {
51  return mesh().boundaryMesh().patchSet(linerPatches_);
52 }
53 
54 
55 Foam::labelHashSet Foam::fvMeshMovers::multiValveEngine::findSlidingPatchSet()
56 {
57  return mesh().boundaryMesh().patchSet(slidingPatches_);
58 }
59 
60 
61 Foam::labelHashSet Foam::fvMeshMovers::multiValveEngine::findStaticPatchSet()
62 {
63  labelHashSet movingPatchSet(piston_.patchSet_);
64 
65  forAll(valves_, valvei)
66  {
67  movingPatchSet += valves_[valvei].patchSet_;
68  }
69 
70  labelHashSet staticPatchSet;
71 
72  forAll(mesh().boundaryMesh(), patchi)
73  {
74  const polyPatch& pp = mesh().boundaryMesh()[patchi];
75 
76  // Exclude non-static patches
77  if
78  (
79  !polyPatch::constraintType(pp.type())
80  && !slidingPatchSet_.found(pp.index())
81  && !movingPatchSet.found(pp.index())
82  )
83  {
84  staticPatchSet.insert(pp.index());
85  }
86  }
87 
88  return staticPatchSet;
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
95 (
96  fvMesh& mesh,
97  const dictionary& dict
98 )
99 :
100  fvMeshMover(mesh),
101  linerPatches_(dict.lookup("linerPatches")),
102  linerPatchSet_(findLinerPatchSet()),
103  slidingPatches_(dict.lookup("slidingPatches")),
104  slidingPatchSet_(findSlidingPatchSet()),
105  piston_("piston", *this, dict.subDict("piston")),
106  valves_(*this, dict.subOrEmptyDict("valves")),
107  staticPatchSet_(findStaticPatchSet()),
108  frozenPointZones_
109  (
110  dict.lookupOrDefault("frozenZones", wordReList::null())
111  ),
112  linerPatchSet(linerPatchSet_),
113  slidingPatchSet(slidingPatchSet_),
114  piston(piston_),
115  valves(valves_),
116  staticPatchSet(staticPatchSet_)
117 {}
118 
119 
120 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
121 
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 {
130  return mesh().time().userTimeValue();
131 }
132 
133 
135 {
136  return mesh().time().userDeltaTValue();
137 }
138 
139 
141 {
142  // Accumulate point motion from each moving object into newPoints field.
143  pointField newPoints(mesh().points());
144 
145  piston_.updatePoints(newPoints);
146 
147  forAll(valves_, valvei)
148  {
149  valves_[valvei].updatePoints(newPoints);
150  }
151 
152  // Update the mesh according to the newPoints field.
153  mesh().movePoints(newPoints);
154 
155  return true;
156 }
157 
158 
160 {
162 }
163 
164 
166 {
167  slidingPatchSet_ = findSlidingPatchSet();
168  linerPatchSet_ = findLinerPatchSet();
169  staticPatchSet_ = findStaticPatchSet();
170 
171  piston_.mapMesh(map);
172 
173  forAll(valves_, valvei)
174  {
175  valves_[valvei].mapMesh(map);
176  }
177 }
178 
179 
181 (
182  const polyDistributionMap&
183 )
184 {
186 }
187 
188 
189 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:819
scalar userDeltaTValue() const
Return user time step value.
Definition: Time.C:825
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for fvMesh movers.
Definition: fvMeshMover.H:53
fvMesh & mesh()
Return the fvMesh.
Definition: fvMeshMover.H:129
A mesh mover using explicit node translation based on scaled distance functions per moving object....
multiValveEngine(fvMesh &mesh, const dictionary &dict)
Construct from fvMesh.
scalar userTime() const
Return current user-time, CAD, s or ...
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
static word cylinderHeadName
Name of the cylinder head pointZone.
virtual bool update()
Update the mesh for both mesh motion and topology change.
scalar userDeltaT() const
Return current user-time-step, CAD, s, ...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:1244
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the patch set corresponding to the given names.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
label patchi
const pointField & points
defineTypeNameAndDebug(none, 0)
addToRunTimeSelectionTable(fvMeshMover, none, fvMesh)
Namespace for OpenFOAM.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
dictionary dict