cellCellStencil.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::cellCellStencil
28 
29 Description
30  Calculation of interpolation stencils.
31 
32  Looks up zoneID labelIOList to give the zoning. Wrapped in
33  MeshObject as cellCellStencilObject. Kept separate so meshes can
34  implement more clever methods (e.g. solid body motion does not require
35  full recalculation)
36 
37 SourceFiles
38  cellCellStencil.C
39  cellCellStencilObject.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef cellCellStencil_H
44 #define cellCellStencil_H
45 
46 #include "scalarList.H"
47 #include "mapDistribute.H"
48 #include "pointList.H"
49 #include "volFields.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 class mapDistribute;
57 
58 /*---------------------------------------------------------------------------*\
59  Class cellCellStencil Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 class cellCellStencil
63 {
64 public:
65 
66  enum patchCellType
67  {
68  OTHER = 0, // not on special patch
69  PATCH = 1, // next to (non-coupled) boundary
70  OVERSET = 2 // next to 'overset' boundary
71  };
72 
73  enum cellType
74  {
75  CALCULATED = 0, // normal operation
76  INTERPOLATED = 1, // interpolated
77  HOLE = 2 // hole
78  };
79 
80 
81 protected:
82 
83  // Protected data
84 
85  //- Mode type names
86  static const Enum<cellType> cellTypeNames_;
87 
88  //- Reference to the mesh
89  const fvMesh& mesh_;
90 
91  //- Set of fields that should not be interpolated
93 
94 
95  // Protected Member Functions
96 
97  //- Count occurrences (in parallel)
98  static labelList count(const label size, const labelUList& lst);
99 
100  //- Helper: create volScalarField for postprocessing.
101  template<class Type>
103  (
104  const fvMesh& mesh,
105  const word& name,
106  const UList<Type>&
107  );
108 
109 
110 private:
111 
112  // Private Member Functions
113 
114  //- No copy construct
115  cellCellStencil(const cellCellStencil&) = delete;
116 
117  //- No copy assignment
118  void operator=(const cellCellStencil&) = delete;
119 
120 
121 public:
122 
123  //- Runtime type information
124  TypeName("cellCellStencil");
125 
126 
127  // Declare run-time constructor selection tables
128 
130  (
131  autoPtr,
133  mesh,
134  (
135  const fvMesh& mesh,
136  const dictionary& dict,
137  const bool update
138  ),
139  (mesh, dict, update)
140  );
141 
142 
143  // Constructors
144 
145  //- Construct from fvMesh
146  cellCellStencil(const fvMesh&);
147 
148  //- New function which constructs and returns pointer to a
149  // cellCellStencil
151  (
152  const fvMesh&,
153  const dictionary& dict,
154  const bool update = true
155  );
156 
157 
158  //- Destructor
159  virtual ~cellCellStencil();
160 
161 
162  // Member Functions
163 
164  //- Update stencils. Return false if nothing changed.
165  virtual bool update() = 0;
166 
167  //- Return the cell type list
168  virtual const labelUList& cellTypes() const = 0;
169 
170  //- Indices of interpolated cells
171  virtual const labelUList& interpolationCells() const = 0;
172 
173  //- Return a communication schedule
174  virtual const mapDistribute& cellInterpolationMap() const = 0;
175 
176  //- Per interpolated cell the neighbour cells (in terms of slots as
177  // constructed by above cellInterpolationMap) to interpolate
178  virtual const labelListList& cellStencil() const = 0;
179 
180  //- Weights for cellStencil
181  virtual const List<scalarList>& cellInterpolationWeights() const = 0;
182 
183  //- Per interpolated cell the interpolation factor. (0 = use
184  // calculated, 1 = use interpolated)
185  virtual const scalarList& cellInterpolationWeight() const = 0;
186 
187  //- Calculate weights for a single acceptor
188  virtual void stencilWeights
189  (
190  const point& sample,
191  const pointList& donorCcs,
192  scalarList& weights
193  ) const = 0;
194 
195  //- Return the names of any (stencil or mesh specific) fields that
196  // should not be interpolated
197  virtual const wordHashSet& nonInterpolatedFields() const;
198 
199  //- Return non-const non-interpolating fields
201 
202  //- Helper: is stencil fully local
203  bool localStencil(const labelUList&) const;
204 
205  //- Helper: get reference to registered zoneID. Loads volScalarField
206  // if not registered.
207  static const labelIOList& zoneID(const fvMesh&);
208 
209  //- Helper: get reference to registered zoneID. Loads volScalarField
210  // if not registered.
211  const labelIOList& zoneID() const
212  {
213  return zoneID(mesh_);
214  }
215 
216  //- Helper: create cell-cell addressing in global numbering
217  static void globalCellCells
218  (
219  const globalIndex& gi,
220  const polyMesh& mesh,
221  const boolList& isValidDonor,
222  const labelList& selectedCells,
223  labelListList& cellCells,
224  pointListList& cellCellCentres
225  );
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #ifdef NoRepository
236 # include "cellCellStencilTemplates.C"
237 #endif
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
Foam::cellCellStencil::OVERSET
Definition: cellCellStencil.H:69
volFields.H
Foam::Enum< cellType >
Foam::cellCellStencil::HOLE
Definition: cellCellStencil.H:76
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::cellCellStencil::nonInterpolatedFields_
wordHashSet nonInterpolatedFields_
Set of fields that should not be interpolated.
Definition: cellCellStencil.H:91
Foam::cellCellStencil::OTHER
Definition: cellCellStencil.H:67
Foam::cellCellStencil::localStencil
bool localStencil(const labelUList &) const
Helper: is stencil fully local.
Definition: cellCellStencil.C:170
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::cellCellStencil::count
static labelList count(const label size, const labelUList &lst)
Count occurrences (in parallel)
Definition: cellCellStencil.C:143
Foam::cellCellStencil::nonInterpolatedFields
virtual const wordHashSet & nonInterpolatedFields() const
Return the names of any (stencil or mesh specific) fields that.
Definition: cellCellStencil.C:158
Foam::cellCellStencil::patchCellType
patchCellType
Definition: cellCellStencil.H:65
Foam::HashSet< word, Hash< word > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
scalarList.H
Foam::cellCellStencil::TypeName
TypeName("cellCellStencil")
Runtime type information.
Foam::cellCellStencil::cellStencil
virtual const labelListList & cellStencil() const =0
Per interpolated cell the neighbour cells (in terms of slots as.
Foam::cellCellStencil::cellTypes
virtual const labelUList & cellTypes() const =0
Return the cell type list.
Foam::cellCellStencil::globalCellCells
static void globalCellCells(const globalIndex &gi, const polyMesh &mesh, const boolList &isValidDonor, const labelList &selectedCells, labelListList &cellCells, pointListList &cellCellCentres)
Helper: create cell-cell addressing in global numbering.
Definition: cellCellStencil.C:184
Foam::cellCellStencil::INTERPOLATED
Definition: cellCellStencil.H:75
Foam::cellCellStencil::zoneID
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
Definition: cellCellStencil.H:210
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
pointList.H
Foam::cellCellStencil::~cellCellStencil
virtual ~cellCellStencil()
Destructor.
Definition: cellCellStencil.C:93
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::cellCellStencil::stencilWeights
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const =0
Calculate weights for a single acceptor.
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::cellCellStencil::cellInterpolationWeight
virtual const scalarList & cellInterpolationWeight() const =0
Per interpolated cell the interpolation factor. (0 = use.
cellCellStencilTemplates.C
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cellCellStencil::cellInterpolationMap
virtual const mapDistribute & cellInterpolationMap() const =0
Return a communication schedule.
Foam::cellCellStencil::interpolationCells
virtual const labelUList & interpolationCells() const =0
Indices of interpolated cells.
mapDistribute.H
Foam::cellCellStencil
Calculation of interpolation stencils.
Definition: cellCellStencil.H:61
Foam::Vector< scalar >
Foam::List< label >
Foam::UList< label >
Foam::cellCellStencil::createField
static tmp< volScalarField > createField(const fvMesh &mesh, const word &name, const UList< Type > &)
Helper: create volScalarField for postprocessing.
Foam::PDRpatchDef
Bookkeeping for patch definitions.
Definition: PDRpatchDef.H:53
Foam::IOList< label >
Foam::cellCellStencil::cellInterpolationWeights
virtual const List< scalarList > & cellInterpolationWeights() const =0
Weights for cellStencil.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::cellCellStencil::mesh_
const fvMesh & mesh_
Reference to the mesh.
Definition: cellCellStencil.H:88
Foam::cellCellStencil::CALCULATED
Definition: cellCellStencil.H:74
Foam::cellCellStencil::update
virtual bool update()=0
Update stencils. Return false if nothing changed.
Foam::cellCellStencil::cellType
cellType
Definition: cellCellStencil.H:72
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::cellCellStencil::cellTypeNames_
static const Enum< cellType > cellTypeNames_
Mode type names.
Definition: cellCellStencil.H:85
Foam::cellCellStencil::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, cellCellStencil, mesh,(const fvMesh &mesh, const dictionary &dict, const bool update),(mesh, dict, update))
sample
Minimal example by using system/controlDict.functions: