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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fvPatchMapper.H"
29 #include "fvPatch.H"
30 #include "fvBoundaryMesh.H"
31 #include "fvMesh.H"
32 #include "mapPolyMesh.H"
33 #include "faceMapper.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 void Foam::fvPatchMapper::calcAddressing() const
38 {
39  if
40  (
41  directAddrPtr_
42  || interpolationAddrPtr_
43  || weightsPtr_
44  )
45  {
47  << "Addressing already calculated"
48  << abort(FatalError);
49  }
50 
51  // Mapping
52  const label oldPatchStart =
53  faceMap_.oldPatchStarts()[patch_.index()];
54 
55  const label oldPatchEnd =
56  oldPatchStart + faceMap_.oldPatchSizes()[patch_.index()];
57 
58  hasUnmapped_ = false;
59 
60  // Assemble the maps: slice to patch
61  if (direct())
62  {
63  // Direct mapping - slice to size
64  directAddrPtr_ = new labelList
65  (
66  patch_.patchSlice
67  (
68  static_cast<const labelList&>(faceMap_.directAddressing())
69  )
70  );
71  labelList& addr = *directAddrPtr_;
72 
73  // Adjust mapping to manage hits into other patches and into
74  // internal
75  forAll(addr, facei)
76  {
77  if
78  (
79  addr[facei] >= oldPatchStart
80  && addr[facei] < oldPatchEnd
81  )
82  {
83  addr[facei] -= oldPatchStart;
84  }
85  else
86  {
87  //addr[facei] = 0;
88  addr[facei] = -1;
89  hasUnmapped_ = true;
90  }
91  }
92 
93  if (fvMesh::debug)
94  {
95  if (min(addr) < 0)
96  {
98  << "Unmapped entry in patch mapping for patch "
99  << patch_.index() << " named " << patch_.name()
100  << endl;
101  }
102  }
103  }
104  else
105  {
106  // Interpolative mapping
107  interpolationAddrPtr_ =
108  new labelListList
109  (
110  patch_.patchSlice(faceMap_.addressing())
111  );
112  labelListList& addr = *interpolationAddrPtr_;
113 
114  weightsPtr_ =
115  new scalarListList
116  (
117  patch_.patchSlice(faceMap_.weights())
118  );
119  scalarListList& w = *weightsPtr_;
120 
121  // Adjust mapping to manage hits into other patches and into
122  // internal
123  forAll(addr, facei)
124  {
125  labelList& curAddr = addr[facei];
126  scalarList& curW = w[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  labelList newAddr(curAddr.size(), false);
144  scalarField newWeights(curAddr.size());
145  label nActive = 0;
146 
147  forAll(curAddr, lfI)
148  {
149  if
150  (
151  curAddr[lfI] >= oldPatchStart
152  && curAddr[lfI] < oldPatchEnd
153  )
154  {
155  newAddr[nActive] = curAddr[lfI] - oldPatchStart;
156  newWeights[nActive] = curW[lfI];
157  nActive++;
158  }
159  }
160 
161  newAddr.setSize(nActive);
162  newWeights.setSize(nActive);
163 
164  // Re-scale the weights
165  if (nActive > 0)
166  {
167  newWeights /= sum(newWeights);
168  }
169  else
170  {
171  hasUnmapped_ = true;
172  }
173 
174  // Reset addressing and weights
175  curAddr = newAddr;
176  curW = newWeights;
177  }
178  }
179 
180  if (fvMesh::debug)
181  {
182  forAll(addr, i)
183  {
184  if (min(addr[i]) < 0)
185  {
187  << "Error in patch mapping for patch "
188  << patch_.index() << " named " << patch_.name()
189  << abort(FatalError);
190  }
191  }
192  }
193  }
194 }
195 
196 
197 void Foam::fvPatchMapper::clearOut()
198 {
199  deleteDemandDrivenData(directAddrPtr_);
200  deleteDemandDrivenData(interpolationAddrPtr_);
201  deleteDemandDrivenData(weightsPtr_);
202  hasUnmapped_ = false;
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
207 
208 Foam::fvPatchMapper::fvPatchMapper
209 (
210  const fvPatch& patch,
211  const faceMapper& faceMap
212 )
213 :
214  patch_(patch),
215  faceMap_(faceMap),
216  sizeBeforeMapping_(faceMap.oldPatchSizes()[patch_.index()]),
217  hasUnmapped_(false),
218  directAddrPtr_(nullptr),
219  interpolationAddrPtr_(nullptr),
220  weightsPtr_(nullptr)
221 {}
222 
223 
224 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
225 
227 {
228  clearOut();
229 }
230 
231 
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233 
235 {
236  if (!direct())
237  {
239  << "Requested direct addressing for an interpolative mapper."
240  << abort(FatalError);
241  }
242 
243  if (!directAddrPtr_)
244  {
245  calcAddressing();
246  }
247 
248  return *directAddrPtr_;
249 }
250 
251 
253 {
254  if (direct())
255  {
257  << "Requested interpolative addressing for a direct mapper."
258  << abort(FatalError);
259  }
260 
261  if (!interpolationAddrPtr_)
262  {
263  calcAddressing();
264  }
265 
266  return *interpolationAddrPtr_;
267 }
268 
269 
271 {
272  if (direct())
273  {
275  << "Requested interpolative weights for a direct mapper."
276  << abort(FatalError);
277  }
278 
279  if (!weightsPtr_)
280  {
281  calcAddressing();
282  }
283 
284  return *weightsPtr_;
285 }
286 
287 
288 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
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: blockMeshMergeFast.C:94
Foam::fvPatchMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: fvPatchMapper.C:234
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:337
Foam::fvPatch::name
virtual const word & name() const
Return name.
Definition: fvPatch.H:151
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:252
Foam::fvPatch::patchSlice
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:194
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
fvPatchMapper.H
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::fvSchemes::debug
static int debug
Debug switch.
Definition: fvSchemes.H:105
Foam::label
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:62
Foam::fvPatchMapper::weights
virtual const scalarListList & weights() const
Return interpolation weights.
Definition: fvPatchMapper.C:270
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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:226
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:181
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:137
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:355
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:294