Merge branch 'master' of ssh://opencfd@repo.or.cz/srv/git/OpenFOAM-1.5.x
[OpenFOAM-1.5.x.git] / applications / test / testPointEdgeWave / testPointEdgeWave.C
blobed68787171b7831ae40b59289a68551344439aaa
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Description
26     Test pointEdgeWave.
28 \*---------------------------------------------------------------------------*/
30 #include "argList.H"
31 #include "Time.H"
32 #include "polyMesh.H"
33 #include "pointMesh.H"
34 #include "OSspecific.H"
35 #include "IFstream.H"
36 #include "pointEdgePoint.H"
37 #include "PointEdgeWave.H"
39 using namespace Foam;
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 int main(int argc, char *argv[])
45     argList::validArgs.append("patch");
47 #   include "setRootCase.H"
48 #   include "createTime.H"
49 #   include "createPolyMesh.H"
51     pointMesh pMesh(mesh);
53     const polyBoundaryMesh& patches = mesh.boundaryMesh();
55     // Get name of patch
56     word patchName(args.additionalArgs()[0]);
57     
58     // Find the label in patches by name.
59     label patchI = patches.findPatchID(patchName);
61     {
62         // Test whether any processor has patch
63         label maxPatchI = patchI;
65         reduce(maxPatchI, maxOp<label>());
67         if (maxPatchI == -1)
68         {
69             FatalErrorIn(args.executable())
70                 << "Cannot find patch named " << patchName << exit(FatalError);
71         }
72     }
75     // Set initial changed points to all the patch points(if patch present)
76     List<pointEdgePoint> wallInfo;
77     labelList wallPoints;
79     if (patchI != -1)
80     {
81         // Retrieve the patch now we have its index in patches.
82         const polyPatch& pp = mesh.boundaryMesh()[patchI];
84         wallPoints = pp.meshPoints();
86         wallInfo.setSize(pp.nPoints());
88         forAll(pp.localPoints(), ppI)
89         {
90             wallInfo[ppI] = pointEdgePoint(pp.localPoints()[ppI], 0.0);
91         }
92     }
94     // Current info on points
95     List<pointEdgePoint> allPointInfo(mesh.nPoints());
97     // Current info on edges
98     List<pointEdgePoint> allEdgeInfo(mesh.nEdges());
100     PointEdgeWave<pointEdgePoint> wallCalc
101     (
102         pMesh,
103         wallPoints,
104         wallInfo,
106         allPointInfo,
107         allEdgeInfo,
108         mesh.nPoints()  // max iterations
109     );
112     pointScalarField psf
113     (
114         IOobject
115         (
116             "wallDist",
117             runTime.timeName(),
118             mesh,
119             IOobject::NO_READ,
120             IOobject::AUTO_WRITE
121         ),
122         pMesh,
123         dimensionedScalar("wallDist", dimLength, 0.0)
124     );
126     forAll(allPointInfo, pointI)
127     {
128         psf[pointI] = Foam::sqrt(allPointInfo[pointI].distSqr());
129     }
131     Info<< "Writing wallDist pointScalarField to " << runTime.value()
132         << endl;
134     psf.write();
136     Info << nl << "End" << endl;
138     return 0;
142 // ************************************************************************* //