patchDataWave.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "patchDataWave.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 // Set initial set of changed faces (= all wall faces)
31 template<class TransferType>
33 (
34  const labelHashSet& patchIDs,
35  labelList& changedFaces,
36  List<TransferType>& faceDist
37 ) const
38 {
39  const polyMesh& mesh = cellDistFuncs::mesh();
40 
41  label nChangedFaces = 0;
42 
43  forAll(mesh.boundaryMesh(), patchi)
44  {
45  if (patchIDs.found(patchi))
46  {
47  const polyPatch& patch = mesh.boundaryMesh()[patchi];
48 
49  const Field<Type>& patchField = initialPatchValuePtrs_[patchi];
50 
51  forAll(patch.faceCentres(), patchFacei)
52  {
53  label meshFacei = patch.start() + patchFacei;
54 
55  changedFaces[nChangedFaces] = meshFacei;
56 
57  faceDist[nChangedFaces] =
58  TransferType
59  (
60  patch.faceCentres()[patchFacei],
61  patchField[patchFacei],
62  0.0
63  );
64 
65  nChangedFaces++;
66  }
67  }
68  }
69 }
70 
71 
72 // Copy from MeshWave data into *this (distance) and field_ (transported data)
73 template<class TransferType>
75 (
76  const MeshWave<TransferType>& waveInfo
77 )
78 {
79  const polyMesh& mesh = cellDistFuncs::mesh();
80 
81  const List<TransferType>& cellInfo = waveInfo.allCellInfo();
82  const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
83 
84  label nIllegal = 0;
85 
86  // Copy cell values
87  distance_.setSize(cellInfo.size());
88 
89  forAll(cellInfo, celli)
90  {
91  const TransferType & wpn = cellInfo[celli];
92 
93  scalar dist = wpn.distSqr();
94 
95  if (cellInfo[celli].valid(waveInfo.data()))
96  {
97  distance_[celli] = Foam::sqrt(dist);
98 
99  cellData_[celli] = cellInfo[celli].data();
100  }
101  else
102  {
103  // Illegal/unset value. What to do with data?
104  // Note: mag for now. Should maybe be member of TransferType?
105 
106  distance_[celli] = mag(dist);
107 
108  // cellData_[celli] = point::max;
109  cellData_[celli] = cellInfo[celli].data();
110 
111  nIllegal++;
112  }
113  }
114 
115  // Copy boundary values
116  forAll(patchDistance_, patchi)
117  {
118  const polyPatch& patch = mesh.boundaryMesh()[patchi];
119 
120  // Allocate storage for patchDistance
121  scalarField* patchFieldPtr = new scalarField(patch.size());
122 
123  patchDistance_.set(patchi, patchFieldPtr);
124 
125  scalarField& patchField = *patchFieldPtr;
126 
127  // Allocate storage for patchData
128  Field<Type>* patchDataFieldPtr = new Field<Type>(patch.size());
129 
130  patchData_.set(patchi, patchDataFieldPtr);
131 
132  Field<Type>& patchDataField = *patchDataFieldPtr;
133 
134  // Copy distance and data
135  forAll(patchField, patchFacei)
136  {
137  label meshFacei = patch.start() + patchFacei;
138 
139  scalar dist = faceInfo[meshFacei].distSqr();
140 
141  if (faceInfo[meshFacei].valid(waveInfo.data()))
142  {
143  // Adding small to avoid problems with /0 in the turbulence
144  // models
145  patchField[patchFacei] = Foam::sqrt(dist) + small;
146 
147  patchDataField[patchFacei] = faceInfo[meshFacei].data();
148  }
149  else
150  {
151  // Illegal/unset value. What to do with data?
152 
153  patchField[patchFacei] = mag(dist);
154 
155  // patchDataField[patchFacei] = point::max;
156  patchDataField[patchFacei] = faceInfo[meshFacei].data();
157 
158  nIllegal++;
159  }
160  }
161  }
162 
163  return nIllegal;
164 }
165 
166 
167 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168 
169 // Construct from components
170 template<class TransferType>
172 (
173  const polyMesh& mesh,
174  const labelHashSet& patchIDs,
175  const UPtrList<Field<Type>>& initialPatchValuePtrs,
176  const bool correctWalls
177 )
178 :
179  cellDistFuncs(mesh),
180  patchIDs_(patchIDs),
181  initialPatchValuePtrs_(initialPatchValuePtrs),
182  correctWalls_(correctWalls),
183  nUnset_(0),
184  distance_(mesh.nCells()),
185  patchDistance_(mesh.boundaryMesh().size()),
186  cellData_(mesh.nCells()),
187  patchData_(mesh.boundaryMesh().size())
188 {
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
194 
195 template<class TransferType>
197 {}
198 
199 
200 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201 
202 // Correct for mesh geom/topo changes
203 template<class TransferType>
205 {
206  //
207  // Set initial changed faces: set TransferType for wall faces
208  // to wall centre.
209  //
210 
211  // Count walls
212  label nWalls = sumPatchSize(patchIDs_);
213 
214  List<TransferType> faceDist(nWalls);
215  labelList changedFaces(nWalls);
216 
217  setChangedFaces(patchIDs_, changedFaces, faceDist);
218 
219  //
220  // Do calculate wall distance by 'growing' from faces.
221  //
222 
223  MeshWave<TransferType> waveInfo
224  (
225  mesh(),
226  changedFaces,
227  faceDist,
228  mesh().globalData().nTotalCells()+1 // max iterations
229  );
230 
231 
232  //
233  // Copy distance into return field
234  //
235 
236  nUnset_ = getValues(waveInfo);
237 
238  //
239  // Correct wall cells for true distance
240  //
241 
242  if (correctWalls_)
243  {
244  Map<label> nearestFace(2 * nWalls);
245 
246  // Get distance and indices of nearest face
247  correctBoundaryFaceCells
248  (
249  patchIDs_,
250  distance_,
251  nearestFace
252  );
253 
254  correctBoundaryPointCells
255  (
256  patchIDs_,
257  distance_,
258  nearestFace
259  );
260 
261  // Transfer data from nearest face to cell
262  const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
263 
264  const labelList wallCells(nearestFace.toc());
265 
266  forAll(wallCells, wallCelli)
267  {
268  label celli = wallCells[wallCelli];
269 
270  label facei = nearestFace[celli];
271 
272  cellData_[celli] = faceInfo[facei].data();
273  }
274  }
275 }
276 
277 
278 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
label nCells() const
dimensionedScalar sqrt(const dimensionedScalar &ds)
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:63
patchDataWave(const polyMesh &mesh, const labelHashSet &patchIDs, const UPtrList< Field< Type >> &initialPatchValuePtrs, bool correctWalls=true)
Construct from mesh, information on patches to initialize and flag.
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:61
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
const List< Type > & allFaceInfo() const
Get allFaceInfo.
Definition: MeshWave.H:119
virtual void correct()
Correct for mesh geom/topo changes.
dynamicFvMesh & mesh
Pre-declare SubField and related Field type.
Definition: Field.H:57
FaceCellWave plus data.
Definition: MeshWave.H:56
List< label > labelList
A List of labels.
Definition: labelList.H:56
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
volScalarField scalarField(fieldObject, mesh)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
T * data()
Return a pointer to the first data element,.
Definition: UListI.H:149
label patchi
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
virtual ~patchDataWave()
Destructor.