ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / OpenFOAM / meshes / polyMesh / mapPolyMesh / faceMapper / faceMapper.C
blobdb90ba0dd8c8d6b48886bd8cfc8a702a3fd3d99d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "faceMapper.H"
27 #include "demandDrivenData.H"
28 #include "polyMesh.H"
29 #include "mapPolyMesh.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
36 void Foam::faceMapper::calcAddressing() const
38     if
39     (
40         directAddrPtr_
41      || interpolationAddrPtr_
42      || weightsPtr_
43      || insertedFaceLabelsPtr_
44     )
45     {
46         FatalErrorIn("void faceMapper::calcAddressing() const")
47             << "Addressing already calculated."
48             << abort(FatalError);
49     }
51     if (direct())
52     {
53         // Direct addressing, no weights
55         directAddrPtr_ = new labelList(mpm_.faceMap());
56         labelList& directAddr = *directAddrPtr_;
58         // Reset the size of addressing list to contain only live faces
59         directAddr.setSize(mesh_.nFaces());
61         insertedFaceLabelsPtr_ = new labelList(mesh_.nFaces());
62         labelList& insertedFaces = *insertedFaceLabelsPtr_;
64         label nInsertedFaces = 0;
66         forAll(directAddr, faceI)
67         {
68             if (directAddr[faceI] < 0)
69             {
70                 // Found inserted face
71                 directAddr[faceI] = 0;
72                 insertedFaces[nInsertedFaces] = faceI;
73                 nInsertedFaces++;
74             }
75         }
77         insertedFaces.setSize(nInsertedFaces);
78     }
79     else
80     {
81         // Interpolative addressing
83         interpolationAddrPtr_ = new labelListList(mesh_.nFaces());
84         labelListList& addr = *interpolationAddrPtr_;
86         weightsPtr_ = new scalarListList(mesh_.nFaces());
87         scalarListList& w = *weightsPtr_;
89         const List<objectMap>& ffp = mpm_.facesFromPointsMap();
91         forAll(ffp, ffpI)
92         {
93             // Get addressing
94             const labelList& mo = ffp[ffpI].masterObjects();
96             label faceI = ffp[ffpI].index();
98             if (addr[faceI].size())
99             {
100                 FatalErrorIn("void faceMapper::calcAddressing() const")
101                     << "Master face " << faceI
102                     << " mapped from point faces " << mo
103                     << " already destination of mapping." << abort(FatalError);
104             }
106             // Map from masters, uniform weights
107             addr[faceI] = mo;
108             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
109         }
111         const List<objectMap>& ffe = mpm_.facesFromEdgesMap();
113         forAll(ffe, ffeI)
114         {
115             // Get addressing
116             const labelList& mo = ffe[ffeI].masterObjects();
118             label faceI = ffe[ffeI].index();
120             if (addr[faceI].size())
121             {
122                 FatalErrorIn("void faceMapper::calcAddressing() const")
123                     << "Master face " << faceI
124                     << " mapped from edge faces " << mo
125                     << " already destination of mapping." << abort(FatalError);
126             }
128             // Map from masters, uniform weights
129             addr[faceI] = mo;
130             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
131         }
133         const List<objectMap>& fff = mpm_.facesFromFacesMap();
135         forAll(fff, fffI)
136         {
137             // Get addressing
138             const labelList& mo = fff[fffI].masterObjects();
140             label faceI = fff[fffI].index();
142             if (addr[faceI].size())
143             {
144                 FatalErrorIn("void faceMapper::calcAddressing() const")
145                     << "Master face " << faceI
146                     << " mapped from face faces " << mo
147                     << " already destination of mapping." << abort(FatalError);
148             }
150             // Map from masters, uniform weights
151             addr[faceI] = mo;
152             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
153         }
156         // Do mapped faces. Note that can already be set from facesFromFaces
157         // so check if addressing size still zero.
158         const labelList& fm = mpm_.faceMap();
160         forAll(fm, faceI)
161         {
162             if (fm[faceI] > -1 && addr[faceI].empty())
163             {
164                 // Mapped from a single face
165                 addr[faceI] = labelList(1, fm[faceI]);
166                 w[faceI] = scalarList(1, 1.0);
167             }
168         }
171         // Grab inserted points (for them the size of addressing is still zero)
173         insertedFaceLabelsPtr_ = new labelList(mesh_.nFaces());
174         labelList& insertedFaces = *insertedFaceLabelsPtr_;
176         label nInsertedFaces = 0;
178         forAll(addr, faceI)
179         {
180             if (addr[faceI].empty())
181             {
182                 // Mapped from a dummy face
183                 addr[faceI] = labelList(1, 0);
184                 w[faceI] = scalarList(1, 1.0);
186                 insertedFaces[nInsertedFaces] = faceI;
187                 nInsertedFaces++;
188             }
189         }
191         insertedFaces.setSize(nInsertedFaces);
192     }
196 void Foam::faceMapper::clearOut()
198     deleteDemandDrivenData(directAddrPtr_);
199     deleteDemandDrivenData(interpolationAddrPtr_);
200     deleteDemandDrivenData(weightsPtr_);
201     deleteDemandDrivenData(insertedFaceLabelsPtr_);
205 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
207 // Construct from components
208 Foam::faceMapper::faceMapper(const mapPolyMesh& mpm)
210     mesh_(mpm.mesh()),
211     mpm_(mpm),
212     insertedFaces_(true),
213     direct_(false),
214     directAddrPtr_(NULL),
215     interpolationAddrPtr_(NULL),
216     weightsPtr_(NULL),
217     insertedFaceLabelsPtr_(NULL)
219     // Check for possibility of direct mapping
220     if
221     (
222         mpm_.facesFromPointsMap().empty()
223      && mpm_.facesFromEdgesMap().empty()
224      && mpm_.facesFromFacesMap().empty()
225     )
226     {
227         direct_ = true;
228     }
229     else
230     {
231         direct_ = false;
232     }
234     // Check for inserted faces
235     if (direct_ && (mpm_.faceMap().empty() || min(mpm_.faceMap()) > -1))
236     {
237         insertedFaces_ = false;
238     }
239     else
240     {
241         // Need to check all 3 lists to see if there are inserted faces
242         // with no owner
244         // Make a copy of the face map, add the entries for faces from points
245         // and faces from edges and check for left-overs
246         labelList fm(mesh_.nFaces(), -1);
248         const List<objectMap>& ffp = mpm_.facesFromPointsMap();
250         forAll(ffp, ffpI)
251         {
252             fm[ffp[ffpI].index()] = 0;
253         }
255         const List<objectMap>& ffe = mpm_.facesFromEdgesMap();
257         forAll(ffe, ffeI)
258         {
259             fm[ffe[ffeI].index()] = 0;
260         }
262         const List<objectMap>& fff = mpm_.facesFromFacesMap();
264         forAll(fff, fffI)
265         {
266             fm[fff[fffI].index()] = 0;
267         }
269         if (min(fm) < 0)
270         {
271             insertedFaces_ = true;
272         }
273     }
277 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
279 Foam::faceMapper::~faceMapper()
281     clearOut();
285 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
287 Foam::label Foam::faceMapper::size() const
289     return mesh_.nFaces();
293 Foam::label Foam::faceMapper::sizeBeforeMapping() const
295     return mpm_.nOldFaces();
299 Foam::label Foam::faceMapper::internalSizeBeforeMapping() const
301     return mpm_.nOldInternalFaces();
305 const Foam::labelUList& Foam::faceMapper::directAddressing() const
307     if (!direct())
308     {
309         FatalErrorIn
310         (
311             "const labelUList& faceMapper::directAddressing() const"
312         )   << "Requested direct addressing for an interpolative mapper."
313             << abort(FatalError);
314     }
316     if (!insertedObjects())
317     {
318         // No inserted faces.  Re-use faceMap
319         return mpm_.faceMap();
320     }
321     else
322     {
323         if (!directAddrPtr_)
324         {
325             calcAddressing();
326         }
328         return *directAddrPtr_;
329     }
333 const Foam::labelListList& Foam::faceMapper::addressing() const
335     if (direct())
336     {
337         FatalErrorIn
338         (
339             "const labelListList& faceMapper::addressing() const"
340         )   << "Requested interpolative addressing for a direct mapper."
341             << abort(FatalError);
342     }
344     if (!interpolationAddrPtr_)
345     {
346         calcAddressing();
347     }
349     return *interpolationAddrPtr_;
353 const Foam::scalarListList& Foam::faceMapper::weights() const
355     if (direct())
356     {
357         FatalErrorIn
358         (
359             "const scalarListList& faceMapper::weights() const"
360         )   << "Requested interpolative weights for a direct mapper."
361             << abort(FatalError);
362     }
364     if (!weightsPtr_)
365     {
366         calcAddressing();
367     }
369     return *weightsPtr_;
373 const Foam::labelList& Foam::faceMapper::insertedObjectLabels() const
375     if (!insertedFaceLabelsPtr_)
376     {
377         if (!insertedObjects())
378         {
379             // There are no inserted faces
380             insertedFaceLabelsPtr_ = new labelList(0);
381         }
382         else
383         {
384             calcAddressing();
385         }
386     }
388     return *insertedFaceLabelsPtr_;
392 const Foam::labelHashSet& Foam::faceMapper::flipFaceFlux() const
394     return mpm_.flipFaceFlux();
398 Foam::label Foam::faceMapper::nOldInternalFaces() const
400     return mpm_.nOldInternalFaces();
404 const Foam::labelList& Foam::faceMapper::oldPatchStarts() const
406     return mpm_.oldPatchStarts();
410 const Foam::labelList& Foam::faceMapper::oldPatchSizes() const
412     return mpm_.oldPatchSizes();
416 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
419 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
422 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
425 // ************************************************************************* //