projectEdge.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) 2016 OpenFOAM Foundation
9 Copyright (C) 2021 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 "projectEdge.H"
30#include "pointConstraint.H"
32#include "OBJstream.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace blockEdges
41{
44}
45}
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::blockEdges::projectEdge::findNearest
50(
51 const point& pt,
52 point& near,
53 pointConstraint& constraint
54) const
55{
56 if (surfaces_.size())
57 {
58 const scalar distSqr = Foam::magSqr(lastPoint()-firstPoint());
59
60 pointField boundaryNear(1);
61 List<pointConstraint> boundaryConstraint(1);
62
64 (
65 geometry_,
66 surfaces_,
67 pointField(1, pt),
68 scalarField(1, distSqr),
69 boundaryNear,
70 boundaryConstraint
71 );
72 near = boundaryNear[0];
73 constraint = boundaryConstraint[0];
74 }
75 else
76 {
77 near = pt;
78 constraint = pointConstraint();
79 }
80}
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
86(
87 const dictionary& dict,
88 const label index,
89 const searchableSurfaces& geometry,
90 const pointField& points,
91 Istream& is
92)
93:
94 blockEdge(dict, index, points, is),
95 geometry_(geometry)
96{
97 wordList names(is);
98 surfaces_.resize(names.size());
99 forAll(names, i)
100 {
101 surfaces_[i] = geometry_.findSurfaceID(names[i]);
102
103 if (surfaces_[i] == -1)
104 {
106 << "Cannot find surface " << names[i] << " in geometry"
107 << exit(FatalIOError);
108 }
109 }
110}
111
112
113// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114
116{
117 // Initial guess
119
120 point near(start);
121
122 if (lambda >= SMALL && lambda < 1.0-SMALL)
123 {
124 pointConstraint constraint;
125 findNearest(start, near, constraint);
126 }
127
128 return near;
129}
130
131
134{
135 // For debugging to tag the output
136 static label eIter = 0;
137
138 autoPtr<OBJstream> debugStr;
139 if (debug)
140 {
141 debugStr.reset
142 (
143 new OBJstream("projectEdge_" + Foam::name(eIter++) + ".obj")
144 );
145 Info<< "Writing lines from straight-line start points"
146 << " to projected points to " << debugStr().name() << endl;
147 }
148
149
150 auto tpoints = tmp<pointField>::New(lambdas.size());
151 auto& points = tpoints.ref();
152
153 const scalar distSqr = Foam::magSqr(lastPoint()-firstPoint());
154
155 // Initial guess
156 forAll(lambdas, i)
157 {
158 points[i] = blockEdge::linearPosition(lambdas[i]);
159 }
160
161
162 // Upper limit for number of iterations
163 constexpr label maxIter = 10;
164
165 // Residual tolerance
166 constexpr scalar relTol = 0.1;
167 constexpr scalar absTol = 1e-4;
168
169 scalar initialResidual = 0.0;
170
171 for (label iter = 0; iter < maxIter; iter++)
172 {
173 // Do projection
174 {
175 List<pointConstraint> constraints(lambdas.size());
176 pointField start(points);
178 (
179 geometry_,
180 surfaces_,
181 start,
182 scalarField(start.size(), distSqr),
183 points,
184 constraints
185 );
186
187 // Reset start and end point
188 if (lambdas[0] < SMALL)
189 {
190 points[0] = firstPoint();
191 }
192 if (lambdas.last() > 1.0-SMALL)
193 {
194 points.last() = lastPoint();
195 }
196
197 if (debugStr)
198 {
199 forAll(points, i)
200 {
201 debugStr().write(linePointRef(start[i], points[i]));
202 }
203 }
204 }
205
206 // Calculate lambdas (normalised coordinate along edge)
207 scalarField projLambdas(points.size());
208 {
209 projLambdas[0] = 0.0;
210 for (label i = 1; i < points.size(); i++)
211 {
212 projLambdas[i] = projLambdas[i-1] + mag(points[i]-points[i-1]);
213 }
214 projLambdas /= projLambdas.last();
215 }
216 linearInterpolationWeights interpolator(projLambdas);
217
218 // Compare actual distances and move points (along straight line;
219 // not along surface)
220 vectorField residual(points.size(), Zero);
221 labelList indices;
222 scalarField weights;
223 for (label i = 1; i < points.size() - 1; i++)
224 {
225 interpolator.valueWeights(lambdas[i], indices, weights);
226
227 point predicted(Zero);
228 forAll(indices, indexi)
229 {
230 predicted += weights[indexi]*points[indices[indexi]];
231 }
232 residual[i] = predicted-points[i];
233 }
234
235 scalar scalarResidual = sum(mag(residual));
236
237 if (debug)
238 {
239 Pout<< "Iter:" << iter << " initialResidual:" << initialResidual
240 << " residual:" << scalarResidual << endl;
241 }
242
243 if (scalarResidual < absTol*0.5*lambdas.size())
244 {
245 break;
246 }
247 else if (iter == 0)
248 {
249 initialResidual = scalarResidual;
250 }
251 else if (scalarResidual/initialResidual < relTol)
252 {
253 break;
254 }
255
256
257 if (debugStr)
258 {
259 forAll(points, i)
260 {
261 const point predicted(points[i] + residual[i]);
262 debugStr().write(linePointRef(points[i], predicted));
263 }
264 }
265
266 points += residual;
267 }
268
269 return tpoints;
270}
271
272
274{
276 return 1;
277}
278
279
280// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
OFstream that keeps track of vertices.
Definition: OBJstream.H:61
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:216
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Define a curved edge that is parameterized for 0<lambda<1 between the start/end points.
Definition: blockEdge.H:64
point linearPosition(const scalar lambda) const
The point position in the straight line.
Definition: blockEdgeI.H:88
const point & lastPoint() const
The location of the last point.
Definition: blockEdgeI.H:55
const point & firstPoint() const
The location of the first point.
Definition: blockEdgeI.H:49
Defines the edge from the projection onto a surface (single surface) or intersection of two surfaces.
Definition: projectEdge.H:61
virtual scalar length() const
The length of the edge.
Definition: projectEdge.C:273
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool valueWeights(const scalar t, labelList &indices, scalarField &weights) const
Calculate weights and indices to calculate t from samples.
Accumulates point constraints through successive applications of the applyConstraint function.
static void findNearest(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &)
Find nearest. Return -1 (and a miss()) or surface and nearest.
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
label findSurfaceID(const word &name) const
Find index of surface. Return -1 if not found.
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
const pointField & points
Namespace for OpenFOAM.
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
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
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dictionary dict
volScalarField & e
Definition: createFields.H:11
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333