1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Description
Test pointEdgeWave.
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "Time.H"
#include "polyMesh.H"
#include "pointMesh.H"
#include "OSspecific.H"
#include "IFstream.H"
#include "pointEdgePoint.H"
#include "PointEdgeWave.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
argList::validArgs.append("(patches)");
#include "setRootCase.H"
#include "createTime.H"
#include "createPolyMesh.H"
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
labelList patchIDs
(
pbm.patchSet(wordReList(IStringStream(args[1])())).sortedToc()
);
Info<< "Starting walk from patches "
<< UIndirectList<word>(pbm.names(), patchIDs)
<< nl
<< endl;
label nPoints = 0;
forAll(patchIDs, i)
{
nPoints += pbm[patchIDs[i]].nPoints();
}
Info<< "Seeding " << returnReduce(nPoints, sumOp<label>())
<< " patch points" << nl << endl;
// Set initial changed points to all the patch points(if patch present)
List<pointEdgePoint> wallInfo(nPoints);
labelList wallPoints(nPoints);
nPoints = 0;
forAll(patchIDs, i)
{
// Retrieve the patch now we have its index in patches.
const polyPatch& pp = pbm[patchIDs[i]];
forAll(pp.meshPoints(), ppI)
{
label meshPointi = pp.meshPoints()[ppI];
wallPoints[nPoints] = meshPointi;
wallInfo[nPoints] = pointEdgePoint(mesh.points()[meshPointi], 0.0);
nPoints++;
}
}
// Current info on points
List<pointEdgePoint> allPointInfo(mesh.nPoints());
// Current info on edges
List<pointEdgePoint> allEdgeInfo(mesh.nEdges());
PointEdgeWave<pointEdgePoint> wallCalc
(
mesh,
wallPoints,
wallInfo,
allPointInfo,
allEdgeInfo,
mesh.nPoints() // max iterations
);
pointScalarField psf
(
IOobject
(
"wallDist",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
pointMesh::New(mesh),
dimensionedScalar("wallDist", dimLength, 0.0)
);
forAll(allPointInfo, pointi)
{
psf[pointi] = Foam::sqrt(allPointInfo[pointi].distSqr());
}
Info<< "Writing wallDist pointScalarField to " << runTime.value()
<< endl;
psf.write();
Info<< "\nEnd\n" << endl;
return 0;
}
// ************************************************************************* //
|