faPatch.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-2017 Wikki Ltd
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::faPatch
29 
30 Description
31  Finite area patch class. Used for 2-D non-Euclidian finite area method.
32 
33 Author
34  Zeljko Tukovic, FMENA
35  Hrvoje Jasak, Wikki Ltd.
36 
37 SourceFiles
38  faPatch.C
39  faPatchNew.C
40  faPatchTemplates.C
41  faPatchFaMeshTemplates.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef Foam_faPatch_H
46 #define Foam_faPatch_H
47 
48 #include "patchIdentifier.H"
49 #include "labelList.H"
50 #include "pointField.H"
51 #include "typeInfo.H"
52 #include "PtrList.H"
53 #include "faPatchFieldsFwd.H"
54 #include "autoPtr.H"
55 #include "runTimeSelectionTables.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 namespace Foam
60 {
61 
62 // Forward Declarations
63 class faBoundaryMesh;
64 class faPatch;
65 
66 //- Store lists of faPatch as a PtrList
68 
70 
71 /*---------------------------------------------------------------------------*\
72  Class faPatch Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class faPatch
76 :
77  public patchIdentifier,
78  public labelList
79 {
80  // Private Data
81 
82  //- Neighbour polyPatch index
83  const label nbrPolyPatchId_;
84 
85  //- Reference to boundary mesh
86  const faBoundaryMesh& boundaryMesh_;
87 
88  //- Demand-driven: edge-face addressing
89  mutable labelList::subList* edgeFacesPtr_;
90 
91  //- Demand-driven: local points labels
92  mutable labelList* pointLabelsPtr_;
93 
94  //- Demand-driven: point-edge addressing
95  mutable labelListList* pointEdgesPtr_;
96 
97 
98  // Private Member Functions
99 
100  //- No copy construct
101  faPatch(const faPatch&) = delete;
102 
103  //- No copy assignment
104  void operator=(const faPatch&) = delete;
105 
106  //- Clear out topological patch data
107  void clearOut();
108 
109 
110 protected:
111 
112  // Protected Member Functions
113 
114  //- The faPatch geometry initialisation is called by faBoundaryMesh
115  friend class faBoundaryMesh;
116 
117  //- Calculate patch point labels
118  void calcPointLabels() const;
119 
120  //- Calculate patch point-edge addressing
121  void calcPointEdges() const;
122 
123  //- Initialise the calculation of the patch geometry
124  virtual void initGeometry(PstreamBuffers&)
125  {}
126 
127  //- Calculate the patch geometry
128  virtual void calcGeometry(PstreamBuffers&)
129  {}
131  //- Initialise the patches for moving points
132  virtual void initMovePoints(PstreamBuffers&, const pointField&)
133  {}
134 
135  //- Correct patch after moving points
136  virtual void movePoints(PstreamBuffers&, const pointField&);
137 
138  //- Initialise the update of the patch topology
139  virtual void initUpdateMesh(PstreamBuffers&)
140  {}
141 
142  //- Update of the patch topology
143  virtual void updateMesh(PstreamBuffers&)
144  {}
146 
147 public:
148 
149  //- The boundary type associated with the patch
152 
153  //- Runtime type information
154  TypeName("patch");
155 
156  // Declare run-time constructor selection tables
159  (
160  autoPtr,
161  faPatch,
162  dictionary,
163  (
164  const word& name,
165  const dictionary& dict,
166  const label index,
167  const faBoundaryMesh& bm,
168  const word& patchType
169  ),
170  (name, dict, index, bm, patchType)
171  );
172 
173 
174  // Constructors
175 
176  //- Construct from components
177  faPatch
178  (
179  const word& name,
180  const labelUList& edgeLabels,
181  const label index,
182  const faBoundaryMesh& bm,
183  const label nbrPolyPatchi,
184  const word& patchType
185  );
186 
187  //- Construct from dictionary
188  faPatch
189  (
190  const word& name,
191  const dictionary& dict,
192  const label index,
193  const faBoundaryMesh& bm,
194  const word& patchType
195  );
196 
197  //- Copy construct, resetting the boundary mesh
198  faPatch(const faPatch& p, const faBoundaryMesh& bm);
199 
200  //- Copy construct, resetting boundary mesh and addressing
201  faPatch
202  (
203  const faPatch& p,
204  const faBoundaryMesh& bm,
205  const label index,
206  const labelUList& edgeLabels,
207  const label nbrPolyPatchi
208  );
209 
210 
211  //- Construct and return a clone, resetting the boundary mesh
212  virtual autoPtr<faPatch> clone(const faBoundaryMesh& bm) const
213  {
214  return autoPtr<faPatch>::New(*this, bm);
215  }
216 
217  //- Construct and return a clone, resetting the edge list
218  //- and boundary mesh
219  virtual autoPtr<faPatch> clone
220  (
221  const faBoundaryMesh& bm,
222  const labelUList& edgeLabels,
223  const label index,
224  const label nbrPolyPatchi
225  ) const
226  {
227  return autoPtr<faPatch>::New
228  (
229  *this,
230  bm,
231  index,
232  edgeLabels,
233  nbrPolyPatchi
234  );
235  }
236 
237 
238  // Selectors
239 
240  //- Return pointer to a new patch created on freestore from dictionary
241  static autoPtr<faPatch> New
242  (
243  const word& name,
244  const dictionary& dict,
245  const label index,
246  const faBoundaryMesh& bm
247  );
248 
249  //- Return pointer to a new patch created on freestore from dictionary
250  static autoPtr<faPatch> New
251  (
252  const word& patchType,
253  const word& name,
254  const dictionary& dict,
255  const label index,
256  const faBoundaryMesh& bm
257  );
258 
259 
260  //- Destructor
261  virtual ~faPatch();
262 
263 
264  // Static Member Functions
265 
266  //- Return true if the given type is a constraint type
267  static bool constraintType(const word& patchType);
268 
269  //- Return a list of all the constraint patch types
270  static wordList constraintTypes();
271 
272 
273  // Member Functions
274 
275  //- Return the list of edges
276  const labelList& edgeLabels() const noexcept
277  {
278  return static_cast<const labelList&>(*this);
279  }
280 
281  //- Number of patch points
282  label nPoints() const
283  {
284  return pointLabels().size();
285  }
286 
287  //- Number of edge labels (boundary edges) addressed by this patch
288  label nEdges() const noexcept
289  {
290  return labelList::size();
291  }
292 
293  //- The neighbour polyPatch index
294  label ngbPolyPatchIndex() const noexcept
295  {
296  return nbrPolyPatchId_;
297  }
298 
299  //- Return boundaryMesh reference
300  const faBoundaryMesh& boundaryMesh() const noexcept;
301 
302  //- Return true if this patch is coupled
303  virtual bool coupled() const
304  {
305  return false;
306  }
307 
308  //- Patch start in edge list
309  label start() const;
310 
311  //- Patch size is the number of edge labels, but can be overloaded
312  virtual label size() const
313  {
314  return labelList::size();
315  }
316 
317  //- Return label of edge in patch from global edge label
318  label whichEdge(const label edgei) const
319  {
320  return edgei - start();
321  }
322 
323  //- Slice List to patch, using the virtual patch size
324  template<class T>
325  typename List<T>::subList patchSlice(const List<T>& l) const
326  {
327  return typename List<T>::subList(l, size(), start());
328  }
329 
330  //- Slice List to patch, using the number of patch edges
331  template<class T>
332  typename List<T>::subList patchRawSlice(const List<T>& l) const
333  {
334  return typename List<T>::subList(l, nEdges(), start());
335  }
336 
337 
338  //- Write
339  virtual void write(Ostream&) const;
340 
341 
342  // Topological information
344  //- List of proc/face for the boundary edge neighbours
345  //- in locally reordered edge numbering.
347 
348  //- Boundary edge neighbour processors
349  //- (does not include own proc)
350  labelList boundaryProcs() const;
352  //- List of proc/size for the boundary edge neighbour processors
353  //- (does not include own proc)
355 
356 
357  // Access functions for geometrical data
358 
359  //- Return patch point labels
360  const labelList& pointLabels() const;
361 
362  //- Return patch point-edge addressing
363  const labelListList& pointEdges() const;
364 
365  //- Return normals of neighbour polyPatch faces
366  // Same as faMesh::haloFaceNormals()
368 
369  //- Return normals of neighbour polyPatch joined points
371 
372  //- Return edge-face addressing
373  const labelUList& edgeFaces() const;
374 
375  //- Return edge centres
376  const vectorField& edgeCentres() const;
377 
378  //- Return edge length vectors
379  const vectorField& edgeLengths() const;
380 
381  //- Return edge length magnitudes
382  const scalarField& magEdgeLengths() const;
383 
384  //- Return edge normals
386 
387  //- Return neighbour face centres
389 
390  //- Return cell-centre to face-centre vector
391  // except for coupled patches for which the cell-centre
392  // to coupled-cell-centre vector is returned
393  virtual tmp<vectorField> delta() const;
394 
395 
396  // Access functions for demand driven data
397 
398  //- Make patch weighting factors
399  virtual void makeWeights(scalarField&) const;
400 
401  //- Return patch weighting factors
402  const scalarField& weights() const;
403 
404  //- Make patch edge - neighbour face distances
405  virtual void makeDeltaCoeffs(scalarField&) const;
406 
407  void makeCorrectionVectors(vectorField&) const;
408 
409  //- Return patch edge - neighbour face distances
410  const scalarField& deltaCoeffs() const;
412 
413  // Topological changes
414 
415  //- Reset the list of edges (use with caution)
416  void resetEdges(const labelUList& newEdges);
417 
418  //- Reset the list of edges (use with caution)
419  void resetEdges(labelList&& newEdges);
420 
421 
422  // Evaluation
423 
424  //- Extract internal field next to patch using edgeFaces mapping
425  template<class Type>
426  inline void patchInternalField
427  (
428  const UList<Type>& f,
429  const labelUList& edgeFaces,
430  Field<Type>& pfld
431  ) const;
432 
433  //- Return given internal field next to patch as patch field
434  template<class Type>
436 
437  //- Return given internal field next to patch as patch field
438  //- using provided addressing
439  template<class Type>
441  (
442  const UList<Type>& f,
443  const labelUList& edgeFaces
444  ) const;
445 
446  //- Extract internal field next to patch as patch field
447  template<class Type>
448  void patchInternalField(const UList<Type>&, Field<Type>&) const;
449 
450  //- Return the patch field of the GeometricField
451  //- corresponding to this patch.
452  template<class GeometricField, class AnyType = bool>
453  const typename GeometricField::Patch& patchField
454  (
455  const GeometricField& gf
456  ) const;
457 
458  //- Lookup the named field from the local registry and
459  //- return the patch field corresponding to this patch.
460  // N.B. The dummy pointer arguments are used if this function is
461  // instantiated within a templated function to avoid a bug in gcc.
462  template<class GeometricField, class AnyType = bool>
464  (
465  const word& name,
466  const GeometricField* = nullptr,
467  const AnyType* = nullptr
468  ) const;
469 
470 
471  // Ostream Operator
472 
473  friend Ostream& operator<<(Ostream&, const faPatch&);
474 };
475 
476 
477 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
478 
479 } // End namespace Foam
480 
481 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
482 
483 #ifdef NoRepository
484  #include "faPatchTemplates.C"
485 #endif
486 
487 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
488 
489 #endif
490 
491 // ************************************************************************* //
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition: faPatch.H:59
static autoPtr< faPatch > New(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm)
Return pointer to a new patch created on freestore from dictionary.
Definition: faPatchNew.C:28
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
const GeometricField::Patch & patchField(const GeometricField &gf) const
Return the patch field of the GeometricField corresponding to this patch.
friend class faBoundaryMesh
The faPatch geometry initialisation is called by faBoundaryMesh.
Definition: faPatch.H:130
label nPoints() const
Number of patch points.
Definition: faPatch.H:343
const scalarField & weights() const
Return patch weighting factors.
Definition: faPatch.C:518
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:470
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
const scalarField & deltaCoeffs() const
Return patch edge - neighbour face distances.
Definition: faPatch.C:506
Identifies a patch by name and index, with optional physical type and group information.
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: faPatch.C:46
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition: faPatch.H:157
const labelList & pointLabels() const
Return patch point labels.
Definition: faPatch.C:266
virtual label size() const
Patch size is the number of edge labels, but can be overloaded.
Definition: faPatch.H:385
virtual void makeDeltaCoeffs(scalarField &) const
Make patch edge - neighbour face distances.
Definition: faPatch.C:490
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
const scalarField & magEdgeLengths() const
Return edge length magnitudes.
Definition: faPatch.C:448
SubList< T > subList
Declare type of subList.
Definition: List.H:122
friend Ostream & operator<<(Ostream &, const faPatch &)
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: faPatch.H:151
label nEdges() const noexcept
Number of edge labels (boundary edges) addressed by this patch.
Definition: faPatch.H:351
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: areaFieldsFwd.H:56
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const AnyType *=nullptr) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
const labelList & edgeLabels() const noexcept
Return the list of edges.
Definition: faPatch.H:335
label whichEdge(const label edgei) const
Return label of edge in patch from global edge label.
Definition: faPatch.H:393
friend class List< T >
Declare friendship with the List class.
Definition: UList.H:207
List< T >::subList patchSlice(const List< T > &l) const
Slice List to patch, using the virtual patch size.
Definition: faPatch.H:402
List< labelPair > boundaryProcSizes() const
List of proc/size for the boundary edge neighbour processors (does not include own proc) ...
Definition: faPatch.C:236
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: faPatch.C:55
void makeCorrectionVectors(vectorField &) const
Definition: faPatch.C:496
A List obtained as a section of another List.
Definition: SubList.H:50
virtual bool coupled() const
Return true if this patch is coupled.
Definition: faPatch.H:372
const vectorField & edgeCentres() const
Return edge centres.
Definition: faPatch.C:436
tmp< vectorField > edgeFaceCentres() const
Return neighbour face centres.
Definition: faPatch.C:464
void calcPointLabels() const
Calculate patch point labels.
Definition: faPatch.C:288
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void write(Ostream &) const
Write.
Definition: faPatch.C:542
List< T >::subList patchRawSlice(const List< T > &l) const
Slice List to patch, using the number of patch edges.
Definition: faPatch.H:411
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition: faPatch.C:422
tmp< vectorField > edgeNormals() const
Return edge normals.
Definition: faPatch.C:454
const word & name() const noexcept
The patch name.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
const direction noexcept
Definition: Scalar.H:258
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: faPatch.H:174
labelList boundaryProcs() const
Boundary edge neighbour processors (does not include own proc)
Definition: faPatch.C:216
List< labelPair > boundaryConnections() const
List of proc/face for the boundary edge neighbours in locally reordered edge numbering.
Definition: faPatch.C:196
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:72
labelList f(nPoints)
label start() const
Patch start in edge list.
Definition: faPatch.C:190
Buffers for inter-processor communications streams (UOPstream, UIPstream).
virtual ~faPatch()
Destructor.
Definition: faPatch.C:176
TypeName("patch")
Runtime type information.
void patchInternalField(const UList< Type > &f, const labelUList &edgeFaces, Field< Type > &pfld) const
Extract internal field next to patch using edgeFaces mapping.
tmp< vectorField > ngbPolyPatchPointNormals() const
Return normals of neighbour polyPatch joined points.
Definition: faPatch.C:387
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:76
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: faPatch.C:512
faBoundaryMesh BoundaryMesh
The boundary type associated with the patch.
Definition: faPatch.H:183
const vectorField & edgeLengths() const
Return edge length vectors.
Definition: faPatch.C:442
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Finite area boundary mesh.
autoPtr< List< label > > clone() const
Clone.
Definition: ListI.H:93
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition: faPatch.C:277
label size() const noexcept
The number of elements in the UList.
Definition: UListI.H:413
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Macros to ease declaration of run-time selection tables.
label index() const noexcept
The index of this patch in the boundaryMesh.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
declareRunTimeSelectionTable(autoPtr, faPatch, dictionary,(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm, const word &patchType),(name, dict, index, bm, patchType))
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
tmp< vectorField > ngbPolyPatchFaceNormals() const
Return normals of neighbour polyPatch faces.
Definition: faPatch.C:376
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
Definition: faPatch.C:524
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: faPatch.H:168
void calcPointEdges() const
Calculate patch point-edge addressing.
Definition: faPatch.C:345
void resetEdges(const labelUList &newEdges)
Reset the list of edges (use with caution)
Definition: faPatch.C:528
const faBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: faPatch.C:184
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: faPatch.H:145
Namespace for OpenFOAM.
label ngbPolyPatchIndex() const noexcept
The neighbour polyPatch index.
Definition: faPatch.H:359