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-2024 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 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 void Foam::fvPatchMapper::calcAddressing() const
39 {
40  if
41  (
42  directAddrPtr_
43  || interpAddrPtr_
44  || weightsPtr_
45  )
46  {
48  << "Addressing already calculated"
49  << abort(FatalError);
50  }
51 
52  // Mapping
53  const label oldPatchStart =
54  faceMap_.oldPatchStarts()[patch_.index()];
55 
56  const label oldPatchEnd =
57  oldPatchStart + faceMap_.oldPatchSizes()[patch_.index()];
58 
59  hasUnmapped_ = false;
60 
61  // Assemble the maps: slice to patch
62  if (direct())
63  {
64  // Direct mapping - slice to size
65  directAddrPtr_ = std::make_unique<labelList>
66  (
67  patch_.patchSlice
68  (
69  static_cast<const labelList&>(faceMap_.directAddressing())
70  )
71  );
72  auto& addr = *directAddrPtr_;
73 
74  // Adjust mapping to manage hits into other patches and into
75  // internal
76  forAll(addr, facei)
77  {
78  if
79  (
80  addr[facei] >= oldPatchStart
81  && addr[facei] < oldPatchEnd
82  )
83  {
84  addr[facei] -= oldPatchStart;
85  }
86  else
87  {
88  //addr[facei] = 0;
89  addr[facei] = -1;
90  hasUnmapped_ = true;
91  }
92  }
93 
94  if (fvMesh::debug)
95  {
96  if (min(addr) < 0)
97  {
99  << "Unmapped entry in patch mapping for patch "
100  << patch_.index() << " named " << patch_.name()
101  << endl;
102  }
103  }
104  }
105  else
106  {
107  // Interpolative mapping
108  interpAddrPtr_ = std::make_unique<labelListList>
109  (
110  patch_.patchSlice(faceMap_.addressing())
111  );
112  auto& addr = *interpAddrPtr_;
113 
114  weightsPtr_ = std::make_unique<scalarListList>
115  (
116  patch_.patchSlice(faceMap_.weights())
117  );
118  auto& wght = *weightsPtr_;
119 
120  // Adjust mapping to manage hits into other patches and into
121  // internal
122 
123  forAll(addr, facei)
124  {
125  auto& curAddr = addr[facei];
126  auto& curWght = wght[facei];
127 
128  if
129  (
130  min(curAddr) >= oldPatchStart
131  && max(curAddr) < oldPatchEnd
132  )
133  {
134  // No adjustment of weights, just subtract patch start
135  forAll(curAddr, i)
136  {
137  curAddr[i] -= oldPatchStart;
138  }
139  }
140  else
141  {
142  // Need to recalculate weights to exclude hits into internal
143 
144  label nActive = 0;
145  scalar sumWeight = 0;
146 
147  forAll(curAddr, i)
148  {
149  if
150  (
151  curAddr[i] >= oldPatchStart
152  && curAddr[i] < oldPatchEnd
153  )
154  {
155  curAddr[nActive] = curAddr[i] - oldPatchStart;
156  curWght[nActive] = curWght[i];
157 
158  sumWeight += curWght[i];
159  ++nActive;
160  }
161  }
162 
163  // Reset addressing and weights
164  curAddr.resize(nActive);
165  curWght.resize(nActive);
166 
167  // Re-scale the weights
168  if (nActive)
169  {
170  for (auto& w : curWght)
171  {
172  w /= sumWeight;
173  }
174  }
175  else
176  {
177  hasUnmapped_ = true;
178  }
179  }
180  }
181 
182  if (fvMesh::debug)
183  {
184  forAll(addr, i)
185  {
186  if (min(addr[i]) < 0)
187  {
189  << "Error in patch mapping for patch "
190  << patch_.index() << " named " << patch_.name()
191  << abort(FatalError);
192  }
193  }
194  }
195  }
196 }
197 
198 
199 // void Foam::fvPatchMapper::clearOut()
200 // {
201 // directAddrPtr_.reset(nullptr);
202 // interpAddrPtr_.reset(nullptr);
203 // weightsPtr_.reset(nullptr);
204 // hasUnmapped_ = false;
205 // }
206 
207 
208 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
209 
210 Foam::fvPatchMapper::fvPatchMapper
211 (
212  const fvPatch& patch,
213  const faceMapper& faceMap
214 )
215 :
216  patch_(patch),
217  faceMap_(faceMap),
218  sizeBeforeMapping_(faceMap.oldPatchSizes()[patch_.index()]),
219  hasUnmapped_(false)
220 {}
221 
222 
223 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
226 {}
227 
228 
229 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
230 
232 {
233  if (!direct())
234  {
236  << "Requested direct addressing for an interpolative mapper."
237  << abort(FatalError);
238  }
239 
240  if (!directAddrPtr_)
241  {
242  calcAddressing();
243  }
244 
245  return *directAddrPtr_;
246 }
247 
248 
250 {
251  if (direct())
252  {
254  << "Requested interpolative addressing for a direct mapper."
255  << abort(FatalError);
256  }
257 
258  if (!interpAddrPtr_)
259  {
260  calcAddressing();
261  }
262 
263  return *interpAddrPtr_;
264 }
265 
266 
268 {
269  if (direct())
270  {
272  << "Requested interpolative weights for a direct mapper."
273  << abort(FatalError);
274  }
275 
276  if (!weightsPtr_)
277  {
278  calcAddressing();
279  }
280 
281  return *weightsPtr_;
282 }
283 
284 
285 // ************************************************************************* //
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:608
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.
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:337
virtual const labelListList & addressing() const
Return interpolated addressing.
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:355
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
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:402
const std::string patch
OpenFOAM patch number as a std::string.
virtual const word & name() const
Return name.
Definition: fvPatch.H:210
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:311
static int debug
Debug switch.
virtual const labelList & oldPatchSizes() const
Return old patch sizes.
Definition: faceMapper.C:408