Loading...
Searching...
No Matches
cuttingPlaneCuts.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) 2018-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "cuttingPlane.H"
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
32
33namespace Foam
34{
35 // Check for face/plane intersection based on crossings
36 // Took (-1,0,+1) from plane::sign and packed as (0,1,2).
37 // Now use for left shift to obtain (1,2,4).
38 //
39 // Test accumulated value for an intersection with the plane.
40 inline bool intersectsFace
41 (
42 const PackedList<2>& sides,
43 const labelUList& indices
44 )
45 {
46 unsigned accum = 0u;
47
48 for (const label pointi : indices)
49 {
50 accum |= (1u << sides[pointi]);
51 }
52
53 // Accumulated value 3,5,6,7 indicates an intersection
54 return (accum == 3 || accum >= 5);
55 }
56
57} // End namespace Foam
58
59
60// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61
62Foam::label Foam::cuttingPlane::calcCellCuts
63(
64 const primitiveMesh& mesh,
65 PackedList<2>& sides,
66 bitSet& cellCuts
67)
68{
69 const faceList& faces = mesh.faces();
70 const pointField& pts = mesh.points();
71
72 const label nCells = mesh.nCells();
73 const label nFaces = mesh.nFaces();
74 const label nInternFaces = mesh.nInternalFaces();
75
76
77 // Classify sides of plane (0=BACK, 1=ONPLANE, 2=FRONT) for each point
78 {
79 const plane& pln = *this;
80 const label len = pts.size();
81
82 sides.resize(len);
83
84 // From signed (-1,0,+1) to (0,1,2) for PackedList
85 for (label i=0; i < len; ++i)
86 {
87 sides.set(i, unsigned(1 + pln.sign(pts[i], SMALL)));
88 }
89 }
90
91
92 // Detect cells cuts from the face cuts
93
94 label nFaceCuts = 0;
95
96 // 1st face cut of cell
97 bitSet hasCut1(nCells);
98
99 // 2nd face cut of cell
100 bitSet hasCut2(nCells);
101
102 for (label facei = 0; facei < nInternFaces; ++facei)
103 {
104 if (intersectsFace(sides, faces[facei]))
105 {
106 const label own = mesh.faceOwner()[facei];
107 const label nei = mesh.faceNeighbour()[facei];
108
109 ++nFaceCuts;
110
111 if (!hasCut1.set(own))
112 {
113 hasCut2.set(own);
114 }
115 if (!hasCut1.set(nei))
116 {
117 hasCut2.set(nei);
118 }
119 }
120 }
121
122 for (label facei = nInternFaces; facei < nFaces; ++facei)
123 {
124 if (intersectsFace(sides, faces[facei]))
125 {
126 const label own = mesh.faceOwner()[facei];
127
128 ++nFaceCuts;
129
130 if (!hasCut1.set(own))
131 {
132 hasCut2.set(own);
133 }
134 }
135 }
136
137 hasCut1.clearStorage(); // Not needed now
138
139 if (cellCuts.size())
140 {
141 cellCuts.resize(nCells); // safety
142 cellCuts &= hasCut2; // restrict to cell subset
143
144 if (debug)
145 {
146 Pout<<"detected " << cellCuts.count() << "/" << nCells
147 << " cells cut, subsetted from "
148 << hasCut2.count() << "/" << nCells << " cells." << endl;
149 }
150 }
151 else
152 {
153 cellCuts = std::move(hasCut2);
154
155 if (debug)
156 {
157 Pout<<"detected " << cellCuts.count() << "/" << nCells
158 << " cells cut." << endl;
159 }
160 }
161
162
163 const fvMesh* fvMeshPtr = nullptr;
164 if (debug && (fvMeshPtr = isA<fvMesh>(mesh)) != nullptr)
165 {
166 auto tcellCutsDebug = volScalarField::New
167 (
168 "cuttingPlane.cellCuts",
170 *fvMeshPtr,
171 dimensionedScalar(dimless, Foam::zero{})
172 );
173 auto& cellCutsDebug = tcellCutsDebug.ref();
174
175 auto& fld = cellCutsDebug.primitiveFieldRef();
176
177 for (const label celli : cellCuts)
178 {
179 fld[celli] = 1;
180 }
181
182 Pout<< "Writing cut types:" << cellCutsDebug.objectPath() << endl;
183 cellCutsDebug.write();
184 }
185
186
187 return nFaceCuts;
188}
189
190
191// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition PackedList.H:146
static int debug
Debug information.
plane()
Construct zero-initialised.
Definition planeI.H:23
dynamicFvMesh & mesh
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
bool intersectsFace(const PackedList< 2 > &sides, const labelUList &indices)
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const pointField & pts