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 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
33 namespace 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 
62 Foam::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  if (debug && isA<fvMesh>(mesh))
164  {
165  const fvMesh& fvm = dynamicCast<const fvMesh&>(mesh);
166 
167  volScalarField cCuts
168  (
169  IOobject
170  (
171  "cuttingPlane.cellCuts",
172  fvm.time().timeName(),
173  fvm.time(),
176  false
177  ),
178  fvm,
180  );
181 
182  auto& cCutsFld = cCuts.primitiveFieldRef();
183 
184  for (const label celli : cellCuts)
185  {
186  cCutsFld[celli] = 1;
187  }
188 
189  Pout<< "Writing cut types:" << cCuts.objectPath() << endl;
190  cCuts.write();
191  }
192 
193 
194  return nFaceCuts;
195 }
196 
197 
198 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
cuttingPlane.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::intersectsFace
bool intersectsFace(const PackedList< 2 > &sides, const labelUList &indices)
Definition: cuttingPlaneCuts.C:41
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::PackedList< 2 >
Foam::cuttingSurfaceBase::debug
static int debug
Debug information.
Definition: cuttingSurfaceBase.H:181
Foam::UList< label >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::plane::plane
plane()
Construct zero-initialised.
Definition: planeI.H:30