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

faceMapper.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 "faceMapper.H"
00027 #include <OpenFOAM/demandDrivenData.H>
00028 #include <OpenFOAM/polyMesh.H>
00029 #include <OpenFOAM/mapPolyMesh.H>
00030 
00031 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00032 
00033 
00034 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00035 
00036 void Foam::faceMapper::calcAddressing() const
00037 {
00038     if
00039     (
00040         directAddrPtr_
00041      || interpolationAddrPtr_
00042      || weightsPtr_
00043      || insertedFaceLabelsPtr_
00044     )
00045     {
00046         FatalErrorIn("void faceMapper::calcAddressing() const")
00047             << "Addressing already calculated."
00048             << abort(FatalError);
00049     }
00050 
00051     if (direct())
00052     {
00053         // Direct addressing, no weights
00054 
00055         directAddrPtr_ = new labelList(mpm_.faceMap());
00056         labelList& directAddr = *directAddrPtr_;
00057 
00058         // Reset the size of addressing list to contain only live faces
00059         directAddr.setSize(mesh_.nFaces());
00060 
00061         insertedFaceLabelsPtr_ = new labelList(mesh_.nFaces());
00062         labelList& insertedFaces = *insertedFaceLabelsPtr_;
00063 
00064         label nInsertedFaces = 0;
00065 
00066         forAll (directAddr, faceI)
00067         {
00068             if (directAddr[faceI] < 0)
00069             {
00070                 // Found inserted face
00071                 directAddr[faceI] = 0;
00072                 insertedFaces[nInsertedFaces] = faceI;
00073                 nInsertedFaces++;
00074             }
00075         }
00076 
00077         insertedFaces.setSize(nInsertedFaces);
00078     }
00079     else
00080     {
00081         // Interpolative addressing
00082 
00083         interpolationAddrPtr_ = new labelListList(mesh_.nFaces());
00084         labelListList& addr = *interpolationAddrPtr_;
00085 
00086         weightsPtr_ = new scalarListList(mesh_.nFaces());
00087         scalarListList& w = *weightsPtr_;
00088 
00089         const List<objectMap>& ffp = mpm_.facesFromPointsMap();
00090 
00091         forAll (ffp, ffpI)
00092         {
00093             // Get addressing
00094             const labelList& mo = ffp[ffpI].masterObjects();
00095 
00096             label faceI = ffp[ffpI].index();
00097 
00098             if (addr[faceI].size())
00099             {
00100                 FatalErrorIn("void faceMapper::calcAddressing() const")
00101                     << "Master face " << faceI
00102                     << " mapped from point faces " << mo
00103                     << " already destination of mapping." << abort(FatalError);
00104             }
00105 
00106             // Map from masters, uniform weights
00107             addr[faceI] = mo;
00108             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
00109         }
00110 
00111         const List<objectMap>& ffe = mpm_.facesFromEdgesMap();
00112 
00113         forAll (ffe, ffeI)
00114         {
00115             // Get addressing
00116             const labelList& mo = ffe[ffeI].masterObjects();
00117 
00118             label faceI = ffe[ffeI].index();
00119 
00120             if (addr[faceI].size())
00121             {
00122                 FatalErrorIn("void faceMapper::calcAddressing() const")
00123                     << "Master face " << faceI
00124                     << " mapped from edge faces " << mo
00125                     << " already destination of mapping." << abort(FatalError);
00126             }
00127 
00128             // Map from masters, uniform weights
00129             addr[faceI] = mo;
00130             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
00131         }
00132 
00133         const List<objectMap>& fff = mpm_.facesFromFacesMap();
00134 
00135         forAll (fff, fffI)
00136         {
00137             // Get addressing
00138             const labelList& mo = fff[fffI].masterObjects();
00139 
00140             label faceI = fff[fffI].index();
00141 
00142             if (addr[faceI].size())
00143             {
00144                 FatalErrorIn("void faceMapper::calcAddressing() const")
00145                     << "Master face " << faceI
00146                     << " mapped from face faces " << mo
00147                     << " already destination of mapping." << abort(FatalError);
00148             }
00149 
00150             // Map from masters, uniform weights
00151             addr[faceI] = mo;
00152             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
00153         }
00154 
00155 
00156         // Do mapped faces. Note that can already be set from facesFromFaces
00157         // so check if addressing size still zero.
00158         const labelList& fm = mpm_.faceMap();
00159 
00160         forAll (fm, faceI)
00161         {
00162             if (fm[faceI] > -1 && addr[faceI].empty())
00163             {
00164                 // Mapped from a single face
00165                 addr[faceI] = labelList(1, fm[faceI]);
00166                 w[faceI] = scalarList(1, 1.0);
00167             }
00168         }
00169 
00170 
00171         // Grab inserted points (for them the size of addressing is still zero)
00172 
00173         insertedFaceLabelsPtr_ = new labelList(mesh_.nFaces());
00174         labelList& insertedFaces = *insertedFaceLabelsPtr_;
00175 
00176         label nInsertedFaces = 0;
00177 
00178         forAll (addr, faceI)
00179         {
00180             if (addr[faceI].empty())
00181             {
00182                 // Mapped from a dummy face
00183                 addr[faceI] = labelList(1, 0);
00184                 w[faceI] = scalarList(1, 1.0);
00185 
00186                 insertedFaces[nInsertedFaces] = faceI;
00187                 nInsertedFaces++;
00188             }
00189         }
00190 
00191         insertedFaces.setSize(nInsertedFaces);
00192     }
00193 }
00194 
00195 
00196 void Foam::faceMapper::clearOut()
00197 {
00198     deleteDemandDrivenData(directAddrPtr_);
00199     deleteDemandDrivenData(interpolationAddrPtr_);
00200     deleteDemandDrivenData(weightsPtr_);
00201     deleteDemandDrivenData(insertedFaceLabelsPtr_);
00202 }
00203 
00204 
00205 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00206 
00207 // Construct from components
00208 Foam::faceMapper::faceMapper(const mapPolyMesh& mpm)
00209 :
00210     mesh_(mpm.mesh()),
00211     mpm_(mpm),
00212     insertedFaces_(true),
00213     direct_(false),
00214     directAddrPtr_(NULL),
00215     interpolationAddrPtr_(NULL),
00216     weightsPtr_(NULL),
00217     insertedFaceLabelsPtr_(NULL)
00218 {
00219     // Check for possibility of direct mapping
00220     if
00221     (
00222         mpm_.facesFromPointsMap().empty() 
00223      && mpm_.facesFromEdgesMap().empty()
00224      && mpm_.facesFromFacesMap().empty()
00225     )
00226     {
00227         direct_ = true;
00228     }
00229     else
00230     {
00231         direct_ = false;
00232     }
00233 
00234     // Check for inserted faces
00235     if (direct_ && (mpm_.faceMap().empty() || min(mpm_.faceMap()) > -1))
00236     {
00237         insertedFaces_ = false;
00238     }
00239     else
00240     {
00241         // Need to check all 3 lists to see if there are inserted faces
00242         // with no owner
00243 
00244         // Make a copy of the face map, add the entries for faces from points
00245         // and faces from edges and check for left-overs
00246         labelList fm(mesh_.nFaces(), -1);
00247 
00248         const List<objectMap>& ffp = mpm_.facesFromPointsMap();
00249 
00250         forAll (ffp, ffpI)
00251         {
00252             fm[ffp[ffpI].index()] = 0;
00253         }
00254 
00255         const List<objectMap>& ffe = mpm_.facesFromEdgesMap();
00256 
00257         forAll (ffe, ffeI)
00258         {
00259             fm[ffe[ffeI].index()] = 0;
00260         }
00261 
00262         const List<objectMap>& fff = mpm_.facesFromFacesMap();
00263 
00264         forAll (fff, fffI)
00265         {
00266             fm[fff[fffI].index()] = 0;
00267         }
00268 
00269         if (min(fm) < 0)
00270         {
00271             insertedFaces_ = true;
00272         }
00273     }
00274 }
00275 
00276 
00277 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00278 
00279 Foam::faceMapper::~faceMapper()
00280 {
00281     clearOut();
00282 }
00283 
00284 
00285 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00286 
00287 Foam::label Foam::faceMapper::size() const
00288 {
00289     return mesh_.nFaces();
00290 }
00291 
00292 
00293 Foam::label Foam::faceMapper::sizeBeforeMapping() const
00294 {
00295     return mpm_.nOldFaces();
00296 }
00297 
00298 
00299 Foam::label Foam::faceMapper::internalSizeBeforeMapping() const
00300 {
00301     return mpm_.nOldInternalFaces();
00302 }
00303 
00304 
00305 const Foam::unallocLabelList& Foam::faceMapper::directAddressing() const
00306 {
00307     if (!direct())
00308     {
00309         FatalErrorIn
00310         (
00311             "const unallocLabelList& faceMapper::directAddressing() const"
00312         )   << "Requested direct addressing for an interpolative mapper."
00313             << abort(FatalError);
00314     }
00315 
00316     if (!insertedObjects())
00317     {
00318         // No inserted faces.  Re-use faceMap
00319         return mpm_.faceMap();
00320     }
00321     else
00322     {
00323         if (!directAddrPtr_)
00324         {
00325             calcAddressing();
00326         }
00327 
00328         return *directAddrPtr_;
00329     }
00330 }
00331 
00332 
00333 const Foam::labelListList& Foam::faceMapper::addressing() const
00334 {
00335     if (direct())
00336     {
00337         FatalErrorIn
00338         (
00339             "const labelListList& faceMapper::addressing() const"
00340         )   << "Requested interpolative addressing for a direct mapper."
00341             << abort(FatalError);
00342     }
00343 
00344     if (!interpolationAddrPtr_)
00345     {
00346         calcAddressing();
00347     }
00348 
00349     return *interpolationAddrPtr_;
00350 }
00351 
00352 
00353 const Foam::scalarListList& Foam::faceMapper::weights() const
00354 {
00355     if (direct())
00356     {
00357         FatalErrorIn
00358         (
00359             "const scalarListList& faceMapper::weights() const"
00360         )   << "Requested interpolative weights for a direct mapper."
00361             << abort(FatalError);
00362     }
00363 
00364     if (!weightsPtr_)
00365     {
00366         calcAddressing();
00367     }
00368 
00369     return *weightsPtr_;
00370 }
00371 
00372 
00373 const Foam::labelList& Foam::faceMapper::insertedObjectLabels() const
00374 {
00375     if (!insertedFaceLabelsPtr_)
00376     {
00377         if (!insertedObjects())
00378         {
00379             // There are no inserted faces
00380             insertedFaceLabelsPtr_ = new labelList(0);
00381         }
00382         else
00383         {
00384             calcAddressing();
00385         }
00386     }
00387 
00388     return *insertedFaceLabelsPtr_;
00389 }
00390 
00391 
00392 const Foam::labelHashSet& Foam::faceMapper::flipFaceFlux() const
00393 {
00394     return mpm_.flipFaceFlux();
00395 }
00396 
00397 
00398 Foam::label Foam::faceMapper::nOldInternalFaces() const
00399 {
00400     return mpm_.nOldInternalFaces();
00401 }
00402 
00403 
00404 const Foam::labelList& Foam::faceMapper::oldPatchStarts() const
00405 {
00406     return mpm_.oldPatchStarts();
00407 }
00408 
00409 
00410 const Foam::labelList& Foam::faceMapper::oldPatchSizes() const
00411 {
00412     return mpm_.oldPatchSizes();
00413 }
00414 
00415 
00416 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
00417 
00418 
00419 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
00420 
00421 
00422 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
00423 
00424 
00425 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines