ShapesEvaluated Class Reference

ShapesEvaluated.h DG++. More...

#include <ShapesEvaluated.h>

Inheritance diagram for ShapesEvaluated:
BasisFunctions ShapesEvaluated__< ShapeObj, QuadObj >

List of all members.

Public Member Functions

const std::vector< double > & getShapes () const
 Shape functions at quadrature points getShapes()[q*getBasisDimension()+a] gives the value of shape function a at quadrature point q.
const std::vector< double > & getDShapes () const
 Derivatives of shape functions at quadrature points getDShapes()[q*getBasisDimension()*getNumberOfDerivativesPerFunction()+ +a*getNumberOfDerivativesPerFunction()+i] gives the derivative in the i-th direction of degree of freedom a at quadrature point q.
const std::vector< double > & getIntegrationWeights () const
const std::vector< double > & getQuadraturePointCoordinates () const
 Coordinates of quadrature points in the real configuration getQuadraturePointCoordinates() q*ElementGeometrygetEmbeddingDimension()+i] returns the i-th coordinate in real space of quadrature point q.
size_t getBasisDimension () const
 returns the number of shape functions provided
size_t getNumberOfDerivativesPerFunction () const
 returns the number of directional derivative for each shape function
size_t getSpatialDimensions () const
 returns the number of number of coordinates for each Gauss point

Protected Member Functions

 ShapesEvaluated ()
virtual ~ShapesEvaluated ()
 ShapesEvaluated (const ShapesEvaluated &SEI)
virtual const ShapeaccessShape () const =0
 Returns the specific Shape objects in derived classes.
virtual const QuadratureaccessQuadrature () const =0
 Quadrature type Returns the specific Quadrature objects in derived classes.
void createObject (const ElementGeometry &eg)
 Since it is not possible to have a virtual constructor, one is emulated below only accessible from derived classes The virtual aspect are the calls to accessShape and accessQuadrature.

Private Attributes

std::vector< double > LocalShapes
std::vector< double > LocalDShapes
std::vector< double > LocalWeights
std::vector< double > LocalCoordinates

Detailed Description

ShapesEvaluated.h DG++.

Created by Adrian Lew on 9/7/06.

Copyright (c) 2006 Adrian Lew

Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE. ShapesEvaluated: Evaluation of the shape functions and derivatives at the quadrature points. Abstract class.

ShapesEvaluated is the class that takes the element geometry, a quadrature rule and shape functions on the reference element, and computes the values of the shape functions and their derivatives at the quadrature points. As opposed to Shapes, there will be one or more ShapesEvaluated objects per element in the mesh. Each different interpolation needed in an element will have a different ShapesEvaluated object.

This class provides all functionality for the templated derived classes ShapesEvaluated__. Only two abstract functions are left to be defined by derived classes, accessShape and accessQuadrature.

Objects in this class should only be constructed by derived classes.

ShapesEvaluated evaluates ElementGeometry::(D)map at the quadrature point (map)Coordinates using Quadrature::getQuadraturePoint().

ShapesEvaluated evaluates Shape functions at the quadrature point (Shape)Coordinates using Quadrature::getQuadraturePointShape().

Warning:
: The type of coordinates used to evaluate functions in Shape and ElementGeometry::(D)map (barycentric, Cartesian, etc.) should be consistent with those provided in Quadrature::getQuadraturePoint. In other words, if the Quadrature object returns Cartesian coordinates, the Shape and ElementGeometry objects should evaluate functions taking Cartesian coordinates as arguments.
Todo:
It would be nice to provide an interface of iterators to navigate the vectors, so that it is not necessary to remember in which order the data is stored in the vector. In the interest of simiplicity, this is for the moment skipped.
Todo:
It would useful to have the option of not computing the derivatives of the shape functions if not needed. Right now, it is not possible.
Todo:
Make coordinate types so that it is not necessary to check whether one is using the right convention between the three related classes.
Todo:
The computation of derivatives of shape functions with respect to the coordinates of the embedding space can only be performed if the ElementGeometry::map has domain and ranges of the same dimension. Otherwise, the derivatives should be computed with respect to some system of coordinates on the manifold. This will likely need to be revisited in the case of shells.
Todo:
Need to fix the Lapack interface... it is too particular to Mac here, and even not the best way to do it in Mac...

Constructor & Destructor Documentation

ShapesEvaluated::ShapesEvaluated (  )  [inline, protected]
virtual ShapesEvaluated::~ShapesEvaluated (  )  [inline, protected, virtual]
ShapesEvaluated::ShapesEvaluated ( const ShapesEvaluated SEI  )  [protected]

Member Function Documentation

virtual const Quadrature& ShapesEvaluated::accessQuadrature (  )  const [protected, pure virtual]

Quadrature type Returns the specific Quadrature objects in derived classes.

Implemented in ShapesEvaluated__< ShapeObj, QuadObj >.

virtual const Shape& ShapesEvaluated::accessShape (  )  const [protected, pure virtual]

Returns the specific Shape objects in derived classes.

Implemented in ShapesEvaluated__< ShapeObj, QuadObj >.

void ShapesEvaluated::createObject ( const ElementGeometry eg  )  [protected]

Since it is not possible to have a virtual constructor, one is emulated below only accessible from derived classes The virtual aspect are the calls to accessShape and accessQuadrature.

size_t ShapesEvaluated::getBasisDimension (  )  const [inline, virtual]

returns the number of shape functions provided

Implements BasisFunctions.

const std::vector<double>& ShapesEvaluated::getDShapes (  )  const [inline, virtual]

Derivatives of shape functions at quadrature points getDShapes()[q*getBasisDimension()*getNumberOfDerivativesPerFunction()+ +a*getNumberOfDerivativesPerFunction()+i] gives the derivative in the i-th direction of degree of freedom a at quadrature point q.

getDShapes returns an empty vector if no derivatives are available

Implements BasisFunctions.

const std::vector<double>& ShapesEvaluated::getIntegrationWeights (  )  const [inline, virtual]
Returns:
vector of integration weights Integration weights

Implements BasisFunctions.

size_t ShapesEvaluated::getNumberOfDerivativesPerFunction (  )  const [inline, virtual]

returns the number of directional derivative for each shape function

Implements BasisFunctions.

const std::vector<double>& ShapesEvaluated::getQuadraturePointCoordinates (  )  const [inline, virtual]

Coordinates of quadrature points in the real configuration getQuadraturePointCoordinates() q*ElementGeometrygetEmbeddingDimension()+i] returns the i-th coordinate in real space of quadrature point q.

Implements BasisFunctions.

const std::vector<double>& ShapesEvaluated::getShapes (  )  const [inline, virtual]

Shape functions at quadrature points getShapes()[q*getBasisDimension()+a] gives the value of shape function a at quadrature point q.

getShapes returns an empty vector if no shape functions are available

Implements BasisFunctions.

size_t ShapesEvaluated::getSpatialDimensions (  )  const [inline, virtual]

returns the number of number of coordinates for each Gauss point

Implements BasisFunctions.


Member Data Documentation

std::vector<double> ShapesEvaluated::LocalCoordinates [private]
std::vector<double> ShapesEvaluated::LocalDShapes [private]
std::vector<double> ShapesEvaluated::LocalShapes [private]
std::vector<double> ShapesEvaluated::LocalWeights [private]

The documentation for this class was generated from the following files:
Generated on Tue Aug 2 11:51:28 2011 for Galois by  doxygen 1.6.3