sampledCutPlane.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledCutPlane.H"
27 #include "emptyFvPatchFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace sampledSurfaces
35 {
38 
39  // Backwards compatible lookup as "cuttingPlane"
41  (
43  cutPlane,
44  word,
45  cuttingPlane
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
54 Foam::sampledSurfaces::cutPlane::calcIsoSurf() const
55 {
56  // Compute the distance from the mesh points to the plane
57  scalarField pointDistance(mesh().nPoints());
58  {
59  forAll(mesh().points(), pointi)
60  {
61  pointDistance[pointi] =
62  plane_.signedDistance(mesh().points()[pointi]);
63  }
64  }
65 
66  // Construct an iso-surface at the given distance
67  return autoPtr<cutPolyIsoSurface>
68  (
69  new cutPolyIsoSurface(mesh(), pointDistance, 0, zoneIDs())
70  );
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
77 (
78  const word& name,
79  const polyMesh& mesh,
80  const dictionary& dict
81 )
82 :
84  plane_(dict)
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 {
98  return timeIndex() == -1;
99 }
100 
101 
103 {
104  os << "cutPlane: " << name() << " :"
105  << " plane:" << plane_
106  << " faces:" << faces().size()
107  << " points:" << points().size();
108 }
109 
110 
111 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
scalar signedDistance(const point &p) const
Return signed distance from the given point to the plane.
Definition: plane.C:322
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
An abstract class for surfaces with sampling.
const polyMesh & mesh() const
Access to the underlying mesh.
A sampledSurface defined by a plane.
virtual ~cutPlane()
Destructor.
virtual bool needsUpdate() const
Does the surface need an update?
virtual void print(Ostream &) const
Write.
cutPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
A base class for sampled surfaces constructed from iso-surfaces.
virtual const pointField & points() const
Points of surface.
const labelList & zoneIDs() const
Access the zone indices.
A class for handling words, derived from string.
Definition: word.H:62
const pointField & points
label nPoints
addNamedToRunTimeSelectionTable(sampledSurface, cutPlane, word, cuttingPlane)
defineTypeNameAndDebug(cutPlane, 0)
addToRunTimeSelectionTable(sampledSurface, cutPlane, word)
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict