Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / foam / meshes / polyMesh / mapPolyMesh / faceMapper / faceMapper.C
blob1cb63e024a5d75dfca21a61720c1bd38b40037a1
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  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 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
33 void Foam::faceMapper::calcAddressing() const
35     if
36     (
37         directAddrPtr_
38      || interpolationAddrPtr_
39      || weightsPtr_
40      || insertedFaceLabelsPtr_
41     )
42     {
43         FatalErrorIn("void faceMapper::calcAddressing() const")
44             << "Addressing already calculated."
45             << abort(FatalError);
46     }
48     if (direct())
49     {
50         // Direct addressing, no weights
52         directAddrPtr_ = new labelList(mpm_.faceMap());
53         labelList& directAddr = *directAddrPtr_;
55         // Reset the size of addressing list to contain only live faces
56         directAddr.setSize(mesh_.nFaces());
58         insertedFaceLabelsPtr_ = new labelList(mesh_.nFaces());
59         labelList& insertedFaces = *insertedFaceLabelsPtr_;
61         label nInsertedFaces = 0;
63         forAll (directAddr, faceI)
64         {
65             if (directAddr[faceI] < 0)
66             {
67                 // Found inserted face
68                 directAddr[faceI] = 0;
69                 insertedFaces[nInsertedFaces] = faceI;
70                 nInsertedFaces++;
71             }
72         }
74         insertedFaces.setSize(nInsertedFaces);
75     }
76     else
77     {
78         // Interpolative addressing
80         interpolationAddrPtr_ = new labelListList(mesh_.nFaces());
81         labelListList& addr = *interpolationAddrPtr_;
83         weightsPtr_ = new scalarListList(mesh_.nFaces());
84         scalarListList& w = *weightsPtr_;
86         const List<objectMap>& ffp = mpm_.facesFromPointsMap();
88         forAll (ffp, ffpI)
89         {
90             // Get addressing
91             const labelList& mo = ffp[ffpI].masterObjects();
93             label faceI = ffp[ffpI].index();
95             if (addr[faceI].size())
96             {
97                 FatalErrorIn("void faceMapper::calcAddressing() const")
98                     << "Master face " << faceI
99                     << " mapped from point faces " << mo
100                     << " already destination of mapping." << abort(FatalError);
101             }
103             // Map from masters, uniform weights
104             addr[faceI] = mo;
105             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
106         }
108         const List<objectMap>& ffe = mpm_.facesFromEdgesMap();
110         forAll (ffe, ffeI)
111         {
112             // Get addressing
113             const labelList& mo = ffe[ffeI].masterObjects();
115             label faceI = ffe[ffeI].index();
117             if (addr[faceI].size())
118             {
119                 FatalErrorIn("void faceMapper::calcAddressing() const")
120                     << "Master face " << faceI
121                     << " mapped from edge faces " << mo
122                     << " already destination of mapping." << abort(FatalError);
123             }
125             // Map from masters, uniform weights
126             addr[faceI] = mo;
127             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
128         }
130         const List<objectMap>& fff = mpm_.facesFromFacesMap();
132         forAll (fff, fffI)
133         {
134             // Get addressing
135             const labelList& mo = fff[fffI].masterObjects();
137             label faceI = fff[fffI].index();
139             if (addr[faceI].size())
140             {
141                 FatalErrorIn("void faceMapper::calcAddressing() const")
142                     << "Master face " << faceI
143                     << " mapped from face faces " << mo
144                     << " already destination of mapping." << abort(FatalError);
145             }
147             // Map from masters, uniform weights
148             addr[faceI] = mo;
149             w[faceI] = scalarList(mo.size(), 1.0/mo.size());
150         }
153         // Do mapped faces. Note that can already be set from facesFromFaces
154         // so check if addressing size still zero.
155         const labelList& fm = mpm_.faceMap();
157         forAll (fm, faceI)
158         {
159             if (fm[faceI] > -1 && addr[faceI].empty())
160             {
161                 // Mapped from a single face
162                 addr[faceI] = labelList(1, fm[faceI]);
163                 w[faceI] = scalarList(1, 1.0);
164             }
165         }
168         // Grab inserted faces (for them the size of addressing is still zero)
170         insertedFaceLabelsPtr_ = new labelList(mesh_.nFaces());
171         labelList& insertedFaces = *insertedFaceLabelsPtr_;
173         label nInsertedFaces = 0;
175         forAll (addr, faceI)
176         {
177             if (addr[faceI].empty())
178             {
179                 // Mapped from a dummy face
180                 addr[faceI] = labelList(1, 0);
181                 w[faceI] = scalarList(1, 1.0);
183                 insertedFaces[nInsertedFaces] = faceI;
184                 nInsertedFaces++;
185             }
186         }
188         insertedFaces.setSize(nInsertedFaces);
189     }
193 void Foam::faceMapper::clearOut()
195     deleteDemandDrivenData(directAddrPtr_);
196     deleteDemandDrivenData(interpolationAddrPtr_);
197     deleteDemandDrivenData(weightsPtr_);
198     deleteDemandDrivenData(insertedFaceLabelsPtr_);
202 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
204 // Construct from components
205 Foam::faceMapper::faceMapper(const mapPolyMesh& mpm)
207     mesh_(mpm.mesh()),
208     mpm_(mpm),
209     insertedFaces_(true),
210     direct_(false),
211     directAddrPtr_(NULL),
212     interpolationAddrPtr_(NULL),
213     weightsPtr_(NULL),
214     insertedFaceLabelsPtr_(NULL)
216     // Check for possibility of direct mapping
217     if
218     (
219         mpm_.facesFromPointsMap().empty()
220      && mpm_.facesFromEdgesMap().empty()
221      && mpm_.facesFromFacesMap().empty()
222     )
223     {
224         direct_ = true;
225     }
226     else
227     {
228         direct_ = false;
229     }
231     // Check for inserted faces
232     // Bug fix: mapping with insertion or truncation.  HJ, 5/Sep/2007
233     if
234     (
235         direct_
236      && (mpm_.faceMap().size() || min(mpm_.faceMap()) > -1)
237      && mpm_.faceMap().size() == mesh_.nFaces()
238     )
239     {
240         insertedFaces_ = false;
241     }
242     else
243     {
244         // Need to check all 3 lists to see if there are inserted faces
245         // with no owner
247         // Make a copy of the face map, add the entries for faces from points
248         // and faces from edges and check for left-overs
249         labelList fm(mesh_.nFaces(), -1);
251         const List<objectMap>& ffp = mpm_.facesFromPointsMap();
253         forAll (ffp, ffpI)
254         {
255             fm[ffp[ffpI].index()] = 0;
256         }
258         const List<objectMap>& ffe = mpm_.facesFromEdgesMap();
260         forAll (ffe, ffeI)
261         {
262             fm[ffe[ffeI].index()] = 0;
263         }
265         const List<objectMap>& fff = mpm_.facesFromFacesMap();
267         forAll (fff, fffI)
268         {
269             fm[fff[fffI].index()] = 0;
270         }
272         if (min(fm) < 0)
273         {
274             insertedFaces_ = true;
275         }
276     }
280 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
282 Foam::faceMapper::~faceMapper()
284     clearOut();
288 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
290 Foam::label Foam::faceMapper::size() const
292     return mesh_.nFaces();
296 Foam::label Foam::faceMapper::sizeBeforeMapping() const
298     return mpm_.nOldFaces();
302 Foam::label Foam::faceMapper::internalSizeBeforeMapping() const
304     return mpm_.nOldInternalFaces();
308 const Foam::unallocLabelList& Foam::faceMapper::directAddressing() const
310     if (!direct())
311     {
312         FatalErrorIn
313         (
314             "const unallocLabelList& faceMapper::directAddressing() const"
315         )   << "Requested direct addressing for an interpolative mapper."
316             << abort(FatalError);
317     }
319     if (!insertedObjects())
320     {
321         // No inserted faces.  Re-use faceMap
322         return mpm_.faceMap();
323     }
324     else
325     {
326         if (!directAddrPtr_)
327         {
328             calcAddressing();
329         }
331         return *directAddrPtr_;
332     }
336 const Foam::labelListList& Foam::faceMapper::addressing() const
338     if (direct())
339     {
340         FatalErrorIn
341         (
342             "const labelListList& faceMapper::addressing() const"
343         )   << "Requested interpolative addressing for a direct mapper."
344             << abort(FatalError);
345     }
347     if (!interpolationAddrPtr_)
348     {
349         calcAddressing();
350     }
352     return *interpolationAddrPtr_;
356 const Foam::scalarListList& Foam::faceMapper::weights() const
358     if (direct())
359     {
360         FatalErrorIn
361         (
362             "const scalarListList& faceMapper::weights() const"
363         )   << "Requested interpolative weights for a direct mapper."
364             << abort(FatalError);
365     }
367     if (!weightsPtr_)
368     {
369         calcAddressing();
370     }
372     return *weightsPtr_;
376 const Foam::labelList& Foam::faceMapper::insertedObjectLabels() const
378     if (!insertedFaceLabelsPtr_)
379     {
380         if (!insertedObjects())
381         {
382             // There are no inserted faces
383             insertedFaceLabelsPtr_ = new labelList(0);
384         }
385         else
386         {
387             calcAddressing();
388         }
389     }
391     return *insertedFaceLabelsPtr_;
395 const Foam::labelHashSet& Foam::faceMapper::flipFaceFlux() const
397     return mpm_.flipFaceFlux();
401 Foam::label Foam::faceMapper::nOldInternalFaces() const
403     return mpm_.nOldInternalFaces();
407 const Foam::labelList& Foam::faceMapper::oldPatchStarts() const
409     return mpm_.oldPatchStarts();
413 const Foam::labelList& Foam::faceMapper::oldPatchSizes() const
415     return mpm_.oldPatchSizes();
419 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
422 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
425 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
428 // ************************************************************************* //