smoothDelta.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::LESModels::smoothDelta
26 
27 Description
28  Smoothed delta which takes a given simple geometric delta and applies
29  smoothing to it such that the ratio of deltas between two cells is no
30  larger than a specified amount, typically 1.15.
31 
32 SourceFiles
33  smoothDelta.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef smoothDelta_H
38 #define smoothDelta_H
39 
40 #include "LESdelta.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace LESModels
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class smoothDelta Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 class smoothDelta
54 :
55  public LESdelta
56 {
57 public:
58 
59  //- Public member class used by mesh-wave to propagate the delta-ratio
60  class deltaData
61  {
62  // Private Member Data
63 
64  //-
65  scalar delta_;
66 
67 
68  // Private Member Functions
69 
70  //- Update. Gets information from neighbouring face/cell and
71  // uses this to update itself (if necessary) and return true.
72  template<class TrackingData>
73  inline bool update
74  (
75  const deltaData& w2,
76  const scalar scale,
77  const scalar tol,
78  TrackingData& td
79  );
80 
81 
82  public:
83 
84  // Constructors
85 
86  //- Construct null
87  inline deltaData();
88 
89  //- Construct from delta value
90  inline deltaData(const scalar delta);
91 
92 
93  // Member Functions
94 
95  // Access
96 
97  //-
98  scalar delta() const
99  {
100  return delta_;
101  }
102 
103 
104  // Needed by FvFaceCellWave
105 
106  //- Check whether origin has been changed at all or
107  // still contains original (invalid) value.
108  template<class TrackingData>
109  inline bool valid(TrackingData& td) const;
110 
111  //- Check for identical geometrical data.
112  // Used for cyclics checking.
113  template<class TrackingData>
114  inline bool sameGeometry
115  (
116  const fvMesh&,
117  const deltaData&,
118  const scalar,
119  TrackingData& td
120  ) const;
121 
122  //- Transform across an interface
123  template<class TrackingData>
124  inline void transform
125  (
126  const fvPatch& patch,
127  const label patchFacei,
128  const transformer& transform,
129  TrackingData& td
130  );
131 
132  //- Influence of neighbouring face.
133  template<class TrackingData>
134  inline bool updateCell
135  (
136  const fvMesh&,
137  const label thisCelli,
138  const labelPair& neighbourPatchAndFacei,
139  const deltaData& neighbourInfo,
140  const scalar tol,
141  TrackingData& td
142  );
143 
144  //- Influence of neighbouring cell.
145  template<class TrackingData>
146  inline bool updateFace
147  (
148  const fvMesh&,
149  const labelPair& thisPatchAndFacei,
150  const label neighbourCelli,
151  const deltaData& neighbourInfo,
152  const scalar tol,
153  TrackingData& td
154  );
155 
156  //- Influence of different value on same face.
157  template<class TrackingData>
158  inline bool updateFace
159  (
160  const fvMesh&,
161  const labelPair& thisPatchAndFacei,
162  const deltaData& neighbourInfo,
163  const scalar tol,
164  TrackingData& td
165  );
166 
167  //- Same (like operator==)
168  template<class TrackingData>
169  inline bool equal(const deltaData&, TrackingData& td) const;
170 
171  // Member Operators
172 
173  inline bool operator==(const deltaData&) const;
174 
175  inline bool operator!=(const deltaData&) const;
176 
177 
178  // IOstream Operators
179 
180  friend Ostream& operator<<(Ostream& os, const deltaData& wDist)
181  {
182  return os << wDist.delta_;
183  }
184 
185  friend Istream& operator>>(Istream& is, deltaData& wDist)
186  {
187  return is >> wDist.delta_;
188  }
189  };
190 
191 
192 private:
193 
194  // Private Data
195 
196  autoPtr<LESdelta> geometricDelta_;
197  scalar maxDeltaRatio_;
198 
199 
200  // Private Member Functions
201 
202  // Calculate the delta values
203  void calcDelta();
204 
205  //- Fill changedPatchAndFaces (with patch and face labels) and
206  // changedFacesInfo (with delta).
207  // This is the initial set of faces from which to start the waves.
208  // Since there might be lots of places with delta jumps we can follow
209  // various strategies for this initial 'seed'.
210  // - start from single cell/face and let FvFaceCellWave pick up all
211  // others from there. might be quite a few waves before everything
212  // settles.
213  // - start from all faces. Lots of initial transfers.
214  // We do something in between:
215  // - start from all faces where there is a jump. Since we cannot easily
216  // determine this across coupled patches (cyclic, processor)
217  // introduce all faces of these and let FvFaceCellWave sort it out.
218  void setChangedFaces
219  (
220  const fvMesh& mesh,
221  const volScalarField& delta,
222  DynamicList<labelPair>& changedPatchAndFaces,
223  DynamicList<deltaData>& changedFacesInfo
224  );
225 
226 
227 public:
228 
229  //- Runtime type information
230  TypeName("smooth");
231 
232 
233  // Constructors
234 
235  //- Construct from name, momentumTransportModel and dictionary
237  (
238  const word& name,
240  const dictionary&
241  );
242 
243  //- Disallow default bitwise copy construction
244  smoothDelta(const smoothDelta&) = delete;
245 
246 
247  //- Destructor
248  virtual ~smoothDelta()
249  {}
250 
251 
252  // Member Functions
253 
254  //- Read the LESdelta dictionary
255  virtual void read(const dictionary&);
256 
257  // Correct values
258  virtual void correct();
259 
260 
261  // Member Operators
262 
263  //- Disallow default bitwise assignment
264  void operator=(const smoothDelta&) = delete;
265 };
266 
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 } // End namespace LESModels
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 //- Data associated with deltaData type are contiguous
275 template<>
276 inline bool contiguous<LESModels::smoothDelta::deltaData>()
277 {
278  return true;
279 }
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 } // End namespace Foam
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 #include "smoothDeltaDeltaDataI.H"
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 #endif
292 
293 // ************************************************************************* //
scalar delta
#define w2
Definition: blockCreate.C:32
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Public member class used by mesh-wave to propagate the delta-ratio.
Definition: smoothDelta.H:60
bool operator==(const deltaData &) const
friend Istream & operator>>(Istream &is, deltaData &wDist)
Definition: smoothDelta.H:184
bool equal(const deltaData &, TrackingData &td) const
Same (like operator==)
bool updateCell(const fvMesh &, const label thisCelli, const labelPair &neighbourPatchAndFacei, const deltaData &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
bool sameGeometry(const fvMesh &, const deltaData &, const scalar, TrackingData &td) const
Check for identical geometrical data.
friend Ostream & operator<<(Ostream &os, const deltaData &wDist)
Definition: smoothDelta.H:179
bool operator!=(const deltaData &) const
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
bool updateFace(const fvMesh &, const labelPair &thisPatchAndFacei, const label neighbourCelli, const deltaData &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
void transform(const fvPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Smoothed delta which takes a given simple geometric delta and applies smoothing to it such that the r...
Definition: smoothDelta.H:55
void operator=(const smoothDelta &)=delete
Disallow default bitwise assignment.
TypeName("smooth")
Runtime type information.
smoothDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
Definition: smoothDelta.C:147
virtual ~smoothDelta()
Destructor.
Definition: smoothDelta.H:247
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: smoothDelta.C:174
Abstract base class for LES deltas.
Definition: LESdelta.H:51
const momentumTransportModel & turbulence() const
Return momentumTransportModel reference.
Definition: LESdelta.H:124
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for turbulence models (RAS, LES and laminar).
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47