Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
surfaceOptimizer Class Reference
Collaboration diagram for surfaceOptimizer:
Collaboration graph
[legend]

Public Member Functions

 surfaceOptimizer (DynList< point > &pts, const DynList< triFace > &trias)
 Construct from transformed points and triangles forming a simplex. More...
 
 ~surfaceOptimizer ()
 
point optimizePoint (const scalar tol=0.1)
 optimizes position of a central point in the simplex More...
 

Private Member Functions

scalar evaluateStabilisationFactor () const
 evaluate stabilisation factor More...
 
scalar evaluateFunc (const scalar &K) const
 evaluate the functional More...
 
void evaluateGradients (const scalar &, vector &, tensor &) const
 evaluate gradients needed for optimisation More...
 
scalar optimiseDivideAndConquer (const scalar tol)
 optimise point position using the divide and conquer technique More...
 
scalar optimiseSteepestDescent (const scalar tol)
 optimise point position via the steepest descent method More...
 
 surfaceOptimizer (const surfaceOptimizer &)
 Disallow default bitwise copy construct. More...
 
void operator= (const surfaceOptimizer &)
 Disallow default bitwise assignment. More...
 

Private Attributes

DynList< point > & pts_
 reference to the simplex points More...
 
const DynList< triFace > & trias_
 reference to the triangles forming a simplex More...
 
point pMin_
 min position of the bnd box More...
 
point pMax_
 max position of the bnd box More...
 

Static Private Attributes

static const vector dirVecs [4]
 direction vectors for divide and conquer algorithm More...
 

Detailed Description

Definition at line 50 of file surfaceOptimizer.H.

Constructor & Destructor Documentation

◆ surfaceOptimizer() [1/2]

surfaceOptimizer ( const surfaceOptimizer )
private

Disallow default bitwise copy construct.

◆ surfaceOptimizer() [2/2]

surfaceOptimizer ( DynList< point > &  pts,
const DynList< triFace > &  trias 
)

Construct from transformed points and triangles forming a simplex.

Definition at line 352 of file surfaceOptimizer.C.

References forAll, Foam::max(), Foam::min(), and tf.

Here is the call graph for this function:

◆ ~surfaceOptimizer()

Definition at line 377 of file surfaceOptimizer.C.

Member Function Documentation

◆ evaluateStabilisationFactor()

scalar evaluateStabilisationFactor ( ) const
private

evaluate stabilisation factor

find the minimum area

K is greater than zero in case the stabilisation is needed

Definition at line 51 of file surfaceOptimizer.C.

References Foam::endl(), forAll, Foam::Info, Foam::magSqr(), Foam::max(), Foam::min(), surfaceOptimizer::pts_, surfaceOptimizer::trias_, Vector< Cmpt >::x(), and Vector< Cmpt >::y().

Referenced by surfaceOptimizer::optimiseDivideAndConquer(), and surfaceOptimizer::optimiseSteepestDescent().

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

◆ evaluateFunc()

scalar evaluateFunc ( const scalar &  K) const
private

evaluate the functional

Definition at line 97 of file surfaceOptimizer.C.

References Foam::endl(), forAll, Foam::func(), Foam::Info, Foam::magSqr(), Foam::max(), surfaceOptimizer::pts_, Foam::sqr(), Foam::sqrt(), surfaceOptimizer::trias_, Vector< Cmpt >::x(), and Vector< Cmpt >::y().

Referenced by surfaceOptimizer::optimiseDivideAndConquer(), and surfaceOptimizer::optimiseSteepestDescent().

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

◆ evaluateGradients()

void evaluateGradients ( const scalar &  K,
vector gradF,
tensor gradGradF 
) const
private

evaluate gradients needed for optimisation

calculate the gradient

calculate the second gradient

stabilise diagonal terms

Definition at line 135 of file surfaceOptimizer.C.

References forAll, Foam::mag(), Foam::magSqr(), Foam::max(), Foam::pow(), Foam::sqr(), Foam::sqrt(), Foam::twoSymm(), Vector< Cmpt >::x(), Tensor::xx(), Vector< Cmpt >::y(), Tensor::yy(), Tensor::zero, and Vector< scalar >::zero.

Referenced by surfaceOptimizer::optimiseSteepestDescent().

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

◆ optimiseDivideAndConquer()

scalar optimiseDivideAndConquer ( const scalar  tol)
private

optimise point position using the divide and conquer technique

find the value of the functional in the centre of the bnd box

set the centre with the minimum value as the centre for future search

halve the search range

calculate the tolerence

Definition at line 208 of file surfaceOptimizer.C.

References surfaceOptimizer::dirVecs, Foam::endl(), surfaceOptimizer::evaluateFunc(), surfaceOptimizer::evaluateStabilisationFactor(), Foam::func(), Foam::Info, Foam::mag(), surfaceOptimizer::pMax_, surfaceOptimizer::pMin_, surfaceOptimizer::pts_, surfaceOptimizer::trias_, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< scalar >::zero.

Referenced by surfaceOptimizer::optimizePoint().

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

◆ optimiseSteepestDescent()

scalar optimiseSteepestDescent ( const scalar  tol)
private

optimise point position via the steepest descent method

find the bounding box

find the minimum value on the 5 x 5 raster

start with steepest descent optimisation

store data into a matrix

calculate the determinant

Definition at line 271 of file surfaceOptimizer.C.

References Foam::det(), matrix2D::determinant(), Foam::endl(), surfaceOptimizer::evaluateFunc(), surfaceOptimizer::evaluateGradients(), surfaceOptimizer::evaluateStabilisationFactor(), Foam::Info, Foam::mag(), surfaceOptimizer::pMax_, surfaceOptimizer::pMin_, surfaceOptimizer::pts_, matrix2D::solveFirst(), matrix2D::solveSecond(), surfaceOptimizer::trias_, Vector< Cmpt >::x(), Tensor::xx(), Tensor::xy(), Vector< Cmpt >::y(), Tensor::yx(), Tensor::yy(), Vector< Cmpt >::z(), and Vector< scalar >::zero.

Referenced by surfaceOptimizer::optimizePoint().

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

◆ operator=()

void operator= ( const surfaceOptimizer )
private

Disallow default bitwise assignment.

◆ optimizePoint()

point optimizePoint ( const scalar  tol = 0.1)

optimizes position of a central point in the simplex

Definition at line 382 of file surfaceOptimizer.C.

References forAll, Foam::mag(), surfaceOptimizer::optimiseDivideAndConquer(), surfaceOptimizer::optimiseSteepestDescent(), surfaceOptimizer::pMax_, surfaceOptimizer::pMin_, surfaceOptimizer::pts_, and surfaceOptimizer::trias_.

Referenced by meshSurfaceOptimizer::newPositionSurfaceOptimizer(), and meshSurfaceOptimizer::nodeDisplacementSurfaceOptimizer().

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

Field Documentation

◆ dirVecs

const vector dirVecs
staticprivate
Initial value:
=
{
vector(-1.0, -1.0, 0.0),
vector(1.0, -1.0, 0.0),
vector(-1.0, 1.0, 0.0),
vector(1.0, 1.0, 0.0)
}

direction vectors for divide and conquer algorithm

Definition at line 54 of file surfaceOptimizer.H.

Referenced by surfaceOptimizer::optimiseDivideAndConquer().

◆ pts_

DynList<point>& pts_
private

◆ trias_

const DynList<triFace>& trias_
private

◆ pMin_

point pMin_
private

◆ pMax_

point pMax_
private

The documentation for this class was generated from the following files:
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49