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-2023 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  //- Helper: populate nonInterpolatedFields_ with motion solver
121  // fields
122  void suppressMotionFields();
123 
124 
125 private:
126 
127  // Private Member Functions
128 
129  //- No copy construct
130  cellCellStencil(const cellCellStencil&) = delete;
131 
132  //- No copy assignment
133  void operator=(const cellCellStencil&) = delete;
134 
135 
136 public:
137 
138  //- Runtime type information
139  TypeName("cellCellStencil");
140 
141 
142  // Declare run-time constructor selection tables
143 
145  (
146  autoPtr,
148  mesh,
149  (
150  const fvMesh& mesh,
151  const dictionary& dict,
152  const bool update
153  ),
154  (mesh, dict, update)
155  );
156 
157 
158  // Constructors
159 
160  //- Construct from fvMesh
161  cellCellStencil(const fvMesh&);
162 
163  //- New function which constructs and returns pointer to a
164  // cellCellStencil
166  (
167  const fvMesh&,
168  const dictionary& dict,
169  const bool update = true
170  );
171 
172 
173  //- Destructor
174  virtual ~cellCellStencil();
175 
176 
177  // Member Functions
178 
179  //- Update stencils. Return false if nothing changed.
180  virtual bool update() = 0;
181 
182  //- Return the cell type list
183  virtual const labelUList& cellTypes() const = 0;
184 
185  //- Indices of interpolated cells
186  virtual const labelUList& interpolationCells() const = 0;
187 
188  //- Return a communication schedule
189  virtual const mapDistribute& cellInterpolationMap() const = 0;
190 
191  //- Per interpolated cell the neighbour cells (in terms of slots as
192  // constructed by above cellInterpolationMap) to interpolate
193  virtual const labelListList& cellStencil() const = 0;
194 
195  //- Weights for cellStencil
196  virtual const List<scalarList>& cellInterpolationWeights() const = 0;
197 
198  //- Per interpolated cell the interpolation factor. (0 = use
199  // calculated, 1 = use interpolated)
200  virtual const scalarList& cellInterpolationWeight() const = 0;
201 
202  //- Calculate weights for a single acceptor
203  virtual void stencilWeights
204  (
205  const point& sample,
206  const pointList& donorCcs,
207  scalarList& weights
208  ) const = 0;
209 
210  //- Return the names of any (stencil or mesh specific) fields that
211  // should not be interpolated
212  virtual const wordHashSet& nonInterpolatedFields() const;
213 
214  //- Return non-const non-interpolating fields
216 
217  //- Helper: is stencil fully local
218  bool localStencil(const labelUList&) const;
219 
220  //- Helper: get reference to registered zoneID. Loads volScalarField
221  // if not registered.
222  static const labelIOList& zoneID(const fvMesh&);
223 
224  //- Helper: get reference to registered zoneID. Loads volScalarField
225  // if not registered.
226  const labelIOList& zoneID() const
227  {
228  return zoneID(mesh_);
229  }
230 
231  //- Count occurrences (in parallel)
232  static labelList count(const label size, const labelUList& lst);
233 
234  //- Helper: create cell-cell addressing in global numbering
235  static void globalCellCells
236  (
237  const globalIndex& gi,
238  const polyMesh& mesh,
239  const boolList& isValidDonor,
240  const labelList& selectedCells,
241  labelListList& cellCells,
242  pointListList& cellCellCentres
243  );
244 
245  //- Interpolation of acceptor cells from donor cells
246  template<class T>
247  static void interpolate
248  (
249  Field<T>& psi,
250  const fvMesh& mesh,
251  const cellCellStencil& overlap,
252  const List<scalarList>& wghts
253  );
254 
255  //- Explicit interpolation of acceptor cells from donor cells. Looks
256  //- up cellCellStencil.
257  template<class T>
258  void interpolate(const fvMesh& mesh, Field<T>& psi) const;
259 
260  //- Explicit interpolation of acceptor cells from donor cells with
261  //- boundary condition handling
262  template<class GeoField>
263  void interpolate(GeoField& psi) const;
264 
265  //- Explicit interpolation of all registered fields. No boundary
266  //- conditions. Excludes selected fields (and their old-time fields)
267  template<class GeoField>
268  void interpolate
269  (
270  const fvMesh& mesh,
271  const wordHashSet& suppressed
272  ) const;
273 
274  //- Version of correctBoundaryConditions that excludes 'overset' bcs
275  template<class GeoField, class SuppressBC>
276  static void correctBoundaryConditions(GeoField& psi);
277 
278  //- Surround holes with layer(s) of interpolated cells
279  void walkFront
280  (
281  const globalIndex& globalCells,
282  const scalar layerRelax,
283  const labelListList& allStencil,
284  labelList& allCellTypes,
285  scalarField& allWeight,
286  const scalarList& compactCellVol,
287  const labelListList& compactStencil,
288  const labelList& zoneID,
289  const label holeLayers,
290  const label useLayer
291  ) const;
292 
293  //- Set up front using allCellTypes
294  void setUpFront
295  (
296  const labelList& allCellTypes,
297  bitSet& isFront
298 
299  ) const;
300 
301  //- Set up front on overset patches
303  (
304  const labelList& allCellTypes,
305  bitSet& isFront
306  ) const;
307 
308  //- Seed faces of cell with wantedFraction (if higher than current)
309  void seedCell
310  (
311  const label cellI,
312  const scalar wantedFraction,
313  bitSet& isFront,
314  scalarField& fraction
315  ) const;
316 
317 
318  // Write
319 
320  //- Return info proxy,
321  //- used to print stencil information to a stream
322  InfoProxy<cellCellStencil> info() const noexcept { return *this; }
323 };
324 
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 } // End namespace Foam
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
332 #ifdef NoRepository
333 # include "cellCellStencilTemplates.C"
334 #endif
335 
336 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 
338 #endif
339 
340 // ************************************************************************* //
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:129
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:61
static labelList count(const label size, const labelUList &lst)
Count occurrences (in parallel)
void suppressMotionFields()
Helper: populate nonInterpolatedFields_ with motion solver.
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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:78
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:74
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.