moab
moab::Element::SpectralHex Class Reference

#include <ElemUtil.hpp>

Inheritance diagram for moab::Element::SpectralHex:
moab::Element::Map

List of all members.

Public Member Functions

 SpectralHex (const std::vector< CartVect > &vertices)
 SpectralHex (int order, double *x, double *y, double *z)
 SpectralHex (int order)
 SpectralHex ()
virtual ~SpectralHex ()
void set_gl_points (double *x, double *y, double *z)
virtual CartVect evaluate (const CartVect &xi) const
 Evaluate the map on (calculate $ x = F($ )$ )
virtual CartVect ievaluate (const CartVect &x) const
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 ()

Protected Attributes

real_xyz [3]

Static Protected Attributes

static int _n
static real_z [3]
static lagrange_data _ld [3]
static opt_data_3 _data
static real_odwork
static bool _init = false

Detailed Description

Definition at line 194 of file ElemUtil.hpp.


Constructor & Destructor Documentation

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

Definition at line 196 of file ElemUtil.hpp.

: Map(vertices){};
moab::Element::SpectralHex::SpectralHex ( int  order,
double *  x,
double *  y,
double *  z 
)

Definition at line 812 of file ElemUtil.cpp.

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

Definition at line 817 of file ElemUtil.cpp.

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

Definition at line 808 of file ElemUtil.cpp.

                           : Map(0)
  {
  }

Definition at line 822 of file ElemUtil.cpp.

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

Member Function Documentation

CartVect moab::Element::SpectralHex::evaluate ( const CartVect xi) const [virtual]

Evaluate the map on (calculate $ x = F($ )$ )

Implements moab::Element::Map.

Definition at line 57 of file SpectralHex.cpp.

{
    //piece that we shouldn't want to cache
  int d=0;
  for(d=0; d<3; d++){
    lagrange_0(&_ld[d], params[d]);
  }
  CartVect result;
  for (d=0; d<3; d++)
  {
    result[d] = tensor_i3(_ld[0].J,_ld[0].n,
                          _ld[1].J,_ld[1].n,
                          _ld[2].J,_ld[2].n,
                          _xyz[d],   // this is the "field"
                          _odwork);
  }
  return result;
}
double moab::Element::SpectralHex::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 929 of file ElemUtil.cpp.

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

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

Definition at line 41 of file SpectralHex.cpp.

{
  for(int d=0; d<3; d++){
    free(_z[d]);
    lagrange_free(&_ld[d]);
  }
  opt_free_3(&_data);
  free(_odwork);
}

Definition at line 886 of file ElemUtil.cpp.

  {
    CartVect result(0.);

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

    real r[3] = {0, 0, 0 }; // initial guess for parametric coords
    unsigned c = opt_no_constraints_3;
    real dist = opt_findpt_3(&_data, (const real **)_xyz, x_star, r, &c);
    // if it did not converge, get out with throw...
    if (dist > 0.9e+30)
      throw Map::EvaluationError();
    //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
    result = r;

    return  result;
  }
void moab::Element::SpectralHex::Init ( int  order)

Definition at line 828 of file ElemUtil.cpp.

  {
    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;
    //triplicates! n is the same in all directions !!!
    for(int d=0; d<3; d++){
      _z[d] = tmalloc(real, _n);
      lobatto_nodes(_z[d], _n);
      lagrange_setup(&_ld[d], _z[d], _n);
    }
    opt_alloc_3(&_data, _ld);

    unsigned int nf = _n*_n, ne = _n, nw = 2*_n*_n + 3*_n;
    _odwork = tmalloc(real, 6*nf + 9*ne + nw);
  }
bool moab::Element::SpectralHex::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 993 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) &&
           ( xi[2]>=-1.-tol) && (xi[2]<=1.+tol);
  }
double moab::Element::SpectralHex::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 944 of file ElemUtil.cpp.

  {
    // set the position of GL points
    // set the positions of GL nodes, before evaluations
    _data.elx[0]=_xyz[0];
    _data.elx[1]=_xyz[1];
    _data.elx[2]=_xyz[2];
    double xi[3];
    //triple loop; the most inner loop is in r direction, then s, then t
    double integral = 0.;
    //double volume = 0;
    int index=0; // used fr the inner loop
    for (int k=0; k<_n; k++ )
    {
      xi[2]=_ld[2].z[k];
      //double wk= _ld[2].w[k];
      for (int j=0; j<_n; j++)
      {
        xi[1]=_ld[1].z[j];
        //double wj= _ld[1].w[j];
        for (int i=0; i<_n; i++)
        {
          xi[0]=_ld[0].z[i];
          //double wi= _ld[0].w[i];
          opt_vol_set_intp_3(&_data,xi);
          double wk= _ld[2].J[k];
          double wj= _ld[1].J[j];
          double wi= _ld[0].J[i];
          Matrix3 J(0.);
          // it is organized differently
          J(0,0) = _data.jac[0]; // dx/dr
          J(0,1) = _data.jac[1]; // dx/ds
          J(0,2) = _data.jac[2]; // dx/dt
          J(1,0) = _data.jac[3]; // dy/dr
          J(1,1) = _data.jac[4]; // dy/ds
          J(1,2) = _data.jac[5]; // dy/dt
          J(2,0) = _data.jac[6]; // dz/dr
          J(2,1) = _data.jac[7]; // dz/ds
          J(2,2) = _data.jac[8]; // dz/dt
          double bm = wk*wj*wi* J.determinant();
          integral+= bm*field_vertex_values[index++];
          //volume +=bm;
        }
      }
    }
    //std::cout << "volume: " << volume << "\n";
    return integral;
  }
Matrix3 moab::Element::SpectralHex::jacobian ( const CartVect xi) const [virtual]

Evaluate the map's Jacobi matrix.

Implements moab::Element::Map.

Definition at line 98 of file SpectralHex.cpp.

{
  double x_i[3];
  params.get(x_i);
    // set the positions of GL nodes, before evaluations
  _data.elx[0]=_xyz[0];
  _data.elx[1]=_xyz[1];
  _data.elx[2]=_xyz[2];
  opt_vol_set_intp_3(&_data,x_i);
  Matrix3 J(0.);
    // it is organized differently
  J(0,0) = _data.jac[0]; // dx/dr
  J(0,1) = _data.jac[1]; // dx/ds
  J(0,2) = _data.jac[2]; // dx/dt
  J(1,0) = _data.jac[3]; // dy/dr
  J(1,1) = _data.jac[4]; // dy/ds
  J(1,2) = _data.jac[5]; // dy/dt
  J(2,0) = _data.jac[6]; // dz/dr
  J(2,1) = _data.jac[7]; // dz/ds
  J(2,2) = _data.jac[8]; // dz/dt
  return J;
}
void moab::Element::SpectralHex::set_gl_points ( double *  x,
double *  y,
double *  z 
)

Definition at line 51 of file SpectralHex.cpp.

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

Member Data Documentation

Definition at line 218 of file ElemUtil.hpp.

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

Definition at line 222 of file ElemUtil.hpp.

Definition at line 217 of file ElemUtil.hpp.

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

Definition at line 215 of file ElemUtil.hpp.

Definition at line 219 of file ElemUtil.hpp.

Definition at line 224 of file ElemUtil.hpp.

real * moab::Element::SpectralHex::_z [static, protected]

Definition at line 216 of file ElemUtil.hpp.


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