Test-PointEdgeWave.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 | Copyright (C) 2011-2015 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 Description
25  Test pointEdgeWave.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "argList.H"
30 #include "Time.H"
31 #include "polyMesh.H"
32 #include "pointMesh.H"
33 #include "OSspecific.H"
34 #include "IFstream.H"
35 #include "pointEdgePoint.H"
36 #include "PointEdgeWave.H"
37 
38 using namespace Foam;
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 int main(int argc, char *argv[])
43 {
44  argList::validArgs.append("(patches)");
45 
46  #include "setRootCase.H"
47  #include "createTime.H"
48  #include "createPolyMesh.H"
49 
50  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
51 
52  labelList patchIDs
53  (
54  pbm.patchSet(wordReList(IStringStream(args[1])())).sortedToc()
55  );
56 
57  Info<< "Starting walk from patches "
58  << UIndirectList<word>(pbm.names(), patchIDs)
59  << nl
60  << endl;
61 
62  label nPoints = 0;
63  forAll(patchIDs, i)
64  {
65  nPoints += pbm[patchIDs[i]].nPoints();
66  }
67 
68  Info<< "Seeding " << returnReduce(nPoints, sumOp<label>())
69  << " patch points" << nl << endl;
70 
71 
72  // Set initial changed points to all the patch points(if patch present)
73  List<pointEdgePoint> wallInfo(nPoints);
74  labelList wallPoints(nPoints);
75  nPoints = 0;
76 
77  forAll(patchIDs, i)
78  {
79  // Retrieve the patch now we have its index in patches.
80  const polyPatch& pp = pbm[patchIDs[i]];
81 
82  forAll(pp.meshPoints(), ppI)
83  {
84  label meshPointI = pp.meshPoints()[ppI];
85  wallPoints[nPoints] = meshPointI;
86  wallInfo[nPoints] = pointEdgePoint(mesh.points()[meshPointI], 0.0);
87  nPoints++;
88  }
89  }
90 
91  // Current info on points
92  List<pointEdgePoint> allPointInfo(mesh.nPoints());
93 
94  // Current info on edges
95  List<pointEdgePoint> allEdgeInfo(mesh.nEdges());
96 
98  (
99  mesh,
100  wallPoints,
101  wallInfo,
102 
103  allPointInfo,
104  allEdgeInfo,
105  mesh.nPoints() // max iterations
106  );
107 
108 
109  pointScalarField psf
110  (
111  IOobject
112  (
113  "wallDist",
114  runTime.timeName(),
115  mesh,
118  ),
120  dimensionedScalar("wallDist", dimLength, 0.0)
121  );
122 
123  forAll(allPointInfo, pointI)
124  {
125  psf[pointI] = Foam::sqrt(allPointInfo[pointI].distSqr());
126  }
127 
128  Info<< "Writing wallDist pointScalarField to " << runTime.value()
129  << endl;
130 
131  psf.write();
132 
133  Info<< "\nEnd\n" << endl;
134  return 0;
135 }
136 
137 
138 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::PointEdgeWave
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:85
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyMesh.H
Foam::primitiveMesh::nEdges
label nEdges() const
Definition: primitiveMeshI.H:41
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:527
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
pointEdgePoint.H
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:59
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
PointEdgeWave.H
argList.H
main
int main(int argc, char *argv[])
Definition: Test-PointEdgeWave.C:42
IFstream.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointEdgePoint
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Definition: pointEdgePoint.H:59
setRootCase.H
Foam::sumOp
Definition: ops.H:162
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:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
createTime.H
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:750
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
createPolyMesh.H
pointMesh.H