All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
twoDPointCorrector.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-2023 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::twoDPointCorrector
26 
27 Description
28  Class applies a two-dimensional correction to mesh motion point field.
29 
30  The correction guarantees that the mesh does not get twisted during motion
31  and thus introduce a third dimension into a 2-D problem.
32 
33  The operation is performed by looping through all edges approximately
34  normal to the plane and enforcing their orthogonality onto the plane by
35  adjusting points on their ends.
36 
37 SourceFiles
38  twoDPointCorrector.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef twoDPointCorrector_H
43 #define twoDPointCorrector_H
44 
45 #include "DemandDrivenMeshObject.H"
46 #include "pointField.H"
47 #include "labelList.H"
48 #include "vector.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 // Forward class declarations
56 class polyMesh;
57 
58 /*---------------------------------------------------------------------------*\
59  Class twoDPointCorrector Declaration
60 \*---------------------------------------------------------------------------*/
61 
63 :
65  <
66  polyMesh,
67  TopoChangeableMeshObject,
68  twoDPointCorrector
69  >
70 {
71  // Private Data
72 
73  //- Is 2D correction required, i.e. is the mesh
74  bool required_;
75 
76  //- 2-D plane unit normal
77  mutable vector* planeNormalPtr_;
78 
79  //- Indices of edges normal to plane
80  mutable labelList* normalEdgeIndicesPtr_;
81 
82  //- Flag to indicate a wedge geometry
83  mutable bool isWedge_;
84 
85  //- Wedge axis (if wedge geometry)
86  mutable vector wedgeAxis_;
87 
88  //- Wedge angle (if wedge geometry)
89  mutable scalar wedgeAngle_;
90 
91 
92  // Private Member Functions
93 
94  //- Calculate addressing
95  void calcAddressing() const;
96 
97  //- Clear addressing
98  void clearAddressing() const;
99 
100  //- Snap a point to the wedge patch(es)
101  void snapToWedge(const vector& n, const point& A, point& p) const;
102 
103 
104  // Static Data Members
105 
106  //- Edge orthogonality tolerance
107  static const scalar edgeOrthogonalityTol;
108 
109 
110 protected:
111 
112  friend class DemandDrivenMeshObject
113  <
114  polyMesh,
117  >;
118 
119  // Protected Constructors
120 
121  //- Construct from mesh
122  explicit twoDPointCorrector(const polyMesh& mesh);
123 
124 
125 public:
126 
127  // Declare name of the class and its debug switch
128  ClassName("twoDPointCorrector");
129 
130 
131  // Constructors
132 
133  //- Disallow default bitwise copy construction
134  twoDPointCorrector(const twoDPointCorrector&) = delete;
135 
136 
137  //- Destructor
139 
140 
141  // Member Functions
142 
143  //- Is 2D correction required, i.e. is the mesh a wedge or slab
144  bool required() const
145  {
146  return required_;
147  }
148 
149  //- Return plane normal
150  const vector& planeNormal() const;
151 
152  //- Return indices of normal edges.
153  const labelList& normalEdgeIndices() const;
154 
155  //- Return direction normal to plane
156  direction normalDir() const;
157 
158  //- Correct motion points
159  void correctPoints(pointField& p) const;
160 
161  //- Correct motion displacements
162  void correctDisplacement(const pointField& p, vectorField& disp) const;
163 
164  //- Update topology
165  virtual void topoChange(const polyTopoChangeMap&);
166 
167  //- Update from another mesh using the given map
168  virtual void mapMesh(const polyMeshMap&);
169 
170  //- Redistribute or update using the given distribution map
171  virtual void distribute(const polyDistributionMap&);
172 
173  //- Correct weighting factors for moving mesh.
174  virtual bool movePoints();
175 
176 
177  // Member Operators
178 
179  //- Disallow default bitwise assignment
180  void operator=(const twoDPointCorrector&) = delete;
181 };
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace Foam
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 #endif
191 
192 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
label n
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Class applies a two-dimensional correction to mesh motion point field.
virtual bool movePoints()
Correct weighting factors for moving mesh.
void operator=(const twoDPointCorrector &)=delete
Disallow default bitwise assignment.
virtual void topoChange(const polyTopoChangeMap &)
Update topology.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
void correctDisplacement(const pointField &p, vectorField &disp) const
Correct motion displacements.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
twoDPointCorrector(const polyMesh &mesh)
Construct from mesh.
direction normalDir() const
Return direction normal to plane.
bool required() const
Is 2D correction required, i.e. is the mesh a wedge or slab.
void correctPoints(pointField &p) const
Correct motion points.
ClassName("twoDPointCorrector")
const labelList & normalEdgeIndices() const
Return indices of normal edges.
const vector & planeNormal() const
Return plane normal.
Namespace for OpenFOAM.
uint8_t direction
Definition: direction.H:45
volScalarField & p