MOAB: Mesh Oriented datABase  (version 5.4.1)
MeshTransform.cpp
Go to the documentation of this file.
00001 /* *****************************************************************
00002     MESQUITE -- The Mesh Quality Improvement Toolkit
00003 
00004     Copyright 2004 Sandia Corporation and Argonne National
00005     Laboratory.  Under the terms of Contract DE-AC04-94AL85000
00006     with Sandia Corporation, the U.S. Government retains certain
00007     rights in this software.
00008 
00009     This library is free software; you can redistribute it and/or
00010     modify it under the terms of the GNU Lesser General Public
00011     License as published by the Free Software Foundation; either
00012     version 2.1 of the License, or (at your option) any later version.
00013 
00014     This library is distributed in the hope that it will be useful,
00015     but WITHOUT ANY WARRANTY; without even the implied warranty of
00016     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00017     Lesser General Public License for more details.
00018 
00019     You should have received a copy of the GNU Lesser General Public License
00020     (lgpl.txt) along with this library; if not, write to the Free Software
00021     Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00022 
00023     [email protected], [email protected], [email protected],
00024     [email protected], [email protected], [email protected]
00025 
00026   ***************************************************************** */
00027 /*!
00028   \file   MeshTransform.cpp
00029   \brief
00030 
00031   The MeshTransform Class is the base class for all the smoothing algorythms
00032 
00033   \author Michael Brewer
00034   \date   2004-11-06
00035 */
00036 
00037 #include "MeshTransform.hpp"
00038 #include "MeshInterface.hpp"
00039 #include "MsqVertex.hpp"
00040 #include "MsqError.hpp"
00041 
00042 namespace MBMesquite
00043 {
00044 
00045 MeshTransform::~MeshTransform() {}
00046 
00047 /*!
00048   Actually apply the affine transformation
00049   */
00050 double MeshTransform::loop_over_mesh( MeshDomainAssoc* mesh_and_domain, const Settings*, MsqError& err )
00051 {
00052     Mesh* mesh = mesh_and_domain->get_mesh();
00053 
00054     std::vector< Mesh::VertexHandle > handle_list;
00055     mesh->get_all_vertices( handle_list, err );
00056     if( MSQ_CHKERR( err ) ) return 1.0;
00057 
00058     std::vector< bool > fixed( 1 );
00059     MsqVertex vertex;
00060     std::vector< Mesh::VertexHandle >::const_iterator iter;
00061     for( iter = handle_list.begin(); iter != handle_list.end(); ++iter )
00062     {
00063         mesh->vertices_get_coordinates( &*iter, &vertex, 1, err );
00064         if( MSQ_CHKERR( err ) ) return 1.0;
00065 
00066         if( skipFixed )
00067         {
00068             mesh->vertices_get_fixed_flag( &*iter, fixed, 1, err );
00069             if( MSQ_CHKERR( err ) ) return 1.0;
00070             if( fixed.front() ) continue;
00071         }
00072 
00073         vertex = mMat * vertex + mVec;
00074 
00075         mesh->vertex_set_coordinates( *iter, vertex, err );
00076         if( MSQ_CHKERR( err ) ) return 1.0;
00077     }
00078 
00079     return 0.0;
00080 }
00081 
00082 void MeshTransform::add_translation( const Vector3D& offset )
00083 {
00084     mVec += offset;
00085 }
00086 
00087 void MeshTransform::add_rotation( const Vector3D& axis, double radians )
00088 {
00089     const double c   = cos( radians );
00090     const double s   = sin( radians );
00091     const Vector3D a = axis / axis.length();
00092     const Matrix3D m1( c, -a[2] * s, a[1] * s, a[2] * s, c, -a[0] * s, -a[1] * s, a[0] * s, c );
00093     Matrix3D m2;
00094     m2.outer_product( a, a );
00095     Matrix3D rot = m1 + ( 1.0 - c ) * m2;
00096     mMat         = rot * mMat;
00097     mVec         = rot * mVec;
00098 }
00099 
00100 void MeshTransform::add_scale( double factor )
00101 {
00102     add_scale( Vector3D( factor ) );
00103 }
00104 
00105 void MeshTransform::add_scale( const Vector3D& f )
00106 {
00107     for( int i = 0; i < 3; ++i )
00108     {
00109         mVec[i] *= f[i];
00110         mMat[i][0] *= f[i];
00111         mMat[i][1] *= f[i];
00112         mMat[i][2] *= f[i];
00113     }
00114 }
00115 
00116 void MeshTransform::initialize_queue( MeshDomainAssoc*, const Settings*, MsqError& ) {}
00117 
00118 }  // namespace MBMesquite
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines