projectCurveEdge.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 
30 #include "projectCurveEdge.H"
31 #include "unitConversion.H"
33 #include "pointConstraint.H"
34 #include "OBJstream.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace blockEdges
43 {
44  defineTypeNameAndDebug(projectCurveEdge, 0);
45  addToRunTimeSelectionTable(blockEdge, projectCurveEdge, Istream);
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 Foam::blockEdges::projectCurveEdge::projectCurveEdge
53 (
54  const dictionary& dict,
55  const label index,
56  const searchableSurfaces& geometry,
57  const pointField& points,
58  Istream& is
59 )
60 :
61  blockEdge(dict, index, points, is),
62  geometry_(geometry)
63 {
64  wordList names(is);
65  surfaces_.resize(names.size());
66  forAll(names, i)
67  {
68  surfaces_[i] = geometry_.findSurfaceID(names[i]);
69 
70  if (surfaces_[i] == -1)
71  {
73  << "Cannot find surface " << names[i] << " in geometry"
74  << exit(FatalIOError);
75  }
76 
77  if (isA<searchableExtrudedCircle>(geometry_[surfaces_[i]]))
78  {
79  Info<< type() << " : Using curved surface "
80  << geometry_[surfaces_[i]].name()
81  << " to predict starting points." << endl;
82  }
83  }
84 }
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
91 {
93  return point::max;
94 }
95 
96 
99 {
100  // For debugging to tag the output
101  static label eIter = 0;
102 
103  autoPtr<OBJstream> debugStr;
104  if (debug)
105  {
106  debugStr.reset
107  (
108  new OBJstream("projectCurveEdge_" + Foam::name(eIter++) + ".obj")
109  );
110  Info<< "Writing lines from straight-line start points"
111  << " to projected points to " << debugStr().name() << endl;
112  }
113 
114 
115  auto tpoints = tmp<pointField>::New(lambdas.size());
116  auto& points = tpoints.ref();
117 
118  const point& startPt = points_[start_];
119  const point& endPt = points_[end_];
120  const vector d = endPt-startPt;
121 
122  // Initial guess
123  forAll(lambdas, i)
124  {
125  points[i] = startPt+lambdas[i]*d;
126  }
127 
128  // Use special interpolation to keep initial guess on same position on
129  // surface
130  forAll(surfaces_, i)
131  {
132  if (isA<searchableExtrudedCircle>(geometry_[surfaces_[i]]))
133  {
134  const searchableExtrudedCircle& s =
135  refCast<const searchableExtrudedCircle>
136  (
137  geometry_[surfaces_[i]]
138  );
139  List<pointIndexHit> nearInfo;
140  s.findParametricNearest
141  (
142  points[0],
143  points.last(),
144  scalarField(lambdas),
145  scalarField(points.size(), magSqr(d)),
146  nearInfo
147  );
148  forAll(nearInfo, i)
149  {
150  if (nearInfo[i].hit())
151  {
152  points[i] = nearInfo[i].hitPoint();
153  }
154  }
155 
156  break;
157  }
158  }
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(), 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.0-SMALL)
193  {
194  points.last() = endPt;
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 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
searchableExtrudedCircle.H
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
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
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::searchableExtrudedCircle
Searching on edgeMesh with constant radius.
Definition: searchableExtrudedCircle.H:92
unitConversion.H
Unit conversion functions.
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::projectCurveEdge::position
virtual point position(const scalar) const
The point position corresponding to the curve parameter.
Definition: projectCurveEdge.C:90
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::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::blockEdges::projectCurveEdge::length
virtual scalar length() const
The length of the curve.
Definition: projectCurveEdge.C:273
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
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
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::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::max
static const Vector< Cmpt > max
Definition: VectorSpace.H:117
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)
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::roots::type
type
Types of root.
Definition: Roots.H:54
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:401
projectCurveEdge.H
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