fvPatchMapper.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvPatchMapper.H"
30 #include "fvPatch.H"
31 #include "fvBoundaryMesh.H"
32 #include "fvMesh.H"
33 #include "mapPolyMesh.H"
34 #include "faceMapper.H"
35 #include "demandDrivenData.H"
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 void Foam::fvPatchMapper::calcAddressing() const
40 {
41  if
42  (
43  directAddrPtr_
44  || interpolationAddrPtr_
45  || weightsPtr_
46  )
47  {
49  << "Addressing already calculated"
50  << abort(FatalError);
51  }
52 
53  // Mapping
54  const label oldPatchStart =
55  faceMap_.oldPatchStarts()[patch_.index()];
56 
57  const label oldPatchEnd =
58  oldPatchStart + faceMap_.oldPatchSizes()[patch_.index()];
59 
60  hasUnmapped_ = false;
61 
62  // Assemble the maps: slice to patch
63  if (direct())
64  {
65  // Direct mapping - slice to size
66  directAddrPtr_ = new labelList
67  (
68  patch_.patchSlice
69  (
70  static_cast<const labelList&>(faceMap_.directAddressing())
71  )
72  );
73  labelList& addr = *directAddrPtr_;
74 
75  // Adjust mapping to manage hits into other patches and into
76  // internal
77  forAll(addr, facei)
78  {
79  if
80  (
81  addr[facei] >= oldPatchStart
82  && addr[facei] < oldPatchEnd
83  )
84  {
85  addr[facei] -= oldPatchStart;
86  }
87  else
88  {
89  //addr[facei] = 0;
90  addr[facei] = -1;
91  hasUnmapped_ = true;
92  }
93  }
94 
95  if (fvMesh::debug)
96  {
97  if (min(addr) < 0)
98  {
100  << "Unmapped entry in patch mapping for patch "
101  << patch_.index() << " named " << patch_.name()
102  << endl;
103  }
104  }
105  }
106  else
107  {
108  // Interpolative mapping
109  interpolationAddrPtr_ =
110  new labelListList
111  (
112  patch_.patchSlice(faceMap_.addressing())
113  );
114  labelListList& addr = *interpolationAddrPtr_;
115 
116  weightsPtr_ =
117  new scalarListList
118  (
119  patch_.patchSlice(faceMap_.weights())
120  );
121  scalarListList& w = *weightsPtr_;
122 
123  // Adjust mapping to manage hits into other patches and into
124  // internal
125  forAll(addr, facei)
126  {
127  labelList& curAddr = addr[facei];
128  scalarList& curW = w[facei];
129 
130  if
131  (
132  min(curAddr) >= oldPatchStart
133  && max(curAddr) < oldPatchEnd
134  )
135  {
136  // No adjustment of weights, just subtract patch start
137  forAll(curAddr, i)
138  {
139  curAddr[i] -= oldPatchStart;
140  }
141  }
142  else
143  {
144  // Need to recalculate weights to exclude hits into internal
145  labelList newAddr(curAddr.size());
146  scalarField newWeights(curAddr.size());
147 
148  label nActive = 0;
149  scalar sumWeight = 0;
150 
151  forAll(curAddr, lfI)
152  {
153  if
154  (
155  curAddr[lfI] >= oldPatchStart
156  && curAddr[lfI] < oldPatchEnd
157  )
158  {
159  newAddr[nActive] = curAddr[lfI] - oldPatchStart;
160  newWeights[nActive] = curW[lfI];
161 
162  sumWeight += curW[lfI];
163  ++nActive;
164  }
165  }
166 
167  newAddr.resize(nActive);
168  newWeights.resize(nActive);
169 
170  // Re-scale the weights
171  if (nActive > 0)
172  {
173  newWeights /= sumWeight;
174  }
175  else
176  {
177  hasUnmapped_ = true;
178  }
179 
180  // Reset addressing and weights
181  curAddr = newAddr;
182  curW = newWeights;
183  }
184  }
185 
186  if (fvMesh::debug)
187  {
188  forAll(addr, i)
189  {
190  if (min(addr[i]) < 0)
191  {
193  << "Error in patch mapping for patch "
194  << patch_.index() << " named " << patch_.name()
195  << abort(FatalError);
196  }
197  }
198  }
199  }
200 }
201 
202 
203 void Foam::fvPatchMapper::clearOut()
204 {
205  deleteDemandDrivenData(directAddrPtr_);
206  deleteDemandDrivenData(interpolationAddrPtr_);
207  deleteDemandDrivenData(weightsPtr_);
208  hasUnmapped_ = false;
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
213 
214 Foam::fvPatchMapper::fvPatchMapper
215 (
216  const fvPatch& patch,
217  const faceMapper& faceMap
218 )
219 :
220  patch_(patch),
221  faceMap_(faceMap),
222  sizeBeforeMapping_(faceMap.oldPatchSizes()[patch_.index()]),
223  hasUnmapped_(false),
224  directAddrPtr_(nullptr),
225  interpolationAddrPtr_(nullptr),
226  weightsPtr_(nullptr)
227 {}
228 
229 
230 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
231 
233 {
234  clearOut();
235 }
236 
237 
238 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
239 
241 {
242  if (!direct())
243  {
245  << "Requested direct addressing for an interpolative mapper."
246  << abort(FatalError);
247  }
248 
249  if (!directAddrPtr_)
250  {
251  calcAddressing();
252  }
253 
254  return *directAddrPtr_;
255 }
256 
257 
259 {
260  if (direct())
261  {
263  << "Requested interpolative addressing for a direct mapper."
264  << abort(FatalError);
265  }
266 
267  if (!interpolationAddrPtr_)
268  {
269  calcAddressing();
270  }
271 
272  return *interpolationAddrPtr_;
273 }
274 
275 
277 {
278  if (direct())
279  {
281  << "Requested interpolative weights for a direct mapper."
282  << abort(FatalError);
283  }
284 
285  if (!weightsPtr_)
286  {
287  calcAddressing();
288  }
289 
290  return *weightsPtr_;
291 }
292 
293 
294 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
const List< T >::subList patchSlice(const List< T > &values) const
This patch slice from the complete list, which has size mesh::nFaces(), using the virtual patch size...
Definition: fvPatch.H:271
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition: faceMapper.H:53
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual ~fvPatchMapper()
Destructor.
virtual bool direct() const
Is the mapping direct.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual const labelUList & directAddressing() const
Return direct addressing.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:322
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
virtual const labelListList & addressing() const
Return interpolated addressing.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:340
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Template functions to aid in the implementation of demand driven data.
virtual const scalarListList & weights() const
Return interpolation weights.
#define WarningInFunction
Report a warning using Foam::Warning.
virtual const labelList & oldPatchStarts() const
Return old patch starts.
Definition: faceMapper.C:389
const std::string patch
OpenFOAM patch number as a std::string.
virtual const word & name() const
Return name.
Definition: fvPatch.H:210
List< label > labelList
A List of labels.
Definition: List.H:62
void deleteDemandDrivenData(DataPtr &dataPtr)
label index() const noexcept
The index of this patch in the boundary mesh.
Definition: fvPatch.H:218
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:296
static int debug
Debug switch.
virtual const labelList & oldPatchSizes() const
Return old patch sizes.
Definition: faceMapper.C:395