Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42namespace blockEdges
43{
46}
47}
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52Foam::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"
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
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 scalar distSqr = Foam::magSqr(lastPoint()-firstPoint());
119
120 // Initial guess
121 forAll(lambdas, i)
122 {
123 points[i] = blockEdge::linearPosition(lambdas[i]);
124 }
125
126 // Use special interpolation to keep initial guess on same position on
127 // surface
128 forAll(surfaces_, i)
129 {
130 if (isA<searchableExtrudedCircle>(geometry_[surfaces_[i]]))
131 {
132 const searchableExtrudedCircle& s =
134 (
135 geometry_[surfaces_[i]]
136 );
137 List<pointIndexHit> nearInfo;
138 s.findParametricNearest
139 (
140 points[0],
141 points.last(),
142 scalarField(lambdas),
143 scalarField(points.size(), distSqr),
144 nearInfo
145 );
146 forAll(nearInfo, i)
147 {
148 if (nearInfo[i].hit())
149 {
150 points[i] = nearInfo[i].point();
151 }
152 }
153
154 break;
155 }
156 }
157
158
159
160 // Upper limit for number of iterations
161 constexpr label maxIter = 10;
162
163 // Residual tolerance
164 constexpr scalar relTol = 0.1;
165 constexpr scalar absTol = 1e-4;
166
167 scalar initialResidual = 0.0;
168
169 for (label iter = 0; iter < maxIter; iter++)
170 {
171 // Do projection
172 {
173 List<pointConstraint> constraints(lambdas.size());
174 pointField start(points);
176 (
177 geometry_,
178 surfaces_,
179 start,
180 scalarField(start.size(), distSqr),
181 points,
182 constraints
183 );
184
185 // Reset start and end point
186 if (lambdas[0] < SMALL)
187 {
188 points[0] = firstPoint();
189 }
190 if (lambdas.last() > 1.0-SMALL)
191 {
192 points.last() = lastPoint();
193 }
194
195 if (debugStr)
196 {
197 forAll(points, i)
198 {
199 debugStr().writeLine(start[i], points[i]);
200 }
201 }
202 }
203
204 // Calculate lambdas (normalised coordinate along edge)
205 scalarField projLambdas(points.size());
206 {
207 projLambdas[0] = 0.0;
208 for (label i = 1; i < points.size(); i++)
209 {
210 projLambdas[i] = projLambdas[i-1] + mag(points[i]-points[i-1]);
211 }
212 projLambdas /= projLambdas.last();
213 }
214 linearInterpolationWeights interpolator(projLambdas);
215
216 // Compare actual distances and move points (along straight line;
217 // not along surface)
218 vectorField residual(points.size(), Zero);
219 labelList indices;
220 scalarField weights;
221 for (label i = 1; i < points.size() - 1; i++)
222 {
223 interpolator.valueWeights(lambdas[i], indices, weights);
224
225 point predicted(Zero);
226 forAll(indices, indexi)
227 {
228 predicted += weights[indexi]*points[indices[indexi]];
229 }
230 residual[i] = predicted-points[i];
231 }
232
233 scalar scalarResidual = sum(mag(residual));
234
235 if (debug)
236 {
237 Pout<< "Iter:" << iter << " initialResidual:" << initialResidual
238 << " residual:" << scalarResidual << endl;
239 }
240
241 if (scalarResidual < absTol*0.5*lambdas.size())
242 {
243 break;
244 }
245 else if (iter == 0)
246 {
247 initialResidual = scalarResidual;
248 }
249 else if (scalarResidual/initialResidual < relTol)
250 {
251 break;
252 }
253
254
255 if (debugStr)
256 {
257 forAll(points, i)
258 {
259 const point predicted(points[i] + residual[i]);
260 debugStr().writeLine(points[i], predicted);
261 }
262 }
263
264 points += residual;
265 }
266
267 return tpoints;
268}
269
270
272{
274 return 1;
275}
276
277
278// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
T & last()
Access last element of the list, position [size()-1].
Definition UList.H:971
static const Form max
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
Define a curved edge that is parameterized for 0<lambda<1 between the start/end points.
Definition blockEdge.H:61
blockEdge(const pointField &points, const label from, const label to)
Construct from components.
Definition blockEdge.H:106
point linearPosition(const scalar lambda) const
The point position in the straight line.
Definition blockEdgeI.H:81
const point & lastPoint() const
The location of the last point.
Definition blockEdgeI.H:48
label start() const noexcept
Index of start (first) point.
Definition blockEdgeI.H:30
const point & firstPoint() const
The location of the first point.
Definition blockEdgeI.H:42
Defines the edge from the projection onto a surface (single surface) or intersection of two surfaces.
virtual scalar length() const
The length of the curve.
virtual point position(const scalar) const
The point position corresponding to the curve parameter.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual bool valueWeights(const scalar t, labelList &indices, scalarField &weights) const
Calculate weights and indices to calculate t from samples.
Searching on edgeMesh with constant radius.
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. The collection is specified as a dictionary. For example,...
label findSurfaceID(const word &name) const
Find index of surface. Return -1 if not found.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
auto & names
const pointField & points
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))
A namespace for various blockEdge types.
Definition arcEdge.C:31
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.