insideSurface.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 "insideSurface.H"
27 #include "Time.H"
28 #include "polyMesh.H"
29 #include "searchableSurface.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace zoneGenerators
37  {
40  }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 template<class ZoneType, class UnaryOp>
48 (
49  const insideSurface& zoneGen,
50  const vectorField& pts,
51  const UnaryOp& uop
52 ) const
53 {
54  labelList indices(pts.size());
55 
56  List<volumeType> volType;
57  surfacePtr_->getVolumeType(pts, volType);
58 
59  label nInZone = 0;
60  forAll(volType, i)
61  {
62  if (uop(volType[i] == volumeType::inside))
63  {
64  indices[nInZone++] = i;
65  }
66  }
67 
68  indices.setSize(nInZone);
69 
70  return indices;
71 }
72 
73 
74 template<class ZoneType, class UnaryOp>
76 (
77  const insideSurface& zoneGen,
78  const zoneGeneratorList& zoneGenerators,
79  const vectorField& pts,
80  const UnaryOp& uop
81 ) const
82 {
83  if (!zoneGenerators.size())
84  {
85  return zoneGen.select<ZoneType>(zoneGen, pts, uop);
86  }
87 
88  boolList selected(pts.size(), false);
89 
90  forAll(zoneGenerators, zgi)
91  {
92  zoneSet zs(zoneGenerators[zgi].generate());
93 
94  const labelList& zsIndices(zs.zone<ZoneType>());
95 
96  const pointField zonePoints(pts, zsIndices);
97  List<volumeType> volType;
98  surfacePtr_->getVolumeType(zonePoints, volType);
99 
100  forAll(zsIndices, i)
101  {
102  if (uop(volType[i] == volumeType::inside))
103  {
104  selected[zsIndices[i]] = true;
105  }
106  }
107  }
108 
109  return indices(selected);
110 }
111 
112 
113 template<class UnaryOp>
115 (
116  const insideSurface& zoneGen,
117  const zoneGeneratorList& zoneGenerators,
118  const vectorField& pts,
119  boolList& flipMap,
120  const UnaryOp& uop
121 ) const
122 {
123  if (!zoneGenerators.size())
124  {
125  return select<faceZone>(zoneGen, pts, uop);
126  }
127 
128  bool oriented = true;
129  boolList selected(pts.size(), false);
130 
131  forAll(zoneGenerators, zgi)
132  {
133  zoneSet zs(zoneGenerators[zgi].generate());
134 
135  const faceZone& fZone = zs.fZone();
136  const labelList& zsIndices(fZone);
137 
138  const pointField zonePoints(pts, zsIndices);
139  List<volumeType> volType;
140  surfacePtr_->getVolumeType(zonePoints, volType);
141 
142  if (oriented && fZone.oriented())
143  {
144  flipMap.setSize(mesh_.nFaces(), false);
145  const boolList& zsFlipMap(fZone.flipMap());
146 
147  forAll(zsIndices, i)
148  {
149  if (uop(volType[i] == volumeType::inside))
150  {
151  selected[zsIndices[i]] = true;
152  flipMap[zsIndices[i]] = zsFlipMap[i];
153  }
154  }
155  }
156  else
157  {
158  oriented = false;
159  flipMap.clear();
160 
161  forAll(zsIndices, i)
162  {
163  if (uop(volType[i] == volumeType::inside))
164  {
165  selected[zsIndices[i]] = true;
166  }
167  }
168  }
169  }
170 
171  return indices(selected);
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
176 
178 (
179  const word& name,
180  const polyMesh& mesh,
181  const dictionary& dict
182 )
183 :
184  volume(name, mesh, dict),
185  surfacePtr_
186  (
188  (
189  word(dict.lookup("surface")),
190  IOobject
191  (
192  dict.lookupOrDefault
193  (
194  "surfaceName",
195  mesh.objectRegistry::db().name()
196  ),
197  mesh.time().constant(),
198  searchableSurface::geometryDir(mesh.time()),
199  mesh.objectRegistry::db(),
200  IOobject::MUST_READ,
201  IOobject::NO_WRITE
202  ),
203  dict
204  )
205  )
206 {}
207 
208 
209 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
210 
212 {}
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
218 {
219  return volume::generate(*this);
220 }
221 
222 
223 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
void setSize(const label)
Reset size of List.
Definition: List.C:281
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Named list of face indices representing a sub-set of the mesh faces.
Definition: faceZone.H:66
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:262
bool oriented() const
Return true if the faceZone is oriented, i.e. the flipMap is set.
Definition: faceZone.H:256
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
A class for handling words, derived from string.
Definition: word.H:62
List of zoneGenerators.
Abstract base class for all zoneGenerators, providing runtime selection.
Definition: zoneGenerator.H:57
static labelList indices(const boolList &selected)
Return the list of selected indices.
A zoneGenerator which selects points, cells or faces with centres either inside or outside a surface.
virtual zoneSet generate() const
Generate and return the zoneSet.
labelList select(const insideSurface &zoneGen, const vectorField &pts, const UnaryOp &uop) const
virtual ~insideSurface()
Destructor.
insideSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
A zoneGenerator which looks-up and returns a zoneSet containing point, and/or cell and/or faces zones...
Definition: lookup.H:139
Abstract zoneGenerator which selects points, cells or faces with centres either inside a volume.
Definition: volume.H:114
virtual zoneSet generate() const =0
Generate and return the zoneSet.
Zone container returned by zoneGenerator::generate.
Definition: zoneSet.H:94
const faceZone & fZone() const
Return a reference to the faceZone if allocated.
Definition: zoneSetI.H:230
const ZoneType & zone() const
Return a reference to the zone type.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dictionary dict