offsetSurface.C
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) 2014 OpenFOAM Foundation
9 Copyright (C) 2015-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#include "offsetSurface.H"
31#include "triSurface.H"
32#include "triSurfaceSearch.H"
33#include "barycentric2D.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace extrudeModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
45
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52:
53 extrudeModel(typeName, dict),
54 project_(coeffDict_.getOrDefault("project", false))
55{
56 // Read surface
57 fileName baseName(coeffDict_.get<fileName>("baseSurface").expand());
58 baseSurfPtr_.reset(new triSurface(baseName));
59
60 // Read offsetted surface
61 fileName offsetName(coeffDict_.get<fileName>("offsetSurface").expand());
62 offsetSurfPtr_.reset(new triSurface(offsetName));
63
64
65 const triSurface& b = *baseSurfPtr_;
66 const triSurface& o = *offsetSurfPtr_;
67
68 if
69 (
70 b.size() != o.size()
71 || b.nPoints() != o.nPoints()
72 || b.nEdges() != o.nEdges()
73 )
74 {
76 << "offsetSurface:\n " << offsetName
77 << " has different topology than the baseSurface:\n "
78 << baseName << endl
80 }
81
82 // Construct search engine
83 baseSearchPtr_.reset(new triSurfaceSearch(b));
84
85 // Construct search engine
86 offsetSearchPtr_.reset(new triSurfaceSearch(o));
87}
88
89
90// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91
93{}
94
95
96// * * * * * * * * * * * * * * * * Operators * * * * * * * * * * * * * * * * //
97
99(
100 const point& surfacePoint,
101 const vector& surfaceNormal,
102 const label layer
103) const
104{
105 if (layer == 0)
106 {
107 return surfacePoint;
108 }
109 else
110 {
111 pointField samples(1, surfacePoint);
112 scalarField nearestDistSqr(1, GREAT);
114 baseSearchPtr_().findNearest(samples, nearestDistSqr, info);
115
116 label triI = info[0].index();
117
118
119 const triSurface& base = baseSurfPtr_();
120 const triPointRef baseTri(base[triI].tri(base.points()));
121 const barycentric2D bary = baseTri.pointToBarycentric(surfacePoint);
122
123 const triSurface& offset = offsetSurfPtr_();
124 const triPointRef offsetTri(offset[triI].tri(offset.points()));
125
126 const point offsetPoint
127 (
128 bary[0]*offsetTri.a()
129 + bary[1]*offsetTri.b()
130 + bary[2]*offsetTri.c()
131 );
132
133 point interpolatedPoint
134 (
135 surfacePoint + sumThickness(layer)*(offsetPoint-surfacePoint)
136 );
137
138
139 // Either return interpolatedPoint or re-project onto surface (since
140 // snapping might not have do so exactly)
141
142 if (project_)
143 {
144 // Re-project onto surface
145 offsetSearchPtr_().findNearest
146 (
147 pointField(1, interpolatedPoint),
148 scalarField(1, GREAT),
149 info
150 );
151 return info[0].hitPoint();
152 }
153 else
154 {
155 return interpolatedPoint;
156 }
157 }
158}
159
160
161// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162
163} // End namespace extrudeModels
164} // End namespace Foam
165
166// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant.
Definition: Barycentric2D.H:57
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
label nEdges() const
Number of edges in patch.
label nPoints() const
Number of points supporting patch faces.
const Field< point_type > & points() const noexcept
Return reference to global points.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Top level extrusion model class.
Definition: extrudeModel.H:77
const dictionary & coeffDict_
Definition: extrudeModel.H:86
Extrudes by interpolating points from one surface to the other. Surfaces have to be topologically ide...
Definition: offsetSurface.H:88
virtual ~offsetSurface()
Destructor.
Definition: offsetSurface.C:92
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
A class for handling file names.
Definition: fileName.H:76
string & expand(const bool allowEmpty=false)
Definition: string.C:173
Helper class to search on triSurface.
Triangulated surface description with patch information.
Definition: triSurface.H:79
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:80
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:272
const Point & a() const
Return first vertex.
Definition: triangleI.H:86
const Point & c() const
Return third vertex.
Definition: triangleI.H:98
const Point & b() const
Return second vertex.
Definition: triangleI.H:92
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Namespace for OpenFOAM.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
volScalarField & b
Definition: createFields.H:27
scalarField samples(nIntervals, Zero)