cellCellStencilObject.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) 2017-2019 OpenCFD Ltd.
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 Class
27  Foam::cellCellStencilObject
28 
29 Description
30 
31 SourceFiles
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef cellCellStencilObject_H
36 #define cellCellStencilObject_H
37 
38 #include "cellCellStencil.H"
39 #include "MeshObject.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
47 typedef MeshObject
48 <
49  fvMesh,
52 > Stencil;
53 
54 /*---------------------------------------------------------------------------*\
55  Class cellCellStencilObject Declaration
56 \*---------------------------------------------------------------------------*/
57 
59 :
60  public Stencil,
61  public cellCellStencil
62 {
63  // Private data
64 
65  autoPtr<cellCellStencil> stencilPtr_;
66 
67 
68 public:
69 
70  TypeName("cellCellStencilObject");
71 
72  // Constructors
73 
74  //- Construct with mesh
76  (
77  const fvMesh& mesh,
78  const bool update = true
79  )
80  :
82  <
83  fvMesh,
86  >(mesh),
88  stencilPtr_
89  (
91  (
92  mesh,
93  mesh.schemesDict().subDict
94  (
95  "oversetInterpolation"
96  ),
97  update
98  )
99  )
100  {}
101 
102 
103  //- Destructor
104  virtual ~cellCellStencilObject() = default;
105 
106 
107  // Member Functions
108 
109  //- Callback for geometry motion
110  virtual bool movePoints()
111  {
112  return stencilPtr_().update();
113  }
114 
115  //- Update stencils. Return false if nothing changed.
116  virtual bool update()
117  {
118  return stencilPtr_().update();
119  }
120 
121  //- Return the cell type list
122  virtual const labelUList& cellTypes() const
123  {
124  return stencilPtr_().cellTypes();
125  }
126 
127  //- Indices of interpolated cells
128  virtual const labelUList& interpolationCells() const
129  {
130  return stencilPtr_().interpolationCells();
131  }
132 
133  //- Return a communication schedule
134  virtual const mapDistribute& cellInterpolationMap() const
135  {
136  return stencilPtr_().cellInterpolationMap();
137  }
138 
139  //- Per interpolated cell the neighbour cells (in terms of slots as
140  // constructed by above cellInterpolationMap) to interpolate
141  virtual const labelListList& cellStencil() const
142  {
143  return stencilPtr_().cellStencil();
144  }
145 
146  //- Weights for cellStencil
147  virtual const List<scalarList>& cellInterpolationWeights() const
148  {
149  return stencilPtr_().cellInterpolationWeights();
150  }
151 
152  //- Per interpolated cell the interpolation factor. (0 = use
153  // calculated, 1 = use interpolated)
154  virtual const scalarList& cellInterpolationWeight() const
155  {
156  return stencilPtr_().cellInterpolationWeight();
157  }
158 
159  //- Calculate weights for a single acceptor
160  virtual void stencilWeights
161  (
162  const point& sample,
163  const pointList& donorCcs,
164  scalarList& weights
165  ) const
166  {
167  stencilPtr_().stencilWeights(sample, donorCcs, weights);
168  }
169 
170  //- Return the names of any (stencil or mesh specific) fields that
171  // should not be interpolated
172  virtual const wordHashSet& nonInterpolatedFields() const
173  {
174  return stencilPtr_().nonInterpolatedFields();
175  }
176 };
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace Foam
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 #endif
186 
187 // ************************************************************************* //
Foam::cellCellStencilObject::nonInterpolatedFields
virtual const wordHashSet & nonInterpolatedFields() const
Return the names of any (stencil or mesh specific) fields that.
Definition: cellCellStencilObject.H:171
Foam::cellCellStencilObject::cellTypes
virtual const labelUList & cellTypes() const
Return the cell type list.
Definition: cellCellStencilObject.H:121
Foam::MeshObject::New
static const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::MoveableMeshObject
Definition: MeshObject.H:220
Foam::cellCellStencilObject::cellInterpolationWeight
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
Definition: cellCellStencilObject.H:153
Foam::HashSet< word, Hash< word > >
Foam::cellCellStencilObject::TypeName
TypeName("cellCellStencilObject")
Foam::cellCellStencilObject::update
virtual bool update()
Update stencils. Return false if nothing changed.
Definition: cellCellStencilObject.H:115
Foam::MeshObject::mesh
const Mesh & mesh() const
Definition: MeshObject.H:122
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::cellCellStencilObject::cellInterpolationMap
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
Definition: cellCellStencilObject.H:133
Foam::cellCellStencilObject::cellCellStencilObject
cellCellStencilObject(const fvMesh &mesh, const bool update=true)
Construct with mesh.
Definition: cellCellStencilObject.H:75
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cellCellStencilObject::cellStencil
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
Definition: cellCellStencilObject.H:140
Foam::cellCellStencilObject::~cellCellStencilObject
virtual ~cellCellStencilObject()=default
Destructor.
Foam::cellCellStencil
Calculation of interpolation stencils.
Definition: cellCellStencil.H:61
Foam::Vector< scalar >
Foam::List< labelList >
Foam::Stencil
MeshObject< fvMesh, Foam::MoveableMeshObject, cellCellStencilObject > Stencil
Definition: cellCellStencilObject.H:45
Foam::cellCellStencilObject::cellInterpolationWeights
virtual const List< scalarList > & cellInterpolationWeights() const
Weights for cellStencil.
Definition: cellCellStencilObject.H:146
Foam::UList< label >
MeshObject.H
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:88
Foam::cellCellStencilObject::interpolationCells
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
Definition: cellCellStencilObject.H:127
cellCellStencil.H
Foam::cellCellStencilObject::movePoints
virtual bool movePoints()
Callback for geometry motion.
Definition: cellCellStencilObject.H:109
Foam::cellCellStencil::New
static autoPtr< cellCellStencil > New(const fvMesh &, const dictionary &dict, const bool update=true)
New function which constructs and returns pointer to a.
Definition: cellCellStencil.C:64
Foam::cellCellStencilObject::stencilWeights
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate weights for a single acceptor.
Definition: cellCellStencilObject.H:160
Foam::cellCellStencilObject
Definition: cellCellStencilObject.H:57
sample
Minimal example by using system/controlDict.functions: