ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / dynamicMesh / meshCut / wallLayerCells / wallLayerCells.C
blob1e2354f080fcd6c4124bf04456a257f3b6cfefcc
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 "wallLayerCells.H"
27 #include "DynamicList.H"
28 #include "MeshWave.H"
29 #include "wallNormalInfo.H"
30 #include "OFstream.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 namespace Foam
38 defineTypeNameAndDebug(wallLayerCells, 0);
43 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
45 bool Foam::wallLayerCells::usesCoupledPatch(const label cellI) const
47     const polyBoundaryMesh& patches = mesh().boundaryMesh();
49     const cell& cFaces = mesh().cells()[cellI];
51     forAll(cFaces, cFaceI)
52     {
53         label faceI = cFaces[cFaceI];
55         label patchID = patches.whichPatch(faceI);
57         if ((patchID >= 0) && (patches[patchID].coupled()))
58         {
59             return true;
60         }
61     }
62     return false;
65 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
67 // Construct from components
68 Foam::wallLayerCells::wallLayerCells
70     const polyMesh& mesh,
71     const List<word>& patchNames,
72     const label nLayers
75     edgeVertex(mesh),
76     List<refineCell>()
78     // Find out cells connected to walls.
80     const polyPatchList& patches = mesh.boundaryMesh();
82     // Make map from name to local patch ID
83     HashTable<label> patchNameToIndex(patches.size());
85     forAll(patches, patchI)
86     {
87         patchNameToIndex.insert(patches[patchI].name(), patchI);
88     }
91     // Count size of walls to set
92     label nWalls = 0;
94     forAll(patchNames, patchNameI)
95     {
96         const word& name = patchNames[patchNameI];
98         if (patchNameToIndex.found(name))
99         {
100             label patchI = patchNameToIndex[name];
102             nWalls += patches[patchI].size();
103         }
104     }
106     // Allocate storage for start of wave on faces
107     List<wallNormalInfo> changedFacesInfo(nWalls);
108     labelList changedFaces(nWalls);
110     // Fill changedFaces info
111     label nChangedFaces = 0;
113     forAll(patchNames, patchNameI)
114     {
115         const word& name = patchNames[patchNameI];
117         if (patchNameToIndex.found(name))
118         {
119             label patchI = patchNameToIndex[name];
121             const polyPatch& pp = patches[patchI];
123             forAll(pp, patchFaceI)
124             {
125                 label meshFaceI = pp.start() + patchFaceI;
127                 changedFaces[nChangedFaces] = meshFaceI;
129                 // Set transported information to the wall normal.
130                 const vector& norm = pp.faceNormals()[patchFaceI];
132                 changedFacesInfo[nChangedFaces] = wallNormalInfo(norm);
134                 nChangedFaces++;
135             }
136         }
137     }
140     // Do a wave of nLayers, transporting the index in patchNames
141     // (cannot use local patchIDs since we might get info from neighbouring
142     //  processor)
144     MeshWave<wallNormalInfo> regionCalc
145     (
146         mesh,
147         changedFaces,
148         changedFacesInfo,
149         0
150     );
152     regionCalc.iterate(nLayers);
155     // Now regionCalc should hold info on faces that are reachable from
156     // changedFaces within nLayers iterations. We use face info since that is
157     // guaranteed to be consistent across processor boundaries.
159     const List<wallNormalInfo>& faceInfo = regionCalc.allFaceInfo();
161     if (debug)
162     {
163         Info<< "wallLayerCells::getRefinement : dumping selected faces to "
164             << "selectedFaces.obj" << endl;
166         OFstream fcStream("selectedFaces.obj");
168         label vertI = 0;
170         forAll(faceInfo, faceI)
171         {
172             const wallNormalInfo& info = faceInfo[faceI];
174             if (info.valid(regionCalc.data()))
175             {
176                 const face& f = mesh.faces()[faceI];
178                 point mid(0.0, 0.0, 0.0);
180                 forAll(f, fp)
181                 {
182                     mid += mesh.points()[f[fp]];
183                 }
184                 mid /= f.size();
186                 fcStream
187                     << "v " << mid.x() << ' ' << mid.y() << ' ' << mid.z()
188                     << endl;
189                 vertI++;
191                 point end(mid + info.normal());
193                 fcStream
194                     << "v " << end.x() << ' ' << end.y() << ' ' << end.z()
195                     << endl;
196                 vertI++;
198                 fcStream << "l " << vertI << ' ' <<vertI-1 << endl;
199             }
200         }
201     }
204     //
205     // Copy meshWave information to List<refineCell>
206     //
208     // Estimate size
210     DynamicList<refineCell> refineCells(3*nWalls);
212     const List<wallNormalInfo>& cellInfo = regionCalc.allCellInfo();
214     forAll(cellInfo, cellI)
215     {
216         const wallNormalInfo& info = cellInfo[cellI];
218         if (info.valid(regionCalc.data()) && !usesCoupledPatch(cellI))
219         {
220             refineCells.append(refineCell(cellI, info.normal()));
221         }
222     }
224     // Transfer refineCells storage to this.
225     transfer(refineCells);
229 // ************************************************************************* //