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) 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 "projectEdge.H"
30 #include "pointConstraint.H"
32 #include "OBJstream.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace blockEdges
41 {
42  defineTypeNameAndDebug(projectEdge, 0);
43  addToRunTimeSelectionTable(blockEdge, projectEdge, Istream);
44 }
45 }
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void 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 = magSqr(points_[end_]-points_[start_]);
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 
85 Foam::blockEdges::projectEdge::projectEdge
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_.setSize(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
118  const point start(points_[start_] + lambda*(points_[end_]-points_[start_]));
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 point& startPt = points_[start_];
154  const point& endPt = points_[end_];
155  const vector d = endPt-startPt;
156 
157  // Initial guess
158  forAll(lambdas, i)
159  {
160  points[i] = startPt+lambdas[i]*d;
161  }
162 
163 
164  // Upper limit for number of iterations
165  constexpr label maxIter = 10;
166 
167  // Residual tolerance
168  constexpr scalar relTol = 0.1;
169  constexpr scalar absTol = 1e-4;
170 
171  scalar initialResidual = 0.0;
172 
173  for (label iter = 0; iter < maxIter; iter++)
174  {
175  // Do projection
176  {
177  List<pointConstraint> constraints(lambdas.size());
178  pointField start(points);
180  (
181  geometry_,
182  surfaces_,
183  start,
184  scalarField(start.size(), magSqr(d)),
185  points,
186  constraints
187  );
188 
189  // Reset start and end point
190  if (lambdas[0] < SMALL)
191  {
192  points[0] = startPt;
193  }
194  if (lambdas.last() > 1.0-SMALL)
195  {
196  points.last() = endPt;
197  }
198 
199  if (debugStr)
200  {
201  forAll(points, i)
202  {
203  debugStr().write(linePointRef(start[i], points[i]));
204  }
205  }
206  }
207 
208  // Calculate lambdas (normalised coordinate along edge)
209  scalarField projLambdas(points.size());
210  {
211  projLambdas[0] = 0.0;
212  for (label i = 1; i < points.size(); i++)
213  {
214  projLambdas[i] = projLambdas[i-1] + mag(points[i]-points[i-1]);
215  }
216  projLambdas /= projLambdas.last();
217  }
218  linearInterpolationWeights interpolator(projLambdas);
219 
220  // Compare actual distances and move points (along straight line;
221  // not along surface)
222  vectorField residual(points.size(), Zero);
223  labelList indices;
224  scalarField weights;
225  for (label i = 1; i < points.size() - 1; i++)
226  {
227  interpolator.valueWeights(lambdas[i], indices, weights);
228 
229  point predicted(Zero);
230  forAll(indices, indexi)
231  {
232  predicted += weights[indexi]*points[indices[indexi]];
233  }
234  residual[i] = predicted-points[i];
235  }
236 
237  scalar scalarResidual = sum(mag(residual));
238 
239  if (debug)
240  {
241  Pout<< "Iter:" << iter << " initialResidual:" << initialResidual
242  << " residual:" << scalarResidual << endl;
243  }
244 
245  if (scalarResidual < absTol*0.5*lambdas.size())
246  {
247  break;
248  }
249  else if (iter == 0)
250  {
251  initialResidual = scalarResidual;
252  }
253  else if (scalarResidual/initialResidual < relTol)
254  {
255  break;
256  }
257 
258 
259  if (debugStr)
260  {
261  forAll(points, i)
262  {
263  const point predicted(points[i] + residual[i]);
264  debugStr().write(linePointRef(points[i], predicted));
265  }
266  }
267 
268  points += residual;
269  }
270 
271  return tpoints;
272 }
273 
274 
276 {
278  return 1;
279 }
280 
281 
282 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::autoPtr::reset
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:109
pointConstraint.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
projectEdge.H
Foam::OBJstream
OFstream that keeps track of vertices.
Definition: OBJstream.H:58
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::pointConstraint
Accumulates point constraints through successive applications of the applyConstraint function.
Definition: pointConstraint.H:60
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::blockEdge
Define a curved edge that is parameterized for 0<lambda<1 between the start/end points.
Definition: blockEdge.H:62
searchableSurfacesQueries.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
linearInterpolationWeights.H
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:445
Foam::blockEdges::defineTypeNameAndDebug
defineTypeNameAndDebug(arcEdge, 0)
Foam::Field< vector >
Foam::blockEdges::addToRunTimeSelectionTable
addToRunTimeSelectionTable(blockEdge, arcEdge, Istream)
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::blockEdge::points_
const pointField & points_
The referenced point field.
Definition: blockEdge.H:69
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::blockEdge::start_
const label start_
Index of the start point.
Definition: blockEdge.H:72
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::linearInterpolationWeights::valueWeights
virtual bool valueWeights(const scalar t, labelList &indices, scalarField &weights) const
Calculate weights and indices to calculate t from samples.
Definition: linearInterpolationWeights.C:89
Foam::linearInterpolationWeights
Definition: linearInterpolationWeights.H:50
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::blockEdge::end_
const label end_
Index of the end point.
Definition: blockEdge.H:75
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::linePointRef
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
Foam::Vector< scalar >
Foam::List< word >
Foam::searchableSurfaces
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Definition: searchableSurfaces.H:92
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::blockEdges::projectEdge::length
virtual scalar length() const
The length of the edge.
Definition: projectEdge.C:275
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::blockEdges::projectEdge::position
virtual point position(const scalar) const
The point position corresponding to the curve parameter.
Definition: projectEdge.C:115
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:401
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::searchableSurfacesQueries::findNearest
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.
Definition: searchableSurfacesQueries.C:349
OBJstream.H