Loading...
Searching...
No Matches
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-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 "fvSurfaceMapper.H"
30#include "fvMesh.H"
31#include "mapPolyMesh.H"
32#include "faceMapper.H"
33
34// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35
36void Foam::fvSurfaceMapper::calcAddressing() const
37{
38 if
39 (
40 directAddrPtr_
41 || interpAddrPtr_
42 || weightsPtr_
43 || insertedObjectsPtr_
44 )
45 {
47 << "Addressing already calculated"
48 << abort(FatalError);
49 }
50
51 // Mapping
52
53 const label oldNInternal = faceMap_.nOldInternalFaces();
54
55 // Assemble the maps
56 if (direct())
57 {
58 // Direct mapping - slice to size
59 directAddrPtr_ = std::make_unique<labelList>
60 (
61 labelList::subList(faceMap_.directAddressing(), size())
62 );
63 auto& addr = *directAddrPtr_;
64
65 // Adjust for creation of an internal face from a boundary face
66 forAll(addr, facei)
67 {
68 if (addr[facei] > oldNInternal)
69 {
70 addr[facei] = 0;
71 }
72 }
73 }
74 else
75 {
76 // Interpolative mapping - slice to size
77 interpAddrPtr_ = std::make_unique<labelListList>
78 (
79 labelListList::subList(faceMap_.addressing(), size())
80 );
81 auto& addr = *interpAddrPtr_;
82
83 weightsPtr_ = std::make_unique<scalarListList>
84 (
85 scalarListList::subList(faceMap_.weights(), size())
86 );
87 auto& wght = *weightsPtr_;
88
89 // Adjust for creation of an internal face from a boundary face
90 forAll(addr, facei)
91 {
92 if (max(addr[facei]) >= oldNInternal)
93 {
94 addr[facei] = labelList(1, Foam::zero{});
95 wght[facei] = scalarList(1, scalar(1));
96 }
97 }
98 }
99
100 // Inserted objects
101
102 insertedObjectsPtr_ = std::make_unique<labelList>();
103 auto& inserted = *insertedObjectsPtr_;
104
105 // If there are, assemble the labels
106 if (faceMap_.insertedObjects())
107 {
108 const labelList& insFaces = faceMap_.insertedObjectLabels();
109
110 inserted.resize(insFaces.size());
111
112 label count = 0;
113 for (const label facei : insFaces)
114 {
115 // If the face is internal, keep it here
116 if (facei < size())
117 {
118 inserted[count] = facei;
119 ++count;
120 }
121 }
122
123 inserted.resize(count);
124 }
125}
126
127
128// void Foam::fvSurfaceMapper::clearOut()
129// {
130// directAddrPtr_.reset(nullptr);
131// interpAddrPtr_.reset(nullptr);
132// weightsPtr_.reset(nullptr);
133// insertedObjectsPtr_.reset(nullptr);
134// }
135
136
137// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
138
139Foam::fvSurfaceMapper::fvSurfaceMapper
140(
141 const fvMesh& mesh,
142 const faceMapper& mapper
143)
144:
145 mesh_(mesh),
146 faceMap_(mapper)
147{}
148
149
150// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153{}
154
155
156// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157
159{
160 if (!direct())
161 {
163 << "Requested direct addressing for an interpolative mapper."
164 << abort(FatalError);
165 }
166
167 if (!directAddrPtr_)
168 {
169 calcAddressing();
170 }
171
172 return *directAddrPtr_;
173}
174
175
177{
178 if (direct())
179 {
181 << "Requested interpolative addressing for a direct mapper."
182 << abort(FatalError);
183 }
184
185 if (!interpAddrPtr_)
186 {
187 calcAddressing();
188 }
189
190 return *interpAddrPtr_;
191}
192
193
195{
196 if (direct())
197 {
199 << "Requested interpolative weights for a direct mapper."
200 << abort(FatalError);
201 }
202
203 if (!weightsPtr_)
204 {
205 calcAddressing();
206 }
207
208 return *weightsPtr_;
209}
210
211
213{
214 if (!insertedObjectsPtr_)
215 {
216 calcAddressing();
217 }
218
219 return *insertedObjectsPtr_;
220}
221
222
223// ************************************************************************* //
SubList< label > subList
Definition List.H:129
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition faceMapper.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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 direct() const
Is the mapping direct.
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
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
List< label > labelList
A List of labels.
Definition List.H:62
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
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299