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 | 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 Description
25  FV surface mapper.
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 
37 {
38  if
39  (
42  || weightsPtr_
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
60  new labelList
61  (
63  );
64  labelList& addr = *directAddrPtr_;
65 
66  // Adjust for creation of an internal face from a boundary face
67  forAll(addr, faceI)
68  {
69  if (addr[faceI] > oldNInternal)
70  {
71  addr[faceI] = 0;
72  }
73  }
74  }
75  else
76  {
77  // Interpolative mapping - slice to size
79  new labelListList
80  (
82  );
84 
85  weightsPtr_ =
86  new scalarListList
87  (
89  );
91 
92  // Adjust for creation of an internal face from a boundary face
93  forAll(addr, faceI)
94  {
95  if (max(addr[faceI]) >= oldNInternal)
96  {
97  addr[faceI] = labelList(1, label(0));
98  w[faceI] = scalarList(1, 1.0);
99  }
100  }
101  }
102 
103  // Inserted objects
104 
105  // If there are, assemble the labels
106  if (insertedObjects())
107  {
108  const labelList& insFaces = faceMap_.insertedObjectLabels();
109 
110  insertedObjectLabelsPtr_ = new labelList(insFaces.size());
112 
113  label nIns = 0;
114 
115  forAll(insFaces, faceI)
116  {
117  // If the face is internal, keep it here
118  if (insFaces[faceI] < size())
119  {
120  ins[nIns] = insFaces[faceI];
121  nIns++;
122  }
123  }
124 
125  ins.setSize(nIns);
126  }
127  else
128  {
129  // No inserted objects
131  }
132 }
133 
134 
136 {
137  deleteDemandDrivenData(directAddrPtr_);
138  deleteDemandDrivenData(interpolationAddrPtr_);
139  deleteDemandDrivenData(weightsPtr_);
140 
141  deleteDemandDrivenData(insertedObjectLabelsPtr_);
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
146 
147 // Construct from components
149 (
150  const fvMesh& mesh,
151  const faceMapper& fMapper
152 )
153 :
154  mesh_(mesh),
155  faceMap_(fMapper),
156  directAddrPtr_(NULL),
157  interpolationAddrPtr_(NULL),
158  weightsPtr_(NULL),
159  insertedObjectLabelsPtr_(NULL)
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 // ************************************************************************* //
Foam::fvSurfaceMapper::faceMap_
const faceMapper & faceMap_
Reference to face mapper.
Definition: fvSurfaceMapper.H:64
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Foam::faceMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:302
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::fvSurfaceMapper::weightsPtr_
scalarListList * weightsPtr_
Interpolation weights.
Definition: fvSurfaceMapper.H:76
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
mapPolyMesh.H
Foam::fvSurfaceMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: fvSurfaceMapper.C:191
Foam::fvSurfaceMapper::fvSurfaceMapper
fvSurfaceMapper(const fvSurfaceMapper &)
Disallow default bitwise copy construct.
Foam::fvSurfaceMapper::directAddrPtr_
labelList * directAddrPtr_
Direct addressing (only one for of addressing is used)
Definition: fvSurfaceMapper.H:70
Foam::fvSurfaceMapper::clearOut
void clearOut()
Clear out local storage.
Definition: fvSurfaceMapper.C:135
Foam::fvSurfaceMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: fvSurfaceMapper.C:173
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
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::fvSurfaceMapper::size
virtual label size() const
Return size.
Definition: fvSurfaceMapper.H:118
Foam::fvSurfaceMapper::direct
virtual bool direct() const
Is the mapping direct.
Definition: fvSurfaceMapper.H:130
Foam::fvSurfaceMapper::~fvSurfaceMapper
virtual ~fvSurfaceMapper()
Destructor.
Definition: fvSurfaceMapper.C:165
faceMapper.H
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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
Foam::fvSurfaceMapper::interpolationAddrPtr_
labelListList * interpolationAddrPtr_
Interpolated addressing (only one for of addressing is used)
Definition: fvSurfaceMapper.H:73
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fvSurfaceMapper::calcAddressing
void calcAddressing() const
Calculate addressing.
Definition: fvSurfaceMapper.C:36
Foam::faceMapper::insertedObjectLabels
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: faceMapper.C:364
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
fvSurfaceMapper.H
Foam::fvSurfaceMapper::insertedObjectLabels
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: fvSurfaceMapper.C:227
Foam::fvSurfaceMapper::insertedObjects
virtual bool insertedObjects() const
Are there any inserted faces.
Definition: fvSurfaceMapper.H:151
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::faceMapper::nOldInternalFaces
virtual label nOldInternalFaces() const
Return number of old internalFaces.
Definition: faceMapper.C:389
Foam::faceMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:328
Foam::fvSurfaceMapper::weights
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: fvSurfaceMapper.C:209
Foam::fvSurfaceMapper::insertedObjectLabelsPtr_
labelList * insertedObjectLabelsPtr_
Inserted faces.
Definition: fvSurfaceMapper.H:79