Jive reference manual
List of all members | Public Member Functions | Static Public Attributes | Protected Member Functions
jive::geom::ParametricLine Class Reference

#include <jive/geom/ParametricLine.h>

Inheritance diagram for jive::geom::ParametricLine:
Inheritance graph

Public Member Functions

 ParametricLine ()
 
 ParametricLine (const String &name, const Matrix &ischeme, const Ref< StdShape > &xshape, const Ref< StdShape > &sshape=NIL)
 
 ParametricLine (const String &name, const Matrix &ischeme, const Boundary &boundary, const Ref< StdShape > &xshape, const Ref< StdShape > &sshape=NIL)
 
virtual void getBoundaryCoords (const Matrix &x, idx_t ibnd, const Matrix &c) const
 
virtual void getGlobalIntegrationPoints (const Matrix &x, const Matrix &c) const
 Computes the integration points in the global coordinate system. More...
 
virtual void getIntegrationWeights (const Vector &w, const Matrix &c) const
 Computes the integration weights in the global coordinate system. More...
 
virtual void getGlobalPoint (const Vector &x, const Vector &u, const Matrix &c) const
 Converts a local coordinate vector to a global coordinate vector. More...
 
virtual bool getLocalPoint (const Vector &u, const Vector &x, double eps, const Matrix &c) const
 Converts a global coordinate vector to a local coordinate vector. More...
 
virtual void evalShapeGradients (const Matrix &g, const Vector &u, const Matrix &c) const
 Computes the gradients of the shape functions in a specific point within this shape. More...
 
virtual void evalShapeGradGrads (const Matrix &g2, const Matrix &g, const Vector &u, const Matrix &c) const
 
virtual void calcGradients (const Cubix &g, const Vector *w, const Matrix &c, const PointSet &points) const
 
virtual void calcGradGrads (const Cubix &g2, const Cubix &g, const Vector *w, const Matrix &c, const PointSet &points) const
 
- Public Member Functions inherited from jive::geom::ParametricShape
 ParametricShape ()
 
 ParametricShape (const String &name, const Matrix &ischeme, const Ref< StdShape > &xshape, const Ref< StdShape > &sshape=NIL)
 
 ParametricShape (const String &name, const Matrix &ischeme, const Boundary &boundary, const Ref< StdShape > &xshape, const Ref< StdShape > &sshape=NIL)
 
virtual void readFrom (ObjectInput &in)
 
virtual void writeTo (ObjectOutput &out) const
 
virtual idx_t localRank () const
 Returns the local rank of this shape. More...
 
virtual idx_t nodeCount () const
 Returns the number of nodes. More...
 
virtual idx_t vertexCount () const
 
virtual idx_t ipointCount () const
 
virtual idx_t shapeFuncCount () const
 Returns the number of shape functions. More...
 
virtual idx_t boundaryCount () const
 Returns the number of boundaries. More...
 
virtual String getGeometry () const
 
virtual BoundaryShapegetBoundaryShape (idx_t ibnd) const
 
virtual IdxVector getBoundaryNodes (idx_t ibnd) const
 
virtual Topology getBoundaryTopology () const
 Returns the boundary topology of this internal shape. More...
 
virtual void mapBoundaryPoint (const Vector &v, idx_t ibnd, const Vector &u) const
 
Vector getLocalWeights () const
 
virtual Matrix getVertexCoords () const
 
virtual Matrix getIntegrationScheme () const
 
virtual void getIntegrationWeights (const Vector &w, const Matrix &c) const
 
virtual void getGlobalPoint (const Vector &x, const Vector &u, const Matrix &c) const
 
virtual bool getLocalPoint (const Vector &u, const Vector &x, double eps, const Matrix &c) const
 
virtual bool containsLocalPoint (const Vector &u) const
 
virtual Matrix getShapeFunctions () const
 Returns the shape functions evaluated in the integration points. More...
 
virtual Matrix getVertexFunctions () const
 
virtual void evalShapeFunctions (const Vector &h, const Vector &u) const
 
virtual void evalShapeGradients (const Matrix &g, const Vector &u, const Matrix &c) const
 
virtual void getShapeGradients (const Cubix &g, const Vector &w, const Matrix &c) const
 
virtual void getVertexGradients (const Cubix &g, const Matrix &c) const
 
virtual void getBoundaryGradients (const Cubix &g, idx_t ibn, const Matrix &c) const
 
virtual void evalShapeGradGrads (const Matrix &g2, const Matrix &g, const Vector &u, const Matrix &c) const
 
virtual void getShapeGradGrads (const Cubix &g2, const Cubix &g, const Vector &w, const Matrix &c) const
 
virtual void getVertexGradGrads (const Cubix &g2, const Cubix &g, const Matrix &c) const
 
virtual void getBoundaryGradGrads (const Cubix &g2, const Cubix &g, idx_t ibn, const Matrix &c) const
 
virtual void calcGradients (const Cubix &g, const Vector *w, const Matrix &c, const PointSet &points) const
 
virtual void calcGradGrads (const Cubix &g2, const Cubix &g, const Vector *w, const Matrix &c, const PointSet &points) const
 
virtual void * getExtByID (ExtensionID extID) const
 
StdShapegetXShape () const
 
StdShapegetSShape () const
 
- Public Member Functions inherited from jive::geom::InternalShape
virtual int globalRank () const
 Returns the global rank of this shape. More...
 
virtual BoundaryShapegetBoundaryShape (int ibnd) const =0
 Returns a boundary shape object given a boundary index. More...
 
virtual IntVector getBoundaryNodes (int ibnd) const =0
 Returns the indices of the nodes on one of the boundaries of this shape. More...
 
virtual void getBoundaryCoords (const Matrix &x, int ibnd, const Matrix &c) const
 Collects the coordinates of the nodes on a specific boundary of this shape. More...
 
virtual void getShapeGradients (const Cubix &g, const Vector &w, const Matrix &c) const
 Computes the gradients of the shape functions in the integration points of this shape. More...
 
virtual void getNodeGradients (const Cubix &g, const Matrix &c) const
 Computes the gradients of the shape functions in the nodes of this shape. More...
 
virtual void getBoundaryGradients (const Cubix &g, int ibnd, const Matrix &c) const =0
 Computes the gradients of the shape functions in the integration points on a boundary of this shape. More...
 
virtual jem::ClassgetClass () const
 Returns the Class instance representing the runtime class of this object. More...
 
- Public Member Functions inherited from jive::geom::Shape
virtual int integrationPointCount () const =0
 Returns the number of integration points. More...
 
virtual Matrix getLocalNodeCoords () const =0
 Returns the coordinates of the nodes in the local coordinate system. More...
 
virtual Matrix getIntegrationPoints () const =0
 Returns the integration points in the local coordinate system. More...
 
virtual void evalShapeFunctions (const Vector &h, const Vector &u) const =0
 Computes the shape functions in a given point. More...
 
virtual bool containsLocalPoint (const Vector &u) const =0
 Tests whether a point lies within this shape. More...
 
- Public Member Functions inherited from jem::Object
virtual String toString () const
 Returns a short textual description of this object. More...
 
virtual long hashValue () const
 Returns a hash value for this object. More...
 
virtual bool equals (const Ref< Object > &obj) const
 Tests whether two objects are equal. More...
 
Ref< Objectclone () const
 Returns a copy of this object. More...
 
- Public Member Functions inherited from jive::geom::Grad2Extension
virtual void evalShapeGradGrads (const Matrix &g2, const Matrix &g, const Vector &u, const Matrix &c) const =0
 
virtual void getShapeGradGrads (const Cubix &g2, const Cubix &g, const Vector &w, const Matrix &c) const =0
 
virtual void getVertexGradGrads (const Cubix &g2, const Cubix &g, const Matrix &c) const =0
 
virtual void getBoundaryGradGrads (const Cubix &g2, const Cubix &g, idx_t ibnd, const Matrix &c) const =0
 

Static Public Attributes

static const int RANK = 1
 
- Static Public Attributes inherited from jive::geom::Grad2Extension
static const jem::byte ID [1]
 

Protected Member Functions

virtual ~ParametricLine ()
 
- Protected Member Functions inherited from jive::geom::ParametricShape
virtual ~ParametricShape ()
 
void updateFuncs_ (const PointSet_ &points) const
 
void updateGrads_ (const PointSet_ &points) const
 
void updateGrads2_ (const PointSet_ &points) const
 
- Protected Member Functions inherited from jive::geom::InternalShape
virtual ~InternalShape ()
 
- Protected Member Functions inherited from jem::Collectable
 Collectable ()
 Creates an empty Collectable. More...
 
 ~Collectable ()
 Frees resources. More...
 
- Protected Member Functions inherited from jive::geom::Grad2Extension
virtual ~Grad2Extension ()
 
- Protected Member Functions inherited from jem::Interface
virtual ~Interface ()
 Empty destructor. More...
 
- Protected Member Functions inherited from jem::io::Serializable
virtual ~Serializable ()
 
virtual void emitVtableFunc_ ()
 

Additional Inherited Members

- Public Types inherited from jive::geom::ParametricShape
typedef ShapeBoundary Boundary
 
- Public Types inherited from jive::geom::InternalShape
typedef jem::numeric::SparseStructure Topology
 A type representing the boundary topology of an internal shape. More...
 
typedef util::IntVector IntVector
 A type representing a one-dimensional integer array. More...
 
typedef util::Cubix Cubix
 A three-dimensional array type. More...
 
- Public Types inherited from jive::geom::Shape
typedef util::Vector Vector
 A vector type. More...
 
typedef util::Matrix Matrix
 A matrix type. More...
 
- Static Public Member Functions inherited from jive::geom::InternalShape
static jem::ClassgetType ()
 Returns a pointer to a Class object representing this class. More...
 
- Static Public Member Functions inherited from jive::geom::Shape
static jem::ClassgetType ()
 Returns a pointer to a Class object representing this class. More...
 
- Static Public Member Functions inherited from jem::Object
static ClassgetType ()
 Returns the Class instance representing the Object class. More...
 
- Protected Attributes inherited from jive::geom::ParametricShape
Ref< StdShapexshape_
 
Ref< StdShapesshape_
 
Vector center_
 
Matrix ischeme_
 
Vector iweights_
 
Boundary boundary_
 
idx_t rank_
 
idx_t rank2_
 
idx_t ipCount_
 
idx_t nodeCount_
 
idx_t funcCount_
 
idx_t vertexCount_
 
PointSet_ ipoints_
 
jem::Array< PointSet_bpoints_
 
PointSet_ vertices_
 
Scratch_ scratch_
 
- Static Protected Attributes inherited from jive::geom::ParametricShape
static const int MAX_NR_ITER_
 
static const double NR_ALPHA_
 
static const double NR_BETA_
 
static const int U_NONE_
 
static const int U_FUNCS_
 
static const int U_GRADS_
 
static const int U_GRADS2_
 
static const int U_XFUNCS_
 
static const int U_XGRADS_
 
static const int U_XGRADS2_
 
static const int U_SFUNCS_
 
static const int U_SGRADS_
 
static const int U_SGRADS2_
 

Constructor & Destructor Documentation

jive::geom::ParametricLine::ParametricLine ( )
jive::geom::ParametricLine::ParametricLine ( const String name,
const Matrix ischeme,
const Ref< StdShape > &  xshape,
const Ref< StdShape > &  sshape = NIL 
)
jive::geom::ParametricLine::ParametricLine ( const String name,
const Matrix ischeme,
const Boundary boundary,
const Ref< StdShape > &  xshape,
const Ref< StdShape > &  sshape = NIL 
)
virtual jive::geom::ParametricLine::~ParametricLine ( )
protectedvirtual

Member Function Documentation

virtual void jive::geom::ParametricLine::getBoundaryCoords ( const Matrix x,
idx_t  ibnd,
const Matrix c 
) const
virtual

Reimplemented from jive::geom::ParametricShape.

virtual void jive::geom::ParametricLine::getGlobalIntegrationPoints ( const Matrix x,
const Matrix c 
) const
virtual

Fills the matrix x with the global coordinates of the integration points encapsulated by this shape: x(i,j) is set to the i-th coordinate of the j-th integration point. The matrix c should contain the global coordinates of the nodes attached to this shape: c(i,j) should be equal to the i-th coordinate of the j-th node.

The default implementation of this function is given by:

int i;
for ( i = 0; i < u.size(1); i++ )
{
this->getGlobalPoint ( x(slice(ALL),i), u(slice(ALL),i), c );
}
Parameters
x- a matrix that will be filled with the global coordinates of the integration points.
c- a matrix containing the global node coordinates.
Precondition
x.size(0) == this->globalRank() &&
x.size(1) == this->integrationPointCount() &&
c.size(0) == this->globalRank() &&
c.size(1) == this->nodeCount()

Reimplemented from jive::geom::ParametricShape.

virtual void jive::geom::ParametricLine::getIntegrationWeights ( const Vector w,
const Matrix c 
) const
virtual

Fills the vector w with the integration weights in the global coordinate system: w[i] is set to the weight of the i-th integration point. The matrix c should contain the global node coordinates: c(i,j) should be equal to the i-th coordinate of the j-th node.

The following code fragment demonstrates how the integration weights can be used to evaluate the integral of a function over the domain of a shape.

// Evaluates the integral of the function
//
// f(x,y,z,...) = x + y + z + ...
//
// over the domain of a given shape.
double evalIntegral ( const Shape& shape, const Matrix& c )
{
const int rank = shape.globalRank ();
const int ipCount = shape.integrationPointCount ();
Matrix x ( rank, ipCount );
Vector w ( ipCount );
double r;
double f;
int i;
shape.getGlobalIntegrationPoints ( x, c );
shape.getIntegrationWeights ( w, c );
r = 0.0;
for ( i = 0; i < ipCount; i++ )
{
f = sum( x(slice(ALL),i) );
r += w[i] * f;
}
return r;
}
Parameters
w- a vector that will be filled with the global integration weights.
c- a matrix containing the global node coordinates.
Precondition
w.size() == this->integrationPointCount() &&
c.size(0) == this->globalRank() &&
c.size(1) == this->nodeCount()

Implements jive::geom::Shape.

virtual void jive::geom::ParametricLine::getGlobalPoint ( const Vector x,
const Vector u,
const Matrix c 
) const
virtual

Fills the vector x with the global coordinates of a specific point within this shape: x[i] is set to the i-th global coordinate. The vector u should contain the local coordinates of that point: u[i] should be equal to the i-th local coordinate. The matrix c should contain the global node coordinates: c(i,j) should be equal to the i-th global coordinate of the j-th node.

Parameters
x- a vector that will be set to the global coordinates of the specified point.
u- a vector containing the local coordinates of that point.
c- a matrix containing the global node coordinates.
Precondition
x.size() == this->globalRank() &&
u.size() == this->localRank() &&
c.size(0) == this->globalRank() &&
c.size(1) == this->nodeCount() &&
this->containsLocalPoint ( u )

Implements jive::geom::Shape.

virtual bool jive::geom::ParametricLine::getLocalPoint ( const Vector u,
const Vector x,
double  eps,
const Matrix c 
) const
virtual

This function tries to convert a global coordinate vector to a local coordinate vector. It returns true if the conversion succeeds and the global coordinate vector points to a location `near' this shape. Otherwise it returns false.

If the return value equals true, the vector u will contain the local coordinate vector: u[i] will be set to the i-th local coordinate. The vector x should contain the global coordinate vector to be converted: x[i] should be equal to the i-th global coordinate. The double eps specifies the maximum allowed distance – in the global coordinate system – between this shape and the location pointed to by the vector x. The matrix c should contain the global node coordinates: c(i,j) should be equal to the i-th coordinate of the j-th node.

Note that a shape class does not have to support this function. In fact, the default implementation of this function simply throws a jem::IllegalOperationException.

Parameters
u- a vector that will be filled with the local coordinates of the point specified by the global coordinate vector x.
x- a global coordinate vector.
eps- the maximum allowed distance between the specified point and this shape.
c- a matrix containing the global node coordinates.
Precondition
u.size() == this->localRank() &&
x.size() == this->globalRank() &&
eps > 0.0 &&
c.size(0) == this->globalRank() &&
c.size(1) == this->nodeCount()
Returns
true if the conversion succeeds, and false otherwise.

Reimplemented from jive::geom::Shape.

virtual void jive::geom::ParametricLine::evalShapeGradients ( const Matrix g,
const Vector u,
const Matrix c 
) const
virtual

Fills the matrix g with the spatial derivatives of the shape functions in a specific point in this shape: g(i,j) is set to the derivative with respect of the i-th coordinate of the j-th shape function in that point. The vector u should contain the local coordinates of that point: u[i] should be equal to the i-th local coordinate. The matrix c should contain the global node coordinates of this shape; c(i,j) should be equal to the i-th coordinate of the j-th node.

Note that the derivatives are computed in the global coordinate system.

Parameters
g- a matrix that will be filled with the gradients of the shape functions.
u- a local coordinate vector.
c- a matrix containing the global node coordinates.
Precondition
g.size(0) == this->globalRank() &&
g.size(1) == this->nodeCount() &&
u.size() == this->localRank() &&
c.size(0) == this->globalRank() &&
c.size(1) == this->nodeCount()

Implements jive::geom::InternalShape.

virtual void jive::geom::ParametricLine::evalShapeGradGrads ( const Matrix g2,
const Matrix g,
const Vector u,
const Matrix c 
) const
virtual
virtual void jive::geom::ParametricLine::calcGradients ( const Cubix g,
const Vector w,
const Matrix c,
const PointSet points 
) const
virtual
virtual void jive::geom::ParametricLine::calcGradGrads ( const Cubix g2,
const Cubix g,
const Vector w,
const Matrix c,
const PointSet points 
) const
virtual

Member Data Documentation

const int jive::geom::ParametricLine::RANK = 1
static