Loading...
Searching...
No Matches
marchingCells.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2020-2023 PCOpt/NTUA
9 Copyright (C) 2020-2023 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "marchingCells.H"
31#include "cellSet.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
37 defineTypeNameAndDebug(marchingCells, 1);
38}
39
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
44{
45 label nChangedFaces(0);
46 labelList changedFaces(mesh_.nFaces(), -1);
47
48 // Gather initial faces from the seeding patches
49 for (const label patchI : seedPatches_)
50 {
51 const fvPatch& patch = mesh_.boundary()[patchI];
52 const label start = patch.start();
53 forAll(patch, faceI)
54 {
55 changedFaces[nChangedFaces++] = start + faceI;
56 }
57 }
58 // Gather the boundary faces from the seeding cell zones
59 for (label cellZoneID : seedCellZoneIDs_)
60 {
61 const labelList& zoneCells = mesh_.cellZones()[cellZoneID];
62 // Create sub mesh from the given zone
63 fvMeshSubset subSetMesh(mesh_, zoneCells);
64 const fvMesh& subMesh = subSetMesh.subMesh();
65 // Get the exposed faces of the subMesh and use them as seeds to
66 // patchWave
67 const word patchName(fvMeshSubset::exposedPatchName);
68 const label patchID(subMesh.boundaryMesh().findPatchID(patchName));
69 const polyPatch& subMeshPatch = subMesh.boundaryMesh()[patchID];
70 // Map to faces of the original mesh
71 const labelList& faceMap = subSetMesh.faceMap();
72 const label start = subMeshPatch.start();
73 // Go through a primitivePatch field (faceCentres) since the type
74 // of the patch containing the exposed faces is empty and a zero size
75 // is returned for the fields accessed through fvPatch
76 forAll(subMeshPatch.faceCentres(), faceI)
77 {
78 changedFaces[nChangedFaces++] = faceMap[start + faceI];
79 }
80 }
81 // Gather faces from the seeding face zones
82 for (label faceZoneID : seedFaceZoneIDs_)
83 {
84 const labelList& zoneFaces = mesh_.faceZones()[faceZoneID];
85 for (label faceI : zoneFaces)
86 {
87 changedFaces[nChangedFaces++] = faceI;
88 }
89 }
90 changedFaces.setSize(nChangedFaces);
91 List<wallPointData<bool>> changedFacesInfo(nChangedFaces);
92 const vectorField& faceCentres = mesh_.faceCentres();
93 forAll(changedFaces, fI)
94 {
95 changedFacesInfo[fI] =
96 wallPointData<bool>(faceCentres[changedFaces[fI]], true, 0.0);
97 }
98
99 meshWave_.setFaceInfo(changedFaces, changedFacesInfo);
100
101 // Initialization has been completed
102 initialised_ = true;
103}
104
105
106void Foam::marchingCells::appendSeedCell(const label cellID)
107{
108 if (!isFixedCell_[cellID])
110 isActiveCell_[cellID] = true;
111 activeCells_.append(cellID);
112 }
113}
114
115
117(
118 label nVisited,
119 const label cI,
120 labelList& newlyAddedCells
121)
122{
123 ++nVisited;
124 if (nVisited < marchingStep_ + 1)
125 {
126 const labelList& cellCells = mesh_.cellCells()[cI];
127 for (label neiCellI : cellCells)
128 {
129 if (!isFixedCell_[neiCellI])
130 {
131 if (!isActiveCell_[neiCellI])
132 {
133 isActiveCell_[neiCellI] = true;
134 newlyAddedCells.append(neiCellI);
135 }
136 march(nVisited, neiCellI, newlyAddedCells);
137 }
139 }
140}
141
142
143// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
144
145Foam::marchingCells::marchingCells
146(
147 const fvMesh& mesh,
148 const dictionary& dict
149)
150:
151 mesh_(mesh),
152 seedPatches_
153 (
154 mesh_.boundaryMesh().patchSet
155 (
156 dict.getOrDefault<wordRes>("seedPatches", wordRes(0))
157 )
158 ),
159 seedCellZoneIDs_
160 (
161 mesh_.cellZones().indices
162 (
163 dict.getOrDefault<wordRes>("seedCellZones", wordRes(0))
164 )
165 ),
166 seedFaceZoneIDs_
167 (
168 mesh_.faceZones().indices
169 (
170 dict.getOrDefault<wordRes>("seedFaceZones", wordRes(0))
171 )
172 ),
173 marchingStep_(dict.get<label>("marchingStep")),
174 isActiveCell_(mesh_.nCells(), false),
175 isFixedCell_(mesh_.nCells(), false),
176 activeCells_(0),
177 addedCells_(0),
178 initialised_(false),
179 nIters_(0),
180 allFaceInfo_(mesh_.nFaces()),
183{}
184
185
186// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187
188void Foam::marchingCells::update(const label iters)
189{
190 // If the first seeding has been performed, do it now
191 if (!initialised_)
192 {
193 initialise();
194 }
195
196 // Iterate the meshWave algorithm
197 meshWave_.iterate(iters*marchingStep_);
198
199 // Grab the newly added cells
200 addedCells_ = labelList(mesh_.nCells(), -1);
201 label nAddedCells(0);
202 forAll(allCellInfo_, cI)
203 {
204 if (allCellInfo_[cI].data() && !isActiveCell_[cI] && !isFixedCell_[cI])
205 {
206 addedCells_[nAddedCells++] = cI;
207 isActiveCell_[cI] = true;
208 }
209 }
210 addedCells_.setSize(nAddedCells);
211
212 // Add newly found cells to the list of activeCells
213 activeCells_.append(addedCells_);
214
215 // Write cell set
216 if (debug)
217 {
218 cellSet activeCellList
219 (
220 mesh_,
221 "activeCells" + name(nIters_),
222 activeCells_.size()
223 );
224 for (label cellI : activeCells_)
225 {
226 activeCellList.insert(cellI);
227 }
228 activeCellList.write();
229 }
230
231 nIters_ += iters;
232}
233
234
236(
238 const labelList& fixedCellZoneIDs
239)
240{
241 for (const label cellZoneID : fixedCellZoneIDs)
242 {
243 for (const label cI : cellZoneMesh[cellZoneID])
244 {
245 isFixedCell_[cI] = true;
246 isActiveCell_[cI] = false;
247 }
248 }
249}
250
251
252void Foam::marchingCells::addFixedCells(const labelList& fixedCellIDs)
253{
254 for (const label cI : fixedCellIDs)
255 {
256 isFixedCell_[cI] = true;
257 isActiveCell_[cI] = false;
258 }
259}
260
261
262// ************************************************************************* //
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A collection of cell labels.
Definition cellSet.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
static word exposedPatchName
Name for exposed internal faces (default: oldInternalFaces).
const labelList & faceMap() const
Return face map.
const fvMesh & subMesh() const
Return reference to subset mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
FaceCellWave< wallPointData< bool > > meshWave_
Engine propagating the active cells.
List< wallPointData< bool > > allFaceInfo_
Information for all faces.
void appendSeedCell(const label cellID)
Append cell to seed cells.
List< wallPointData< bool > > allCellInfo_
Information for all cells.
void initialise()
Initialise the active cells from the seeding patches.
labelList addedCells_
Which are the added cells.
const fvMesh & mesh_
bool initialised_
Has the initial seeding been conducted.
label nIters_
Iterations conducted thus far.
boolList isFixedCell_
Should this cell remain incative.
label marchingStep_
Marching step.
boolList isActiveCell_
Whether each cell is curently active or not.
void addFixedCells(const cellZoneMesh &cellZoneMesh, const labelList &fixedCellZoneIDs)
Add fixed cells through cellZone IDs.
labelList seedFaceZoneIDs_
Face zones used as seeds in the marching algorithm.
void update(const label iters=1)
Update active cells.
DynamicList< label > activeCells_
Which are the active cells.
labelHashSet seedPatches_
Patches used as seeds in the marching algorithm.
void march(label nVisited, const label cI, labelList &newlyAddedCells)
labelList seedCellZoneIDs_
Cell zones, the boundary faces of which are used as seeds in the marching algorithm.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const vectorField::subField faceCentres() const
Return face centres.
Definition polyPatch.C:320
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
label nFaces() const noexcept
Number of mesh faces.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Holds information (coordinate and normal) regarding nearest wall point.
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
auto & name
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< label > labelList
A List of labels.
Definition List.H:62
DynamicID< cellZoneMesh > cellZoneID
Foam::cellZoneID.
Definition ZoneIDs.H:38
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
Field< vector > vectorField
Specialisation of Field<T> for vector.
DynamicID< faceZoneMesh > faceZoneID
Foam::faceZoneID.
Definition ZoneIDs.H:43
static bool initialised_(false)
loopControl iters(runTime, aMesh.solutionDict(), "solution")
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299