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-2022 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 Foam_cellCellStencil_H
44 #define Foam_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 // Forward Declarations
57 class cellCellStencil;
58 class mapDistribute;
59 
60 template<>
61 Ostream& operator<<(Ostream&, const InfoProxy<cellCellStencil>&);
62 
63 
64 /*---------------------------------------------------------------------------*\
65  Class cellCellStencil Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class cellCellStencil
69 {
70 public:
71 
72  enum patchCellType
73  {
74  OTHER = 0, // not on special patch
75  PATCH = 1, // next to (non-coupled) boundary
76  OVERSET = 2 // next to 'overset' boundary
77  };
78 
79  enum cellType
80  {
81  CALCULATED = 0, // normal operation
82  INTERPOLATED = 1, // interpolated
83  HOLE = 2, // hole
84  SPECIAL = 3, // hole that is donor
85  POROUS = 4 // tag for special treatment due to lack of donors
86  };
87 
88 
89 protected:
90 
91  // Protected data
92 
93  //- Mode type names
94  static const Enum<cellType> cellTypeNames_;
95 
96  //- Reference to the mesh
97  const fvMesh& mesh_;
98 
99  //- Set of fields that should not be interpolated
102  //- Dictionary of motion control parameters
103  const dictionary dict_;
104 
105 
106  // Protected Member Functions
107 
108  //- Helper: create volScalarField for postprocessing.
109  template<class Type>
111  (
112  const fvMesh& mesh,
113  const word& name,
114  const UList<Type>&
115  );
116 
117  //- Helper: strip off trailing _0
118  static word baseName(const word& name);
119 
120 
121 private:
122 
123  // Private Member Functions
124 
125  //- No copy construct
126  cellCellStencil(const cellCellStencil&) = delete;
127 
128  //- No copy assignment
129  void operator=(const cellCellStencil&) = delete;
130 
131 
132 public:
133 
134  //- Runtime type information
135  TypeName("cellCellStencil");
136 
137 
138  // Declare run-time constructor selection tables
139 
141  (
142  autoPtr,
144  mesh,
145  (
146  const fvMesh& mesh,
147  const dictionary& dict,
148  const bool update
149  ),
150  (mesh, dict, update)
151  );
152 
153 
154  // Constructors
155 
156  //- Construct from fvMesh
157  cellCellStencil(const fvMesh&);
158 
159  //- New function which constructs and returns pointer to a
160  // cellCellStencil
162  (
163  const fvMesh&,
164  const dictionary& dict,
165  const bool update = true
166  );
167 
168 
169  //- Destructor
170  virtual ~cellCellStencil();
171 
172 
173  // Member Functions
174 
175  //- Update stencils. Return false if nothing changed.
176  virtual bool update() = 0;
177 
178  //- Return the cell type list
179  virtual const labelUList& cellTypes() const = 0;
180 
181  //- Indices of interpolated cells
182  virtual const labelUList& interpolationCells() const = 0;
183 
184  //- Return a communication schedule
185  virtual const mapDistribute& cellInterpolationMap() const = 0;
186 
187  //- Per interpolated cell the neighbour cells (in terms of slots as
188  // constructed by above cellInterpolationMap) to interpolate
189  virtual const labelListList& cellStencil() const = 0;
190 
191  //- Weights for cellStencil
192  virtual const List<scalarList>& cellInterpolationWeights() const = 0;
193 
194  //- Per interpolated cell the interpolation factor. (0 = use
195  // calculated, 1 = use interpolated)
196  virtual const scalarList& cellInterpolationWeight() const = 0;
197 
198  //- Calculate weights for a single acceptor
199  virtual void stencilWeights
200  (
201  const point& sample,
202  const pointList& donorCcs,
203  scalarList& weights
204  ) const = 0;
205 
206  //- Return the names of any (stencil or mesh specific) fields that
207  // should not be interpolated
208  virtual const wordHashSet& nonInterpolatedFields() const;
209 
210  //- Return non-const non-interpolating fields
212 
213  //- Helper: is stencil fully local
214  bool localStencil(const labelUList&) const;
215 
216  //- Helper: get reference to registered zoneID. Loads volScalarField
217  // if not registered.
218  static const labelIOList& zoneID(const fvMesh&);
219 
220  //- Helper: get reference to registered zoneID. Loads volScalarField
221  // if not registered.
222  const labelIOList& zoneID() const
223  {
224  return zoneID(mesh_);
225  }
226 
227  //- Count occurrences (in parallel)
228  static labelList count(const label size, const labelUList& lst);
229 
230  //- Helper: create cell-cell addressing in global numbering
231  static void globalCellCells
232  (
233  const globalIndex& gi,
234  const polyMesh& mesh,
235  const boolList& isValidDonor,
236  const labelList& selectedCells,
237  labelListList& cellCells,
238  pointListList& cellCellCentres
239  );
240 
241  //- Interpolation of acceptor cells from donor cells
242  template<class T>
243  static void interpolate
244  (
245  Field<T>& psi,
246  const fvMesh& mesh,
247  const cellCellStencil& overlap,
248  const List<scalarList>& wghts
249  );
250 
251  //- Explicit interpolation of acceptor cells from donor cells. Looks
252  //- up cellCellStencil.
253  template<class T>
254  void interpolate(const fvMesh& mesh, Field<T>& psi) const;
255 
256  //- Explicit interpolation of acceptor cells from donor cells with
257  //- boundary condition handling
258  template<class GeoField>
259  void interpolate(GeoField& psi) const;
260 
261  //- Explicit interpolation of all registered fields. No boundary
262  //- conditions. Excludes selected fields (and their old-time fields)
263  template<class GeoField>
264  void interpolate
265  (
266  const fvMesh& mesh,
267  const wordHashSet& suppressed
268  ) const;
269 
270  //- Version of correctBoundaryConditions that excludes 'overset' bcs
271  template<class GeoField, class SuppressBC>
272  static void correctBoundaryConditions(GeoField& psi);
274  //- Surround holes with layer(s) of interpolated cells
275  void walkFront
276  (
277  const globalIndex& globalCells,
278  const scalar layerRelax,
279  const labelListList& allStencil,
280  labelList& allCellTypes,
281  scalarField& allWeight,
282  const scalarList& compactCellVol,
283  const labelListList& compactStencil,
284  const labelList& zoneID,
285  const label holeLayers,
286  const label useLayer
287  ) const;
288 
289  //- Set up front using allCellTypes
290  void setUpFront
291  (
292  const labelList& allCellTypes,
293  bitSet& isFront
294 
295  ) const;
296 
297  //- Set up front on overset patches
299  (
300  const labelList& allCellTypes,
301  bitSet& isFront
302  ) const;
303 
304  //- Seed faces of cell with wantedFraction (if higher than current)
305  void seedCell
306  (
307  const label cellI,
308  const scalar wantedFraction,
309  bitSet& isFront,
310  scalarField& fraction
311  ) const;
312 
313 
314  // Write
315 
316  //- Return info proxy,
317  //- used to print stencil information to a stream
318  InfoProxy<cellCellStencil> info() const noexcept { return *this; }
319 };
320 
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 } // End namespace Foam
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 #ifdef NoRepository
329 # include "cellCellStencilTemplates.C"
330 #endif
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 #endif
335 
336 // ************************************************************************* //
virtual const wordHashSet & nonInterpolatedFields() const
Return the names of any (stencil or mesh specific) fields that.
dictionary dict
virtual const List< scalarList > & cellInterpolationWeights() const =0
Weights for cellStencil.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
void seedCell(const label cellI, const scalar wantedFraction, bitSet &isFront, scalarField &fraction) const
Seed faces of cell with wantedFraction (if higher than current)
static autoPtr< cellCellStencil > New(const fvMesh &, const dictionary &dict, const bool update=true)
New function which constructs and returns pointer to a.
void setUpFront(const labelList &allCellTypes, bitSet &isFront) const
Set up front using allCellTypes.
virtual ~cellCellStencil()
Destructor.
static void interpolate(Field< T > &psi, const fvMesh &mesh, const cellCellStencil &overlap, const List< scalarList > &wghts)
Interpolation of acceptor cells from donor cells.
wordHashSet nonInterpolatedFields_
Set of fields that should not be interpolated.
bool localStencil(const labelUList &) const
Helper: is stencil fully local.
void setUpFrontOnOversetPatch(const labelList &allCellTypes, bitSet &isFront) const
Set up front on overset patches.
virtual const scalarList & cellInterpolationWeight() const =0
Per interpolated cell the interpolation factor. (0 = use.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:62
static labelList count(const label size, const labelUList &lst)
Count occurrences (in parallel)
dynamicFvMesh & mesh
virtual const mapDistribute & cellInterpolationMap() const =0
Return a communication schedule.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
Calculation of interpolation stencils.
InfoProxy< cellCellStencil > info() const noexcept
Return info proxy, used to print stencil information to a stream.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const cellCellStencilObject & overlap
Definition: correctPhi.H:57
static tmp< volScalarField > createField(const fvMesh &mesh, const word &name, const UList< Type > &)
Helper: create volScalarField for postprocessing.
virtual const labelUList & interpolationCells() const =0
Indices of interpolated cells.
const fvMesh & mesh_
Reference to the mesh.
static word baseName(const word &name)
Helper: strip off trailing _0.
const direction noexcept
Definition: Scalar.H:258
virtual bool update()=0
Update stencils. Return false if nothing changed.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
TypeName("cellCellStencil")
Runtime type information.
virtual const labelListList & cellStencil() const =0
Per interpolated cell the neighbour cells (in terms of slots as.
virtual const labelUList & cellTypes() const =0
Return the cell type list.
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.
declareRunTimeSelectionTable(autoPtr, cellCellStencil, mesh,(const fvMesh &mesh, const dictionary &dict, const bool update),(mesh, dict, update))
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
Class containing processor-to-processor mapping information.
A helper class for outputting values to Ostream.
Definition: ensightCells.H:43
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
static const Enum< cellType > cellTypeNames_
Mode type names.
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const =0
Calculate weights for a single acceptor.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
static void correctBoundaryConditions(GeoField &psi)
Version of correctBoundaryConditions that excludes &#39;overset&#39; bcs.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const volScalarField & psi
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
void walkFront(const globalIndex &globalCells, const scalar layerRelax, const labelListList &allStencil, labelList &allCellTypes, scalarField &allWeight, const scalarList &compactCellVol, const labelListList &compactStencil, const labelList &zoneID, const label holeLayers, const label useLayer) const
Surround holes with layer(s) of interpolated cells.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Bookkeeping for patch definitions.
Definition: PDRpatchDef.H:49
Namespace for OpenFOAM.
const dictionary dict_
Dictionary of motion control parameters.