volPointInterpolationAdjoint.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-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Class
28  Foam::volPointInterpolationAdjoint
29 
30 Description
31  Interpolate from cell centres to points (vertices) using inverse distance
32  weighting
33 
34 SourceFiles
35  volPointInterpolationAdjoint.C
36  volPointInterpolate.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef volPointInterpolationAdjoint_H
41 #define volPointInterpolationAdjoint_H
42 
43 #include "MeshObject.H"
44 #include "scalarList.H"
45 #include "volFields.H"
46 #include "pointFields.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 class fvMesh;
54 class pointMesh;
55 
56 /*---------------------------------------------------------------------------*\
57  Class volPointInterpolationAdjoint Declaration
58 \*---------------------------------------------------------------------------*/
59 
61 :
62  public MeshObject<fvMesh, UpdateableMeshObject, volPointInterpolationAdjoint>
63 {
64 protected:
65 
66  // Protected data
67 
68  //- Boundary addressing
70 
71  //- Per boundary face whether is on non-coupled, non-empty patch
73 
74  //- Per mesh(!) point whether is on non-coupled, non-empty patch (on
75  // any processor)
77 
78  //- Per mesh(!) point whether is on symmetry plane
79  // any processor)
81 
82  //- Per boundary point the weights per pointFaces.
84 
85 
86  // Protected Member Functions
87 
88  //- Construct addressing over all boundary faces
90 
91  //- Make weights for points on uncoupled patches
92  void makeBoundaryWeights(scalarField& sumWeights);
93 
94  //- Construct all point weighting factors
95  void makeWeights();
96 
97  //- Helper: push master point data to collocated points
98  template<class Type>
99  void pushUntransformedData(List<Type>&) const;
100 
101  //- Get boundary field in same order as boundary faces. Field is
102  // zero on all coupled and empty patches
103  template<class Type>
105  (
107  ) const;
108 
109  //- Add separated contributions
110  template<class Type>
111  void addSeparated
112  (
114  ) const;
115 
116  //- No copy construct
118 
119  //- No copy assignment
120  void operator=(const volPointInterpolationAdjoint&) = delete;
121 
122 
123 public:
124 
125  // Declare name of the class and its debug switch
126  ClassName("volPointInterpolationAdjoint");
127 
128 
129  // Constructors
130 
131  //- Constructor given fvMesh and pointMesh.
132  explicit volPointInterpolationAdjoint(const fvMesh&);
133 
134 
135  //- Destructor
137 
138 
139  // Member functions
140 
141  // Edit
142 
143  //- Update mesh topology using the morph engine
144  void updateMesh(const mapPolyMesh&);
145 
146  //- Correct weighting factors for moving mesh.
147  bool movePoints();
148 
149 
150  // Interpolation Functions
151 
152  //- Interpolate sensitivties from points to faces
153  // conditions
154  template<class Type>
156  (
159  const labelHashSet& patchIDs
160  ) const;
161 
162  template<class Type>
164  (
167  ) const;
168 
169  //- Interpolate sensitivties from points to faces
171  (
172  const vectorField& pf,
173  vectorField& vf,
174  const labelHashSet& patchIDs
175  ) const;
176 };
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace Foam
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 #ifdef NoRepository
187 #endif
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 #endif
192 
193 // ************************************************************************* //
const labelList patchIDs(pbm.indices(polyPatchNames, true))
bitSet isPatchPoint_
Per mesh(!) point whether is on non-coupled, non-empty patch (on.
void calcBoundaryAddressing()
Construct addressing over all boundary faces.
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
void operator=(const volPointInterpolationAdjoint &)=delete
No copy assignment.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
bitSet isSymmetryPoint_
Per mesh(!) point whether is on symmetry plane.
autoPtr< primitivePatch > boundaryPtr_
Boundary addressing.
Interpolate from cell centres to points (vertices) using inverse distance weighting.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:85
void addSeparated(GeometricField< Type, pointPatchField, pointMesh > &) const
Add separated contributions.
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
void makeWeights()
Construct all point weighting factors.
bool movePoints()
Correct weighting factors for moving mesh.
bitSet boundaryIsPatchFace_
Per boundary face whether is on non-coupled, non-empty patch.
ClassName("volPointInterpolationAdjoint")
void pushUntransformedData(List< Type > &) const
Helper: push master point data to collocated points.
Generic GeometricBoundaryField class.
Definition: areaFieldsFwd.H:46
void makeBoundaryWeights(scalarField &sumWeights)
Make weights for points on uncoupled patches.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
volPointInterpolationAdjoint(const volPointInterpolationAdjoint &)=delete
No copy construct.
tmp< Field< Type > > flatBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Get boundary field in same order as boundary faces. Field is.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void interpolateSensitivitiesField(const GeometricField< Type, pointPatchField, pointMesh > &pf, typename GeometricField< Type, fvPatchField, volMesh >::Boundary &vf, const labelHashSet &patchIDs) const
Interpolate sensitivties from points to faces.
scalarListList boundaryPointWeights_
Per boundary point the weights per pointFaces.
Namespace for OpenFOAM.