| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

In This Package:

SoUtils Namespace Reference

from HepVis package More...


Classes

class  SbProjector

Functions

StatusCode eigensystem (const Gaudi::SymMatrix3x3 &matrix, Gaudi::Vector3 &eigenvalues, std::vector< Gaudi::Vector3 > &eigenvectors)
 calculate eigenvalues and eigenvectors of symmetric matrix
StatusCode eigenvalues (const Gaudi::SymMatrix3x3 &matrix, Gaudi::Vector3 &eigenvalues)
 calculate eigenvalues of symmetric matrix
StatusCode ellipsoid (const Gaudi::XYZPoint &center, const Gaudi::SymMatrix3x3 &cov, SoEllipsoid *&node)
 create SoEllipsoid Inventor node from the covarinace matrix and position values of error codes from "eigensystem" method
StatusCode ellipticalprism (const Gaudi::XYZPoint &center, const Gaudi::SymMatrix2x2 &cov, const double extent, SoEllipticalPrism *&node)
 create SoEllipticalPrism Inventor node from the covarinace matrix and position values of error codes from "eigensystem" method

Detailed Description

from HepVis package

class SbProjector

To handle special (non-linear) projections (rz, phiz, etc...)


Function Documentation

StatusCode SoUtils::eigensystem ( const Gaudi::SymMatrix3x3 matrix,
Gaudi::Vector3 eigenvalues,
std::vector< Gaudi::Vector3 > &  eigenvectors 
)

calculate eigenvalues and eigenvectors of symmetric matrix

return value of status code

Parameters:
matrix reference to symmetric matrix
eigenvalues (return) vector of eigenvalues
einenvectors (return) vector of eigenvectors
status code

Definition at line 58 of file EigenSystems.cpp.

00061 {
00062   // dimension of matrix 
00063   //  const size_t num = matrix.num_row() ;
00064   const size_t num = 3 ;
00065   // reset the output values 
00066   eigenvalues  = Gaudi::Vector3 (0,0,0);
00067   eigenvectors = std::vector<Gaudi::Vector3>( num , Gaudi::Vector3( 0, 0 , 0 ) ) ;
00068   if( 0 == num ) { return StatusCode::SUCCESS ; }                    // RETURN
00069   
00070   // allocate GSL matrices 
00071   gsl_matrix* mtrx = gsl_matrix_alloc( num , num ) ;
00072   gsl_matrix* evec = gsl_matrix_alloc( num , num ) ;
00073   if ( 0 == mtrx  ) 
00074   { GSL_ERROR_VAL("SoCaloUtils::eigenvalues: the matrix(1) is not allocated", 
00075                   1000 , StatusCode( 1000 ) ) ; }                  // RETURN
00076   if ( 0 == evec  ) 
00077   { GSL_ERROR_VAL("SoCaloUtils::eigenvalues: the matrix(2) is not allocated", 
00078                   1000 , StatusCode( 1000 ) ) ; }                  // RETURN
00079   
00080   // allocate working space 
00081   gsl_eigen_symmv_workspace* wspace = gsl_eigen_symmv_alloc( num ) ;
00082   if ( 0 == wspace  ) 
00083   { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: workspace is not allocated",
00084                   1001 , StatusCode ( 1001 ) ) ; }               // RETURN
00085   
00086   // allocate the vector of eigenvalues 
00087   gsl_vector* eval  = gsl_vector_alloc ( num ) ;
00088   gsl_vector* vcol  = gsl_vector_alloc ( num ) ;
00089   if ( 0 == eval    ) 
00090   { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: vector(1) is not allocated",
00091                   1002 , StatusCode ( 1002 ) ) ; }             // RETURN  
00092   if ( 0 == vcol  ) 
00093   { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: vector(2) is not allocated",
00094                   1002 , StatusCode ( 1002 ) ) ; }             // RETURN  
00095   
00096   // copy  symmetric matrix into GSL matrix 
00097   { for( size_t i = 0 ; i < num ; ++i ) 
00098   { for( size_t j = 0 ; j < num ; ++j ) 
00099   { gsl_matrix_set ( mtrx , i , j , matrix( i + 1 , j + 1 ) ) ; } } }
00100   
00101   int ierror = gsl_eigen_symmv ( mtrx , eval , evec , wspace ) ;
00102   if ( ierror      )
00103   { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: values are not computed",
00104                   1010 + ierror , StatusCode ( 1010 + ierror ) ) ; } 
00105   
00106   // decode  GSL values into Vector
00107   { for( size_t i = 0 ; i < num ; ++i ) 
00108   { eigenvalues( i + 1 ) = gsl_vector_get( eval , i ) ; } }
00109   
00110   
00111   { // get eigenvectors 
00112     for( size_t i = 0 ; i < num ; ++i ) 
00113     {
00114       ierror = gsl_matrix_get_col( vcol , evec , i ) ;
00115       if ( ierror ) 
00116       { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: invalid matrix operation",
00117                       1100 + ierror , StatusCode ( 1100 + ierror ) ) ; } 
00118       // copy eigenvector 
00119       for( size_t j = 0 ; j < num ; ++j ) 
00120       { eigenvectors[i]( j + 1 ) = gsl_vector_get( vcol , j ) ; }
00121     } 
00122   }
00123   
00124   // free everything at the end 
00125   if ( 0 != mtrx   ) { gsl_matrix_free      ( mtrx   ) ; }
00126   if ( 0 != evec   ) { gsl_matrix_free      ( evec   ) ; }
00127   if ( 0 != wspace ) { gsl_eigen_symmv_free ( wspace ) ; }
00128   if ( 0 != eval   ) { gsl_vector_free      ( eval   ) ; }
00129   if ( 0 != vcol   ) { gsl_vector_free      ( vcol   ) ; }
00130   
00131   return StatusCode::SUCCESS ;
00132   
00133 };

StatusCode SoUtils::eigenvalues ( const Gaudi::SymMatrix3x3 matrix,
Gaudi::Vector3 eigenvalues 
)

calculate eigenvalues of symmetric matrix

return value of status code

Parameters:
matrix reference to symmetric matrix
eigenvalues (return) vector of eigenvalues
status code

Definition at line 151 of file EigenSystems.cpp.

00153 {
00154   // dimension of matrix 
00155   //  const size_t num = matrix.num_row() ;
00156   const size_t num = 3 ;
00157   // reset the output values 
00158   eigenvalues = Gaudi::Vector3( 0, 0 , 0 );
00159   if ( 0 == num ) { return StatusCode::SUCCESS ; }             // RETURN
00160   
00161   // alocate GSL matrix 
00162   gsl_matrix* mtrx = gsl_matrix_alloc( num , num ) ;
00163   if ( 0 == mtrx    ) 
00164   { GSL_ERROR_VAL("SoCaloUtils::eigenvalues: the matrix is not allocated", 
00165                   1000 , StatusCode( 1000 ) ) ; }               // RETURN
00166   
00167   // allocate working space 
00168   gsl_eigen_symm_workspace* wspace = gsl_eigen_symm_alloc( num ) ;
00169   if ( 0 == wspace  ) 
00170   { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: workspace is not allocated",
00171                   1001 , StatusCode ( 1001 ) ) ; }             // RETURN  
00172   
00173   // allocate the vector of eigenvalues 
00174   gsl_vector* eval  = gsl_vector_alloc ( num ) ;
00175   if ( 0 == eval    ) 
00176   { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: vector is not allocated",
00177                   1002 , StatusCode ( 1002 ) ) ; }             // RETURN  
00178   
00179   // copy CLHEP symmetric matrix into GSL matrix 
00180   { for( size_t i = 0 ; i < num ; ++i ) 
00181   { for( size_t j = 0 ; j < num ; ++j ) 
00182   { gsl_matrix_set ( mtrx , i , j , matrix( i + 1 , j + 1 ) ) ; } } }
00183   
00184   // use GSL 
00185   int ierror = gsl_eigen_symm( mtrx , eval , wspace ) ;
00186   if ( ierror      ) 
00187   { GSL_ERROR_VAL("SoCaloUtil::eigenvalues: values are not computed",
00188                   1010 + ierror , StatusCode ( 1010 + ierror ) ) ; }
00189   
00190   // decode  GSL values into Vector
00191   { for( size_t i = 0 ; i < num ; ++i ) 
00192   { eigenvalues( i + 1 ) = gsl_vector_get( eval , i ) ; } }
00193   
00194   // free everything at the end 
00195   if ( 0 != mtrx   ) { gsl_matrix_free     ( mtrx   ) ; }
00196   if ( 0 != wspace ) { gsl_eigen_symm_free ( wspace ) ; }
00197   if ( 0 != eval   ) { gsl_vector_free     ( eval   ) ; }
00198   
00199   return StatusCode::SUCCESS ;
00200 };

StatusCode SoUtils::ellipsoid ( const Gaudi::XYZPoint center,
const Gaudi::SymMatrix3x3 cov,
SoEllipsoid *&  node 
)

create SoEllipsoid Inventor node from the covarinace matrix and position values of error codes from "eigensystem" method

Parameters:
center positinoi of ellipsoid center
cov "covarinace" matrix for the ellipsoid
node (return) pointer to created SoEllipsoid Inventor node
Returns:
status code

reset the output

Definition at line 60 of file Ellipsoid.cpp.

00063 {
00065   node = 0 ;
00066   //  if( 3 != cov.num_row()    ) { return StatusCode::FAILURE ; } ///< RETURN !
00068   Gaudi::Vector3              evals  ( 0, 0 , 0 ) ;
00069   typedef std::vector<Gaudi::Vector3> Vectors ;
00070   Vectors evects           ;
00071   StatusCode  sc = SoUtils::eigensystem( cov , evals , evects );
00072   if( sc.isFailure       () ) { return sc                  ; } 
00073   if( 3 != evects.size   () ) { return StatusCode::FAILURE ; } 
00074   //  if( 3 != evals.num_row () ) { return StatusCode::FAILURE ; } ///< RETURN ! 
00076 
00077   SoEllipsoid* ell = new SoEllipsoid();
00078   ell->center.
00079     setValue( (float)center.x() , (float)center.y() , (float)center.z() );
00080   ell->eigenvalues.
00081     setValue( (float)sqrt( evals( 1 )  ) , 
00082               (float)sqrt( evals( 2 )  ) , 
00083               (float)sqrt( evals( 3 )  ) ) ;
00085   Gaudi::Rotation3D rot;
00086   Gaudi::XYZVector vx ( evects[0](1) , evects[0](2) , evects[0](3) ) ;
00087   Gaudi::XYZVector vy ( evects[1](1) , evects[1](2) , evects[1](3) ) ;
00089   //OD  rot.rotateAxes ( vx , vy , vx.Cross( vy ) );
00090 
00091   Gaudi::XYZVector axis  ;
00092   double      angle = 0;
00093   //OD  rot.getAngleAxis( angle , axis );
00094   
00095   ell->rotation.setValue
00096     (SbVec3f((float)axis.x(),(float)axis.y(),(float)axis.z()),(float)angle);
00098   node = ell ;
00100   return StatusCode::SUCCESS ;
00101 };

StatusCode SoUtils::ellipticalprism ( const Gaudi::XYZPoint center,
const Gaudi::SymMatrix2x2 cov,
const double  extent,
SoEllipticalPrism *&  node 
)

create SoEllipticalPrism Inventor node from the covarinace matrix and position values of error codes from "eigensystem" method

Parameters:
center position of elliptical prism
cov "covarinace" matrix (2D) for the ellipsoid
extent extent for elliptical prism
node (return) pointer to created SoEllipsoid Inventor node
Returns:
status code

reset the output

return StatusCode::FAILURE;

Definition at line 61 of file EllipticalPrism.cpp.

00065 {
00067   node = 0 ;
00068   
00069   if( cov == 0  ) { 
00071   } 
00072   Gaudi::Vector3              evals  ( 0 , 0 ,0 ) ;
00073   std::vector<Gaudi::Vector3> evects           ;
00074   //OD StatusCode sc = SoUtils::eigensystem( cov , evals , evects );
00075   //if( sc.isFailure       () ) { return sc                  ; } ///< RETURN !
00076   if( 2 != evects.size   () ) { return StatusCode::FAILURE ; } 
00077   //  if( 2 != evals.num_row () ) { return StatusCode::FAILURE ; } ///< RETURN ! 
00079   SoEllipticalPrism* ell = new SoEllipticalPrism();
00080   ell->center.
00081     setValue( (float)center.x() , (float)center.y() , (float)center.z() );
00082   ell->eigenvalues.
00083     setValue( (float)sqrt( evals( 1 )  ) , 
00084               (float)sqrt( evals( 2 )  ) ) ;
00086 
00087   Gaudi::Rotation3D rot;
00088   Gaudi::XYZVector vx ( evects[0](1) , evects[0](2) , 0 ) ;
00089   Gaudi::XYZVector vy ( evects[1](1) , evects[1](2) , 0 ) ;
00091   //OD  rot.rotateAxes ( vx , vy , vx.cross( vy ) );
00092   Gaudi::XYZVector axis  ;
00093   double      angle = 0;
00094   //OD rot.getAngleAxis( angle , axis );
00095   ell->rotation.setValue
00096     (SbVec3f((float)axis.x(),(float)axis.y(),(float)axis.z()),(float)angle);
00097   ell->extent.setValue( (float)extent );
00099   node = ell ;
00101   return StatusCode::SUCCESS ;
00102 };  

| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

Generated on Mon Apr 11 20:01:49 2011 for SoUtils by doxygen 1.4.7