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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Description
27  FV surface mapper.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "fvSurfaceMapper.H"
32 #include "fvMesh.H"
33 #include "mapPolyMesh.H"
34 #include "faceMapper.H"
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 void Foam::fvSurfaceMapper::calcAddressing() const
39 {
40  if
41  (
42  directAddrPtr_
43  || interpolationAddrPtr_
44  || weightsPtr_
45  || insertedObjectLabelsPtr_
46  )
47  {
49  << "Addressing already calculated"
50  << abort(FatalError);
51  }
52 
53  // Mapping
54 
55  const label oldNInternal = faceMap_.nOldInternalFaces();
56 
57  // Assemble the maps
58  if (direct())
59  {
60  // Direct mapping - slice to size
61  directAddrPtr_ =
62  new labelList
63  (
65  );
66  labelList& addr = *directAddrPtr_;
67 
68  // Adjust for creation of an internal face from a boundary face
69  forAll(addr, facei)
70  {
71  if (addr[facei] > oldNInternal)
72  {
73  addr[facei] = 0;
74  }
75  }
76  }
77  else
78  {
79  // Interpolative mapping - slice to size
80  interpolationAddrPtr_ =
81  new labelListList
82  (
84  );
85  labelListList& addr = *interpolationAddrPtr_;
86 
87  weightsPtr_ =
88  new scalarListList
89  (
90  scalarListList::subList(faceMap_.weights(), size())
91  );
92  scalarListList& w = *weightsPtr_;
93 
94  // Adjust for creation of an internal face from a boundary face
95  forAll(addr, facei)
96  {
97  if (max(addr[facei]) >= oldNInternal)
98  {
99  addr[facei] = labelList(1, Zero);
100  w[facei] = scalarList(1, scalar(1));
101  }
102  }
103  }
104 
105  // Inserted objects
106 
107  // If there are, assemble the labels
108  if (insertedObjects())
109  {
110  const labelList& insFaces = faceMap_.insertedObjectLabels();
111 
112  insertedObjectLabelsPtr_ = new labelList(insFaces.size());
113  labelList& ins = *insertedObjectLabelsPtr_;
114 
115  label nIns = 0;
116 
117  forAll(insFaces, facei)
118  {
119  // If the face is internal, keep it here
120  if (insFaces[facei] < size())
121  {
122  ins[nIns] = insFaces[facei];
123  nIns++;
124  }
125  }
126 
127  ins.setSize(nIns);
128  }
129  else
130  {
131  // No inserted objects
132  insertedObjectLabelsPtr_ = new labelList(0);
133  }
134 }
135 
136 
137 void Foam::fvSurfaceMapper::clearOut()
138 {
139  deleteDemandDrivenData(directAddrPtr_);
140  deleteDemandDrivenData(interpolationAddrPtr_);
141  deleteDemandDrivenData(weightsPtr_);
142 
143  deleteDemandDrivenData(insertedObjectLabelsPtr_);
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 
149 Foam::fvSurfaceMapper::fvSurfaceMapper
150 (
151  const fvMesh& mesh,
152  const faceMapper& fMapper
153 )
154 :
155  mesh_(mesh),
156  faceMap_(fMapper),
157  directAddrPtr_(nullptr),
158  interpolationAddrPtr_(nullptr),
159  weightsPtr_(nullptr),
160  insertedObjectLabelsPtr_(nullptr)
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
165 
167 {
168  clearOut();
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173 
175 {
176  if (!direct())
177  {
179  << "Requested direct addressing for an interpolative mapper."
180  << abort(FatalError);
181  }
182 
183  if (!directAddrPtr_)
184  {
185  calcAddressing();
186  }
187 
188  return *directAddrPtr_;
189 }
190 
191 
193 {
194  if (direct())
195  {
197  << "Requested interpolative addressing for a direct mapper."
198  << abort(FatalError);
199  }
200 
201  if (!interpolationAddrPtr_)
202  {
203  calcAddressing();
204  }
205 
206  return *interpolationAddrPtr_;
207 }
208 
209 
211 {
212  if (direct())
213  {
215  << "Requested interpolative weights for a direct mapper."
216  << abort(FatalError);
217  }
218 
219  if (!weightsPtr_)
220  {
221  calcAddressing();
222  }
223 
224  return *weightsPtr_;
225 }
226 
227 
229 {
230  if (!insertedObjectLabelsPtr_)
231  {
232  calcAddressing();
233  }
234 
235  return *insertedObjectLabelsPtr_;
236 }
237 
238 
239 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
240 
241 
242 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
243 
244 
245 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
246 
247 
248 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::scalarListList
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::faceMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:303
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
mapPolyMesh.H
Foam::fvSurfaceMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: fvSurfaceMapper.C:192
Foam::fvSurfaceMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: fvSurfaceMapper.C:174
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
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:62
Foam::fvSurfaceMapper::size
virtual label size() const
Return size.
Definition: fvSurfaceMapper.H:120
Foam::fvSurfaceMapper::direct
virtual bool direct() const
Is the mapping direct.
Definition: fvSurfaceMapper.H:132
Foam::fvSurfaceMapper::~fvSurfaceMapper
virtual ~fvSurfaceMapper()
Destructor.
Definition: fvSurfaceMapper.C:166
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
faceMapper.H
Foam::List< label >::subList
SubList< label > subList
Declare type of subList.
Definition: List.H:117
Foam::FatalError
error FatalError
Foam::faceMapper::weights
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:347
Foam::faceMapper
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition: faceMapper.H:57
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::faceMapper::insertedObjectLabels
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: faceMapper.C:365
Foam::List< labelList >
Foam::UList< label >
fvSurfaceMapper.H
Foam::fvSurfaceMapper::insertedObjectLabels
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: fvSurfaceMapper.C:228
Foam::fvSurfaceMapper::insertedObjects
virtual bool insertedObjects() const
Are there any inserted faces.
Definition: fvSurfaceMapper.H:153
Foam::faceMapper::nOldInternalFaces
virtual label nOldInternalFaces() const
Return number of old internalFaces.
Definition: faceMapper.C:390
Foam::faceMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:329
Foam::fvSurfaceMapper::weights
virtual const scalarListList & weights() const
Return interpolation weights.
Definition: fvSurfaceMapper.C:210