projectEdge.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "projectEdge.H"
29 #include "pointConstraint.H"
30 #include "OBJstream.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::projectEdge::findNearest
45 (
46  const point& pt,
47  point& near,
48  pointConstraint& constraint
49 ) const
50 {
51  if (surfaces_.size())
52  {
53  const scalar distSqr = magSqr(points_[end_] - points_[start_]);
54 
55  pointField boundaryNear(1);
56  List<pointConstraint> boundaryConstraint(1);
57 
59  (
60  geometry_,
61  surfaces_,
62  pointField(1, pt),
63  scalarField(1, distSqr),
64  boundaryNear,
65  boundaryConstraint
66  );
67  near = boundaryNear[0];
68  constraint = boundaryConstraint[0];
69  }
70  else
71  {
72  near = pt;
73  constraint = pointConstraint();
74  }
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  const dictionary& dict,
83  const label index,
84  const searchableSurfaces& geometry,
85  const pointField& points,
86  Istream& is
87 )
88 :
89  blockEdge(dict, index, points, is),
90  geometry_(geometry)
91 {
92  wordList names(is);
93  surfaces_.setSize(names.size());
94  forAll(names, i)
95  {
96  surfaces_[i] = geometry_.findSurfaceID(names[i]);
97 
98  if (surfaces_[i] == -1)
99  {
101  << "Cannot find surface " << names[i] << " in geometry"
102  << exit(FatalIOError);
103  }
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
111 {
112  // Initial guess
113  const point start
114  (
115  points_[start_] + lambda*(points_[end_] - points_[start_])
116  );
117 
118  point near(start);
119 
120  if (lambda >= small && lambda < 1 - small)
121  {
122  pointConstraint constraint;
123  findNearest(start, near, constraint);
124  }
125 
126  return near;
127 }
128 
129 
132 {
133  // For debugging to tag the output
134  static label eIter = 0;
135 
136  autoPtr<OBJstream> debugStr;
137  if (debug)
138  {
139  debugStr.reset
140  (
141  new OBJstream("projectEdge_" + Foam::name(eIter++) + ".obj")
142  );
143  Info<< "Writing lines from straight-line start points"
144  << " to projected points to " << debugStr().name() << endl;
145  }
146 
147 
148  tmp<pointField> tpoints(new pointField(lambdas.size()));
149  pointField& points = tpoints.ref();
150 
151  const point& startPt = points_[start_];
152  const point& endPt = points_[end_];
153  const vector d = endPt-startPt;
154 
155  // Initial guess
156  forAll(lambdas, i)
157  {
158  points[i] = startPt+lambdas[i]*d;
159  }
160 
161 
162  // Upper limit for number of iterations
163  const label maxIter = 10;
164 
165  // Residual tolerance
166  const scalar relTol = 0.1;
167  const scalar absTol = 1e-4;
168 
169  scalar initialResidual = 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(), magSqr(d)),
183  points,
184  constraints
185  );
186 
187  // Reset start and end point
188  if (lambdas[0] < small)
189  {
190  points[0] = startPt;
191  }
192  if (lambdas.last() > 1 - small)
193  {
194  points.last() = endPt;
195  }
196 
197  if (debugStr.valid())
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;
210  for (label i = 1; i < points.size(); i++)
211  {
212  projLambdas[i] =
213  projLambdas[i-1] + mag(points[i] - points[i-1]);
214  }
215  projLambdas /= projLambdas.last();
216  }
217  linearInterpolationWeights interpolator(projLambdas);
218 
219  // Compare actual distances and move points (along straight line;
220  // not along surface)
221  vectorField residual(points.size(), vector::zero);
222  labelList indices;
223  scalarField weights;
224  for (label i = 1; i < points.size() - 1; i++)
225  {
226  interpolator.valueWeights(lambdas[i], indices, weights);
227 
228  point predicted = vector::zero;
229  forAll(indices, indexi)
230  {
231  predicted += weights[indexi]*points[indices[indexi]];
232  }
233  residual[i] = predicted - points[i];
234  }
235 
236  const scalar scalarResidual = sum(mag(residual));
237 
238  if (debug)
239  {
240  Pout<< "Iter:" << iter << " initialResidual:" << initialResidual
241  << " residual:" << scalarResidual << endl;
242  }
243 
244  if (scalarResidual < absTol*0.5*lambdas.size())
245  {
246  break;
247  }
248  else if (iter == 0)
249  {
250  initialResidual = scalarResidual;
251  }
252  else if (scalarResidual/initialResidual < relTol)
253  {
254  break;
255  }
256 
257  if (debugStr.valid())
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 
273 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
OFstream which keeps track of vertices.
Definition: OBJstream.H:56
T & last()
Return the last element of the list.
Definition: UListI.H:128
static const Form zero
Definition: VectorSpace.H:118
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Define a curved edge that is parameterised for 0<lambda<1 between the start and end point.
Definition: blockEdge.H:57
const pointField & points_
Definition: blockEdge.H:62
const label end_
Definition: blockEdge.H:65
const label start_
Definition: blockEdge.H:64
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
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.
Defines the edge from the projection onto a surface (single surface) or intersection of two surfaces.
Definition: projectEdge.H:54
projectEdge(const dictionary &dict, const label index, const searchableSurfaces &geometry, const pointField &points, Istream &)
Construct from Istream setting pointsList.
Definition: projectEdge.C:81
virtual point position(const scalar) const
Return the point positions corresponding to the curve parameters.
Definition: projectEdge.C:110
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.
label findSurfaceID(const word &name) const
Find index of surface. Return -1 if not found.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const pointField & points
dimensionedScalar lambda(viscosity->lookup("lambda"))
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
IOerror FatalIOError
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict