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

primitiveMeshFaceCentresAndAreas.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     Calulate the face centres and areas.
00026 
00027     Calculate the centre by breaking the face into triangles using the face
00028     centre and area-weighted averaging their centres.  This method copes with
00029     small face-concavity.
00030 
00031 \*---------------------------------------------------------------------------*/
00032 
00033 #include "primitiveMesh.H"
00034 
00035 
00036 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00037 
00038 void Foam::primitiveMesh::calcFaceCentresAndAreas() const
00039 {
00040     if (debug)
00041     {
00042         Pout<< "primitiveMesh::calcFaceCentresAndAreas() : "
00043             << "Calculating face centres and face areas"
00044             << endl;
00045     }
00046 
00047     // It is an error to attempt to recalculate faceCentres
00048     // if the pointer is already set
00049     if (faceCentresPtr_ || faceAreasPtr_)
00050     {
00051         FatalErrorIn("primitiveMesh::calcFaceCentresAndAreas() const")
00052             << "Face centres or face areas already calculated"
00053             << abort(FatalError);
00054     }
00055 
00056     faceCentresPtr_ = new vectorField(nFaces());
00057     vectorField& fCtrs = *faceCentresPtr_;
00058 
00059     faceAreasPtr_ = new vectorField(nFaces());
00060     vectorField& fAreas = *faceAreasPtr_;
00061 
00062     makeFaceCentresAndAreas(points(), fCtrs, fAreas);
00063 
00064     if (debug)
00065     {
00066         Pout<< "primitiveMesh::calcFaceCentresAndAreas() : "
00067             << "Finished calculating face centres and face areas"
00068             << endl;
00069     }
00070 }
00071 
00072 
00073 void Foam::primitiveMesh::makeFaceCentresAndAreas
00074 (
00075     const pointField& p,
00076     vectorField& fCtrs,
00077     vectorField& fAreas
00078 ) const
00079 {
00080     const faceList& fs = faces();
00081 
00082     forAll (fs, facei)
00083     {
00084         const labelList& f = fs[facei];
00085         label nPoints = f.size();
00086 
00087         // If the face is a triangle, do a direct calculation for efficiency
00088         // and to avoid round-off error-related problems
00089         if (nPoints == 3)
00090         {
00091             fCtrs[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
00092             fAreas[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
00093         }
00094         else
00095         {
00096             vector sumN = vector::zero;
00097             scalar sumA = 0.0;
00098             vector sumAc = vector::zero;
00099 
00100             point fCentre = p[f[0]];
00101             for (label pi = 1; pi < nPoints; pi++)
00102             {
00103                 fCentre += p[f[pi]];
00104             }
00105 
00106             fCentre /= nPoints;
00107 
00108             for (label pi = 0; pi < nPoints; pi++)
00109             {
00110                 const point& nextPoint = p[f[(pi + 1) % nPoints]];
00111 
00112                 vector c = p[f[pi]] + nextPoint + fCentre;
00113                 vector n = (nextPoint - p[f[pi]])^(fCentre - p[f[pi]]);
00114                 scalar a = mag(n);
00115 
00116                 sumN += n;
00117                 sumA += a;
00118                 sumAc += a*c;
00119             }
00120 
00121             fCtrs[facei] = (1.0/3.0)*sumAc/(sumA + VSMALL);
00122             fAreas[facei] = 0.5*sumN;
00123         }
00124     }
00125 }
00126 
00127 
00128 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00129 
00130 const Foam::vectorField& Foam::primitiveMesh::faceCentres() const
00131 {
00132     if (!faceCentresPtr_)
00133     {
00134         calcFaceCentresAndAreas();
00135     }
00136 
00137     return *faceCentresPtr_;
00138 }
00139 
00140 
00141 const Foam::vectorField& Foam::primitiveMesh::faceAreas() const
00142 {
00143     if (!faceAreasPtr_)
00144     {
00145         calcFaceCentresAndAreas();
00146     }
00147 
00148     return *faceAreasPtr_;
00149 }
00150 
00151 
00152 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines