Test-PatchEdgeFaceWave.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 PatchEdgeFaceWave.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "argList.H"
30 #include "Time.H"
31 #include "fvMesh.H"
32 #include "volFields.H"
33 #include "PatchEdgeFaceWave.H"
34 #include "patchEdgeFaceInfo.H"
35 #include "patchPatchDist.H"
36 
37 using namespace Foam;
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 int main(int argc, char *argv[])
42 {
43  argList::validArgs.append("patch");
44 
45  #include "setRootCase.H"
46  #include "createTime.H"
47  #include "createMesh.H"
48 
50 
51  // Get name of patch
52  const word patchName = args[1];
53  const polyPatch& patch = patches[patchName];
54 
55  // 1. Walk from a single edge
56  {
57  // Data on all edges and faces
58  List<patchEdgeFaceInfo> allEdgeInfo(patch.nEdges());
59  List<patchEdgeFaceInfo> allFaceInfo(patch.size());
60 
61  // Initial seed
62  DynamicList<label> initialEdges;
63  DynamicList<patchEdgeFaceInfo> initialEdgesInfo;
64 
65 
66  // Just set an edge on the master
67  if (Pstream::master())
68  {
69  label edgeI = 0;
70  Info<< "Starting walk on edge " << edgeI << endl;
71 
72  initialEdges.append(edgeI);
73  const edge& e = patch.edges()[edgeI];
74  initialEdgesInfo.append
75  (
77  (
78  e.centre(patch.localPoints()),
79  0.0
80  )
81  );
82  }
83 
84 
85  // Walk
87  <
90  > calc
91  (
92  mesh,
93  patch,
94  initialEdges,
95  initialEdgesInfo,
96  allEdgeInfo,
97  allFaceInfo,
98  returnReduce(patch.nEdges(), sumOp<label>())
99  );
100 
101 
102  // Extract as patchField
103  volScalarField vsf
104  (
105  IOobject
106  (
107  "patchDist",
108  runTime.timeName(),
109  mesh,
112  ),
113  mesh,
114  dimensionedScalar("patchDist", dimLength, 0.0)
115  );
116  scalarField pf(vsf.boundaryField()[patch.index()].size());
117  forAll(pf, faceI)
118  {
119  pf[faceI] = Foam::sqrt(allFaceInfo[faceI].distSqr());
120  }
121  vsf.boundaryField()[patch.index()] = pf;
122 
123  Info<< "Writing patchDist volScalarField to " << runTime.value()
124  << endl;
125 
126  vsf.write();
127  }
128 
129 
130  // 2. Use a wrapper to walk from all boundary edges on selected patches
131  {
132  labelHashSet otherPatchIDs(identity(mesh.boundaryMesh().size()));
133  otherPatchIDs.erase(patch.index());
134 
135  Info<< "Walking on patch " << patch.index()
136  << " from edges shared with patches " << otherPatchIDs
137  << endl;
138 
139  patchPatchDist pwd(patch, otherPatchIDs);
140 
141  // Extract as patchField
142  volScalarField vsf
143  (
144  IOobject
145  (
146  "otherPatchDist",
147  runTime.timeName(),
148  mesh,
151  ),
152  mesh,
153  dimensionedScalar("otherPatchDist", dimLength, 0.0)
154  );
155  vsf.boundaryField()[patch.index()] = pwd;
156 
157  Info<< "Writing otherPatchDist volScalarField to " << runTime.value()
158  << endl;
159 
160  vsf.write();
161  }
162 
163  Info<< "\nEnd\n" << endl;
164  return 0;
165 }
166 
167 
168 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
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::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
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::DynamicList< label >
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::PatchEdgeFaceWave
Wave propagation of information along patch. Every iteration information goes through one layer of fa...
Definition: PatchEdgeFaceWave.H:69
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::calc
void calc(const argList &args, const Time &runTime, const fvMesh &mesh)
Definition: Lambda2.C:38
patchPatchDist.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::primitivePatch
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Foam::primitivePatch.
Definition: primitivePatch.H:45
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
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
Foam::HashSet< label, Hash< label > >
PatchEdgeFaceWave.H
Foam::patchPatchDist
Like wallDist but calculates on a patch the distance to nearest neighbouring patches....
Definition: patchPatchDist.H:53
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
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
argList.H
Foam::HashTable::erase
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
patchEdgeFaceInfo.H
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
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
createMesh.H
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
createTime.H
Foam::patchEdgeFaceInfo
Definition: patchEdgeFaceInfo.H:57
main
int main(int argc, char *argv[])
Definition: Test-PatchEdgeFaceWave.C:41
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
Foam::patchIdentifier::index
label index() const
Return the index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133