1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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 "patchDataWave.H"
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 // Set initial set of changed faces (= all wall faces)
32 template<class TransferType>
33 void Foam::patchDataWave<TransferType>::setChangedFaces
35 const labelHashSet& patchIDs,
36 labelList& changedFaces,
37 List<TransferType>& faceDist
40 const polyMesh& mesh = cellDistFuncs::mesh();
42 label nChangedFaces = 0;
44 forAll(mesh.boundaryMesh(), patchI)
46 if (patchIDs.found(patchI))
48 const polyPatch& patch = mesh.boundaryMesh()[patchI];
50 const Field<Type>& patchField = initialPatchValuePtrs_[patchI];
52 forAll(patch.faceCentres(), patchFaceI)
54 label meshFaceI = patch.start() + patchFaceI;
56 changedFaces[nChangedFaces] = meshFaceI;
58 faceDist[nChangedFaces] =
61 patch.faceCentres()[patchFaceI],
62 patchField[patchFaceI],
73 // Copy from MeshWave data into *this (distance) and field_ (transported data)
74 template<class TransferType>
75 Foam::label Foam::patchDataWave<TransferType>::getValues
77 const MeshWave<TransferType>& waveInfo
80 const polyMesh& mesh = cellDistFuncs::mesh();
82 const List<TransferType>& cellInfo = waveInfo.allCellInfo();
83 const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
88 distance_.setSize(cellInfo.size());
90 forAll(cellInfo, cellI)
92 const TransferType & wpn = cellInfo[cellI];
94 scalar dist = wpn.distSqr();
96 if (cellInfo[cellI].valid())
98 distance_[cellI] = Foam::sqrt(dist);
100 cellData_[cellI] = cellInfo[cellI].data();
104 // Illegal/unset value. What to do with data?
106 distance_[cellI] = dist;
108 //cellData_[cellI] = wallPoint::greatPoint;
109 cellData_[cellI] = cellInfo[cellI].data();
115 // Copy boundary values
116 forAll(patchDistance_, patchI)
118 const polyPatch& patch = mesh.boundaryMesh()[patchI];
120 // Allocate storage for patchDistance
121 scalarField* patchFieldPtr = new scalarField(patch.size());
123 patchDistance_.set(patchI, patchFieldPtr);
125 scalarField& patchField = *patchFieldPtr;
127 // Allocate storage for patchData
128 Field<Type>* patchDataFieldPtr = new Field<Type>(patch.size());
130 patchData_.set(patchI, patchDataFieldPtr);
132 Field<Type>& patchDataField = *patchDataFieldPtr;
134 // Copy distance and data
135 forAll(patchField, patchFaceI)
137 label meshFaceI = patch.start() + patchFaceI;
139 scalar dist = faceInfo[meshFaceI].distSqr();
141 if (faceInfo[meshFaceI].valid())
143 // Adding SMALL to avoid problems with /0 in the turbulence
145 patchField[patchFaceI] = Foam::sqrt(dist) + SMALL;
147 patchDataField[patchFaceI] = faceInfo[meshFaceI].data();
151 // Illegal/unset value. What to do with data?
153 patchField[patchFaceI] = dist;
155 //patchDataField[patchFaceI] = wallPoint::greatPoint;
156 patchDataField[patchFaceI] = faceInfo[meshFaceI].data();
167 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
169 // Construct from components
170 template<class TransferType>
171 Foam::patchDataWave<TransferType>::patchDataWave
173 const polyMesh& mesh,
174 const labelHashSet& patchIDs,
175 const UPtrList<Field<Type> >& initialPatchValuePtrs,
176 const bool correctWalls
181 initialPatchValuePtrs_(initialPatchValuePtrs),
182 correctWalls_(correctWalls),
184 distance_(mesh.nCells()),
185 patchDistance_(mesh.boundaryMesh().size()),
186 cellData_(mesh.nCells()),
187 patchData_(mesh.boundaryMesh().size())
189 patchDataWave<TransferType>::correct();
193 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
195 template<class TransferType>
196 Foam::patchDataWave<TransferType>::~patchDataWave()
200 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
202 // Correct for mesh geom/topo changes
203 template<class TransferType>
204 void Foam::patchDataWave<TransferType>::correct()
207 // Set initial changed faces: set TransferType for wall faces
212 label nWalls = sumPatchSize(patchIDs_);
214 List<TransferType> faceDist(nWalls);
215 labelList changedFaces(nWalls);
217 setChangedFaces(patchIDs_, changedFaces, faceDist);
220 // Do calculate wall distance by 'growing' from faces.
223 MeshWave<TransferType> waveInfo
228 mesh().globalData().nTotalCells() // max iterations
233 // Copy distance into return field
236 nUnset_ = getValues(waveInfo);
239 // Correct wall cells for true distance
244 Map<label> nearestFace(2 * nWalls);
246 // Get distance and indices of nearest face
247 correctBoundaryFaceCells
254 correctBoundaryPointCells
261 // Transfer data from nearest face to cell
262 const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
264 const labelList wallCells(nearestFace.toc());
266 forAll(wallCells, wallCellI)
268 label cellI = wallCells[wallCellI];
270 label faceI = nearestFace[cellI];
272 cellData_[cellI] = faceInfo[faceI].data();
278 // ************************************************************************* //