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 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
28 
29 #include "offsetSurface.H"
31 #include "triSurface.H"
32 #include "triSurfaceSearch.H"
33 #include "barycentric2D.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace 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
79  << exit(FatalIOError);
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 
98 point offsetSurface::operator()
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);
113  List<pointIndexHit> info;
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 // ************************************************************************* //
label nPoints() const
Number of points supporting patch faces.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual ~offsetSurface()
Destructor.
Definition: offsetSurface.C:85
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
A class for handling file names.
Definition: fileName.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const Cmpt & b() const noexcept
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
scalarField samples(nIntervals, Zero)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const dictionary & coeffDict_
Definition: extrudeModel.H:82
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
defineTypeNameAndDebug(cyclicSector, 0)
Macros for easy insertion into run-time selection tables.
Top level extrusion model class.
Definition: extrudeModel.H:72
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant...
Definition: Barycentric2D.H:50
Extrudes by interpolating points from one surface to the other. Surfaces have to be topologically ide...
Definition: offsetSurface.H:79
Helper class to search on triSurface.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Field< point_type > & points() const noexcept
Return reference to global points.
label nEdges() const
Number of edges in patch.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
addToRunTimeSelectionTable(extrudeModel, cyclicSector, dictionary)
const Cmpt & c() const noexcept
string & expand(const bool allowEmpty=false)
Inplace expand initial tags, tildes, and all occurrences of environment variables as per stringOps::e...
Definition: string.C:166
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
offsetSurface(const dictionary &dict)
Construct from dictionary.
Definition: offsetSurface.C:44
Triangulated surface description with patch information.
Definition: triSurface.H:71
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...