cuttingSurfaceCuts.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-------------------------------------------------------------------------------
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 "fvMesh.H"
30#include "volFields.H"
31#include "pointFields.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35void Foam::cuttingSurface::calcCellCuts
36(
37 const fvMesh& fvm,
38 scalarField& pointDist,
39 bitSet& cellCuts
40)
41{
42 const pointField& cc = fvm.C();
43 const pointField& pts = fvm.points();
44
45 const label nCells = fvm.nCells();
46 const label nPoints = fvm.nPoints();
47
48 // The point distances
49
50 List<pointIndexHit> nearest;
51 surfPtr_().findNearest
52 (
53 pts,
54 scalarField(nPoints, GREAT),
55 nearest
56 );
57
58 vectorField norms;
59 surfPtr_().getNormal(nearest, norms);
60
61 pointDist.resize(nPoints);
62
63 for (label i=0; i < nPoints; ++i)
64 {
65 const point diff(pts[i] - nearest[i].hitPoint());
66
67 // Normal distance
68 pointDist[i] = (diff & norms[i]);
69 }
70
71
72 if (cellCuts.size())
73 {
74 cellCuts.resize(nCells); // safety
75 }
76 else
77 {
78 cellCuts.resize(nCells, true);
79 }
80
81
82 // Check based on cell distance
83
84 surfPtr_().findNearest
85 (
86 cc,
87 scalarField(cc.size(), GREAT),
88 nearest
89 );
90
91
92 boundBox cellBb;
93
94 for (const label celli : cellCuts)
95 {
96 cellBb.clear();
97 cellBb.add(pts, fvm.cellPoints(celli));
98
99 if (!cellBb.contains(nearest[celli].hitPoint()))
100 {
101 cellCuts.unset(celli);
102 }
103 }
104
105
106 if (debug)
107 {
108 volScalarField cCuts
109 (
110 IOobject
111 (
112 "cuttingSurface.cellCuts",
113 fvm.time().timeName(),
114 fvm.time(),
117 false
118 ),
119 fvm,
121 );
122
123 for (const label celli : cellCuts)
124 {
125 cCuts[celli] = 1;
126 }
127
128 Pout<< "Writing cellCuts:" << cCuts.objectPath() << endl;
129 cCuts.write();
130
131 pointScalarField pDist
132 (
133 IOobject
134 (
135 "cuttingSurface.pointDistance",
136 fvm.time().timeName(),
137 fvm.time(),
140 false
141 ),
142 pointMesh::New(fvm),
144 );
145 pDist.primitiveFieldRef() = pointDist;
146
147 Pout<< "Writing point distance:" << pDist.objectPath() << endl;
148 pDist.write();
149 }
150}
151
152
153// ************************************************************************* //
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label nPoints() const
Number of points supporting patch faces.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static int debug
Debug information.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.