fvSurfaceMapper.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 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 "fvSurfaceMapper.H"
30#include "fvMesh.H"
31#include "mapPolyMesh.H"
32#include "faceMapper.H"
33#include "demandDrivenData.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37void Foam::fvSurfaceMapper::calcAddressing() const
38{
39 if
40 (
41 directAddrPtr_
42 || interpolationAddrPtr_
43 || weightsPtr_
44 || insertedObjectLabelsPtr_
45 )
46 {
48 << "Addressing already calculated"
49 << abort(FatalError);
50 }
51
52 // Mapping
53
54 const label oldNInternal = faceMap_.nOldInternalFaces();
55
56 // Assemble the maps
57 if (direct())
58 {
59 // Direct mapping - slice to size
60 directAddrPtr_ =
61 new labelList
62 (
64 );
65 labelList& addr = *directAddrPtr_;
66
67 // Adjust for creation of an internal face from a boundary face
68 forAll(addr, facei)
69 {
70 if (addr[facei] > oldNInternal)
71 {
72 addr[facei] = 0;
73 }
74 }
75 }
76 else
77 {
78 // Interpolative mapping - slice to size
79 interpolationAddrPtr_ =
81 (
83 );
84 labelListList& addr = *interpolationAddrPtr_;
85
86 weightsPtr_ =
88 (
90 );
91 scalarListList& w = *weightsPtr_;
92
93 // Adjust for creation of an internal face from a boundary face
94 forAll(addr, facei)
95 {
96 if (max(addr[facei]) >= oldNInternal)
97 {
98 addr[facei] = labelList(1, Zero);
99 w[facei] = scalarList(1, scalar(1));
100 }
101 }
102 }
103
104 // Inserted objects
105
106 // If there are, assemble the labels
107 if (insertedObjects())
108 {
109 const labelList& insFaces = faceMap_.insertedObjectLabels();
110
111 insertedObjectLabelsPtr_ = new labelList(insFaces.size());
112 labelList& ins = *insertedObjectLabelsPtr_;
113
114 label nIns = 0;
115
116 forAll(insFaces, facei)
117 {
118 // If the face is internal, keep it here
119 if (insFaces[facei] < size())
120 {
121 ins[nIns] = insFaces[facei];
122 nIns++;
123 }
124 }
125
126 ins.setSize(nIns);
127 }
128 else
129 {
130 // No inserted objects
131 insertedObjectLabelsPtr_ = new labelList(0);
132 }
133}
134
135
136void Foam::fvSurfaceMapper::clearOut()
137{
138 deleteDemandDrivenData(directAddrPtr_);
139 deleteDemandDrivenData(interpolationAddrPtr_);
140 deleteDemandDrivenData(weightsPtr_);
141
142 deleteDemandDrivenData(insertedObjectLabelsPtr_);
143}
144
145
146// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147
149(
150 const fvMesh& mesh,
151 const faceMapper& fMapper
152)
153:
154 mesh_(mesh),
155 faceMap_(fMapper),
156 directAddrPtr_(nullptr),
157 interpolationAddrPtr_(nullptr),
158 weightsPtr_(nullptr),
159 insertedObjectLabelsPtr_(nullptr)
160{}
161
162
163// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
164
166{
167 clearOut();
168}
169
170
171// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172
174{
175 if (!direct())
176 {
178 << "Requested direct addressing for an interpolative mapper."
179 << abort(FatalError);
180 }
181
182 if (!directAddrPtr_)
183 {
184 calcAddressing();
185 }
186
187 return *directAddrPtr_;
188}
189
190
192{
193 if (direct())
194 {
196 << "Requested interpolative addressing for a direct mapper."
197 << abort(FatalError);
198 }
199
200 if (!interpolationAddrPtr_)
201 {
202 calcAddressing();
203 }
204
205 return *interpolationAddrPtr_;
206}
207
208
210{
211 if (direct())
212 {
214 << "Requested interpolative weights for a direct mapper."
215 << abort(FatalError);
216 }
217
218 if (!weightsPtr_)
219 {
220 calcAddressing();
221 }
222
223 return *weightsPtr_;
224}
225
226
228{
229 if (!insertedObjectLabelsPtr_)
230 {
231 calcAddressing();
232 }
233
234 return *insertedObjectLabelsPtr_;
235}
236
237
238// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
239
240
241// * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
242
243
244// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
245
246
247// ************************************************************************* //
SubList< label > subList
Declare type of subList.
Definition: List.H:111
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition: faceMapper.H:60
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:329
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:347
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:303
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: faceMapper.C:365
virtual label nOldInternalFaces() const
Return number of old internalFaces.
Definition: faceMapper.C:390
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
FV surface mapper.
virtual label size() const
Return size.
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 const labelList & insertedObjectLabels() const
Return list of inserted faces.
virtual ~fvSurfaceMapper()
Destructor.
virtual bool insertedObjects() const
Are there any inserted faces.
virtual bool direct() const
Is the mapping direct.
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
List< label > labelList
A List of labels.
Definition: List.H:66
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
void deleteDemandDrivenData(DataPtr &dataPtr)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333