fvPatchMapper.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 "fvPatchMapper.H"
27 #include "fvPatch.H"
28 #include "fvBoundaryMesh.H"
29 #include "fvMesh.H"
30 #include "mapPolyMesh.H"
31 #include "faceMapper.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::fvPatchMapper::calcAddressing() const
36 {
37  if
38  (
39  directAddrPtr_
40  || interpolationAddrPtr_
41  || weightsPtr_
42  )
43  {
45  << "Addressing already calculated"
46  << abort(FatalError);
47  }
48 
49  // Mapping
50  const label oldPatchStart =
51  faceMap_.oldPatchStarts()[patch_.index()];
52 
53  const label oldPatchEnd =
54  oldPatchStart + faceMap_.oldPatchSizes()[patch_.index()];
55 
56  hasUnmapped_ = false;
57 
58  // Assemble the maps: slice to patch
59  if (direct())
60  {
61  // Direct mapping - slice to size
62  directAddrPtr_ = new labelList
63  (
64  patch_.patchSlice
65  (
66  static_cast<const labelList&>(faceMap_.directAddressing())
67  )
68  );
69  labelList& addr = *directAddrPtr_;
70 
71  // Adjust mapping to manage hits into other patches and into
72  // internal
73  forAll(addr, facei)
74  {
75  if
76  (
77  addr[facei] >= oldPatchStart
78  && addr[facei] < oldPatchEnd
79  )
80  {
81  addr[facei] -= oldPatchStart;
82  }
83  else
84  {
85  // addr[facei] = 0;
86  addr[facei] = -1;
87  hasUnmapped_ = true;
88  }
89  }
90 
91  if (fvMesh::debug)
92  {
93  if (min(addr) < 0)
94  {
96  << "Unmapped entry in patch mapping for patch "
97  << patch_.index() << " named " << patch_.name()
98  << endl;
99  }
100  }
101  }
102  else
103  {
104  // Interpolative mapping
105  interpolationAddrPtr_ =
106  new labelListList
107  (
108  patch_.patchSlice(faceMap_.addressing())
109  );
110  labelListList& addr = *interpolationAddrPtr_;
111 
112  weightsPtr_ =
113  new scalarListList
114  (
115  patch_.patchSlice(faceMap_.weights())
116  );
117  scalarListList& w = *weightsPtr_;
118 
119  // Adjust mapping to manage hits into other patches and into
120  // internal
121  forAll(addr, facei)
122  {
123  labelList& curAddr = addr[facei];
124  scalarList& curW = w[facei];
125 
126  if
127  (
128  min(curAddr) >= oldPatchStart
129  && max(curAddr) < oldPatchEnd
130  )
131  {
132  // No adjustment of weights, just subtract patch start
133  forAll(curAddr, i)
134  {
135  curAddr[i] -= oldPatchStart;
136  }
137  }
138  else
139  {
140  // Need to recalculate weights to exclude hits into internal
141  labelList newAddr(curAddr.size(), false);
142  scalarField newWeights(curAddr.size());
143  label nActive = 0;
144 
145  forAll(curAddr, lfI)
146  {
147  if
148  (
149  curAddr[lfI] >= oldPatchStart
150  && curAddr[lfI] < oldPatchEnd
151  )
152  {
153  newAddr[nActive] = curAddr[lfI] - oldPatchStart;
154  newWeights[nActive] = curW[lfI];
155  nActive++;
156  }
157  }
158 
159  newAddr.setSize(nActive);
160  newWeights.setSize(nActive);
161 
162  // Re-scale the weights
163  if (nActive > 0)
164  {
165  newWeights /= sum(newWeights);
166  }
167  else
168  {
169  hasUnmapped_ = true;
170  }
171 
172  // Reset addressing and weights
173  curAddr = newAddr;
174  curW = newWeights;
175  }
176  }
177 
178  if (fvMesh::debug)
179  {
180  forAll(addr, i)
181  {
182  if (min(addr[i]) < 0)
183  {
185  << "Error in patch mapping for patch "
186  << patch_.index() << " named " << patch_.name()
187  << abort(FatalError);
188  }
189  }
190  }
191  }
192 }
193 
194 
195 void Foam::fvPatchMapper::clearOut()
196 {
197  deleteDemandDrivenData(directAddrPtr_);
198  deleteDemandDrivenData(interpolationAddrPtr_);
199  deleteDemandDrivenData(weightsPtr_);
200  hasUnmapped_ = false;
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
205 
207 (
208  const fvPatch& patch,
209  const faceMapper& faceMap
210 )
211 :
212  patch_(patch),
213  faceMap_(faceMap),
214  sizeBeforeMapping_(faceMap.oldPatchSizes()[patch_.index()]),
215  hasUnmapped_(false),
216  directAddrPtr_(nullptr),
217  interpolationAddrPtr_(nullptr),
218  weightsPtr_(nullptr)
219 {}
220 
221 
222 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
223 
225 {
226  clearOut();
227 }
228 
229 
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
231 
233 {
234  if (!direct())
235  {
237  << "Requested direct addressing for an interpolative mapper."
238  << abort(FatalError);
239  }
240 
241  if (!directAddrPtr_)
242  {
243  calcAddressing();
244  }
245 
246  return *directAddrPtr_;
247 }
248 
249 
251 {
252  if (direct())
253  {
255  << "Requested interpolative addressing for a direct mapper."
256  << abort(FatalError);
257  }
258 
259  if (!interpolationAddrPtr_)
260  {
261  calcAddressing();
262  }
263 
264  return *interpolationAddrPtr_;
265 }
266 
267 
269 {
270  if (direct())
271  {
273  << "Requested interpolative weights for a direct mapper."
274  << abort(FatalError);
275  }
276 
277  if (!weightsPtr_)
278  {
279  calcAddressing();
280  }
281 
282  return *weightsPtr_;
283 }
284 
285 
286 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:174
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition: faceMapper.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
virtual ~fvPatchMapper()
Destructor.
virtual bool direct() const
Is the mapping direct.
virtual const labelUList & directAddressing() const
Return direct addressing.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:322
List< scalarList > scalarListList
Definition: scalarList.H:51
virtual const labelListList & addressing() const
Return interpolated addressing.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
virtual const scalarListList & weights() const
Return interpolation weights.
Definition: faceMapper.C:340
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:190
fvPatchMapper(const fvPatch &patch, const faceMapper &faceMap)
Construct from mappers.
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
virtual const word & name() const
Return name.
Definition: fvPatch.H:144
void deleteDemandDrivenData(DataPtr &dataPtr)
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:296
virtual const labelList & oldPatchSizes() const
Return old patch sizes.
Definition: faceMapper.C:395