initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / fvMesh / wallDist / nearWallDistNoSearch.C
blob1e17f1d7c71a51da5dd3005f49d53c3fd5ce782e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 \*---------------------------------------------------------------------------*/
27 #include "nearWallDistNoSearch.H"
28 #include "fvMesh.H"
29 #include "wallPoint.H"
30 #include "wallFvPatch.H"
31 #include "surfaceFields.H"
33 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
35 void Foam::nearWallDistNoSearch::doAll()
37     const volVectorField& cellCentres = mesh_.C();
38     const fvPatchList& patches = mesh_.boundary();
40     forAll(patches, patchI)
41     {
42         fvPatchScalarField& ypatch = operator[](patchI);
44         if (patches[patchI].type() == wallFvPatch::typeName)
45         {
46             const unallocLabelList& faceCells = patches[patchI].faceCells();
48             const fvPatchVectorField& patchCentres
49                 = cellCentres.boundaryField()[patchI];
51             const fvsPatchVectorField& Apatch
52                 = mesh_.Sf().boundaryField()[patchI];
54             const fvsPatchScalarField& magApatch
55                 = mesh_.magSf().boundaryField()[patchI];
57             forAll(patchCentres, facei)
58             {
59                 ypatch[facei] =
60                 (
61                     Apatch[facei] &
62                     (
63                         patchCentres[facei]
64                       - cellCentres[faceCells[facei]]
65                     )
66                 )/magApatch[facei];
67             }
68         }
69         else
70         {
71             ypatch = 0.0;
72         }
73     }
77 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
79 Foam::nearWallDistNoSearch::nearWallDistNoSearch(const Foam::fvMesh& mesh)
81     volScalarField::GeometricBoundaryField
82     (
83         mesh.boundary(),
84         mesh.V(),           // Dummy internal field
85         calculatedFvPatchScalarField::typeName
86     ),
87     mesh_(mesh)
89     doAll();
93 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
95 Foam::nearWallDistNoSearch::~nearWallDistNoSearch()
99 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
101 void Foam::nearWallDistNoSearch::correct()
103     if (mesh_.changing())
104     {
105         // Update size of GeometricBoundaryField
106         forAll(mesh_.boundary(), patchI)
107         {
108             operator[](patchI).setSize(mesh_.boundary()[patchI].size());
109         }
110     }
112     doAll();
116 // ************************************************************************* //