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-2018 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 
94  // Protected Member Functions
95 
96  void walkFront
97  (
98  const scalar layerRelax,
99  labelList& allCellTypes,
100  scalarField& allWeight
101  ) const;
102 
103  //- Find cells next to cells of type PATCH
104  void findHoles
105  (
106  const globalIndex& globalCells,
107  const fvMesh& mesh,
108  const labelList& zoneID,
109  const labelListList& stencil,
111  ) const;
112 
113  //- according to additionalDocumentation/MEJ_oversetMesh.txt
114  void markPatchCells
115  (
116  const fvMesh& mesh,
117  const labelList& cellMap,
118  labelList& patchCellTypes
119  ) const;
120 
121  void combineCellTypes
122  (
123  const label subZoneID,
124  const fvMesh& subMesh,
125  const labelList& subCellMap,
126 
127  const label donorZoneID,
128  const labelListList& toOtherCells,
129  const List<scalarList>& weights,
130  const labelList& otherCells,
131  const labelList& interpolatedOtherPatchTypes,
132 
133  labelListList& allStencil,
134  scalarListList& allWeights,
135  labelList& allCellTypes,
136  labelList& allDonorID
137  ) const;
138 
139  //- interpolate (= combine) patch types
141  (
142  const labelListList& addressing,
143  const labelList& patchTypes,
144  labelList& result
145  ) const;
146 
147  //- interpolate (= combine) patch types
149  (
150  const autoPtr<mapDistribute>& mapPtr,
151  const labelListList& addressing,
152  const labelList& patchTypes,
153  labelList& result
154  ) const;
155 
156 
157 private:
158 
159  // Private Member Functions
160 
161  //- No copy construct
162  cellVolumeWeight(const cellVolumeWeight&) = delete;
163 
164  //- No copy assignment
165  void operator=(const cellVolumeWeight&) = delete;
166 
167 
168 public:
169 
170  //- Runtime type information
171  TypeName("cellVolumeWeight");
172 
173 
174  // Constructors
175 
176  //- Construct from fvMesh
177  cellVolumeWeight(const fvMesh&, const dictionary&, const bool doUpdate);
178 
179 
180  //- Destructor
181  virtual ~cellVolumeWeight();
182 
183 
184  // Member Functions
185 
186  //- Access to volume overlap tolerance
187  scalar overlapTolerance() const
188  {
189  return overlapTolerance_;
190  }
191 
192  //- Update stencils. Return false if nothing changed.
193  virtual bool update();
194 
195  //- Return the cell type list
196  virtual const labelUList& cellTypes() const
197  {
198  return cellTypes_;
199  }
200 
201  //- Indices of interpolated cells
202  virtual const labelUList& interpolationCells() const
203  {
204  return interpolationCells_;
205  }
206 
207  //- Return a communication schedule
208  virtual const mapDistribute& cellInterpolationMap() const
209  {
211  {
212  const_cast<cellVolumeWeight&>(*this).update();
213  }
214  return *cellInterpolationMap_;
215  }
216 
217  //- Per interpolated cell the neighbour cells (in terms of slots as
218  // constructed by above cellInterpolationMap) to interpolate
219  virtual const labelListList& cellStencil() const
220  {
221  return cellStencil_;
222  }
223 
224  //- Weights for cellStencil
225  virtual const List<scalarList>& cellInterpolationWeights() const
226  {
228  }
229 
230  //- Per interpolated cell the interpolation factor. (0 = use
231  // calculated, 1 = use interpolated)
232  virtual const scalarList& cellInterpolationWeight() const
233  {
235  }
236 
237  //- Calculate inverse distance weights for a single acceptor. Revert
238  // to inverse distance (so not consistent with volume overlap!)
239  virtual void stencilWeights
240  (
241  const point& sample,
242  const pointList& donorCcs,
243  scalarList& weights
244  ) const;
245 };
246 
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 } // End namespace cellCellStencils
251 } // End namespace Foam
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 #endif
256 
257 // ************************************************************************* //
volFields.H
Foam::cellCellStencils::cellVolumeWeight::cellTypes_
labelList cellTypes_
Per cell the cell type.
Definition: cellVolumeWeightCellCellStencil.H:75
Foam::cellCellStencils::cellVolumeWeight::cellStencil_
labelListList cellStencil_
Interpolation stencil.
Definition: cellVolumeWeightCellCellStencil.H:84
Foam::cellCellStencils::cellVolumeWeight::interpolationCells_
labelList interpolationCells_
Indices of interpolated cells.
Definition: cellVolumeWeightCellCellStencil.H:78
Foam::cellCellStencils::cellVolumeWeight::cellInterpolationMap
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
Definition: cellVolumeWeightCellCellStencil.H:207
Foam::cellCellStencils::cellVolumeWeight::cellInterpolationWeights
virtual const List< scalarList > & cellInterpolationWeights() const
Weights for cellStencil.
Definition: cellVolumeWeightCellCellStencil.H:224
Foam::cellCellStencils::cellVolumeWeight::~cellVolumeWeight
virtual ~cellVolumeWeight()
Destructor.
Definition: cellVolumeWeightCellCellStencil.C:741
Foam::cellCellStencils::cellVolumeWeight::TypeName
TypeName("cellVolumeWeight")
Runtime type information.
Foam::cellCellStencils::cellVolumeWeight::cellInterpolationWeights_
List< scalarList > cellInterpolationWeights_
Interpolation weights.
Definition: cellVolumeWeightCellCellStencil.H:87
Foam::cellCellStencils::cellVolumeWeight::cellTypes
virtual const labelUList & cellTypes() const
Return the cell type list.
Definition: cellVolumeWeightCellCellStencil.H:195
Foam::cellCellStencils::cellVolumeWeight
Volume-weighted interpolation stencil.
Definition: cellVolumeWeightCellCellStencil.H:54
Foam::cellCellStencils::cellVolumeWeight::markPatchCells
void markPatchCells(const fvMesh &mesh, const labelList &cellMap, labelList &patchCellTypes) const
according to additionalDocumentation/MEJ_oversetMesh.txt
Definition: cellVolumeWeightCellCellStencil.C:445
patchTypes
wordList patchTypes(nPatches)
Foam::cellCellStencils::cellVolumeWeight::combineCellTypes
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
Definition: cellVolumeWeightCellCellStencil.C:553
Foam::cellCellStencils::cellVolumeWeight::interpolatePatchTypes
void interpolatePatchTypes(const labelListList &addressing, const labelList &patchTypes, labelList &result) const
interpolate (= combine) patch types
Definition: cellVolumeWeightCellCellStencil.C:485
Foam::Field< scalar >
Foam::cellCellStencils::cellVolumeWeight::dict_
const dictionary dict_
Dictionary of motion control parameters.
Definition: cellVolumeWeightCellCellStencil.H:69
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
Foam::cellCellStencils::cellVolumeWeight::walkFront
void walkFront(const scalar layerRelax, labelList &allCellTypes, scalarField &allWeight) const
Definition: cellVolumeWeightCellCellStencil.C:60
Foam::cellCellStencils::cellVolumeWeight::cellInterpolationMap_
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
Definition: cellVolumeWeightCellCellStencil.H:81
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::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::cellCellStencils::cellVolumeWeight::overlapTolerance
scalar overlapTolerance() const
Access to volume overlap tolerance.
Definition: cellVolumeWeightCellCellStencil.H:186
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cellCellStencils::cellVolumeWeight::cellInterpolationWeight
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
Definition: cellVolumeWeightCellCellStencil.H:231
Foam::cellCellStencils::cellVolumeWeight::findHoles
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.
Definition: cellVolumeWeightCellCellStencil.C:204
Foam::cellCellStencils::cellVolumeWeight::interpolationCells
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
Definition: cellVolumeWeightCellCellStencil.H:201
Foam::cellCellStencil
Calculation of interpolation stencils.
Definition: cellCellStencil.H:61
Foam::cellCellStencils::cellVolumeWeight::defaultOverlapTolerance_
static scalar defaultOverlapTolerance_
Default overlap tolerance. Fraction of volume.
Definition: cellVolumeWeightCellCellStencil.H:63
Foam::Vector< scalar >
Foam::List< label >
Foam::UList< label >
cellCellStencil.H
Foam::cellCellStencils::cellVolumeWeight::stencilWeights
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor. Revert.
Definition: cellVolumeWeightCellCellStencil.C:1279
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::cellCellStencils::cellVolumeWeight::cellInterpolationWeight_
volScalarField cellInterpolationWeight_
Amount of interpolation.
Definition: cellVolumeWeightCellCellStencil.H:90
Foam::cellCellStencils::cellVolumeWeight::overlapTolerance_
scalar overlapTolerance_
Tolerance for volume overlap. Fraction of volume.
Definition: cellVolumeWeightCellCellStencil.H:72
Foam::cellCellStencils::cellVolumeWeight::update
virtual bool update()
Update stencils. Return false if nothing changed.
Definition: cellVolumeWeightCellCellStencil.C:747
sample
Minimal example by using system/controlDict.functions:
Foam::cellCellStencils::cellVolumeWeight::cellStencil
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
Definition: cellVolumeWeightCellCellStencil.H:218