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(), false);
146  scalarField newWeights(curAddr.size());
147  label nActive = 0;
148 
149  forAll(curAddr, lfI)
150  {
151  if
152  (
153  curAddr[lfI] >= oldPatchStart
154  && curAddr[lfI] < oldPatchEnd
155  )
156  {
157  newAddr[nActive] = curAddr[lfI] - oldPatchStart;
158  newWeights[nActive] = curW[lfI];
159  nActive++;
160  }
161  }
162 
163  newAddr.setSize(nActive);
164  newWeights.setSize(nActive);
165 
166  // Re-scale the weights
167  if (nActive > 0)
168  {
169  newWeights /= sum(newWeights);
170  }
171  else
172  {
173  hasUnmapped_ = true;
174  }
175 
176  // Reset addressing and weights
177  curAddr = newAddr;
178  curW = newWeights;
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  deleteDemandDrivenData(directAddrPtr_);
202  deleteDemandDrivenData(interpolationAddrPtr_);
203  deleteDemandDrivenData(weightsPtr_);
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  directAddrPtr_(nullptr),
221  interpolationAddrPtr_(nullptr),
222  weightsPtr_(nullptr)
223 {}
224 
225 
226 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
227 
229 {
230  clearOut();
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 {
238  if (!direct())
239  {
241  << "Requested direct addressing for an interpolative mapper."
242  << abort(FatalError);
243  }
244 
245  if (!directAddrPtr_)
246  {
247  calcAddressing();
248  }
249 
250  return *directAddrPtr_;
251 }
252 
253 
255 {
256  if (direct())
257  {
259  << "Requested interpolative addressing for a direct mapper."
260  << abort(FatalError);
261  }
262 
263  if (!interpolationAddrPtr_)
264  {
265  calcAddressing();
266  }
267 
268  return *interpolationAddrPtr_;
269 }
270 
271 
273 {
274  if (direct())
275  {
277  << "Requested interpolative weights for a direct mapper."
278  << abort(FatalError);
279  }
280 
281  if (!weightsPtr_)
282  {
283  calcAddressing();
284  }
285 
286  return *weightsPtr_;
287 }
288 
289 
290 // ************************************************************************* //
Foam::data::debug
static int debug
Debug switch.
Definition: data.H:77
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarListList
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::faceMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:303
Foam::faceMapper::oldPatchSizes
virtual const labelList & oldPatchSizes() const
Return old patch sizes.
Definition: faceMapper.C:402
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::fvPatchMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: fvPatchMapper.C:236
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::faceMapper::oldPatchStarts
virtual const labelList & oldPatchStarts() const
Return old patch starts.
Definition: faceMapper.C:396
mapPolyMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fvPatch::name
virtual const word & name() const
Return name.
Definition: fvPatch.H:167
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::fvPatchMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: fvPatchMapper.C:254
Foam::fvPatch::patchSlice
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:210
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
fvPatchMapper.H
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::fvPatchMapper::weights
virtual const scalarListList & weights() const
Return interpolation weights.
Definition: fvPatchMapper.C:272
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
faceMapper.H
fvBoundaryMesh.H
Foam::fvPatchMapper::~fvPatchMapper
virtual ~fvPatchMapper()
Destructor.
Definition: fvPatchMapper.C:228
Foam::FatalError
error FatalError
Foam::faceMapper::weights
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:347
Foam::faceMapper
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition: faceMapper.H:57
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
fvMesh.H
Foam::fvPatchMapper::direct
virtual bool direct() const
Is the mapping direct.
Definition: fvPatchMapper.H:136
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< labelList >
Foam::UList< label >
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
fvPatch.H
Foam::faceMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:329
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328