cellVolumeWeightCellCellStencil.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) 2016-2021 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::cellCellStencils::cellVolumeWeight
28 
29 Description
30  Volume-weighted interpolation stencil
31 
32 SourceFiles
33  cellVolumeWeightCellCellStencil.C
34  cellVolumeWeightCellCellStencilTemplates.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef cellCellStencils_cellVolumeWeight_H
39 #define cellCellStencils_cellVolumeWeight_H
40 
41 #include "cellCellStencil.H"
42 #include "volFields.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace cellCellStencils
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class cellVolumeWeight Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class cellVolumeWeight
56 :
57  public cellCellStencil
58 {
59 protected:
60 
61  // Static data members
62 
63  //- Default overlap tolerance. Fraction of volume
64  static scalar defaultOverlapTolerance_;
65 
66 
67  // Protected data
68 
69  //- Dictionary of motion control parameters
70  const dictionary dict_;
71 
72  //- Tolerance for volume overlap. Fraction of volume
73  scalar overlapTolerance_;
74 
75  //- Per cell the cell type
77 
78  //- Indices of interpolated cells
80 
81  //- Fetch interpolated cells
83 
84  //- Interpolation stencil
86 
87  //- Interpolation weights
89 
90  //- Amount of interpolation
92 
93  //- Allow interpolared as donors
95 
96  // Protected Member Functions
97 
98 
99  //- Find cells next to cells of type PATCH
100  void findHoles
101  (
102  const globalIndex& globalCells,
103  const fvMesh& mesh,
105  const labelListList& stencil,
107  ) const;
108 
109  //- according to additionalDocumentation/MEJ_oversetMesh.txt
110  void markPatchCells
111  (
112  const fvMesh& mesh,
113  const labelList& cellMap,
114  labelList& patchCellTypes
115  ) const;
116 
117  void combineCellTypes
118  (
119  const label subZoneID,
120  const fvMesh& subMesh,
121  const labelList& subCellMap,
122 
123  const label donorZoneID,
124  const labelListList& toOtherCells,
125  const List<scalarList>& weights,
126  const labelList& otherCells,
127  const labelList& interpolatedOtherPatchTypes,
128 
129  labelListList& allStencil,
130  scalarListList& allWeights,
131  labelList& allCellTypes,
132  labelList& allDonorID
133  ) const;
134 
135  //- interpolate (= combine) patch types
137  (
138  const labelListList& addressing,
139  const labelList& patchTypes,
140  labelList& result
141  ) const;
142 
143  //- interpolate (= combine) patch types
145  (
146  const autoPtr<mapDistribute>& mapPtr,
147  const labelListList& addressing,
148  const labelList& patchTypes,
149  labelList& result
150  ) const;
151 
152 
153 private:
154 
155  // Private Member Functions
156 
157  //- No copy construct
158  cellVolumeWeight(const cellVolumeWeight&) = delete;
159 
160  //- No copy assignment
161  void operator=(const cellVolumeWeight&) = delete;
162 
163 
164 public:
165 
166  //- Runtime type information
167  TypeName("cellVolumeWeight");
168 
169 
170  // Constructors
171 
172  //- Construct from fvMesh
173  cellVolumeWeight(const fvMesh&, const dictionary&, const bool doUpdate);
174 
175 
176  //- Destructor
177  virtual ~cellVolumeWeight();
178 
179 
180  // Member Functions
181 
182  //- Access to volume overlap tolerance
183  scalar overlapTolerance() const
184  {
185  return overlapTolerance_;
186  }
187 
188  //- Update stencils. Return false if nothing changed.
189  virtual bool update();
190 
191  //- Return the cell type list
192  virtual const labelUList& cellTypes() const
193  {
194  return cellTypes_;
195  }
196 
197  //- Indices of interpolated cells
198  virtual const labelUList& interpolationCells() const
199  {
200  return interpolationCells_;
201  }
202 
203  //- Return a communication schedule
204  virtual const mapDistribute& cellInterpolationMap() const
205  {
207  {
208  const_cast<cellVolumeWeight&>(*this).update();
209  }
210  return *cellInterpolationMap_;
211  }
212 
213  //- Per interpolated cell the neighbour cells (in terms of slots as
214  // constructed by above cellInterpolationMap) to interpolate
215  virtual const labelListList& cellStencil() const
216  {
217  return cellStencil_;
218  }
219 
220  //- Weights for cellStencil
221  virtual const List<scalarList>& cellInterpolationWeights() const
222  {
224  }
225 
226  //- Per interpolated cell the interpolation factor. (0 = use
227  // calculated, 1 = use interpolated)
228  virtual const scalarList& cellInterpolationWeight() const
229  {
231  }
232 
233  //- Calculate inverse distance weights for a single acceptor. Revert
234  // to inverse distance (so not consistent with volume overlap!)
235  virtual void stencilWeights
236  (
237  const point& sample,
238  const pointList& donorCcs,
239  scalarList& weights
240  ) const;
241 };
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace cellCellStencils
247 } // End namespace Foam
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
scalar overlapTolerance_
Tolerance for volume overlap. Fraction of volume.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor. Revert.
void combineCellTypes(const label subZoneID, const fvMesh &subMesh, const labelList &subCellMap, const label donorZoneID, const labelListList &toOtherCells, const List< scalarList > &weights, const labelList &otherCells, const labelList &interpolatedOtherPatchTypes, labelListList &allStencil, scalarListList &allWeights, labelList &allCellTypes, labelList &allDonorID) const
wordList patchTypes(nPatches)
const bool allowInterpolatedDonors_
Allow interpolared as donors.
virtual bool update()
Update stencils. Return false if nothing changed.
void findHoles(const globalIndex &globalCells, const fvMesh &mesh, const labelList &zoneID, const labelListList &stencil, labelList &cellTypes) const
Find cells next to cells of type PATCH.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
const dictionary dict_
Dictionary of motion control parameters.
volScalarField cellInterpolationWeight_
Amount of interpolation.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:61
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
dynamicFvMesh & mesh
labelList interpolationCells_
Indices of interpolated cells.
Calculation of interpolation stencils.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
List< scalarList > cellInterpolationWeights_
Interpolation weights.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
labelListList cellStencil_
Interpolation stencil.
void interpolatePatchTypes(const labelListList &addressing, const labelList &patchTypes, labelList &result) const
interpolate (= combine) patch types
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
Volume-weighted interpolation stencil.
scalar overlapTolerance() const
Access to volume overlap tolerance.
static scalar defaultOverlapTolerance_
Default overlap tolerance. Fraction of volume.
virtual const List< scalarList > & cellInterpolationWeights() const
Weights for cellStencil.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const labelUList & cellTypes() const
Return the cell type list.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
void markPatchCells(const fvMesh &mesh, const labelList &cellMap, labelList &patchCellTypes) const
according to additionalDocumentation/MEJ_oversetMesh.txt
TypeName("cellVolumeWeight")
Runtime type information.
Namespace for OpenFOAM.