Merge branch 'upstream/OpenFOAM' into master
[freefoam.git] / src / finiteVolume / fvMesh / wallDist / wallDist.H
blob8e134bdd32e5fb846263040c8a29b5568ec47898
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 Class
26     Foam::wallDist
28 Description
29     Calculation of distance to nearest wall for all cells and boundary.
30     Uses meshWave to do actual calculation.
32     Distance correction:
34     if correctWalls = true:
35     For each cell with face on wall calculate the true nearest point
36     (by triangle decomposition) on that face and do that same for that face's
37     pointNeighbours. This will find the true nearest distance in almost all
38     cases. Only very skewed cells or cells close to another wall might be
39     missed.
41     For each cell with only point on wall the same is done except now it takes
42     the pointFaces() of the wall point to look for the nearest point.
44 Note
46     correct() : for now does complete recalculation. (which usually is
47     ok since mesh is smoothed). However for topology change where geometry
48     in most of domain does not change you could think of starting from the
49     old cell values. Tried but not done since:
50     - meshWave would have to be called with old cellInfo.
51       This is List\<wallInfo\> of nCells.
52     - cannot construct from distance (y_) only since we don't know a value
53       for origin_. (origin_ = GREAT already used to denote illegal value.)
54     - so we would have to store a List\<wallInfo\> which unfortunately does
55       not get resized/mapped automatically upon mesh changes.
57 SourceFiles
58     wallDist.C
60 \*---------------------------------------------------------------------------*/
62 #ifndef wallDist_H
63 #define wallDist_H
65 #include <finiteVolume/volFields.H>
66 #include <meshTools/cellDistFuncs.H>
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 namespace Foam
74 class fvMesh;
76 /*---------------------------------------------------------------------------*\
77                            Class wallDist Declaration
78 \*---------------------------------------------------------------------------*/
80 class wallDist
82     public volScalarField,
83     public cellDistFuncs
87 private:
89     // Private Member Data
91         //- Do accurate distance calculation for near-wall cells.
92         bool correctWalls_;
94         //- Number of unset cells and faces.
95         label nUnset_;
98     // Private Member Functions
100         //- Disallow default bitwise copy construct
101         wallDist(const wallDist&);
103         //- Disallow default bitwise assignment
104         void operator=(const wallDist&);
107 public:
109     // Constructors
111         //- Construct from mesh and flag whether or not to correct wall.
112         //  Calculate for all cells. correctWalls : correct wall (face&point)
113         //  cells for correct distance, searching neighbours.
114         wallDist(const fvMesh& mesh, bool correctWalls = true);
117     // Destructor
119         virtual ~wallDist();
122     // Member Functions
124         const volScalarField& y() const
125         {
126             return *this;
127         }
129         label nUnset() const
130         {
131             return nUnset_;
132         }
134         //- Correct for mesh geom/topo changes
135         virtual void correct();
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 } // End namespace Foam
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 #endif
147 // ************************ vim: set sw=4 sts=4 et: ************************ //