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-2022 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 "cuttingSurface.H"
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "pointFields.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void 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  for (const label celli : cellCuts)
93  {
94  if (!fvm.cellBb(celli).contains(nearest[celli].point()))
95  {
96  cellCuts.unset(celli);
97  }
98  }
99 
100  if (debug)
101  {
102  volScalarField cCuts
103  (
104  IOobject
105  (
106  "cuttingSurface.cellCuts",
107  fvm.time().timeName(),
108  fvm.time(),
112  ),
113  fvm,
115  );
116 
117  for (const label celli : cellCuts)
118  {
119  cCuts[celli] = 1;
120  }
121 
122  Pout<< "Writing cellCuts:" << cCuts.objectPath() << endl;
123  cCuts.write();
124 
125  pointScalarField pDist
126  (
127  IOobject
128  (
129  "cuttingSurface.pointDistance",
130  fvm.time().timeName(),
131  fvm.time(),
135  ),
136  pointMesh::New(fvm),
138  );
139  pDist.primitiveFieldRef() = pointDist;
140 
141  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
142  pDist.write();
143  }
144 }
145 
146 
147 // ************************************************************************* //
label nPoints() const
Number of points supporting patch faces.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
static int debug
Debug information.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
vector point
Point is a vector.
Definition: point.H:37
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127