Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
36// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37
38void Foam::fvPatchMapper::calcAddressing() const
39{
40 if
41 (
42 directAddrPtr_
43 || interpAddrPtr_
44 || weightsPtr_
45 )
46 {
48 << "Addressing already calculated"
49 << abort(FatalError);
50 }
51
52 // Mapping
53 const label oldPatchStart =
54 faceMap_.oldPatchStarts()[patch_.index()];
55
56 const label oldPatchEnd =
57 oldPatchStart + faceMap_.oldPatchSizes()[patch_.index()];
58
59 hasUnmapped_ = false;
60
61 // Assemble the maps: slice to patch
62 if (direct())
63 {
64 // Direct mapping - slice to size
65 directAddrPtr_ = std::make_unique<labelList>
66 (
67 patch_.patchSlice(faceMap_.directAddressing())
68 );
69 auto& 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 interpAddrPtr_ = std::make_unique<labelListList>
106 (
107 patch_.patchSlice(faceMap_.addressing())
108 );
109 auto& addr = *interpAddrPtr_;
110
111 weightsPtr_ = std::make_unique<scalarListList>
112 (
113 patch_.patchSlice(faceMap_.weights())
114 );
115 auto& wght = *weightsPtr_;
116
117 // Adjust mapping to manage hits into other patches and into
118 // internal
119
120 forAll(addr, facei)
121 {
122 auto& curAddr = addr[facei];
123 auto& curWght = wght[facei];
124
125 if
126 (
127 min(curAddr) >= oldPatchStart
128 && max(curAddr) < oldPatchEnd
129 )
130 {
131 // No adjustment of weights, just subtract patch start
132 forAll(curAddr, i)
133 {
134 curAddr[i] -= oldPatchStart;
135 }
136 }
137 else
138 {
139 // Need to recalculate weights to exclude hits into internal
140
141 label nActive = 0;
142 scalar sumWeight = 0;
143
144 forAll(curAddr, i)
145 {
146 if
147 (
148 curAddr[i] >= oldPatchStart
149 && curAddr[i] < oldPatchEnd
150 )
151 {
152 curAddr[nActive] = curAddr[i] - oldPatchStart;
153 curWght[nActive] = curWght[i];
154
155 sumWeight += curWght[i];
156 ++nActive;
157 }
158 }
159
160 // Reset addressing and weights
161 curAddr.resize(nActive);
162 curWght.resize(nActive);
163
164 // Re-scale the weights
165 if (nActive)
166 {
167 for (auto& w : curWght)
168 {
169 w /= sumWeight;
170 }
171 }
172 else
173 {
174 hasUnmapped_ = true;
175 }
176 }
177 }
178
179 if (fvMesh::debug)
180 {
181 forAll(addr, i)
182 {
183 if (min(addr[i]) < 0)
184 {
186 << "Error in patch mapping for patch "
187 << patch_.index() << " named " << patch_.name()
188 << abort(FatalError);
189 }
190 }
191 }
192 }
193}
194
195
196// void Foam::fvPatchMapper::clearOut()
197// {
198// directAddrPtr_.reset(nullptr);
199// interpAddrPtr_.reset(nullptr);
200// weightsPtr_.reset(nullptr);
201// hasUnmapped_ = false;
202// }
203
204
205// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206
207Foam::fvPatchMapper::fvPatchMapper
208(
209 const fvPatch& patch,
210 const faceMapper& faceMap
211)
212:
213 patch_(patch),
214 faceMap_(faceMap),
215 sizeBeforeMapping_(faceMap.oldPatchSizes()[patch_.index()]),
216 hasUnmapped_(false)
217{}
218
219
220// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
223{}
224
225
226// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
227
229{
230 if (!direct())
231 {
233 << "Requested direct addressing for an interpolative mapper."
234 << abort(FatalError);
235 }
236
237 if (!directAddrPtr_)
238 {
239 calcAddressing();
240 }
241
242 return *directAddrPtr_;
243}
244
245
247{
248 if (direct())
249 {
251 << "Requested interpolative addressing for a direct mapper."
252 << abort(FatalError);
253 }
254
255 if (!interpAddrPtr_)
256 {
257 calcAddressing();
258 }
259
260 return *interpAddrPtr_;
261}
262
263
265{
266 if (direct())
267 {
269 << "Requested interpolative weights for a direct mapper."
270 << abort(FatalError);
271 }
272
273 if (!weightsPtr_)
274 {
275 calcAddressing();
276 }
277
278 return *weightsPtr_;
279}
280
281
282// ************************************************************************* //
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition faceMapper.H:56
virtual const labelListList & addressing() const
Return interpolated addressing.
virtual const scalarListList & weights() const
Return interpolation weights.
virtual const labelUList & directAddressing() const
Return direct addressing.
virtual bool direct() const
Is the mapping direct.
virtual ~fvPatchMapper()
Destructor.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static int debug
Debug switch.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define WarningInFunction
Report a warning using Foam::Warning.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
UList< label > labelUList
A UList of labels.
Definition UList.H:75
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299