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