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 | Copyright (C) 2011-2015 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 
36 {
37  if
38  (
41  || weightsPtr_
42  )
43  {
45  << "Addressing already calculated"
46  << abort(FatalError);
47  }
48 
49  // Mapping
50  const label oldPatchStart =
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
63  (
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
106  new labelListList
107  (
109  );
111 
112  weightsPtr_ =
113  new scalarListList
114  (
116  );
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 
160  //if (nActive == 0)
161  //{
162  // newAddr[nActive] = 0;
163  // newWeights[nActive] = 1;
164  // nActive++;
165  //}
166 
167  newAddr.setSize(nActive);
168  newWeights.setSize(nActive);
169 
170  // Re-scale the weights
171  if (nActive > 0)
172  {
173  newWeights /= sum(newWeights);
174  }
175  else
176  {
177  hasUnmapped_ = true;
178  }
179 
180  // Reset addressing and weights
181  curAddr = newAddr;
182  curW = newWeights;
183  }
184  }
185 
186  if (fvMesh::debug)
187  {
188  forAll(addr, i)
189  {
190  if (min(addr[i]) < 0)
191  {
193  << "Error in patch mapping for patch "
194  << patch_.index() << " named " << patch_.name()
195  << abort(FatalError);
196  }
197  }
198  }
199  }
200 }
201 
202 
204 {
205  deleteDemandDrivenData(directAddrPtr_);
206  deleteDemandDrivenData(interpolationAddrPtr_);
207  deleteDemandDrivenData(weightsPtr_);
208  hasUnmapped_ = false;
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
213 
214 // Construct from components
216 (
217  const fvPatch& patch,
218  const faceMapper& faceMap
219 )
220 :
221  patch_(patch),
222  faceMap_(faceMap),
223  sizeBeforeMapping_(faceMap.oldPatchSizes()[patch_.index()]),
224  hasUnmapped_(false),
225  directAddrPtr_(NULL),
226  interpolationAddrPtr_(NULL),
227  weightsPtr_(NULL)
228 {}
229 
230 
231 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
232 
234 {
235  clearOut();
236 }
237 
238 
239 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
240 
242 {
243  if (!direct())
244  {
246  << "Requested direct addressing for an interpolative mapper."
247  << abort(FatalError);
248  }
249 
250  if (!directAddrPtr_)
251  {
252  calcAddressing();
253  }
254 
255  return *directAddrPtr_;
256 }
257 
258 
260 {
261  if (direct())
262  {
264  << "Requested interpolative addressing for a direct mapper."
265  << abort(FatalError);
266  }
267 
268  if (!interpolationAddrPtr_)
269  {
270  calcAddressing();
271  }
272 
273  return *interpolationAddrPtr_;
274 }
275 
276 
278 {
279  if (direct())
280  {
282  << "Requested interpolative weights for a direct mapper."
283  << abort(FatalError);
284  }
285 
286  if (!weightsPtr_)
287  {
288  calcAddressing();
289  }
290 
291  return *weightsPtr_;
292 }
293 
294 
295 // ************************************************************************* //
Foam::fvPatchMapper::fvPatchMapper
fvPatchMapper(const fvPatchMapper &)
Disallow default bitwise copy construct.
Foam::faceMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:302
Foam::faceMapper::oldPatchSizes
virtual const labelList & oldPatchSizes() const
Return old patch sizes.
Definition: faceMapper.C:401
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::fvPatchMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: fvPatchMapper.C:241
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::faceMapper::oldPatchStarts
virtual const labelList & oldPatchStarts() const
Return old patch starts.
Definition: faceMapper.C:395
mapPolyMesh.H
Foam::fvPatchMapper::hasUnmapped_
bool hasUnmapped_
Definition: fvPatchMapper.H:73
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::fvPatchMapper::patch_
const fvPatch & patch_
Reference to patch.
Definition: fvPatchMapper.H:62
Foam::fvPatchMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: fvPatchMapper.C:259
Foam::fvPatch::patchSlice
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:192
Foam::fvPatchMapper::faceMap_
const faceMapper & faceMap_
Reference to face mapper.
Definition: fvPatchMapper.H:65
fvPatchMapper.H
Foam::fvPatch::name
const word & name() const
Return name.
Definition: fvPatch.H:149
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::fvSchemes::debug
static int debug
Debug switch.
Definition: fvSchemes.H:103
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvPatchMapper::weights
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: fvPatchMapper.C:277
Foam::fvPatchMapper::directAddrPtr_
labelList * directAddrPtr_
Direct addressing (only one for of addressing is used)
Definition: fvPatchMapper.H:76
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::fvPatchMapper::clearOut
void clearOut()
Clear out local storage.
Definition: fvPatchMapper.C:203
faceMapper.H
fvBoundaryMesh.H
Foam::fvPatchMapper::~fvPatchMapper
virtual ~fvPatchMapper()
Destructor.
Definition: fvPatchMapper.C:233
Foam::FatalError
error FatalError
Foam::faceMapper::weights
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:346
Foam::faceMapper
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition: faceMapper.H:55
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
fvMesh.H
Foam::fvPatchMapper::direct
virtual bool direct() const
Is the mapping direct.
Definition: fvPatchMapper.H:137
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::fvPatchMapper::interpolationAddrPtr_
labelListList * interpolationAddrPtr_
Interpolated addressing (only one for of addressing is used)
Definition: fvPatchMapper.H:79
Foam::fvPatchMapper::calcAddressing
void calcAddressing() const
Calculate addressing for mapping with inserted cells.
Definition: fvPatchMapper.C:35
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::UList
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
Foam::fvPatchMapper::weightsPtr_
scalarListList * weightsPtr_
Interpolation weights.
Definition: fvPatchMapper.H:82
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
fvPatch.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::scalarListList
List< scalarList > scalarListList
Definition: scalarList.H:51
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::faceMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:328
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259