changed isA<wallPolyPatch> & isA<wallFvPatch> almost everywhere
[openfoam-extend-OpenFOAM-1.6-ext.git] / src / finiteVolume / fvMesh / wallDist / nearWallDistNoSearch.C
blobf2fa84f75a78e2037aff8fb98f64e7ea26718725
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 \*---------------------------------------------------------------------------*/
27 #include "nearWallDistNoSearch.H"
28 #include "fvMesh.H"
29 #include "wallPoint.H"
30 #include "surfaceFields.H"
32 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
34 void Foam::nearWallDistNoSearch::doAll()
36     const volVectorField& cellCentres = mesh_.C();
37     const fvPatchList& patches = mesh_.boundary();
39     forAll(patches, patchI)
40     {
41         fvPatchScalarField& ypatch = operator[](patchI);
43         if (patches[patchI].isWall())
44         {
45             const unallocLabelList& faceCells = patches[patchI].faceCells();
47             const fvPatchVectorField& patchCentres
48                 = cellCentres.boundaryField()[patchI];
50             const fvsPatchVectorField& Apatch
51                 = mesh_.Sf().boundaryField()[patchI];
53             const fvsPatchScalarField& magApatch
54                 = mesh_.magSf().boundaryField()[patchI];
56             forAll(patchCentres, facei)
57             {
58                 ypatch[facei] =
59                 (
60                     Apatch[facei] &
61                     (
62                         patchCentres[facei]
63                       - cellCentres[faceCells[facei]]
64                     )
65                 )/magApatch[facei];
66             }
67         }
68         else
69         {
70             ypatch = 0.0;
71         }
72     }
76 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
78 Foam::nearWallDistNoSearch::nearWallDistNoSearch(const Foam::fvMesh& mesh)
80     volScalarField::GeometricBoundaryField
81     (
82         mesh.boundary(),
83         mesh.V(),           // Dummy internal field
84         calculatedFvPatchScalarField::typeName
85     ),
86     mesh_(mesh)
88     doAll();
92 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
94 Foam::nearWallDistNoSearch::~nearWallDistNoSearch()
98 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
100 void Foam::nearWallDistNoSearch::correct()
102     if (mesh_.changing())
103     {
104         // Update size of GeometricBoundaryField
105         forAll(mesh_.boundary(), patchI)
106         {
107             operator[](patchI).setSize(mesh_.boundary()[patchI].size());
108         }
109     }
111     doAll();
115 // ************************************************************************* //