Loading...
Searching...
No Matches
pointMapper.H
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-2013 OpenFOAM Foundation
9 Copyright (C) 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
27Class
28 Foam::pointMapper
29
30Description
31 This object provides mapping and fill-in information for point data
32 between the two meshes after the topological change. It is
33 constructed from mapPolyMesh.
34
35SourceFiles
36 pointMapper.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef Foam_pointMapper_H
41#define Foam_pointMapper_H
42
43#include "morphFieldMapper.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49
50// Forward Declarations
51class mapPolyMesh;
52class pointMesh;
54/*---------------------------------------------------------------------------*\
55 Class pointMapper Declaration
56\*---------------------------------------------------------------------------*/
57
58class pointMapper
59:
60 public morphFieldMapper
61{
62 // Private Data
63
64 //- Reference to mapPolyMesh
65 const mapPolyMesh& mpm_;
66
67 //- The size of the mapper = pointMesh::size()
68 const label mapperLen_;
69
70 //- Number of inserted (unmapped) points
71 label nInsertedObjects_;
72
73 //- Is the mapping direct
74 bool direct_;
75
76
77 // Demand-Driven Data
78
79 //- Direct addressing (only one for of addressing is used)
80 mutable std::unique_ptr<labelList> directAddrPtr_;
81
82 //- Interpolated addressing (only one for of addressing is used)
83 mutable std::unique_ptr<labelListList> interpAddrPtr_;
84
85 //- Interpolation weights
86 mutable std::unique_ptr<scalarListList> weightsPtr_;
87
88 //- Inserted points
89 mutable std::unique_ptr<labelList> insertedObjectsPtr_;
90
91
92 // Private Member Functions
93
94 //- Calculate addressing for mapping with inserted points
95 void calcAddressing() const;
96
97public:
98
99 // Generated Methods
100
101 //- No copy construct
102 pointMapper(const pointMapper&) = delete;
103
104 //- No copy assignment
105 void operator=(const pointMapper&) = delete;
106
107
108 // Constructors
109
110 //- Construct from mapPolyMesh
111 pointMapper(const pointMesh&, const mapPolyMesh& mpm);
112
113
114
115 //- Destructor
116 virtual ~pointMapper();
118
119 // Member Functions
120
121 //- The mapper size
122 virtual label size() const;
123
124 //- Return size before mapping
125 virtual label sizeBeforeMapping() const;
126
127 //- Is the mapping direct
128 virtual bool direct() const
129 {
130 return direct_;
131 }
132
133 //- Are there unmapped values? i.e. do all size() elements get value
134 virtual bool hasUnmapped() const
135 {
136 return insertedObjects();
137 }
138
139 //- Return direct addressing
140 virtual const labelUList& directAddressing() const;
141
142 //- Return interpolated addressing
143 virtual const labelListList& addressing() const;
144
145 //- Return interpolation weights
146 virtual const scalarListList& weights() const;
147
148 //- Are there any inserted points
149 bool insertedObjects() const noexcept
150 {
151 return bool(nInsertedObjects_);
152 }
153
154 //- Return list of inserted points
156};
157
158
159// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160
161} // End namespace Foam
162
163// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164
165#endif
166
167// ************************************************************************* //
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
morphFieldMapper()=default
Default construct.
virtual const labelListList & addressing() const
Return interpolated addressing.
virtual const scalarListList & weights() const
Return interpolation weights.
virtual bool hasUnmapped() const
Are there unmapped values? i.e. do all size() elements get value.
virtual const labelUList & directAddressing() const
Return direct addressing.
virtual label size() const
The mapper size.
const labelList & insertedObjectLabels() const
Return list of inserted points.
virtual ~pointMapper()
Destructor.
virtual label sizeBeforeMapping() const
Return size before mapping.
pointMapper(const pointMapper &)=delete
No copy construct.
void operator=(const pointMapper &)=delete
No copy assignment.
bool insertedObjects() const noexcept
Are there any inserted points.
virtual bool direct() const
Is the mapping direct.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
const direction noexcept
Definition scalarImpl.H:265
UList< label > labelUList
A UList of labels.
Definition UList.H:75