MOAB: Mesh Oriented datABase  (version 5.4.1)
moab::Element::SpectralQuad Class Reference

#include <ElemUtil.hpp>

+ Inheritance diagram for moab::Element::SpectralQuad:
+ Collaboration diagram for moab::Element::SpectralQuad:

Public Member Functions

 SpectralQuad (const std::vector< CartVect > &vertices)
 SpectralQuad (int order, double *x, double *y, double *z)
 SpectralQuad (int order)
 SpectralQuad ()
virtual ~SpectralQuad ()
void set_gl_points (double *x, double *y, double *z)
virtual CartVect evaluate (const CartVect &xi) const
 Evaluate the map on \(x_i\) (calculate \(\vec x = F($\vec \xi)\) )
virtual CartVect ievaluate (const CartVect &x, double tol=1e-6, const CartVect &x0=CartVect(0.0)) const
 Evaluate the inverse map (calculate \(\vec \xi = F^-1($\vec x)\) to given tolerance)
virtual Matrix3 jacobian (const CartVect &xi) const
 Evaluate the map's Jacobi matrix.
double evaluate_scalar_field (const CartVect &xi, const double *field_vertex_values) const
 Evaluate a scalar field at a point given field values at the vertices.
double integrate_scalar_field (const double *field_vertex_values) const
 Integrate a scalar field over the element given field values at the vertices.
bool inside_nat_space (const CartVect &xi, double &tol) const
 decide if within the natural param space, with a tolerance
void Init (int order)
void freedata ()
void compute_gl_positions ()
void get_gl_points (double *&x, double *&y, double *&z, int &size)

Protected Attributes

realType_xyz [3]

Static Protected Attributes

static int _n
static realType_z [2]
static lagrange_data _ld [2]
static opt_data_2 _data
static realType_odwork
static bool _init = false
static realType_glpoints

Detailed Description

Definition at line 421 of file ElemUtil.hpp.


Constructor & Destructor Documentation

moab::Element::SpectralQuad::SpectralQuad ( const std::vector< CartVect > &  vertices) [inline]

Definition at line 424 of file ElemUtil.hpp.

References _xyz.

                                                              : Map( vertices )
        {
            _xyz[0] = _xyz[1] = _xyz[2] = NULL;
        };
moab::Element::SpectralQuad::SpectralQuad ( int  order,
double *  x,
double *  y,
double *  z 
)

Definition at line 18 of file SpectralQuad.cpp.

References _xyz, and Init().

                                                                       : Map( 0 )
{
    Init( order );
    _xyz[0] = x;
    _xyz[1] = y;
    _xyz[2] = z;
}

Definition at line 25 of file SpectralQuad.cpp.

References _xyz, and Init().

                                      : Map( 4 )
{
    Init( order );
    _xyz[0] = _xyz[1] = _xyz[2] = NULL;
}

Definition at line 16 of file SpectralQuad.cpp.

: Map( 0 ) {}

Definition at line 30 of file SpectralQuad.cpp.

References _init, and freedata().

{
    if( _init ) freedata();
    _init = false;
}

Member Function Documentation

Definition at line 173 of file SpectralQuad.cpp.

References _glpoints, _n, _xyz, and _z.

Referenced by test_spectral_quad().

{
    // will need to use shape functions on a simple linear quad to compute gl points
    // so we know the position of gl points in parametric space, we will just compute those
    // from the 3d vertex position (corner nodes of the quad), using simple mapping
    assert( this->vertex.size() == 4 );
    static double corner_params[4][2] = { { -1., -1. }, { 1., -1. }, { 1., 1. }, { -1., 1. } };
    // we will use the cached lobatto nodes in parametric space _z[d] (the same in both directions)
    int indexGL = 0;
    int n2      = _n * _n;
    for( int i = 0; i < _n; i++ )
    {
        double csi = _z[0][i];
        for( int j = 0; j < _n; j++ )
        {
            double eta = _z[1][j];  // we could really use the same _z[0] array of lobatto nodes
            CartVect pos( 0.0 );
            for( int k = 0; k < 4; k++ )
            {
                const double N_k = ( 1 + csi * corner_params[k][0] ) * ( 1 + eta * corner_params[k][1] );
                pos += N_k * vertex[k];
            }
            pos *= 0.25;                           // these are x, y, z of gl points; reorder them
            _glpoints[indexGL]          = pos[0];  // x
            _glpoints[indexGL + n2]     = pos[1];  // y
            _glpoints[indexGL + 2 * n2] = pos[2];  // z
            indexGL++;
        }
    }
    // now, we can set the _xyz pointers to internal memory allocated to these!
    _xyz[0] = &( _glpoints[0] );
    _xyz[1] = &( _glpoints[n2] );
    _xyz[2] = &( _glpoints[2 * n2] );
}
CartVect moab::Element::SpectralQuad::evaluate ( const CartVect xi) const [virtual]

Evaluate the map on \(x_i\) (calculate \(\vec x = F($\vec \xi)\) )

Implements moab::Element::Map.

Definition at line 1298 of file ElemUtil.cpp.

References _ld, _odwork, _xyz, lagrange_0(), and tensor_i2().

    {
        // piece that we shouldn't want to cache
        int d = 0;
        for( d = 0; d < 2; d++ )
        {
            lagrange_0( &_ld[d], xi[d] );
        }
        CartVect result;
        for( d = 0; d < 3; d++ )
        {
            result[d] = tensor_i2( _ld[0].J, _ld[0].n, _ld[1].J, _ld[1].n, _xyz[d], _odwork );
        }
        return result;
    }
double moab::Element::SpectralQuad::evaluate_scalar_field ( const CartVect xi,
const double *  field_vertex_values 
) const [virtual]

Evaluate a scalar field at a point given field values at the vertices.

Implements moab::Element::Map.

Definition at line 1343 of file ElemUtil.cpp.

References _ld, _odwork, lagrange_0(), and tensor_i2().

    {
        // piece that we shouldn't want to cache
        int d;
        for( d = 0; d < 2; d++ )
        {
            lagrange_0( &_ld[d], xi[d] );
        }

        double value = tensor_i2( _ld[0].J, _ld[0].n, _ld[1].J, _ld[1].n, field, _odwork );
        return value;
    }

Definition at line 60 of file SpectralQuad.cpp.

References _data, _glpoints, _ld, _odwork, _z, lagrange_free(), and opt_free_2().

Referenced by Init(), and ~SpectralQuad().

{
    for( int d = 0; d < 2; d++ )
    {
        free( _z[d] );
        lagrange_free( &_ld[d] );
    }
    opt_free_2( &_data );
    free( _odwork );
    free( _glpoints );
}
void moab::Element::SpectralQuad::get_gl_points ( double *&  x,
double *&  y,
double *&  z,
int &  size 
)

Definition at line 207 of file SpectralQuad.cpp.

References _n, and _xyz.

Referenced by test_spectral_quad().

{
    x    = (double*)_xyz[0];
    y    = (double*)_xyz[1];
    z    = (double*)_xyz[2];
    size = _n * _n;
}
CartVect moab::Element::SpectralQuad::ievaluate ( const CartVect x,
double  tol = 1e-6,
const CartVect x0 = CartVect( 0.0 ) 
) const [virtual]

Evaluate the inverse map (calculate \(\vec \xi = F^-1($\vec x)\) to given tolerance)

Reimplemented from moab::Element::Map.

Definition at line 1314 of file ElemUtil.cpp.

References _data, _xyz, moab::CartVect::get(), opt_findpt_2(), and opt_no_constraints_3.

    {
        // find nearest point
        realType x_star[3];
        xyz.get( x_star );

        realType r[2] = { 0, 0 };  // initial guess for parametric coords
        unsigned c    = opt_no_constraints_3;
        realType dist = opt_findpt_2( &_data, (const realType**)_xyz, x_star, r, &c );
        // if it did not converge, get out with throw...
        if( dist > 0.9e+30 )
        {
            std::vector< CartVect > dummy;
            throw Map::EvaluationError( xyz, dummy );
        }

        // c tells us if we landed inside the element or exactly on a face, edge, or node
        // also, dist shows the distance to the computed point.
        // copy parametric coords back
        return CartVect( r[0], r[1], 0. );
    }
void moab::Element::SpectralQuad::Init ( int  order)

Definition at line 35 of file SpectralQuad.cpp.

References _data, _glpoints, _init, _ld, _n, _odwork, _z, freedata(), lagrange_setup(), lobatto_nodes(), opt_alloc_2(), and tmalloc.

Referenced by SpectralQuad().

{
    if( _init && _n == order ) return;
    if( _init && _n != order )
    {
        // TODO: free data cached
        freedata();
    }
    // compute stuff that depends only on order
    _init = true;
    _n    = order;
    // duplicates! n is the same in all directions !!!
    for( int d = 0; d < 2; d++ )
    {
        _z[d] = tmalloc( double, _n );
        lobatto_nodes( _z[d], _n );
        lagrange_setup( &_ld[d], _z[d], _n );
    }
    opt_alloc_2( &_data, _ld );

    unsigned int nf = _n * _n, ne = _n, nw = 2 * _n * _n + 3 * _n;
    _odwork   = tmalloc( double, 6 * nf + 9 * ne + nw );
    _glpoints = tmalloc( double, 3 * nf );
}
bool moab::Element::SpectralQuad::inside_nat_space ( const CartVect xi,
double &  tol 
) const [virtual]

decide if within the natural param space, with a tolerance

Implements moab::Element::Map.

Definition at line 1360 of file ElemUtil.cpp.

    {
        // just look at the box+tol here
        return ( xi[0] >= -1. - tol ) && ( xi[0] <= 1. + tol ) && ( xi[1] >= -1. - tol ) && ( xi[1] <= 1. + tol );
    }
double moab::Element::SpectralQuad::integrate_scalar_field ( const double *  field_vertex_values) const [virtual]

Integrate a scalar field over the element given field values at the vertices.

Implements moab::Element::Map.

Definition at line 1355 of file ElemUtil.cpp.

    {
        return 0.;  // not implemented
    }
Matrix3 moab::Element::SpectralQuad::jacobian ( const CartVect xi) const [virtual]

Evaluate the map's Jacobi matrix.

Implements moab::Element::Map.

Definition at line 132 of file SpectralQuad.cpp.

{
    // not implemented
    Matrix3 J( 0. );
    return J;
}
void moab::Element::SpectralQuad::set_gl_points ( double *  x,
double *  y,
double *  z 
)

Definition at line 72 of file SpectralQuad.cpp.

References _xyz.

{
    _xyz[0] = x;
    _xyz[1] = y;
    _xyz[2] = z;
}

Member Data Documentation

Definition at line 461 of file ElemUtil.hpp.

Referenced by compute_gl_positions(), freedata(), and Init().

bool moab::Element::SpectralQuad::_init = false [static, protected]

Definition at line 460 of file ElemUtil.hpp.

Referenced by Init(), and ~SpectralQuad().

int moab::Element::SpectralQuad::_n [static, protected]

Definition at line 453 of file ElemUtil.hpp.

Referenced by compute_gl_positions(), get_gl_points(), and Init().

Definition at line 454 of file ElemUtil.hpp.

Referenced by compute_gl_positions(), freedata(), and Init().

List of all members.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines