blockEdge Class Referenceabstract

Define a curved edge that is parameterized for 0<lambda<1 between the start/end points. More...

Inheritance diagram for blockEdge:
[legend]
Collaboration diagram for blockEdge:
[legend]

Classes

class  iNew
 Class used for the read-construction of. More...
 

Public Member Functions

 TypeName ("blockEdge")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, blockEdge, Istream,(const dictionary &dict, const label index, const searchableSurfaces &geometry, const pointField &points, Istream &is),(dict, index, geometry, points, is))
 
 blockEdge (const pointField &points, const edge &fromTo)
 Construct from components. More...
 
 blockEdge (const dictionary &dict, const label index, const pointField &points, Istream &is)
 Construct from Istream and point field. More...
 
virtual autoPtr< blockEdgeclone () const
 Clone function. More...
 
virtual ~blockEdge ()=default
 Destructor. More...
 
bool valid () const noexcept
 True if first/last indices are unique and non-negative. More...
 
label start () const noexcept
 Index of start (first) point. More...
 
label end () const noexcept
 Index of end (last) point. More...
 
const pointfirstPoint () const
 The location of the first point. More...
 
const pointlastPoint () const
 The location of the last point. More...
 
int compare (const blockEdge &e) const
 Compare the given start/end points with this block edge. More...
 
int compare (const edge &e) const
 Compare the given start/end points with this block edge. More...
 
int compare (const label start, const label end) const
 Compare the given start/end points with this block edge. More...
 
point linearPosition (const scalar lambda) const
 The point position in the straight line. More...
 
virtual point position (const scalar lambda) const =0
 The point position corresponding to the curve parameter. More...
 
virtual tmp< pointFieldposition (const scalarList &lambdas) const
 The point positions corresponding to the curve parameters. More...
 
virtual scalar length () const =0
 The length of the curve. More...
 
void write (Ostream &os, const dictionary &dict) const
 Write edge with variable back-substitution. More...
 

Static Public Member Functions

static autoPtr< blockEdgeNew (const dictionary &dict, const label index, const searchableSurfaces &geometry, const pointField &points, Istream &is)
 New function which constructs and returns pointer to a blockEdge. More...
 

Protected Member Functions

 blockEdge (const pointField &points, const label from, const label to)
 Construct from components. More...
 

Static Protected Member Functions

static pointField appendEndPoints (const pointField &p, const label from, const label to, const pointField &intermediate)
 

Protected Attributes

const pointFieldpoints_
 The referenced point field. More...
 
const label start_
 Index of the first point. More...
 
const label end_
 Index of the last point. More...
 

Friends

Ostreamoperator<< (Ostream &os, const blockEdge &e)
 

Detailed Description

Define a curved edge that is parameterized for 0<lambda<1 between the start/end points.

Source files

Definition at line 63 of file blockEdge.H.

Constructor & Destructor Documentation

◆ blockEdge() [1/3]

blockEdge ( const pointField points,
const label  from,
const label  to 
)
inlineprotected

Construct from components.

Deprecated:
(2021-11) use constructor with edge
Parameters
pointsReferenced point field
fromStart point in point field
toEnd point in point field

Definition at line 96 of file blockEdge.H.

◆ blockEdge() [2/3]

blockEdge ( const pointField points,
const edge fromTo 
)

Construct from components.

Parameters
pointsReferenced point field
fromToStart/end in point field

Definition at line 44 of file blockEdge.C.

◆ blockEdge() [3/3]

blockEdge ( const dictionary dict,
const label  index,
const pointField points,
Istream is 
)

Construct from Istream and point field.

Parameters
pointsReferenced point field

Definition at line 56 of file blockEdge.C.

◆ ~blockEdge()

virtual ~blockEdge ( )
virtualdefault

Destructor.

Member Function Documentation

◆ appendEndPoints()

Foam::pointField appendEndPoints ( const pointField p,
const label  from,
const label  to,
const pointField intermediate 
)
staticprotected

Return a complete point field by appending the start/end points to the given list

Deprecated:
(2020-10) use polyLine::concat
Parameters
pReferenced point field
fromStart point in point field
toEnd point in point field
intermediateIntermediate points (knots)

Definition at line 109 of file blockEdge.C.

References polyLine::concat(), and p.

Here is the call graph for this function:

◆ TypeName()

TypeName ( "blockEdge"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
blockEdge  ,
Istream  ,
(const dictionary &dict, const label index, const searchableSurfaces &geometry, const pointField &points, Istream &is)  ,
(dict, index, geometry, points, is)   
)

◆ clone()

Foam::autoPtr< Foam::blockEdge > clone ( ) const
virtual

Clone function.

Definition at line 70 of file blockEdge.C.

References NotImplemented.

◆ New()

Foam::autoPtr< Foam::blockEdge > New ( const dictionary dict,
const label  index,
const searchableSurfaces geometry,
const pointField points,
Istream is 
)
static

New function which constructs and returns pointer to a blockEdge.

Definition at line 77 of file blockEdge.C.

References Foam::abort(), DebugInFunction, dict, Foam::endl(), Foam::FatalIOError, FatalIOErrorInLookup, and points.

Here is the call graph for this function:

◆ valid()

bool valid ( ) const
inlinenoexcept

True if first/last indices are unique and non-negative.

Definition at line 31 of file blockEdgeI.H.

References blockEdge::end_, and blockEdge::start_.

◆ start()

Foam::label start ( ) const
inlinenoexcept

Index of start (first) point.

Definition at line 37 of file blockEdgeI.H.

◆ end()

Foam::label end ( ) const
inlinenoexcept

Index of end (last) point.

Definition at line 43 of file blockEdgeI.H.

◆ firstPoint()

const Foam::point & firstPoint ( ) const
inline

The location of the first point.

Definition at line 49 of file blockEdgeI.H.

Referenced by arcEdge::arcEdge().

Here is the caller graph for this function:

◆ lastPoint()

const Foam::point & lastPoint ( ) const
inline

The location of the last point.

Definition at line 55 of file blockEdgeI.H.

Referenced by arcEdge::arcEdge().

Here is the caller graph for this function:

◆ compare() [1/3]

int compare ( const blockEdge e) const
inline

Compare the given start/end points with this block edge.

Return:

  • 0: different
  • +1: identical
  • -1: same edge, but different orientation

Definition at line 76 of file blockEdgeI.H.

References blockEdge::compare(), e, and UList< T >::end().

Referenced by blockEdge::compare().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compare() [2/3]

int compare ( const edge e) const
inline

Compare the given start/end points with this block edge.

Return:

  • 0: different
  • +1: identical
  • -1: same edge, but different orientation

Definition at line 82 of file blockEdgeI.H.

References blockEdge::compare(), e, and UList< T >::end().

Here is the call graph for this function:

◆ compare() [3/3]

int compare ( const label  start,
const label  end 
) const
inline

Compare the given start/end points with this block edge.

Return:

  • 0: different
  • +1: identical
  • -1: same edge, but different orientation

Definition at line 61 of file blockEdgeI.H.

◆ linearPosition()

Foam::point linearPosition ( const scalar  lambda) const
inline

The point position in the straight line.

0 <= lambda <= 1

Definition at line 88 of file blockEdgeI.H.

References InfoInFunction, lambda(), and Foam::nl.

Referenced by lineEdge::position(), projectEdge::position(), and projectCurveEdge::position().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ position() [1/2]

virtual point position ( const scalar  lambda) const
pure virtual

The point position corresponding to the curve parameter.

0 <= lambda <= 1

Implemented in arcEdge, bezier, polyLineEdge, BSplineEdge, lineEdge, projectCurveEdge, projectEdge, and splineEdge.

Referenced by lineDivide::lineDivide().

Here is the caller graph for this function:

◆ position() [2/2]

Foam::tmp< Foam::pointField > position ( const scalarList lambdas) const
virtual

The point positions corresponding to the curve parameters.

0 <= lambda <= 1

Reimplemented in projectCurveEdge, and projectEdge.

Definition at line 124 of file blockEdge.C.

References forAll, Time::New(), points, and UList< T >::size().

Here is the call graph for this function:

◆ length()

virtual scalar length ( ) const
pure virtual

The length of the curve.

Implemented in bezier, BSplineEdge, lineEdge, polyLineEdge, projectCurveEdge, projectEdge, splineEdge, and arcEdge.

◆ write()

void write ( Ostream os,
const dictionary dict 
) const

Write edge with variable back-substitution.

Definition at line 137 of file blockEdge.C.

References dict, Foam::endl(), os(), Foam::tab, and ObukhovLength::write().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator<<

Ostream & operator<< ( Ostream os,
const blockEdge e 
)
friend

Member Data Documentation

◆ points_

const pointField& points_
protected

The referenced point field.

Definition at line 70 of file blockEdge.H.

◆ start_

const label start_
protected

Index of the first point.

Definition at line 73 of file blockEdge.H.

Referenced by arcEdge::arcEdge(), and blockEdge::valid().

◆ end_

const label end_
protected

Index of the last point.

Definition at line 76 of file blockEdge.H.

Referenced by arcEdge::arcEdge(), and blockEdge::valid().


The documentation for this class was generated from the following files: