twoDPointCorrector.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) 2011-2015 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::twoDPointCorrector
28
29Description
30 Class applies a two-dimensional correction to mesh motion point field.
31
32 The correction guarantees that the mesh does not get twisted during motion
33 and thus introduce a third dimension into a 2-D problem.
34
35 The operation is performed by looping through all edges approximately
36 normal to the plane and enforcing their orthogonality onto the plane by
37 adjusting points on their ends.
38
39SourceFiles
40 twoDPointCorrector.C
41
42\*---------------------------------------------------------------------------*/
43
44#ifndef twoDPointCorrector_H
45#define twoDPointCorrector_H
46
47#include "MeshObject.H"
48#include "pointField.H"
49#include "labelList.H"
50#include "vector.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57// Forward class declarations
58class polyMesh;
59
60/*---------------------------------------------------------------------------*\
61 Class twoDPointCorrector Declaration
62\*---------------------------------------------------------------------------*/
65:
66 public MeshObject<polyMesh, UpdateableMeshObject, twoDPointCorrector>
67{
68 // Private data
69
70 //- Is 2D correction required, i.e. is the mesh
71 bool required_;
72
73 //- 2-D plane unit normal
74 mutable vector* planeNormalPtr_;
75
76 //- Indices of edges normal to plane
77 mutable labelList* normalEdgeIndicesPtr_;
78
79 //- Flag to indicate a wedge geometry
80 mutable bool isWedge_;
81
82 //- Wedge axis (if wedge geometry)
83 mutable vector wedgeAxis_;
84
85 //- Wedge angle (if wedge geometry)
86 mutable scalar wedgeAngle_;
87
88
89 // Private Member Functions
90
91 //- No copy construct
93
94 //- No copy assignment
95 void operator=(const twoDPointCorrector&) = delete;
96
97
98 //- Calculate addressing
99 void calcAddressing() const;
100
101 //- Clear addressing
102 void clearAddressing() const;
103
104 //- Snap a point to the wedge patch(es)
105 void snapToWedge(const vector& n, const point& A, point& p) const;
106
107
108 // Static data members
109
110 //- Edge orthogonality tolerance
111 static const scalar edgeOrthogonalityTol;
112
113
114public:
115
116 // Declare name of the class and its debug switch
117 ClassName("twoDPointCorrector");
118
119
120 // Constructors
121
122 //- Construct from components
124
125
126 //- Destructor
128
129
130 // Member Functions
131
132 //- Is 2D correction required, i.e. is the mesh a wedge or slab
133 bool required() const
134 {
135 return required_;
136 }
137
138 //- Return plane normal
139 const vector& planeNormal() const;
140
141 //- Return indices of normal edges.
142 const labelList& normalEdgeIndices() const;
143
144 //- Return direction normal to plane
145 direction normalDir() const;
146
147 //- Correct motion points
148 void correctPoints(pointField& p) const;
149
150 //- Correct motion displacements
151 void correctDisplacement(const pointField& p, vectorField& disp) const;
152
153 //- Update topology
154 void updateMesh(const mapPolyMesh&);
155
156 //- Correct weighting factors for moving mesh.
157 bool movePoints();
158};
159
160
161// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162
163} // End namespace Foam
164
165// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166
167#endif
168
169// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
label n
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:91
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Class applies a two-dimensional correction to mesh motion point field.
bool movePoints()
Correct weighting factors for moving mesh.
void correctDisplacement(const pointField &p, vectorField &disp) const
Correct motion displacements.
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")
void updateMesh(const mapPolyMesh &)
Update topology.
const labelList & normalEdgeIndices() const
Return indices of normal edges.
const vector & planeNormal() const
Return plane normal.
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition: className.H:67
volScalarField & p
Namespace for OpenFOAM.
uint8_t direction
Definition: direction.H:56