cuttingSurface.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 "cuttingSurface.H"
29#include "dictionary.H"
30#include "fvMesh.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
35(
36 const polyMesh& mesh,
37 const word& surfaceType,
38 const word& surfaceName
39)
40:
42 surfPtr_
43 (
45 (
46 surfaceType,
48 (
49 surfaceName, // name
50 mesh.time().constant(), // directory
51 "triSurface", // instance
52 mesh.time(), // registry
53 IOobject::MUST_READ,
54 IOobject::NO_WRITE
55 ),
57 )
58 )
59{}
60
61
63(
64 const word& defaultSurfaceName,
65 const polyMesh& mesh,
66 const dictionary& dict
67)
68:
70 surfPtr_
71 (
73 (
74 dict.get<word>("surfaceType"),
76 (
77 dict.getOrDefault("surfaceName", defaultSurfaceName),
78 mesh.time().constant(), // directory
79 "triSurface", // instance
80 mesh.time(), // registry
81 IOobject::MUST_READ,
82 IOobject::NO_WRITE
83 ),
84 dict
85 )
86 )
87{}
88
89
90// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91
93(
94 const primitiveMesh& mesh,
95 const bool triangulate,
96 bitSet&& cellIdLabels
97)
98{
99 const fvMesh& fvm = static_cast<const fvMesh&>(mesh);
100
101 Mesh::clear();
102 meshCells_.clear();
103
104 // Pre-populate with restriction
105 bitSet cellCuts(std::move(cellIdLabels));
106
107 if (cellCuts.size())
108 {
109 cellCuts.resize(mesh.nCells());
110 }
111
112 scalarField pointDist;
113 calcCellCuts(fvm, pointDist, cellCuts);
114
115
116 // Walk cell cuts to create faces
117
118 // Action #1:
119 // - Orient edge so it points in the positive gradient direction.
120 // - Edge intersection when it spans across point-distance == 0.
121 const auto edgeOrientIntersect =
122 [=](edge& e) -> bool
123 {
124 if (pointDist[e.last()] < pointDist[e.first()])
125 {
126 e.flip();
127 }
128
129 const scalar s0 = pointDist[e.first()];
130 const scalar s1 = pointDist[e.last()];
131
132 if
133 (
134 s0 > ROOTVSMALL // Edge is all positive
135 || s1 < ROOTVSMALL // Edge is all negative
136 || Foam::mag(s1-s0) < ROOTVSMALL
137 )
138 {
139 return false;
140 }
141
142 return true;
143 };
144
145 // Action #2:
146 // - The edge intersection alpha for point-distance == 0
147 // Return -1 for error.
148 // This is like the iso-fraction for an iso-surface.
149 const auto edgeAlphaIntersect =
150 [=](const edge& e) -> scalar
151 {
152 const scalar s0 = pointDist[e.first()];
153 const scalar s1 = pointDist[e.last()];
154 const scalar d = s1-s0;
155
156 return Foam::mag(d) < ROOTVSMALL ? -1 : (-s0/d);
157 };
158
159 walkCellCuts
160 (
161 mesh,
162 cellCuts,
163 edgeOrientIntersect,
164 edgeAlphaIntersect,
165 triangulate
166 );
167}
168
169
171{
172 os << " surface:" << surfaceName();
173 if (level)
174 {
175 os << " faces:" << Mesh::surfFaces().size()
176 << " points:" << Mesh::points().size();
177 }
178}
179
180
181// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
T & first()
Return the first element of the list.
Definition: UListI.H:202
T & last()
Return the last element of the list.
Definition: UListI.H:216
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
Description of cuts across cells.
Definition: cellCuts.H:113
Base for creating a MeshedSurface by performing some type of cell cutting/intersection.
Constructs a cutting surface through a mesh.
virtual void performCut(const primitiveMesh &mesh, const bool triangulate, bitSet &&cellIdLabels)
Cut mesh, restricted to a list of cells.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
scalar print()
Print to screen.
constant condensation/saturation model.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
label nCells() const noexcept
Number of mesh cells.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
volScalarField & e
Definition: createFields.H:11