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 -------------------------------------------------------------------------------
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 "projectCurveEdge.H"
30 #include "unitConversion.H"
32 #include "pointConstraint.H"
33 #include "OBJstream.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(projectCurveEdge, 0);
42  addToRunTimeSelectionTable(blockEdge, projectCurveEdge, Istream);
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 Foam::projectCurveEdge::projectCurveEdge
49 (
50  const dictionary& dict,
51  const label index,
52  const searchableSurfaces& geometry,
53  const pointField& points,
54  Istream& is
55 )
56 :
57  blockEdge(dict, index, points, is),
58  geometry_(geometry)
59 {
60  wordList names(is);
61  surfaces_.setSize(names.size());
62  forAll(names, i)
63  {
64  surfaces_[i] = geometry_.findSurfaceID(names[i]);
65 
66  if (surfaces_[i] == -1)
67  {
69  << "Cannot find surface " << names[i] << " in geometry"
70  << exit(FatalIOError);
71  }
72 
73  if (isA<searchableExtrudedCircle>(geometry_[surfaces_[i]]))
74  {
75  Info<< type() << " : Using curved surface "
76  << geometry_[surfaces_[i]].name()
77  << " to predict starting points." << endl;
78  }
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
87 {
88  // For debugging to tag the output
89  static label eIter = 0;
90 
91  autoPtr<OBJstream> debugStr;
92  if (debug)
93  {
94  debugStr.reset
95  (
96  new OBJstream("projectCurveEdge_" + Foam::name(eIter++) + ".obj")
97  );
98  Info<< "Writing lines from straight-line start points"
99  << " to projected points to " << debugStr().name() << endl;
100  }
101 
102 
103  tmp<pointField> tpoints(new pointField(lambdas.size()));
104  pointField& points = tpoints.ref();
105 
106  const point& startPt = points_[start_];
107  const point& endPt = points_[end_];
108  const vector d = endPt-startPt;
109 
110  // Initial guess
111  forAll(lambdas, i)
112  {
113  points[i] = startPt+lambdas[i]*d;
114  }
115 
116  // Use special interpolation to keep initial guess on same position on
117  // surface
118  forAll(surfaces_, i)
119  {
120  if (isA<searchableExtrudedCircle>(geometry_[surfaces_[i]]))
121  {
122  const searchableExtrudedCircle& s =
123  refCast<const searchableExtrudedCircle>
124  (
125  geometry_[surfaces_[i]]
126  );
128  s.findParametricNearest
129  (
130  points[0],
131  points.last(),
132  scalarField(lambdas),
133  scalarField(points.size(), magSqr(d)),
134  nearInfo
135  );
136  forAll(nearInfo, i)
137  {
138  if (nearInfo[i].hit())
139  {
140  points[i] = nearInfo[i].hitPoint();
141  }
142  }
143 
144  break;
145  }
146  }
147 
148 
149 
150  // Upper limit for number of iterations
151  const label maxIter = 10;
152  // Residual tolerance
153  const scalar relTol = 0.1;
154  const scalar absTol = 1e-4;
155 
156  scalar initialResidual = 0.0;
157 
158  for (label iter = 0; iter < maxIter; iter++)
159  {
160  // Do projection
161  {
162  List<pointConstraint> constraints(lambdas.size());
165  (
166  geometry_,
167  surfaces_,
168  start,
169  scalarField(start.size(), magSqr(d)),
170  points,
171  constraints
172  );
173 
174  // Reset start and end point
175  if (lambdas[0] < SMALL)
176  {
177  points[0] = startPt;
178  }
179  if (lambdas.last() > 1.0-SMALL)
180  {
181  points.last() = endPt;
182  }
183 
184  if (debugStr.valid())
185  {
186  forAll(points, i)
187  {
188  debugStr().write(linePointRef(start[i], points[i]));
189  }
190  }
191  }
192 
193  // Calculate lambdas (normalised coordinate along edge)
194  scalarField projLambdas(points.size());
195  {
196  projLambdas[0] = 0.0;
197  for (label i = 1; i < points.size(); i++)
198  {
199  projLambdas[i] = projLambdas[i-1] + mag(points[i]-points[i-1]);
200  }
201  projLambdas /= projLambdas.last();
202  }
203  linearInterpolationWeights interpolator(projLambdas);
204 
205  // Compare actual distances and move points (along straight line;
206  // not along surface)
207  vectorField residual(points.size(), Zero);
208  labelList indices;
209  scalarField weights;
210  for (label i = 1; i < points.size() - 1; i++)
211  {
212  interpolator.valueWeights(lambdas[i], indices, weights);
213 
214  point predicted(Zero);
215  forAll(indices, indexi)
216  {
217  predicted += weights[indexi]*points[indices[indexi]];
218  }
219  residual[i] = predicted-points[i];
220  }
221 
222  scalar scalarResidual = sum(mag(residual));
223 
224  if (debug)
225  {
226  Pout<< "Iter:" << iter << " initialResidual:" << initialResidual
227  << " residual:" << scalarResidual << endl;
228  }
229 
230  if (scalarResidual < absTol*0.5*lambdas.size())
231  {
232  break;
233  }
234  else if (iter == 0)
235  {
236  initialResidual = scalarResidual;
237  }
238  else if (scalarResidual/initialResidual < relTol)
239  {
240  break;
241  }
242 
243 
244  if (debugStr.valid())
245  {
246  forAll(points, i)
247  {
248  const point predicted(points[i] + residual[i]);
249  debugStr().write(linePointRef(points[i], predicted));
250  }
251  }
252 
253  points += residual;
254  }
255 
256  return tpoints;
257 }
258 
259 
260 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
searchableExtrudedCircle.H
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::projectCurveEdge::position
virtual point position(const scalar) const
Return the point positions corresponding to the curve parameters.
Definition: projectCurveEdge.H:102
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::nearInfo
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
Definition: sampledTriSurfaceMesh.C:67
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:56
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::searchableExtrudedCircle
Searching on edgeMesh with constant radius.
Definition: searchableExtrudedCircle.H:87
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::blockEdge::start
label start() const
Return label of start point.
Definition: blockEdgeI.H:30
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
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)
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
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:375
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:57
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
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