FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

patchWave.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "patchWave.H"
00027 #include <OpenFOAM/polyMesh.H>
00028 #include <meshTools/wallPoint.H>
00029 #include <OpenFOAM/MeshWave.H>
00030 #include <OpenFOAM/globalMeshData.H>
00031 
00032 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00033 
00034 void Foam::patchWave::setChangedFaces
00035 (
00036     const labelHashSet& patchIDs,
00037     labelList& changedFaces,
00038     List<wallPoint>& faceDist
00039 ) const
00040 {
00041     const polyMesh& mesh = cellDistFuncs::mesh();
00042 
00043     label nChangedFaces = 0;
00044 
00045     forAll(mesh.boundaryMesh(), patchI)
00046     {
00047         if (patchIDs.found(patchI))
00048         {
00049             const polyPatch& patch = mesh.boundaryMesh()[patchI];
00050 
00051             forAll(patch.faceCentres(), patchFaceI)
00052             {
00053                 label meshFaceI = patch.start() + patchFaceI;
00054 
00055                 changedFaces[nChangedFaces] = meshFaceI;
00056 
00057                 faceDist[nChangedFaces] =
00058                     wallPoint
00059                     (
00060                         patch.faceCentres()[patchFaceI],
00061                         0.0
00062                     );
00063 
00064                 nChangedFaces++;
00065             }
00066         }
00067     }
00068 }
00069 
00070 
00071 Foam::label Foam::patchWave::getValues(const MeshWave<wallPoint>& waveInfo)
00072 {
00073     const List<wallPoint>& cellInfo = waveInfo.allCellInfo();
00074     const List<wallPoint>& faceInfo = waveInfo.allFaceInfo();
00075 
00076     label nIllegal = 0;
00077 
00078     // Copy cell values
00079     distance_.setSize(cellInfo.size());
00080 
00081     forAll(cellInfo, cellI)
00082     {
00083         scalar dist = cellInfo[cellI].distSqr();
00084 
00085         if (cellInfo[cellI].valid())
00086         {
00087             distance_[cellI] = Foam::sqrt(dist);
00088         }
00089         else
00090         {
00091             distance_[cellI] = dist;
00092 
00093             nIllegal++;
00094         }
00095     }
00096 
00097     // Copy boundary values
00098     forAll(patchDistance_, patchI)
00099     {
00100         const polyPatch& patch = mesh().boundaryMesh()[patchI];
00101 
00102         // Allocate storage for patchDistance
00103         scalarField* patchDistPtr = new scalarField(patch.size());
00104 
00105         patchDistance_.set(patchI, patchDistPtr);
00106 
00107         scalarField& patchField = *patchDistPtr;
00108 
00109         forAll(patchField, patchFaceI)
00110         {
00111             label meshFaceI = patch.start() + patchFaceI;
00112 
00113             scalar dist = faceInfo[meshFaceI].distSqr();
00114 
00115             if (faceInfo[meshFaceI].valid())
00116             {
00117                 // Adding SMALL to avoid problems with /0 in the turbulence
00118                 // models
00119                 patchField[patchFaceI] = Foam::sqrt(dist) + SMALL;
00120             }
00121             else
00122             {
00123                 patchField[patchFaceI] = dist;
00124 
00125                 nIllegal++;
00126             }
00127         }
00128     }
00129     return nIllegal;
00130 }
00131 
00132 
00133 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00134 
00135 // Construct from components
00136 Foam::patchWave::patchWave
00137 (
00138     const polyMesh& mesh,
00139     const labelHashSet& patchIDs,
00140     const bool correctWalls
00141 )
00142 :
00143     cellDistFuncs(mesh),
00144     patchIDs_(patchIDs),
00145     correctWalls_(correctWalls),
00146     nUnset_(0),
00147     distance_(mesh.nCells()),
00148     patchDistance_(mesh.boundaryMesh().size())
00149 {
00150     patchWave::correct();
00151 }
00152 
00153 
00154 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00155 
00156 Foam::patchWave::~patchWave()
00157 {}
00158 
00159 
00160 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00161 
00162 // Correct for mesh geom/topo changes. Might be more intelligent in the
00163 // future (if only small topology change)
00164 void Foam::patchWave::correct()
00165 {
00166     //
00167     // Set initial changed faces: set wallPoint for wall faces to wall centre
00168     //
00169 
00170     // Count walls
00171     label nWalls = sumPatchSize(patchIDs_);
00172 
00173     List<wallPoint> faceDist(nWalls);
00174     labelList changedFaces(nWalls);
00175 
00176     // Set to faceDist information to facecentre on walls.
00177     setChangedFaces(patchIDs_, changedFaces, faceDist);
00178 
00179 
00180     //
00181     // Do calculate wall distance by 'growing' from faces.
00182     //
00183 
00184     MeshWave<wallPoint> waveInfo
00185     (
00186         mesh(),
00187         changedFaces,
00188         faceDist,
00189         mesh().globalData().nTotalCells()   // max iterations
00190     );
00191 
00192 
00193     //
00194     // Copy distance into return field
00195     //
00196 
00197     nUnset_ = getValues(waveInfo);
00198 
00199     //
00200     // Correct wall cells for true distance
00201     //
00202 
00203     if (correctWalls_)
00204     {
00205         Map<label> nearestFace(2 * nWalls);
00206 
00207         correctBoundaryFaceCells
00208         (
00209             patchIDs_,
00210             distance_,
00211             nearestFace
00212         );
00213 
00214         correctBoundaryPointCells
00215         (
00216             patchIDs_,
00217             distance_,
00218             nearestFace
00219         );
00220     }
00221 }
00222 
00223 
00224 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines