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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "projectEdge.H"
30 #include "unitConversion.H"
32 #include "pointConstraint.H"
33 #include "OBJstream.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(projectEdge, 0);
41  addToRunTimeSelectionTable(blockEdge, projectEdge, Istream);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::projectEdge::findNearest
48 (
49  const point& pt,
50  point& near,
51  pointConstraint& constraint
52 ) const
53 {
54  if (surfaces_.size())
55  {
56  const scalar distSqr = magSqr(points_[end_]-points_[start_]);
57 
58  pointField boundaryNear(1);
59  List<pointConstraint> boundaryConstraint(1);
60 
62  (
63  geometry_,
64  surfaces_,
65  pointField(1, pt),
66  scalarField(1, distSqr),
67  boundaryNear,
68  boundaryConstraint
69  );
70  near = boundaryNear[0];
71  constraint = boundaryConstraint[0];
72  }
73  else
74  {
75  near = pt;
76  constraint = pointConstraint();
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
83 Foam::projectEdge::projectEdge
84 (
85  const dictionary& dict,
86  const label index,
87  const searchableSurfaces& geometry,
88  const pointField& points,
89  Istream& is
90 )
91 :
92  blockEdge(dict, index, points, is),
93  geometry_(geometry)
94 {
95  wordList names(is);
96  surfaces_.setSize(names.size());
97  forAll(names, i)
98  {
99  surfaces_[i] = geometry_.findSurfaceID(names[i]);
100 
101  if (surfaces_[i] == -1)
102  {
104  << "Cannot find surface " << names[i] << " in geometry"
105  << exit(FatalIOError);
106  }
107  }
108 }
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 {
115  // Initial guess
116  const point start(points_[start_] + lambda*(points_[end_]-points_[start_]));
117 
118  point near(start);
119 
120  if (lambda >= SMALL && lambda < 1.0-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  // Residual tolerance
165  const scalar relTol = 0.1;
166  const scalar absTol = 1e-4;
167 
168  scalar initialResidual = 0.0;
169 
170  for (label iter = 0; iter < maxIter; iter++)
171  {
172  // Do projection
173  {
174  List<pointConstraint> constraints(lambdas.size());
177  (
178  geometry_,
179  surfaces_,
180  start,
181  scalarField(start.size(), magSqr(d)),
182  points,
183  constraints
184  );
185 
186  // Reset start and end point
187  if (lambdas[0] < SMALL)
188  {
189  points[0] = startPt;
190  }
191  if (lambdas.last() > 1.0-SMALL)
192  {
193  points.last() = endPt;
194  }
195 
196  if (debugStr.valid())
197  {
198  forAll(points, i)
199  {
200  debugStr().write(linePointRef(start[i], points[i]));
201  }
202  }
203  }
204 
205  // Calculate lambdas (normalised coordinate along edge)
206  scalarField projLambdas(points.size());
207  {
208  projLambdas[0] = 0.0;
209  for (label i = 1; i < points.size(); i++)
210  {
211  projLambdas[i] = projLambdas[i-1] + mag(points[i]-points[i-1]);
212  }
213  projLambdas /= projLambdas.last();
214  }
215  linearInterpolationWeights interpolator(projLambdas);
216 
217  // Compare actual distances and move points (along straight line;
218  // not along surface)
219  vectorField residual(points.size(), Zero);
220  labelList indices;
221  scalarField weights;
222  for (label i = 1; i < points.size() - 1; i++)
223  {
224  interpolator.valueWeights(lambdas[i], indices, weights);
225 
226  point predicted(Zero);
227  forAll(indices, indexi)
228  {
229  predicted += weights[indexi]*points[indices[indexi]];
230  }
231  residual[i] = predicted-points[i];
232  }
233 
234  scalar scalarResidual = sum(mag(residual));
235 
236  if (debug)
237  {
238  Pout<< "Iter:" << iter << " initialResidual:" << initialResidual
239  << " residual:" << scalarResidual << endl;
240  }
241 
242  if (scalarResidual < absTol*0.5*lambdas.size())
243  {
244  break;
245  }
246  else if (iter == 0)
247  {
248  initialResidual = scalarResidual;
249  }
250  else if (scalarResidual/initialResidual < relTol)
251  {
252  break;
253  }
254 
255 
256  if (debugStr.valid())
257  {
258  forAll(points, i)
259  {
260  const point predicted(points[i] + residual[i]);
261  debugStr().write(linePointRef(points[i], predicted));
262  }
263  }
264 
265  points += residual;
266  }
267 
268  return tpoints;
269 }
270 
271 
272 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
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:158
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:56
Foam::projectEdge::position
virtual point position(const scalar) const
Return the point positions corresponding to the curve parameters.
Definition: projectEdge.C:113
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::pointConstraint
Accumulates point constraints through successive applications of the applyConstraint function.
Definition: pointConstraint.H:60
unitConversion.H
Unit conversion functions.
Foam::autoPtr::valid
bool valid() const noexcept
True if the managed pointer is non-null.
Definition: autoPtrI.H:107
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::blockEdge
Define a curved edge that is parameterized for 0<lambda<1 between the start and end point.
Definition: blockEdge.H:59
searchableSurfacesQueries.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
linearInterpolationWeights.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
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:62
Foam::Field< vector >
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_
Definition: blockEdge.H:65
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_
Definition: blockEdge.H:67
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_
Definition: blockEdge.H:68
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
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::start
label ListType::const_reference const label start
Definition: ListOps.H:408
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
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