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

twoDPointCorrector.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 Description
00025     Class applies a two-dimensional correction to mesh motion point field.
00026     The correction guarantees that the mesh does not get twisted during motion
00027     and thus introduce a third dimension into a 2-D problem.
00028 
00029 \*---------------------------------------------------------------------------*/
00030 
00031 #include "twoDPointCorrector.H"
00032 #include <OpenFOAM/polyMesh.H>
00033 #include <OpenFOAM/wedgePolyPatch.H>
00034 #include <OpenFOAM/emptyPolyPatch.H>
00035 
00036 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00037 
00038 namespace Foam
00039 {
00040 
00041 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00042 
00043 const scalar twoDPointCorrector::edgeOrthogonalityTol = 1.0 - 1e-4;
00044 
00045 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00046 
00047 void twoDPointCorrector::calcAddressing() const
00048 {
00049     // Find geometry normal
00050     planeNormalPtr_ = new vector(0, 0, 0);
00051     vector& pn = *planeNormalPtr_;
00052 
00053     bool isWedge = false;
00054 
00055     // Algorithm:
00056     // Attempt to find wedge patch and work out the normal from it.
00057     // If not found, find an empty patch with faces in it and use the
00058     // normal of the first face.  If neither is found, declare an
00059     // error.
00060 
00061     // Try and find a wedge patch
00062     const polyBoundaryMesh& patches = mesh_.boundaryMesh();
00063 
00064     forAll (patches, patchI)
00065     {
00066         if (isA<wedgePolyPatch>(patches[patchI]))
00067         {
00068             isWedge = true;
00069 
00070             pn = refCast<const wedgePolyPatch>(patches[patchI]).centreNormal();
00071 
00072             if (polyMesh::debug)
00073             {
00074                 Pout << "Found normal from wedge patch " << patchI;
00075             }
00076 
00077             break;
00078         }
00079     }
00080 
00081     // Try to find an empty patch with faces
00082     if (!isWedge)
00083     {
00084         forAll (patches, patchI)
00085         {
00086             if (isA<emptyPolyPatch>(patches[patchI]) && patches[patchI].size())
00087             {
00088                 pn = patches[patchI].faceAreas()[0];
00089 
00090                 if (polyMesh::debug)
00091                 {
00092                     Pout << "Found normal from empty patch " << patchI;
00093                 }
00094 
00095                 break;
00096             }
00097         }
00098     }
00099 
00100 
00101     if (mag(pn) < VSMALL)
00102     {
00103         FatalErrorIn
00104         (
00105             "twoDPointCorrector::twoDPointCorrector(const polyMesh& mesh, "
00106             "const vector& n)"
00107         )   << "Cannot determine normal vector from patches."
00108             << abort(FatalError);
00109     }
00110     else
00111     {
00112         pn /= mag(pn);
00113     }
00114 
00115     if (polyMesh::debug)
00116     {
00117         Pout << " twoDPointCorrector normal: " << pn << endl;
00118     }
00119 
00120     // Select edges to be included in check.
00121     normalEdgeIndicesPtr_ = new labelList(mesh_.nEdges());
00122     labelList& neIndices = *normalEdgeIndicesPtr_;
00123 
00124     const edgeList& meshEdges = mesh_.edges();
00125 
00126     const pointField& meshPoints = mesh_.points();
00127 
00128     label nNormalEdges = 0;
00129 
00130     forAll (meshEdges, edgeI)
00131     {
00132         vector edgeVector =
00133             meshEdges[edgeI].vec(meshPoints)/
00134             (meshEdges[edgeI].mag(meshPoints) + VSMALL);
00135 
00136         if (mag(edgeVector & pn) > edgeOrthogonalityTol)
00137         {
00138             // this edge is normal to the plane. Add it to the list
00139             neIndices[nNormalEdges++] = edgeI;
00140         }
00141     }
00142 
00143     neIndices.setSize(nNormalEdges);
00144 
00145     // Construction check: number of points in a read 2-D or wedge geometry
00146     // should be odd and the number of edges normal to the plane should be
00147     // exactly half the number of points
00148     if (!isWedge)
00149     {
00150         if (meshPoints.size() % 2 != 0)
00151         {
00152             WarningIn
00153             (
00154                 "twoDPointCorrector::twoDPointCorrector("
00155                 "const polyMesh& mesh, const vector& n)"
00156             )   << "the number of vertices in the geometry "
00157                 << "is odd - this should not be the case for a 2-D case. "
00158                 << "Please check the geometry."
00159                 << endl;
00160         }
00161 
00162         if (2*nNormalEdges != meshPoints.size())
00163         {
00164             WarningIn
00165             (
00166                 "twoDPointCorrector::twoDPointCorrector("
00167                 "const polyMesh& mesh, const vector& n)"
00168             )   << "The number of points in the mesh is "
00169                 << "not equal to twice the number of edges normal to the plane "
00170                 << "- this may be OK only for wedge geometries.\n"
00171                 << "    Please check the geometry or adjust "
00172                 << "the orthogonality tolerance.\n" << endl
00173                 << "Number of normal edges: " << nNormalEdges
00174                 << " number of points: " << meshPoints.size()
00175                 << endl;
00176         }
00177     }
00178 }
00179 
00180 
00181 void twoDPointCorrector::clearAddressing() const
00182 {
00183     deleteDemandDrivenData(planeNormalPtr_);
00184     deleteDemandDrivenData(normalEdgeIndicesPtr_);
00185 }
00186 
00187 
00188 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00189 
00190 twoDPointCorrector::twoDPointCorrector(const polyMesh& mesh)
00191 :
00192     mesh_(mesh),
00193     required_(mesh_.nGeometricD() == 2),
00194     planeNormalPtr_(NULL),
00195     normalEdgeIndicesPtr_(NULL)
00196 {}
00197 
00198 
00199 
00200 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00201 
00202 Foam::twoDPointCorrector::~twoDPointCorrector()
00203 {
00204     clearAddressing();
00205 }
00206 
00207 
00208 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00209 
00210 direction twoDPointCorrector::normalDir() const
00211 {
00212     const vector& pn = planeNormal();
00213 
00214     if (mag(pn.x()) >= edgeOrthogonalityTol)
00215     {
00216         return vector::X;
00217     }
00218     else if (mag(pn.y()) >= edgeOrthogonalityTol)
00219     {
00220         return vector::Y;
00221     }
00222     else if (mag(pn.z()) >= edgeOrthogonalityTol)
00223     {
00224         return vector::Z;
00225     }
00226     else
00227     {
00228         FatalErrorIn("direction twoDPointCorrector::normalDir() const")
00229             << "Plane normal not aligned with the coordinate system" << nl
00230             << "    pn = " << pn
00231             << abort(FatalError);
00232 
00233         return vector::Z;
00234     }
00235 }
00236 
00237 
00238 // Return plane normal
00239 const vector& twoDPointCorrector::planeNormal() const
00240 {
00241     if (!planeNormalPtr_)
00242     {
00243         calcAddressing();
00244     }
00245 
00246     return *planeNormalPtr_;
00247 }
00248 
00249 
00250 // Return indices of normal edges.
00251 const labelList& twoDPointCorrector::normalEdgeIndices() const
00252 {
00253     if (!normalEdgeIndicesPtr_)
00254     {
00255         calcAddressing();
00256     }
00257 
00258     return *normalEdgeIndicesPtr_;
00259 }
00260 
00261 
00262 void twoDPointCorrector::correctPoints(pointField& p) const
00263 {
00264     if (!required_) return;
00265 
00266     // Algorithm:
00267     // Loop through all edges. Calculate the average point position A for
00268     // the front and the back. Correct the position of point P (in two planes)
00269     // such that vectors AP and planeNormal are parallel
00270 
00271     // Get reference to edges
00272     const edgeList&  meshEdges = mesh_.edges();
00273 
00274     const labelList& neIndices = normalEdgeIndices();
00275     const vector& pn = planeNormal();
00276 
00277     forAll (neIndices, edgeI)
00278     {
00279         point& pStart = p[meshEdges[neIndices[edgeI]].start()];
00280 
00281         point& pEnd = p[meshEdges[neIndices[edgeI]].end()];
00282 
00283         // calculate average point position
00284         const point A = 0.5*(pStart + pEnd);
00285 
00286         // correct point locations
00287         pStart = A + pn*(pn & (pStart - A));
00288         pEnd = A + pn*(pn & (pEnd - A));
00289     }
00290 }
00291 
00292 
00293 void twoDPointCorrector::updateMesh()
00294 {
00295     clearAddressing();
00296 }
00297 
00298 
00299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00300 
00301 } // End namespace Foam
00302 
00303 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines