MOAB: Mesh Oriented datABase  (version 5.4.1)
MBMesquite::TargetCalculator Class Reference

#include <TargetCalculator.hpp>

+ Inheritance diagram for MBMesquite::TargetCalculator:

Public Member Functions

virtual ~TargetCalculator ()
virtual void initialize_queue (MeshDomainAssoc *mesh_and_domain, const Settings *settings, MsqError &err)
 Called at start of instruction queue processing.
virtual bool get_3D_target (PatchData &pd, size_t element, Sample sample, MsqMatrix< 3, 3 > &W_out, MsqError &err)=0
 Get a target matrix.
virtual bool get_surface_target (PatchData &pd, size_t element, Sample sample, MsqMatrix< 3, 2 > &W_out, MsqError &err)=0
 Get a target matrix.
virtual bool get_2D_target (PatchData &pd, size_t element, Sample sample, MsqMatrix< 2, 2 > &W_out, MsqError &err)=0
 Get a target matrix.
virtual bool have_surface_orient () const =0
 Use 3x2 W for surface elements if true, 2x2 W if false.

Static Public Member Functions

static bool factor_3D (const MsqMatrix< 3, 3 > &W, double &Lambda, MsqMatrix< 3, 3 > &V, MsqMatrix< 3, 3 > &Q, MsqMatrix< 3, 3 > &Delta, MsqError &err)
 Factor some existing target or Jacobian matrix.
static bool factor_surface (const MsqMatrix< 3, 2 > &W, double &Lambda, MsqMatrix< 3, 2 > &V, MsqMatrix< 2, 2 > &Q, MsqMatrix< 2, 2 > &Delta, MsqError &err)
 Factor some existing target or Jacobian matrix.
static bool factor_2D (const MsqMatrix< 2, 2 > &W, double &Lambda, MsqMatrix< 2, 2 > &V, MsqMatrix< 2, 2 > &Q, MsqMatrix< 2, 2 > &Delta, MsqError &err)
 Factor some existing target or Jacobian matrix.
static double size (const MsqMatrix< 3, 3 > &W)
 Get size component of W.
static double size (const MsqMatrix< 3, 2 > &W)
 Get size component of W.
static double size (const MsqMatrix< 2, 2 > &W)
 Get size component of W.
static MsqMatrix< 3, 3 > skew (const MsqMatrix< 3, 3 > &W)
 Get skew component of W.
static MsqMatrix< 2, 2 > skew (const MsqMatrix< 3, 2 > &W)
 Get skew component of W.
static MsqMatrix< 2, 2 > skew (const MsqMatrix< 2, 2 > &W)
 Get skew component of W.
static MsqMatrix< 3, 3 > aspect (const MsqMatrix< 3, 3 > &W)
 Get skew component of W.
static MsqMatrix< 2, 2 > aspect (const MsqMatrix< 3, 2 > &W)
 Get skew component of W.
static MsqMatrix< 2, 2 > aspect (const MsqMatrix< 2, 2 > &W)
 Get skew component of W.
static MsqMatrix< 3, 3 > shape (const MsqMatrix< 3, 3 > &W)
 Get shape (skew and aspect) component of W.
static MsqMatrix< 2, 2 > shape (const MsqMatrix< 3, 2 > &W)
 Get skew component of W.
static MsqMatrix< 2, 2 > shape (const MsqMatrix< 2, 2 > &W)
 Get skew component of W.
static MsqMatrix< 3, 3 > new_orientation_3D (const MsqVector< 3 > &b1, const MsqVector< 3 > &b2)
 Create a new orientation matrix.
static MsqMatrix< 3, 2 > new_orientation_2D (const MsqVector< 3 > &b1, const MsqVector< 3 > &b2)
 Create a new orientation matrix.
static void ideal_skew_3D (EntityTopology element_type, Sample s, const PatchData &pd, MsqMatrix< 3, 3 > &W, MsqError &err)
 Get skew matrix for an ideally shaped element.
static void ideal_skew_2D (EntityTopology element_type, Sample s, const PatchData &pd, MsqMatrix< 2, 2 > &W, MsqError &err)
 Get skew matrix for an ideally shaped element.
static void ideal_shape_3D (EntityTopology element_type, Sample s, const PatchData &pd, MsqMatrix< 3, 3 > &W, MsqError &err)
 Get skew matrix for an ideally shaped element.
static void ideal_shape_2D (EntityTopology element_type, Sample s, const PatchData &pd, MsqMatrix< 2, 2 > &W, MsqError &err)
 Get skew matrix for an ideally shaped element.
static MsqMatrix< 3, 3 > new_aspect_3D (const MsqVector< 3 > &r)
 Create a new aspect ratio matrix.
static MsqMatrix< 2, 2 > new_aspect_2D (const MsqVector< 2 > &r)
 Create a new aspect ratio matrix.
static MsqMatrix< 2, 2 > new_aspect_2D (double rho)
 Create a new aspect ratio matrix.
static void jacobian_3D (PatchData &pd, EntityTopology element_type, int num_nodes, Sample location, const Vector3D *coords, MsqMatrix< 3, 3 > &W_out, MsqError &err)
 Calculate the Jacobian given element vertex coordinates.
static void jacobian_2D (PatchData &pd, EntityTopology element_type, int num_nodes, Sample location, const Vector3D *coords, MsqMatrix< 3, 2 > &W_out, MsqError &err)
 Calculate the Jacobian given element vertex coordinates.
static void get_refmesh_Jacobian_3D (ReferenceMeshInterface *ref_mesh, PatchData &active_mesh, size_t element_no, Sample sample_no, MsqMatrix< 3, 3 > &W_out, MsqError &err)
static void get_refmesh_Jacobian_2D (ReferenceMeshInterface *ref_mesh, PatchData &active_mesh, size_t element_no, Sample sample_no, MsqMatrix< 3, 2 > &W_out, MsqError &err)

Detailed Description

Definition at line 51 of file TargetCalculator.hpp.


Constructor & Destructor Documentation

Definition at line 670 of file TargetCalculator.cpp.

{}

Member Function Documentation

MsqMatrix< 3, 3 > MBMesquite::TargetCalculator::aspect ( const MsqMatrix< 3, 3 > &  W) [static]

Get skew component of W.

Definition at line 117 of file TargetCalculator.cpp.

References MBMesquite::cbrt(), MBMesquite::MsqMatrix< R, C >::column(), and MBMesquite::length().

Referenced by TargetCalculatorTest::test_aspect_2D(), TargetCalculatorTest::test_aspect_3D(), and TargetCalculatorTest::test_aspect_surface().

{
    double a1    = length( W.column( 0 ) );
    double a2    = length( W.column( 1 ) );
    double a3    = length( W.column( 2 ) );
    double coeff = 1.0 / MBMesquite::cbrt( a1 * a2 * a3 );

    MsqMatrix< 3, 3 > result( coeff );
    result( 0, 0 ) *= a1;
    result( 1, 1 ) *= a2;
    result( 2, 2 ) *= a3;
    return result;
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::aspect ( const MsqMatrix< 3, 2 > &  W) [static]

Get skew component of W.

Definition at line 131 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column().

{
    double a1_sqr   = W.column( 0 ) % W.column( 0 );
    double a2_sqr   = W.column( 1 ) % W.column( 1 );
    double sqrt_rho = sqrt( sqrt( a1_sqr / a2_sqr ) );

    MsqMatrix< 2, 2 > result;
    result( 0, 0 ) = sqrt_rho;
    result( 0, 1 ) = 0.0;
    result( 1, 0 ) = 0.0;
    result( 1, 1 ) = 1 / sqrt_rho;
    return result;
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::aspect ( const MsqMatrix< 2, 2 > &  W) [static]

Get skew component of W.

Definition at line 145 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column().

{
    double a1_sqr   = W.column( 0 ) % W.column( 0 );
    double a2_sqr   = W.column( 1 ) % W.column( 1 );
    double sqrt_rho = sqrt( sqrt( a1_sqr / a2_sqr ) );

    MsqMatrix< 2, 2 > result;
    result( 0, 0 ) = sqrt_rho;
    result( 0, 1 ) = 0.0;
    result( 1, 0 ) = 0.0;
    result( 1, 1 ) = 1 / sqrt_rho;
    return result;
}
bool MBMesquite::TargetCalculator::factor_2D ( const MsqMatrix< 2, 2 > &  W,
double &  Lambda,
MsqMatrix< 2, 2 > &  V,
MsqMatrix< 2, 2 > &  Q,
MsqMatrix< 2, 2 > &  Delta,
MsqError err 
) [static]

Factor some existing target or Jacobian matrix.

Utility method for subclasses to use in their implementation of get_2D_target

Parameters:
WThe matrix to factor
LambdaOutput: the size factor of W V Output: the orientation factor of V
QOutput: the skew factor of W
DeltaOutput: the aspect ratio factor of W
Returns:
bool True if W can be factored, false otherwise.

Definition at line 305 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column(), MBMesquite::det(), moab::dot(), MBMesquite::length(), and MBMesquite::MsqMatrix< R, C >::set_column().

Referenced by MBMesquite::LVQDTargetCalculator::evaluate_guide_2D(), TargetCalculatorTest::test_aspect_2D(), TargetCalculatorTest::test_factor_2D(), TargetCalculatorTest::test_factor_2D_zero(), TargetCalculatorTest::test_shape_2D(), TargetCalculatorTest::test_size_2D(), and TargetCalculatorTest::test_skew_2D().

{
    double alpha = det( A );
    Lambda       = sqrt( fabs( alpha ) );
    if( Lambda < DBL_EPSILON ) return false;

    double la1_sqr = A.column( 0 ) % A.column( 0 );
    double la1     = sqrt( la1_sqr );
    double la2     = length( A.column( 1 ) );
    double inv_la1 = 1.0 / la1;
    double dot     = A.column( 0 ) % A.column( 1 );

    V.set_column( 0, A.column( 0 ) * inv_la1 );
    V.set_column( 1, ( la1_sqr * A.column( 1 ) - dot * A.column( 0 ) ) / ( la1 * alpha ) );

    double prod_rt2 = sqrt( la1 * la2 );
    Q( 0, 0 )       = prod_rt2 / Lambda;
    Q( 0, 1 )       = dot / ( prod_rt2 * Lambda );
    Q( 1, 0 )       = 0.0;
    Q( 1, 1 )       = 1.0 / Q( 0, 0 );

    double inv_prod_rt2 = 1.0 / prod_rt2;
    Delta( 0, 0 )       = la1 * inv_prod_rt2;
    Delta( 0, 1 )       = 0.0;
    Delta( 1, 0 )       = 0.0;
    Delta( 1, 1 )       = la2 * inv_prod_rt2;

    return true;
}
bool MBMesquite::TargetCalculator::factor_3D ( const MsqMatrix< 3, 3 > &  W,
double &  Lambda,
MsqMatrix< 3, 3 > &  V,
MsqMatrix< 3, 3 > &  Q,
MsqMatrix< 3, 3 > &  Delta,
MsqError err 
) [static]

Factor some existing target or Jacobian matrix.

Utility method for subclasses to use in their implementation of get_3D_target.

Parameters:
WThe matrix to factor
LambdaOutput: the size factor of W V Output: the orientation factor of V
QOutput: the skew factor of W
DeltaOutput: the aspect ratio factor of W
Returns:
bool True if W can be factored, false otherwise.

Definition at line 216 of file TargetCalculator.cpp.

References MBMesquite::cbrt(), MBMesquite::MsqMatrix< R, C >::column(), MBMesquite::length(), and MBMesquite::MsqMatrix< R, C >::set_column().

Referenced by MBMesquite::LVQDTargetCalculator::get_3D_target(), TargetCalculatorTest::test_aspect_3D(), TargetCalculatorTest::test_factor_3D(), TargetCalculatorTest::test_factor_3D_zero(), TargetCalculatorTest::test_shape_3D(), TargetCalculatorTest::test_size_3D(), and TargetCalculatorTest::test_skew_3D().

{
    MsqVector< 3 > a1xa2 = A.column( 0 ) * A.column( 1 );
    double alpha         = a1xa2 % A.column( 2 );
    Lambda               = MBMesquite::cbrt( fabs( alpha ) );
    if( Lambda < DBL_EPSILON ) return false;

    double la1_sqr = A.column( 0 ) % A.column( 0 );
    double la1     = sqrt( la1_sqr );
    double la2     = length( A.column( 1 ) );
    double la3     = length( A.column( 2 ) );
    double lx      = length( a1xa2 );
    double a1dota2 = A.column( 0 ) % A.column( 1 );

    double inv_la1 = 1.0 / la1;
    double inv_lx  = 1.0 / lx;
    V.set_column( 0, A.column( 0 ) * inv_la1 );
    V.set_column( 1, ( la1_sqr * A.column( 1 ) - a1dota2 * A.column( 0 ) ) * inv_la1 * inv_lx );
    V.set_column( 2, ( alpha / ( fabs( alpha ) * lx ) ) * a1xa2 );

    double inv_la2      = 1.0 / la2;
    double inv_la3      = 1.0 / la3;
    double len_prod_rt3 = MBMesquite::cbrt( la1 * la2 * la3 );
    Q( 0, 0 )           = 1.0;
    Q( 1, 0 ) = Q( 2, 0 ) = 0.0;
    Q( 0, 1 )             = a1dota2 * inv_la1 * inv_la2;
    Q( 1, 1 )             = lx * inv_la1 * inv_la2;
    Q( 2, 1 )             = 0.0;
    Q( 0, 2 )             = ( A.column( 0 ) % A.column( 2 ) ) * inv_la1 * inv_la3;
    Q( 1, 2 )             = ( a1xa2 % ( A.column( 0 ) * A.column( 2 ) ) ) * inv_lx * inv_la1 * inv_la3;
    Q( 2, 2 )             = fabs( alpha ) * inv_lx * inv_la3;
    Q *= len_prod_rt3 / Lambda;

    double inv_prod_rt3 = 1.0 / len_prod_rt3;
    ;
    Delta( 0, 0 ) = la1 * inv_prod_rt3;
    Delta( 0, 1 ) = 0.0;
    Delta( 0, 2 ) = 0.0;
    Delta( 1, 0 ) = 0.0;
    Delta( 1, 1 ) = la2 * inv_prod_rt3;
    Delta( 1, 2 ) = 0.0;
    Delta( 2, 0 ) = 0.0;
    Delta( 2, 1 ) = 0.0;
    Delta( 2, 2 ) = la3 * inv_prod_rt3;

    return true;
}
bool MBMesquite::TargetCalculator::factor_surface ( const MsqMatrix< 3, 2 > &  W,
double &  Lambda,
MsqMatrix< 3, 2 > &  V,
MsqMatrix< 2, 2 > &  Q,
MsqMatrix< 2, 2 > &  Delta,
MsqError err 
) [static]

Factor some existing target or Jacobian matrix.

Utility method for subclasses to use in their implementation of get_2D_target

Parameters:
WThe matrix to factor
LambdaOutput: the size factor of W V Output: the orientation factor of V
QOutput: the skew factor of W
DeltaOutput: the aspect ratio factor of W
Returns:
bool True if W can be factored, false otherwise.

Definition at line 269 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column(), moab::cross(), moab::dot(), MBMesquite::length(), and MBMesquite::MsqMatrix< R, C >::set_column().

Referenced by MBMesquite::LVQDTargetCalculator::evaluate_guide_2D(), MBMesquite::RefMeshTargetCalculator::get_2D_target(), TargetCalculatorTest::test_aspect_surface(), TargetCalculatorTest::test_factor_surface(), TargetCalculatorTest::test_factor_surface_zero(), TargetCalculatorTest::test_shape_surface(), TargetCalculatorTest::test_size_surface(), and TargetCalculatorTest::test_skew_surface().

{
    MsqVector< 3 > cross = A.column( 0 ) * A.column( 1 );
    double alpha         = length( cross );
    Lambda               = sqrt( alpha );
    if( Lambda < DBL_EPSILON ) return false;

    double la1_sqr = A.column( 0 ) % A.column( 0 );
    double la1     = sqrt( la1_sqr );
    double la2     = length( A.column( 1 ) );
    double inv_la1 = 1.0 / la1;
    double dot     = A.column( 0 ) % A.column( 1 );

    V.set_column( 0, A.column( 0 ) * inv_la1 );
    V.set_column( 1, ( la1_sqr * A.column( 1 ) - dot * A.column( 0 ) ) / ( la1 * alpha ) );

    double prod_rt2 = sqrt( la1 * la2 );
    Q( 0, 0 )       = prod_rt2 / Lambda;
    Q( 0, 1 )       = dot / ( prod_rt2 * Lambda );
    Q( 1, 0 )       = 0.0;
    Q( 1, 1 )       = 1.0 / Q( 0, 0 );

    double inv_prod_rt2 = 1.0 / prod_rt2;
    Delta( 0, 0 )       = la1 * inv_prod_rt2;
    Delta( 0, 1 )       = 0.0;
    Delta( 1, 0 )       = 0.0;
    Delta( 1, 1 )       = la2 * inv_prod_rt2;

    return true;
}
virtual bool MBMesquite::TargetCalculator::get_3D_target ( PatchData pd,
size_t  element,
Sample  sample,
MsqMatrix< 3, 3 > &  W_out,
MsqError err 
) [pure virtual]
void MBMesquite::TargetCalculator::get_refmesh_Jacobian_2D ( ReferenceMeshInterface ref_mesh,
PatchData active_mesh,
size_t  element_no,
Sample  sample_no,
MsqMatrix< 3, 2 > &  W_out,
MsqError err 
) [static]

Definition at line 640 of file TargetCalculator.cpp.

References MBMesquite::PatchData::element_by_index(), MBMesquite::MsqMeshEntity::get_element_type(), MBMesquite::ReferenceMeshInterface::get_reference_vertex_coordinates(), MBMesquite::PatchData::get_vertex_handles_array(), MBMesquite::MsqMeshEntity::get_vertex_index_array(), jacobian_2D(), moab::MAX_NODES, MSQ_ERRRTN, n, and MBMesquite::MsqMeshEntity::node_count().

Referenced by MBMesquite::RefMeshTargetCalculator::get_surface_target(), and TargetCalculatorTest::test_get_refmesh_Jacobian_2D().

{
    // get element
    MsqMeshEntity& elem       = pd.element_by_index( element );
    const EntityTopology type = elem.get_element_type();
    const unsigned n          = elem.node_count();

    const unsigned MAX_NODES = 9;
    assert( n <= MAX_NODES );

    // get vertices
    Mesh::VertexHandle elem_verts[MAX_NODES];
    const std::size_t* vtx_idx        = elem.get_vertex_index_array();
    const Mesh::VertexHandle* vtx_hdl = pd.get_vertex_handles_array();
    for( unsigned i = 0; i < n; ++i )
        elem_verts[i] = vtx_hdl[vtx_idx[i]];

    // get vertex coordinates
    Vector3D vert_coords[MAX_NODES];
    ref_mesh->get_reference_vertex_coordinates( elem_verts, n, vert_coords, err );MSQ_ERRRTN( err );

    // calculate Jacobian
    jacobian_2D( pd, type, n, sample, vert_coords, W_out, err );MSQ_ERRRTN( err );
}
void MBMesquite::TargetCalculator::get_refmesh_Jacobian_3D ( ReferenceMeshInterface ref_mesh,
PatchData active_mesh,
size_t  element_no,
Sample  sample_no,
MsqMatrix< 3, 3 > &  W_out,
MsqError err 
) [static]

Definition at line 610 of file TargetCalculator.cpp.

References MBMesquite::PatchData::element_by_index(), MBMesquite::MsqMeshEntity::get_element_type(), MBMesquite::ReferenceMeshInterface::get_reference_vertex_coordinates(), MBMesquite::PatchData::get_vertex_handles_array(), MBMesquite::MsqMeshEntity::get_vertex_index_array(), jacobian_3D(), moab::MAX_NODES, MSQ_ERRRTN, n, and MBMesquite::MsqMeshEntity::node_count().

Referenced by MBMesquite::RefMeshTargetCalculator::get_3D_target(), and TargetCalculatorTest::test_get_refmesh_Jacobian_3D().

{
    // get element
    MsqMeshEntity& elem       = pd.element_by_index( element );
    const EntityTopology type = elem.get_element_type();
    const unsigned n          = elem.node_count();

    const unsigned MAX_NODES = 27;
    assert( n <= MAX_NODES );

    // get vertices
    Mesh::VertexHandle elem_verts[MAX_NODES];
    const std::size_t* vtx_idx        = elem.get_vertex_index_array();
    const Mesh::VertexHandle* vtx_hdl = pd.get_vertex_handles_array();
    for( unsigned i = 0; i < n; ++i )
        elem_verts[i] = vtx_hdl[vtx_idx[i]];

    // get vertex coordinates
    Vector3D vert_coords[MAX_NODES];
    ref_mesh->get_reference_vertex_coordinates( elem_verts, n, vert_coords, err );MSQ_ERRRTN( err );

    // calculate Jacobian
    jacobian_3D( pd, type, n, sample, vert_coords, W_out, err );MSQ_ERRRTN( err );
}
void MBMesquite::TargetCalculator::ideal_shape_2D ( EntityTopology  element_type,
Sample  s,
const PatchData pd,
MsqMatrix< 2, 2 > &  W,
MsqError err 
) [static]

Get skew matrix for an ideally shaped element.

Definition at line 488 of file TargetCalculator.cpp.

References MBMesquite::PatchData::get_mapping_function_2D(), MBMesquite::MappingFunction2D::ideal(), MBMesquite::ideal_constant_skew_I_2D(), MSQ_ERRRTN, MSQ_SETERR, shape(), and MBMesquite::MsqError::UNSUPPORTED_ELEMENT.

Referenced by MBMesquite::IdealShapeTarget::get_2D_target(), TargetCalculatorTest::test_ideal_shape_quad(), and TargetCalculatorTest::test_ideal_shape_tri().

{
    if( !ideal_constant_skew_I_2D( element_type, q ) )
    {
        const MappingFunction2D* map = pd.get_mapping_function_2D( element_type );
        if( !map )
        {
            MSQ_SETERR( err )( MsqError::UNSUPPORTED_ELEMENT );
            return;
        }
        MsqMatrix< 3, 2 > J;
        map->ideal( s, J, err );MSQ_ERRRTN( err );
        q = TargetCalculator::shape( J );
    }
}
void MBMesquite::TargetCalculator::ideal_skew_2D ( EntityTopology  element_type,
Sample  s,
const PatchData pd,
MsqMatrix< 2, 2 > &  W,
MsqError err 
) [static]

Get skew matrix for an ideally shaped element.

Definition at line 449 of file TargetCalculator.cpp.

References MBMesquite::PatchData::get_mapping_function_2D(), MBMesquite::MappingFunction2D::ideal(), MBMesquite::ideal_constant_skew_I_2D(), MSQ_ERRRTN, MSQ_SETERR, skew(), and MBMesquite::MsqError::UNSUPPORTED_ELEMENT.

Referenced by TargetCalculatorTest::test_ideal_shape_tri(), TargetCalculatorTest::test_ideal_skew_quad(), and TargetCalculatorTest::test_ideal_skew_tri().

{
    if( !ideal_constant_skew_I_2D( element_type, q ) )
    {
        const MappingFunction2D* map = pd.get_mapping_function_2D( element_type );
        if( !map )
        {
            MSQ_SETERR( err )( MsqError::UNSUPPORTED_ELEMENT );
            return;
        }
        MsqMatrix< 3, 2 > J;
        map->ideal( s, J, err );MSQ_ERRRTN( err );
        q = TargetCalculator::skew( J );
    }
}
void MBMesquite::TargetCalculator::initialize_queue ( MeshDomainAssoc mesh_and_domain,
const Settings settings,
MsqError err 
) [virtual]

Called at start of instruction queue processing.

Do any preliminary global initialization, consistency checking, etc. Default implementation does nothing.

Definition at line 672 of file TargetCalculator.cpp.

Referenced by MBMesquite::TMPQualityMetric::initialize_queue().

{}
void MBMesquite::TargetCalculator::jacobian_2D ( PatchData pd,
EntityTopology  element_type,
int  num_nodes,
Sample  location,
const Vector3D coords,
MsqMatrix< 3, 2 > &  W_out,
MsqError err 
) [static]

Calculate the Jacobian given element vertex coordinates.

Definition at line 577 of file TargetCalculator.cpp.

References MBMesquite::MappingFunction2D::derivatives(), MBMesquite::PatchData::get_mapping_function_2D(), MBMesquite::get_nodeset(), moab::MAX_NODES, MSQ_ERRRTN, MSQ_SETERR, n, MBMesquite::outer(), and MBMesquite::MsqError::UNSUPPORTED_ELEMENT.

Referenced by get_refmesh_Jacobian_2D(), TargetCalculatorTest::test_get_refmesh_Jacobian_2D(), TargetCalculatorTest::test_ideal_skew_tri(), and TargetCalculatorTest::test_jacobian_2D().

{
    // Get element properties
    NodeSet bits = get_nodeset( type, num_nodes, err );MSQ_ERRRTN( err );
    const MappingFunction2D* mf = pd.get_mapping_function_2D( type );
    if( !mf )
    {
        MSQ_SETERR( err )( MsqError::UNSUPPORTED_ELEMENT );
        return;
    }

    // Get mapping function derivatives
    const int MAX_NODES = 9;
    assert( num_nodes <= MAX_NODES );
    size_t indices[MAX_NODES], n;
    MsqVector< 2 > derivs[MAX_NODES];
    mf->derivatives( location, bits, indices, derivs, n, err );MSQ_ERRRTN( err );

    // calculate Jacobian
    assert( sizeof( Vector3D ) == sizeof( MsqVector< 3 > ) );
    const MsqVector< 3 >* verts = reinterpret_cast< const MsqVector< 3 >* >( coords );
    assert( n > 0 );
    J = outer( verts[indices[0]], derivs[0] );
    for( size_t i = 1; i < n; ++i )
        J += outer( verts[indices[i]], derivs[i] );
}
void MBMesquite::TargetCalculator::jacobian_3D ( PatchData pd,
EntityTopology  element_type,
int  num_nodes,
Sample  location,
const Vector3D coords,
MsqMatrix< 3, 3 > &  W_out,
MsqError err 
) [static]

Calculate the Jacobian given element vertex coordinates.

Definition at line 544 of file TargetCalculator.cpp.

References MBMesquite::MappingFunction3D::derivatives(), MBMesquite::PatchData::get_mapping_function_3D(), MBMesquite::get_nodeset(), moab::MAX_NODES, MSQ_ERRRTN, MSQ_SETERR, n, MBMesquite::outer(), and MBMesquite::MsqError::UNSUPPORTED_ELEMENT.

Referenced by get_refmesh_Jacobian_3D(), TargetCalculatorTest::test_get_refmesh_Jacobian_3D(), TargetCalculatorTest::test_ideal_shape_pyramid(), TargetCalculatorTest::test_ideal_skew_prism(), TargetCalculatorTest::test_ideal_skew_pyramid(), TargetCalculatorTest::test_ideal_skew_tet(), and TargetCalculatorTest::test_jacobian_3D().

{
    // Get element properties
    NodeSet bits = get_nodeset( type, num_nodes, err );MSQ_ERRRTN( err );
    const MappingFunction3D* mf = pd.get_mapping_function_3D( type );
    if( !mf )
    {
        MSQ_SETERR( err )( MsqError::UNSUPPORTED_ELEMENT );
        return;
    }

    // Get mapping function derivatives
    const int MAX_NODES = 27;
    assert( num_nodes <= MAX_NODES );
    size_t indices[MAX_NODES], n;
    MsqVector< 3 > derivs[MAX_NODES];
    mf->derivatives( location, bits, indices, derivs, n, err );MSQ_ERRRTN( err );

    // calculate Jacobian
    assert( sizeof( Vector3D ) == sizeof( MsqVector< 3 > ) );
    const MsqVector< 3 >* verts = reinterpret_cast< const MsqVector< 3 >* >( coords );
    assert( n > 0 );
    J = outer( verts[indices[0]], derivs[0] );
    for( size_t i = 1; i < n; ++i )
        J += outer( verts[indices[i]], derivs[i] );
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::new_aspect_2D ( const MsqVector< 2 > &  r) [static]

Create a new aspect ratio matrix.

Create an aspect ratio matrix such that the ratio of column lengths is proportional to the ratio of the corresponding pair of values in the passed vector.

Definition at line 517 of file TargetCalculator.cpp.

Referenced by TargetCalculatorTest::test_new_aspect_2D().

{
    return new_aspect_2D( r[0] / r[1] );
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::new_aspect_2D ( double  rho) [static]

Create a new aspect ratio matrix.

Create an aspect ratio matrix such that the ratio of column lengths is the passed value.

Definition at line 522 of file TargetCalculator.cpp.

{
    MsqMatrix< 2, 2 > W( sqrt( rho ) );
    W( 1, 1 ) = 1.0 / W( 0, 0 );
    return W;
}
MsqMatrix< 3, 3 > MBMesquite::TargetCalculator::new_aspect_3D ( const MsqVector< 3 > &  r) [static]

Create a new aspect ratio matrix.

Create an aspect ratio matrix such that the ratio of column lengths is proportional to the ratio of the corresponding pair of values in the passed vector.

Definition at line 508 of file TargetCalculator.cpp.

References MBMesquite::cbrt().

Referenced by TargetCalculatorTest::test_new_aspect_3D().

{
    MsqMatrix< 3, 3 > W( 0.0 );
    W( 0, 0 ) = MBMesquite::cbrt( r[0] * r[0] / ( r[1] * r[2] ) );
    W( 1, 1 ) = MBMesquite::cbrt( r[1] * r[1] / ( r[0] * r[2] ) );
    W( 2, 2 ) = MBMesquite::cbrt( r[2] * r[2] / ( r[1] * r[0] ) );
    return W;
}
MsqMatrix< 3, 2 > MBMesquite::TargetCalculator::new_orientation_2D ( const MsqVector< 3 > &  b1,
const MsqVector< 3 > &  b2 
) [static]

Create a new orientation matrix.

Create an orientation matrix such that the first and second Jacobian columns of W are aligned to the passed vectors.

Definition at line 354 of file TargetCalculator.cpp.

References b1, MBMesquite::length(), and MBMesquite::MsqMatrix< R, C >::set_column().

Referenced by TargetCalculatorTest::test_new_orientatin_2D().

{
    double lb1_sqr = b1 % b1;
    double inv_lb1 = 1.0 / sqrt( lb1_sqr );
    MsqMatrix< 3, 2 > V;
    V.set_column( 0, b1 * inv_lb1 );
    V.set_column( 1, ( lb1_sqr * b2 - ( b1 % b2 ) * b1 ) * ( inv_lb1 / length( b1 * b2 ) ) );
    return V;
}
MsqMatrix< 3, 3 > MBMesquite::TargetCalculator::new_orientation_3D ( const MsqVector< 3 > &  b1,
const MsqVector< 3 > &  b2 
) [static]

Create a new orientation matrix.

Create an orientation matrix such that the first and second Jacobian columns of W are aligned to the passed vectors.

Definition at line 340 of file TargetCalculator.cpp.

References b1, b2, moab::cross(), MBMesquite::length(), and MBMesquite::MsqMatrix< R, C >::set_column().

Referenced by TargetCalculatorTest::test_new_orientatin_3D().

{
    double lb1_sqr       = b1 % b1;
    MsqVector< 3 > cross = b1 * b2;
    double lb1           = sqrt( lb1_sqr );
    double inv_lb1       = 1.0 / lb1;
    double inv_lx        = 1.0 / length( cross );
    MsqMatrix< 3, 3 > V;
    V.set_column( 0, inv_lb1 * b1 );
    V.set_column( 1, ( lb1_sqr * b2 - ( b1 % b2 ) * lb1 ) * inv_lb1 * inv_lx );
    V.set_column( 2, cross * inv_lx );
    return V;
}
MsqMatrix< 3, 3 > MBMesquite::TargetCalculator::shape ( const MsqMatrix< 3, 3 > &  W) [static]

Get shape (skew and aspect) component of W.

Definition at line 159 of file TargetCalculator.cpp.

References MBMesquite::cbrt(), MBMesquite::MsqMatrix< R, C >::column(), and MBMesquite::length().

Referenced by ideal_shape_2D(), ideal_shape_3D(), TargetCalculatorTest::test_ideal_shape_pyramid(), TargetCalculatorTest::test_shape_2D(), TargetCalculatorTest::test_shape_3D(), and TargetCalculatorTest::test_shape_surface().

{
    MsqVector< 3 > a1    = W.column( 0 );
    MsqVector< 3 > a2    = W.column( 1 );
    MsqVector< 3 > a3    = W.column( 2 );
    MsqVector< 3 > a1xa2 = a1 * a2;

    double len1  = length( a1 );
    double lenx  = length( a1xa2 );
    double alpha = fabs( a1xa2 % a3 );
    double coeff = MBMesquite::cbrt( 1 / alpha );
    double inv1  = 1.0 / len1;
    double invx  = 1.0 / lenx;

    MsqMatrix< 3, 3 > q;
    q( 0, 0 ) = coeff * len1;
    q( 0, 1 ) = coeff * inv1 * ( a1 % a2 );
    q( 0, 2 ) = coeff * inv1 * ( a1 % a3 );
    q( 1, 0 ) = 0.0;
    q( 1, 1 ) = coeff * inv1 * lenx;
    q( 1, 2 ) = coeff * invx * inv1 * ( a1xa2 % ( a1 * a3 ) );
    q( 2, 0 ) = 0.0;
    q( 2, 1 ) = 0.0;
    q( 2, 2 ) = coeff * alpha * invx;
    return q;
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::shape ( const MsqMatrix< 3, 2 > &  W) [static]

Get skew component of W.

Definition at line 186 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column(), and MBMesquite::length().

{
    double len1       = length( W.column( 0 ) );
    double inv1       = 1.0 / len1;
    double root_alpha = sqrt( length( W.column( 0 ) * W.column( 1 ) ) );
    double coeff      = 1.0 / root_alpha;

    MsqMatrix< 2, 2 > result;
    result( 0, 0 ) = coeff * len1;
    result( 0, 1 ) = coeff * inv1 * ( W.column( 0 ) % W.column( 1 ) );
    result( 1, 0 ) = 0.0;
    result( 1, 1 ) = root_alpha * inv1;
    return result;
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::shape ( const MsqMatrix< 2, 2 > &  W) [static]

Get skew component of W.

Definition at line 201 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column(), MBMesquite::det(), and MBMesquite::length().

{
    double len1       = length( W.column( 0 ) );
    double inv1       = 1.0 / len1;
    double root_alpha = sqrt( fabs( det( W ) ) );
    double coeff      = 1.0 / root_alpha;

    MsqMatrix< 2, 2 > result;
    result( 0, 0 ) = coeff * len1;
    result( 0, 1 ) = coeff * inv1 * ( W.column( 0 ) % W.column( 1 ) );
    result( 1, 0 ) = 0.0;
    result( 1, 1 ) = root_alpha * inv1;
    return result;
}
double MBMesquite::TargetCalculator::size ( const MsqMatrix< 3, 2 > &  W) [static]

Get size component of W.

Definition at line 48 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column(), and MBMesquite::length().

{
    return sqrt( length( W.column( 0 ) * W.column( 1 ) ) );
}
double MBMesquite::TargetCalculator::size ( const MsqMatrix< 2, 2 > &  W) [static]

Get size component of W.

Definition at line 53 of file TargetCalculator.cpp.

References MBMesquite::det().

{
    return sqrt( fabs( det( W ) ) );
}
MsqMatrix< 3, 3 > MBMesquite::TargetCalculator::skew ( const MsqMatrix< 3, 3 > &  W) [static]

Get skew component of W.

Definition at line 58 of file TargetCalculator.cpp.

References MBMesquite::cbrt(), MBMesquite::MsqMatrix< R, C >::column(), moab::dot(), and MBMesquite::length().

Referenced by IdealTargetTest::get_ideal_target(), ideal_skew_2D(), ideal_skew_3D(), TargetCalculatorTest::test_ideal_skew_prism(), TargetCalculatorTest::test_ideal_skew_pyramid(), TargetCalculatorTest::test_ideal_skew_tet(), TargetCalculatorTest::test_ideal_skew_tri(), TargetCalculatorTest::test_skew_2D(), TargetCalculatorTest::test_skew_3D(), and TargetCalculatorTest::test_skew_surface().

{
    MsqVector< 3 > a1    = W.column( 0 ) * ( 1.0 / length( W.column( 0 ) ) );
    MsqVector< 3 > a2    = W.column( 1 ) * ( 1.0 / length( W.column( 1 ) ) );
    MsqVector< 3 > a3    = W.column( 2 ) * ( 1.0 / length( W.column( 2 ) ) );
    MsqVector< 3 > a1xa2 = a1 * a2;

    double lenx  = length( a1xa2 );
    double alpha = fabs( a1xa2 % a3 );
    double coeff = MBMesquite::cbrt( 1 / alpha );
    double dot   = a1xa2 % ( a1 * a3 );

    MsqMatrix< 3, 3 > q;
    q( 0, 0 ) = coeff;
    q( 0, 1 ) = coeff * ( a1 % a2 );
    q( 0, 2 ) = coeff * ( a1 % a3 );
    q( 1, 0 ) = 0.0;
    q( 1, 1 ) = coeff * lenx;
    q( 1, 2 ) = coeff * dot / lenx;
    q( 2, 0 ) = 0.0;
    q( 2, 1 ) = 0.0;
    q( 2, 2 ) = coeff * alpha / lenx;
    return q;
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::skew ( const MsqMatrix< 3, 2 > &  W) [static]

Get skew component of W.

Definition at line 83 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column(), moab::dot(), and MBMesquite::length().

{
    MsqVector< 3 > alpha = W.column( 0 ) * W.column( 1 );
    double a1_sqr        = W.column( 0 ) % W.column( 0 );
    double a2_sqr        = W.column( 1 ) % W.column( 1 );
    double dot           = W.column( 0 ) % W.column( 1 );
    double a1a2          = sqrt( a1_sqr * a2_sqr );
    double coeff         = sqrt( a1a2 / length( alpha ) );

    MsqMatrix< 2, 2 > result;
    result( 0, 0 ) = coeff;
    result( 0, 1 ) = coeff * dot / a1a2;
    result( 1, 0 ) = 0.0;
    result( 1, 1 ) = 1 / coeff;
    return result;
}
MsqMatrix< 2, 2 > MBMesquite::TargetCalculator::skew ( const MsqMatrix< 2, 2 > &  W) [static]

Get skew component of W.

Definition at line 100 of file TargetCalculator.cpp.

References MBMesquite::MsqMatrix< R, C >::column(), MBMesquite::det(), and moab::dot().

{
    double alpha  = fabs( det( W ) );
    double a1_sqr = W.column( 0 ) % W.column( 0 );
    double a2_sqr = W.column( 1 ) % W.column( 1 );
    double dot    = W.column( 0 ) % W.column( 1 );
    double a1a2   = sqrt( a1_sqr * a2_sqr );
    double coeff  = sqrt( a1a2 / alpha );

    MsqMatrix< 2, 2 > result;
    result( 0, 0 ) = coeff;
    result( 0, 1 ) = coeff * dot / a1a2;
    result( 1, 0 ) = 0.0;
    result( 1, 1 ) = 1 / coeff;
    return result;
}

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