pointPatchDist.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "pointPatchDist.H"
30 #include "externalPointEdgePoint.H"
31 #include "pointMesh.H"
32 #include "PointEdgeWave.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 Foam::pointPatchDist::pointPatchDist
37 (
38  const pointMesh& pMesh,
39  const labelHashSet& patchIDs,
40  const pointField& points
41 )
42 :
44  (
45  IOobject
46  (
47  "pointDistance",
48  pMesh.db().time().timeName(),
49  pMesh.db()
50  ),
51  pMesh,
52  dimensionedScalar("y", dimLength, GREAT)
53  ),
54  points_(points),
55  patchIDs_(patchIDs),
56  nUnset_(0)
57 {
58  correct();
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 {
72  const pointBoundaryMesh& pbm = mesh().boundary();
73 
74  label nPoints = 0;
75 
76  for (const label patchi : patchIDs_)
77  {
78  nPoints += pbm[patchi].meshPoints().size();
79  }
80 
82 
83  // Set initial changed points to all the patch points(if patch present)
86  nPoints = 0;
87 
88  for (const label patchi : patchIDs_)
89  {
90  // Retrieve the patch now we have its index in patches.
91 
92  const labelList& mp = pbm[patchi].meshPoints();
93 
94  forAll(mp, ppI)
95  {
96  label meshPointi = mp[ppI];
97  wallPoints[nPoints] = meshPointi;
99  (
100  td.points_[meshPointi],
101  0.0
102  );
103  nPoints++;
104  }
105  }
106 
107  // Current info on points
108  List<externalPointEdgePoint> allPointInfo(mesh()().nPoints());
109 
110  // Current info on edges
111  List<externalPointEdgePoint> allEdgeInfo(mesh()().nEdges());
112 
114  <
117  > wallCalc
118  (
119  mesh()(),
120  wallPoints,
121  wallInfo,
122 
123  allPointInfo,
124  allEdgeInfo,
125  mesh().globalData().nTotalPoints(), // max iterations
126  td
127  );
128 
129  pointScalarField& psf = *this;
130 
131 
132  forAll(allPointInfo, pointi)
133  {
134  if (allPointInfo[pointi].valid(td))
135  {
136  psf[pointi] = Foam::sqrt(allPointInfo[pointi].distSqr());
137  }
138  else
139  {
140  nUnset_++;
141  }
142  }
143 }
144 
145 
146 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
Foam::pointPatchDist::correct
void correct()
Correct for mesh geom/topo changes.
Definition: pointPatchDist.C:70
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::externalPointEdgePoint
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Definition: externalPointEdgePoint.H:63
Foam::PointEdgeWave
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:87
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:186
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::HashSet< label, Hash< label > >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:432
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::externalPointEdgePoint::trackingData
Class used to pass data into container.
Definition: externalPointEdgePoint.H:102
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< vector >
PointEdgeWave.H
externalPointEdgePoint.H
correct
fvOptions correct(rho)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::externalPointEdgePoint::trackingData::points_
const pointField & points_
Definition: externalPointEdgePoint.H:105
Foam::pointPatchDist::~pointPatchDist
virtual ~pointPatchDist()
Destructor.
Definition: pointPatchDist.C:64
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:50
pointPatchDist.H
Foam::pointScalarField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:51
Foam::pointBoundaryMesh
Foam::pointBoundaryMesh.
Definition: pointBoundaryMesh.H:55
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::GeometricField< scalar, pointPatchField, pointMesh >
Foam::wallPoints
For use with FaceCellWave. Determines topological distance to starting faces.
Definition: wallPoints.H:65
pointMesh.H