cuttingPlane.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-2020 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 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 :
34  plane(pln)
35 {}
36 
37 
39 (
40  const plane& pln,
41  const primitiveMesh& mesh,
42  const bool triangulate,
43  const bitSet& cellIdLabels
44 )
45 :
46  plane(pln)
47 {
48  performCut(mesh, triangulate, cellIdLabels);
49 }
50 
51 
53 (
54  const plane& pln,
55  const primitiveMesh& mesh,
56  const bool triangulate,
57  bitSet&& cellIdLabels
58 )
59 :
60  plane(pln)
61 {
62  performCut(mesh, triangulate, cellIdLabels);
63 }
64 
65 
67 (
68  const plane& pln,
69  const primitiveMesh& mesh,
70  const bool triangulate,
71  const labelUList& cellIdLabels
72 )
73 :
74  plane(pln)
75 {
76  performCut(mesh, triangulate, cellIdLabels);
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
83 (
84  const primitiveMesh& mesh,
85  const bool triangulate,
86  bitSet&& cellIdLabels
87 )
88 {
89  const plane& pln = *this;
90  const pointField& pts = mesh.points();
91 
92  Mesh::clear();
93  meshCells_.clear();
94 
95  // Pre-populate with restriction
96  bitSet cellCuts(std::move(cellIdLabels));
97 
98  if (cellCuts.size())
99  {
100  cellCuts.resize(mesh.nCells());
101  }
102 
103  // For each mesh point, the encoded side (0,1,2) of the plane
104  PackedList<2> sides;
105 
106  // Determine cells that are (likely) cut
107  // - some ambiguity when plane is exactly between cells
108  const label nFaceCuts = calcCellCuts(mesh, sides, cellCuts);
109 
110 
111  // Walk cell cuts to create faces
112 
113  // Action #1:
114  // - Orient edge so it points in the positive normal direction.
115  // - Edge/plane intersection when the sign changes
116  const auto edgeOrientIntersect =
117  [=](edge& e) -> bool
118  {
119  if (sides[e.last()] < sides[e.first()])
120  {
121  e.flip();
122  }
123 
124  return sides[e.first()] != sides[e.last()];
125  };
126 
127  // Action #2:
128  // - The edge intersection alpha
129  const auto edgeAlphaIntersect =
130  [=](const edge& e) -> scalar
131  {
132  return pln.lineIntersect(e.line(pts));
133  };
134 
135  walkCellCuts
136  (
137  mesh,
138  cellCuts,
139  edgeOrientIntersect,
140  edgeAlphaIntersect,
141  triangulate,
142  nFaceCuts
143  );
144 }
145 
146 
147 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
148 
150 {
151  if (this == &rhs)
152  {
153  return; // Self-assignment is a no-op
154  }
155 
156  static_cast<Mesh&>(*this) = rhs;
157  static_cast<plane&>(*this) = rhs;
158  meshCells_ = rhs.meshCells();
159 }
160 
161 
162 // ************************************************************************* //
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::cuttingPlane::cuttingPlane
cuttingPlane(const plane &pln)
Construct from a plane description without any cutting.
Definition: cuttingPlane.C:32
Foam::cuttingSurfaceBase::meshCells
const labelList & meshCells() const
The mesh cells cut.
Definition: cuttingSurfaceBase.H:197
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::cuttingPlane::performCut
virtual void performCut(const primitiveMesh &mesh, const bool triangulate, bitSet &&cellIdLabels)
Cut mesh, restricted to a list of cells.
Definition: cuttingPlane.C:83
Foam::cuttingPlane
Constructs cutting plane through a mesh.
Definition: cuttingPlane.H:58
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
cuttingPlane.H
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field< vector >
Foam::plane::lineIntersect
scalar lineIntersect(const line< Point, PointRef > &l) const
Return the cutting point between the plane and.
Definition: plane.H:257
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
clear
patchWriters clear()
Foam::PackedList< 2 >
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::cuttingPlane::operator=
void operator=(const cuttingPlane &rhs)
Copy assignment.
Definition: cuttingPlane.C:149
Foam::MeshedSurface< face >
Foam::cellCuts
Description of cuts across cells.
Definition: cellCuts.H:110
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78