ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / test / PatchEdgeFaceWave / Test-PatchEdgeFaceWave.C
blob8d37f8fe0130d17a177edd665323218b564d6d2c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Description
25     Test PatchEdgeFaceWave.
27 \*---------------------------------------------------------------------------*/
29 #include "argList.H"
30 #include "Time.H"
31 #include "fvMesh.H"
32 #include "volFields.H"
33 #include "PatchEdgeFaceWave.H"
34 #include "patchEdgeFaceInfo.H"
36 using namespace Foam;
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 int main(int argc, char *argv[])
42     argList::validArgs.append("patch");
44 #   include "setRootCase.H"
45 #   include "createTime.H"
46 #   include "createMesh.H"
48     const polyBoundaryMesh& patches = mesh.boundaryMesh();
50     // Get name of patch
51     const word patchName = args[1];
53     const polyPatch& patch = patches[patchName];
55     // Data on all edges and faces    
56     List<patchEdgeFaceInfo> allEdgeInfo(patch.nEdges());
57     List<patchEdgeFaceInfo> allFaceInfo(patch.size());
59     // Initial seed
60     DynamicList<label> initialEdges;
61     DynamicList<patchEdgeFaceInfo> initialEdgesInfo;
64     // Just set an edge on the master
65     if (Pstream::master())
66     {
67         label edgeI = 0;
68         Info<< "Starting walk on edge " << edgeI << endl;
70         initialEdges.append(edgeI);
71         const edge& e = patch.edges()[edgeI];
72         initialEdgesInfo.append
73         (
74             patchEdgeFaceInfo
75             (
76                 e.centre(patch.localPoints()),
77                 0.0
78             )
79         );
80     }
83     // Walk
84     PatchEdgeFaceWave
85     <
86         primitivePatch,
87         patchEdgeFaceInfo
88     > calc
89     (
90         mesh,
91         patch,
92         initialEdges,
93         initialEdgesInfo,
94         allEdgeInfo,
95         allFaceInfo,
96         returnReduce(patch.nEdges(), sumOp<label>())
97     );
100     // Extract as patchField
101     volScalarField vsf
102     (
103         IOobject
104         (
105             "patchDist",
106             runTime.timeName(),
107             mesh,
108             IOobject::NO_READ,
109             IOobject::AUTO_WRITE
110         ),
111         mesh,
112         dimensionedScalar("patchDist", dimLength, 0.0)
113     );
114     scalarField pf(vsf.boundaryField()[patch.index()].size());
115     forAll(pf, faceI)
116     {
117         pf[faceI] = Foam::sqrt(allFaceInfo[faceI].distSqr());
118     }
119     vsf.boundaryField()[patch.index()] = pf;
121     Info<< "Writing patchDist volScalarField to " << runTime.value()
122         << endl;
124     vsf.write();
127     Info<< "\nEnd\n" << endl;
128     return 0;
132 // ************************************************************************* //