pistonBowlPoints.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) 2025 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 "pistonBowlPoints.H"
27 #include "multiValveEngine.H"
28 #include "pistonPointEdgeData.H"
29 #include "PointEdgeWave.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace zoneGenerators
37  {
40  (
44  );
45  }
46 }
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& name,
53  const polyMesh& mesh,
54  const dictionary& dict
55 )
56 :
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
70 {
71  const fvMesh& mesh = refCast<const fvMesh>(mesh_);
72 
74  (
75  refCast<const fvMeshMovers::multiValveEngine>(mesh.mover())
76  );
77 
79 
80  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
81 
82  // Find the maximum axis coordinate of the piston patch-set
83  // Assumes the piston moves in the positive axis direction
84  scalar maxZ = -great;
86  {
87  const label patchi = iter.key();
88  if (pbm[patchi].localPoints().size())
89  {
90  maxZ = max(maxZ, max(piston.axis & pbm[patchi].localPoints()));
91  }
92  }
93  reduce(maxZ, maxOp<scalar>());
94 
95 
96  // Search for points starting at the piston surface and stopping at maxZ
97  // using a PointEdgeWave
98 
99  label nPistonPatchPoints = 0;
100  forAllConstIter(labelHashSet, piston.patchSet, iter)
101  {
102  const label patchi = iter.key();
103  nPistonPatchPoints += pbm[patchi].meshPoints().size();
104  }
105 
106 
107  const pointField& points = mesh.points();
109 
110  // Set initial changed points to all the patch points(if patch present)
111  List<pistonPointEdgeData> pistonPatchPointData(nPistonPatchPoints);
112  labelList pistonPatchPoints(nPistonPatchPoints);
113 
114  // Add the patch points to the pistonPatchPointData
115  nPistonPatchPoints = 0;
116  forAllConstIter(labelHashSet, piston.patchSet, iter)
117  {
118  const label patchi = iter.key();
119  const labelList& mp = pbm[patchi].meshPoints();
120 
121  forAll(mp, ppi)
122  {
123  pistonPatchPoints[nPistonPatchPoints] = mp[ppi];
124  pistonPatchPointData[nPistonPatchPoints] = pistonPointEdgeData
125  (
126  true
127  );
128  nPistonPatchPoints++;
129  }
130  }
131 
132  // Point data for wave
133  List<pistonPointEdgeData> allPointData(mesh.nPoints());
134 
135  // Edge data for wave
136  List<pistonPointEdgeData> allEdgeData(mesh.nEdges());
137 
139  <
142  > patchCalc
143  (
144  mesh,
145  pistonPatchPoints,
146  pistonPatchPointData,
147 
148  allPointData,
149  allEdgeData,
150  mesh.globalData().nTotalPoints(), // max iterations
151  td
152  );
153 
154  // Create a labelHashSet of the point labels in the piston bowl
155  labelList pointIndices(mesh.nPoints());
156  label zpi = 0;
157  forAll(allPointData, pointi)
158  {
159  if (allPointData[pointi].inBowl())
160  {
161  pointIndices[zpi++] = pointi;
162  }
163  }
164  pointIndices.setSize(zpi);
165 
166  return zoneSet
167  (
168  new pointZone
169  (
170  dict_.found("name")
171  ? zoneName_
172  : piston.pistonBowlName,
173  pointIndices,
174  mesh_.pointZones(),
175  moveUpdate_,
176  true
177  )
178  );
179 }
180 
181 
182 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void setSize(const label)
Reset size of List.
Definition: List.C:281
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:89
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const labelHashSet & patchSet
Object patchSet.
static word pistonBowlName
Name of the piston bowl pointZone.
A mesh mover using explicit node translation based on scaled distance functions per moving object....
const pistonObject & piston
Piston object.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvMeshMover & mover() const
Return the mover function class.
Definition: fvMesh.C:1155
label nTotalPoints() const
Return total number of points in decomposed mesh. Not.
Class used to pass data into container.
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Named list of point indices representing a sub-set of the mesh faces.
Definition: pointZone.H:60
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
label nPoints() const
label nEdges() const
A class for handling words, derived from string.
Definition: word.H:62
Abstract base class for all zoneGenerators, providing runtime selection.
Definition: zoneGenerator.H:57
A zoneGenerator which selects the points within the engine piston bowl.
virtual zoneSet generate() const
Generate and return the zoneSet.
pistonBowlPoints(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Zone container returned by zoneGenerator::generate.
Definition: zoneSet.H:94
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
const pointField & points
const dimensionedScalar mp
Proton mass.
defineTypeNameAndDebug(cylinderHeadPoints, 0)
addToRunTimeSelectionTable(zoneGenerator, cylinderHeadPoints, dictionary)
Namespace for OpenFOAM.
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
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary dict