Main Page   Class Hierarchy   Alphabetical List   Compound List   File List   Compound Members   File Members  

Brick.cpp

Go to the documentation of this file.
00001 /* ****************************************************************** **
00002 **    OpenSees - Open System for Earthquake Engineering Simulation    **
00003 **          Pacific Earthquake Engineering Research Center            **
00004 **                                                                    **
00005 **                                                                    **
00006 ** (C) Copyright 1999, The Regents of the University of California    **
00007 ** All Rights Reserved.                                               **
00008 **                                                                    **
00009 ** Commercial use of this program without express permission of the   **
00010 ** University of California, Berkeley, is strictly prohibited.  See   **
00011 ** file 'COPYRIGHT'  in main directory for information on usage and   **
00012 ** redistribution,  and for a DISCLAIMER OF ALL WARRANTIES.           **
00013 **                                                                    **
00014 ** Developed by:                                                      **
00015 **   Frank McKenna (fmckenna@ce.berkeley.edu)                         **
00016 **   Gregory L. Fenves (fenves@ce.berkeley.edu)                       **
00017 **   Filip C. Filippou (filippou@ce.berkeley.edu)                     **
00018 **                                                                    **
00019 ** ****************************************************************** */
00020                                                                         
00021 // $Revision: 1.2 $
00022 // $Date: 2001/07/11 23:15:55 $
00023 // $Source: /usr/local/cvs/OpenSees/SRC/element/brick/Brick.cpp,v $
00024 
00025 // Ed "C++" Love
00026 //
00027 // Eight node Brick element
00028 //
00029 
00030 #include <iostream.h>
00031 #include <stdio.h> 
00032 #include <stdlib.h> 
00033 #include <math.h> 
00034 
00035 #include <ID.h> 
00036 #include <Vector.h>
00037 #include <Matrix.h>
00038 #include <Element.h>
00039 #include <Node.h>
00040 #include <SectionForceDeformation.h>
00041 #include <Domain.h>
00042 #include <ErrorHandler.h>
00043 #include <Brick.h>
00044 
00045 #include <Renderer.h>
00046 
00047 extern void  shp3d( const double ss[3],
00048       double &xsj,
00049       double shp[4][8],
00050       const double xl[3][8]  ) ;
00051 
00052 //static data
00053 double  Brick::xl[3][8] ;
00054 
00055 Matrix  Brick::stiff(24,24) ;
00056 Vector  Brick::resid(24) ;
00057 Matrix  Brick::mass(24,24) ;
00058 Matrix  Brick::damping(24,24) ;
00059 
00060     
00061 //quadrature data
00062 const double  Brick::root3 = sqrt(3.0) ;
00063 const double  Brick::one_over_root3 = 1.0 / root3 ;
00064 
00065 const double  Brick::sg[] = { -one_over_root3,  
00066           one_over_root3  } ;
00067 
00068 const double  Brick::wg[] = { 1.0, 1.0, 1.0, 1.0, 
00069                               1.0, 1.0, 1.0, 1.0  } ;
00070 
00071   
00072 
00073 //null constructor
00074 Brick::Brick( ) :
00075 Element( 0, ELE_TAG_Brick ),
00076 connectedExternalNodes(8) 
00077 { 
00078 
00079 }
00080 
00081 
00082 //*********************************************************************
00083 //full constructor
00084 Brick::Brick(  int tag, 
00085                          int node1,
00086                          int node2,
00087                      int node3,
00088                          int node4,
00089                          int node5,
00090                          int node6,
00091                          int node7,
00092     int node8,
00093     NDMaterial &theMaterial ) :
00094 Element( tag, ELE_TAG_Brick ),
00095 connectedExternalNodes(8) 
00096 {
00097   connectedExternalNodes(0) = node1 ;
00098   connectedExternalNodes(1) = node2 ;
00099   connectedExternalNodes(2) = node3 ;
00100   connectedExternalNodes(3) = node4 ;
00101 
00102   connectedExternalNodes(4) = node5 ;
00103   connectedExternalNodes(5) = node6 ;
00104   connectedExternalNodes(6) = node7 ;
00105   connectedExternalNodes(7) = node8 ;
00106 
00107   int i ;
00108   for ( i=0; i<8; i++ ) {
00109 
00110       materialPointers[i] = theMaterial.getCopy("ThreeDimensional") ;
00111 
00112       if (materialPointers[i] == 0) {
00113 
00114    g3ErrorHandler->fatal("Brick::constructor %s",
00115   "- failed to get a material of type: ShellSection");
00116       } //end if
00117       
00118   } //end for i 
00119 
00120 
00121 }
00122 //******************************************************************
00123 
00124 
00125 //destructor 
00126 Brick::~Brick( )
00127 {
00128   int i ;
00129   for ( i=0 ; i<8; i++ ) {
00130 
00131     delete materialPointers[i] ;
00132     materialPointers[i] = 0 ; 
00133 
00134     nodePointers[i] = 0 ;
00135 
00136   } //end for i
00137 }
00138 
00139 
00140 //set domain
00141 void  Brick::setDomain( Domain *theDomain ) 
00142 {  
00143 
00144   int i ;
00145 
00146   //node pointers
00147   for ( i=0; i<8; i++ ) {
00148      nodePointers[i] = theDomain->getNode( connectedExternalNodes(i) ) ;
00149      if (nodePointers[i] == 0) {
00150        
00151        cerr << "Brick::setDomain - no node " << connectedExternalNodes(i);
00152        cerr << " exists in the domain - segmentation fault likely\n";
00153      }
00154   }
00155 
00156   this->DomainComponent::setDomain(theDomain);
00157 
00158 }
00159 
00160 
00161 //get the number of external nodes
00162 int  Brick::getNumExternalNodes( ) const
00163 {
00164   return 8 ;
00165 } 
00166  
00167 
00168 //return connected external nodes
00169 const ID&  Brick::getExternalNodes( ) 
00170 {
00171   return connectedExternalNodes ;
00172 } 
00173 
00174 
00175 //return number of dofs
00176 int  Brick::getNumDOF( ) 
00177 {
00178   return 24 ;
00179 }
00180 
00181 
00182 //commit state
00183 int  Brick::commitState( )
00184 {
00185   int i ;
00186   int success = 0 ;
00187 
00188   for ( i=0; i<8; i++ ) 
00189     success += materialPointers[i]->commitState( ) ;
00190   
00191   return success ;
00192 }
00193  
00194 
00195 
00196 //revert to last commit 
00197 int  Brick::revertToLastCommit( ) 
00198 {
00199   int i ;
00200   int success = 0 ;
00201 
00202   for ( i=0; i<8; i++ ) 
00203     success += materialPointers[i]->revertToLastCommit( ) ;
00204   
00205   return success ;
00206 }
00207     
00208 
00209 //revert to start 
00210 int  Brick::revertToStart( ) 
00211 {
00212   int i ;
00213   int success = 0 ;
00214 
00215   for ( i=0; i<8; i++ ) 
00216     success += materialPointers[i]->revertToStart( ) ;
00217   
00218   return success ;
00219 }
00220 
00221 //print out element data
00222 void  Brick::Print( ostream &s, int flag )
00223 {
00224   s << '\n' ;
00225   s << "Standard Eight Node Brick  tag: " << this->getTag() << "\n" ;
00226   s << "Node 1 : " << connectedExternalNodes(0) << '\n' ;
00227   s << "Node 2 : " << connectedExternalNodes(1) << '\n' ;
00228   s << "Node 3 : " << connectedExternalNodes(2) << '\n' ;
00229   s << "Node 4 : " << connectedExternalNodes(3) << '\n' ;
00230   s << "Node 5 : " << connectedExternalNodes(4) << '\n' ;
00231   s << "Node 6 : " << connectedExternalNodes(5) << '\n' ;
00232   s << "Node 7 : " << connectedExternalNodes(6) << '\n' ;
00233   s << "Node 8 : " << connectedExternalNodes(7) << '\n' ;
00234 
00235   s << "Material Information : \n " ;
00236   materialPointers[0]->Print( s, flag ) ;
00237 
00238   s << endl ;
00239 }
00240 
00241 //return stiffness matrix 
00242 const Matrix&  Brick::getTangentStiff( ) 
00243 {
00244   int tang_flag = 1 ; //get the tangent 
00245 
00246   //do tangent and residual here
00247   formResidAndTangent( tang_flag ) ;  
00248 
00249   return stiff ;
00250 }    
00251 
00252 
00253 //return secant matrix 
00254 const Matrix&  Brick::getSecantStiff( ) 
00255 {
00256    int tang_flag = 1 ; //get the tangent
00257 
00258   //do tangent and residual here
00259   formResidAndTangent( tang_flag ) ;  
00260     
00261   return stiff ;
00262 }
00263     
00264 
00265 //return damping matrix because frank is a dumb ass 
00266 const Matrix&  Brick::getDamp( ) 
00267 {
00268   //not supported
00269   return damping ;
00270 }    
00271 
00272 
00273 //return mass matrix
00274 const Matrix&  Brick::getMass( ) 
00275 {
00276   int tangFlag = 1 ;
00277 
00278   formInertiaTerms( tangFlag ) ;
00279 
00280   return mass ;
00281 } 
00282 
00283 
00284 //zero the load -- what load?
00285 void  Brick::zeroLoad( )
00286 {
00287   return ;
00288 }
00289 
00290 //add load -- what load?
00291 int  Brick::addLoad( const Vector &addP )
00292 {
00293   return -1 ;
00294 }
00295 
00296 
00297 //get residual
00298 const Vector&  Brick::getResistingForce( ) 
00299 {
00300   int tang_flag = 0 ; //don't get the tangent
00301 
00302   formResidAndTangent( tang_flag ) ;
00303 
00304   return resid ;   
00305 }
00306 
00307 
00308 //get residual with inertia terms
00309 const Vector&  Brick::getResistingForceIncInertia( )
00310 {
00311   int tang_flag = 0 ; //don't get the tangent
00312 
00313   //do tangent and residual here 
00314   formResidAndTangent( tang_flag ) ;
00315 
00316   formInertiaTerms( tang_flag ) ;
00317 
00318   return resid ;
00319 }
00320 
00321 
00322 //*********************************************************************
00323 //form inertia terms
00324 
00325 void   Brick::formInertiaTerms( int tangFlag ) 
00326 {
00327 
00328   static const int ndm = 3 ;
00329 
00330   static const int ndf = 3 ; 
00331 
00332   static const int numberNodes = 8 ;
00333 
00334   static const int numberGauss = 8 ;
00335 
00336   static const int nShape = 4 ;
00337 
00338   static const int massIndex = nShape - 1 ;
00339 
00340   double xsj ;  // determinant jacaobian matrix 
00341 
00342   double dvol[numberGauss] ; //volume element
00343 
00344   static double shp[nShape][numberNodes] ;  //shape functions at a gauss point
00345 
00346   static double Shape[nShape][numberNodes][numberGauss] ; //all the shape functions
00347 
00348   static double gaussPoint[ndm] ;
00349 
00350   static Vector momentum(ndf) ;
00351 
00352   int i, j, k, p, q ;
00353   int jj, kk ;
00354 
00355   double temp, rho, massJK ;
00356 
00357 
00358   //zero mass 
00359   mass.Zero( ) ;
00360 
00361   //compute basis vectors and local nodal coordinates
00362   computeBasis( ) ;
00363 
00364   //gauss loop to compute and save shape functions 
00365 
00366   int count = 0 ;
00367 
00368   for ( i = 0; i < 2; i++ ) {
00369     for ( j = 0; j < 2; j++ ) {
00370       for ( k = 0; k < 2; k++ ) {
00371 
00372         gaussPoint[0] = sg[i] ;        
00373  gaussPoint[1] = sg[j] ;        
00374  gaussPoint[2] = sg[k] ;
00375 
00376  //get shape functions    
00377  shp3d( gaussPoint, xsj, shp, xl ) ;
00378 
00379  //save shape functions
00380  for ( p = 0; p < nShape; p++ ) {
00381    for ( q = 0; q < numberNodes; q++ )
00382      Shape[p][q][count] = shp[p][q] ;
00383  } // end for p
00384 
00385 
00386  //volume element to also be saved
00387  dvol[count] = wg[count] * xsj ;  
00388 
00389  count++ ;
00390 
00391       } //end for k
00392     } //end for j
00393   } // end for i 
00394   
00395 
00396 
00397   //gauss loop 
00398   for ( i = 0; i < numberGauss; i++ ) {
00399 
00400     //extract shape functions from saved array
00401     for ( p = 0; p < nShape; p++ ) {
00402        for ( q = 0; q < numberNodes; q++ )
00403    shp[p][q]  = Shape[p][q][i] ;
00404     } // end for p
00405 
00406 
00407     //node loop to compute acceleration
00408     momentum.Zero( ) ;
00409     for ( j = 0; j < numberNodes; j++ ) 
00410       momentum += shp[massIndex][j] * ( nodePointers[j]->getTrialAccel()  ) ; 
00411 
00412     //density
00413     rho = materialPointers[i]->getRho() ;
00414 
00415 
00416     //multiply acceleration by density to form momentum
00417     momentum *= rho ;
00418 
00419 
00420     //residual and tangent calculations node loops
00421     jj = 0 ;
00422     for ( j = 0; j < numberNodes; j++ ) {
00423 
00424       temp = shp[massIndex][j] * dvol[i] ;
00425 
00426       for ( p = 0; p < ndf; p++ )
00427         resid( jj+p ) += ( temp * momentum(p) )  ;
00428 
00429       
00430       if ( tangFlag == 1 ) {
00431 
00432   //multiply by density
00433   temp *= rho ;
00434 
00435   //node-node mass
00436          kk = 0 ;
00437          for ( k = 0; k < numberNodes; k++ ) {
00438 
00439      massJK = temp * shp[massIndex][k] ;
00440 
00441             for ( p = 0; p < ndf; p++ )  
00442        mass( jj+p, kk+p ) += massJK ;
00443             
00444             kk += ndf ;
00445           } // end for k loop
00446 
00447       } // end if tang_flag 
00448 
00449       jj += ndf ;
00450     } // end for j loop
00451 
00452 
00453   } //end for i gauss loop 
00454 
00455 }
00456 
00457 //*********************************************************************
00458 //form residual and tangent
00459 void  Brick::formResidAndTangent( int tang_flag ) 
00460 {
00461 
00462   //strains ordered : eps11, eps22, eps33, 2*eps12, 2*eps23, 2*eps31 
00463 
00464   static const int ndm = 3 ;
00465 
00466   static const int ndf = 3 ; 
00467 
00468   static const int nstress = 6 ;
00469  
00470   static const int numberNodes = 8 ;
00471 
00472   static const int numberGauss = 8 ;
00473 
00474   static const int nShape = 4 ;
00475 
00476   int i, j, k, p, q ;
00477   int jj, kk ;
00478 
00479   int success ;
00480   
00481   static double volume ;
00482 
00483   static double xsj ;  // determinant jacaobian matrix 
00484 
00485   static double dvol[numberGauss] ; //volume element
00486 
00487   static double gaussPoint[ndm] ;
00488 
00489   static Vector strain(nstress) ;  //strain
00490 
00491   static double shp[nShape][numberNodes] ;  //shape functions at a gauss point
00492 
00493   static double Shape[nShape][numberNodes][numberGauss] ; //all the shape functions
00494 
00495   static Vector residJ(ndf) ; //nodeJ residual 
00496 
00497   static Matrix stiffJK(ndf,ndf) ; //nodeJK stiffness 
00498 
00499   static Vector stress(nstress) ;  //stress
00500 
00501   static Matrix dd(nstress,nstress) ;  //material tangent
00502 
00503 
00504   //---------B-matrices------------------------------------
00505 
00506     static Matrix BJ(nstress,ndf) ;      // B matrix node J
00507 
00508     static Matrix BJtran(ndf,nstress) ;
00509 
00510     static Matrix BK(nstress,ndf) ;      // B matrix node k
00511 
00512     static Matrix BJtranD(ndf,nstress) ;
00513 
00514   //-------------------------------------------------------
00515 
00516   
00517   //zero stiffness and residual 
00518   stiff.Zero( ) ;
00519   resid.Zero( ) ;
00520 
00521   //compute basis vectors and local nodal coordinates
00522   computeBasis( ) ;
00523 
00524   //gauss loop to compute and save shape functions 
00525 
00526   int count = 0 ;
00527   volume = 0.0 ;
00528 
00529   for ( i = 0; i < 2; i++ ) {
00530     for ( j = 0; j < 2; j++ ) {
00531       for ( k = 0; k < 2; k++ ) {
00532 
00533         gaussPoint[0] = sg[i] ;        
00534  gaussPoint[1] = sg[j] ;        
00535  gaussPoint[2] = sg[k] ;
00536 
00537  //get shape functions    
00538  shp3d( gaussPoint, xsj, shp, xl ) ;
00539 
00540  //save shape functions
00541  for ( p = 0; p < nShape; p++ ) {
00542    for ( q = 0; q < numberNodes; q++ )
00543      Shape[p][q][count] = shp[p][q] ;
00544  } // end for p
00545 
00546 
00547  //volume element to also be saved
00548  dvol[count] = wg[count] * xsj ;  
00549 
00550  //volume += dvol[count] ;
00551 
00552  count++ ;
00553 
00554       } //end for k
00555     } //end for j
00556   } // end for i 
00557   
00558 
00559   //gauss loop 
00560   for ( i = 0; i < numberGauss; i++ ) {
00561 
00562     //extract shape functions from saved array
00563     for ( p = 0; p < nShape; p++ ) {
00564        for ( q = 0; q < numberNodes; q++ )
00565    shp[p][q]  = Shape[p][q][i] ;
00566     } // end for p
00567 
00568 
00569     //zero the strains
00570     strain.Zero( ) ;
00571 
00572 
00573     // j-node loop to compute strain 
00574     for ( j = 0; j < numberNodes; j++ )  {
00575 
00576       //compute B matrix 
00577 
00578       BJ = computeB( j, shp ) ;
00579       
00580       //nodal displacements 
00581       const Vector &ul = nodePointers[j]->getTrialDisp( ) ;
00582 
00583       //compute the strain
00584       strain += (BJ*ul) ; 
00585 
00586     } // end for j
00587   
00588 
00589 
00590     //send the strain to the material 
00591     success = materialPointers[i]->setTrialStrain( strain ) ;
00592 
00593     //compute the stress
00594     stress = materialPointers[i]->getStress( ) ;
00595 
00596 
00597     //multiply by volume element
00598     stress  *= dvol[i] ;
00599 
00600     if ( tang_flag == 1 ) {
00601       dd = materialPointers[i]->getTangent( ) ;
00602       dd *= dvol[i] ;
00603     } //end if tang_flag
00604 
00605 
00606     //residual and tangent calculations node loops
00607 
00608     jj = 0 ;
00609     for ( j = 0; j < numberNodes; j++ ) {
00610 
00611       BJ = computeB( j, shp ) ;
00612    
00613       //transpose 
00614       BJtran = transpose( nstress, ndf, BJ ) ;
00615 
00616 
00617       //residual
00618       residJ = BJtran * stress ;
00619 
00620       //residual 
00621       for ( p = 0; p < ndf; p++ )
00622         resid( jj + p ) += residJ(p)  ;
00623 
00624 
00625       if ( tang_flag == 1 ) {
00626 
00627          BJtranD = BJtran * dd ;
00628 
00629          kk = 0 ;
00630          for ( k = 0; k < numberNodes; k++ ) {
00631 
00632             BK = computeB( k, shp ) ;
00633   
00634  
00635             stiffJK =  BJtranD * BK  ;
00636 
00637             for ( p = 0; p < ndf; p++ )  {
00638                for ( q = 0; q < ndf; q++ )
00639                   stiff( jj+p, kk+q ) += stiffJK( p, q ) ;
00640             } //end for p
00641 
00642             kk += ndf ;
00643           } // end for k loop
00644 
00645       } // end if tang_flag 
00646 
00647       jj += ndf ;
00648     } // end for j loop
00649 
00650 
00651   } //end for i gauss loop 
00652 
00653   
00654   return ;
00655 }
00656 
00657 
00658 //************************************************************************
00659 //compute local coordinates and basis
00660 
00661 void   Brick::computeBasis( ) 
00662 {
00663 
00664   //nodal coordinates 
00665 
00666   int i ;
00667   for ( i = 0; i < 8; i++ ) {
00668 
00669        const Vector &coorI = nodePointers[i]->getCrds( ) ;
00670 
00671        xl[0][i] = coorI(0) ;
00672        xl[1][i] = coorI(1) ;
00673        xl[2][i] = coorI(2) ;
00674 
00675   }  //end for i 
00676 
00677 }
00678 
00679 //*************************************************************************
00680 //compute B
00681 
00682 Matrix   Brick::computeB( int node, const double shp[4][8] )
00683 {
00684 
00685   static Matrix B(6,3) ;
00686 
00687 //---B Matrix in standard {1,2,3} mechanics notation---------
00688 //
00689 //                -                   -
00690 //               | N,1      0     0    | 
00691 //   B       =   |   0     N,2    0    |
00692 //               |   0      0     N,3  |   (6x3)
00693 //               | N,2     N,1     0   |
00694 //               |   0     N,3    N,2  |
00695 //               | N,3      0     N,1  |
00696 //                -                   -       
00697 //
00698 //-------------------------------------------------------------------
00699 
00700 
00701   B.Zero( ) ;
00702 
00703   B(0,0) = shp[0][node] ;
00704   B(1,1) = shp[1][node] ;
00705   B(2,2) = shp[2][node] ;
00706 
00707   B(3,0) = shp[1][node] ;
00708   B(3,1) = shp[0][node] ;
00709 
00710   B(4,1) = shp[2][node] ;
00711   B(4,2) = shp[1][node] ;
00712 
00713   B(5,0) = shp[2][node] ;
00714   B(5,2) = shp[0][node] ;
00715 
00716   return B ;
00717 
00718 }
00719 
00720 //***********************************************************************
00721 
00722 Matrix  Brick::transpose( int dim1, 
00723                                        int dim2, 
00724                          const Matrix &M ) 
00725 {
00726   int i ;
00727   int j ;
00728 
00729   Matrix Mtran( dim2, dim1 ) ;
00730 
00731   for ( i = 0; i < dim1; i++ ) {
00732      for ( j = 0; j < dim2; j++ ) 
00733          Mtran(j,i) = M(i,j) ;
00734   } // end for i
00735 
00736   return Mtran ;
00737 }
00738 
00739 //**********************************************************************
00740 
00741 int  Brick::sendSelf (int commitTag, Channel &theChannel)
00742 {
00743     return -1;
00744 }
00745     
00746 int  Brick::recvSelf (int commitTag, 
00747          Channel &theChannel, 
00748          FEM_ObjectBroker &theBroker)
00749 {
00750     return -1;
00751 }
00752 //**************************************************************************
00753 
00754 int
00755 Brick::displaySelf(Renderer &theViewer, int displayMode, float fact)
00756 {
00757     const Vector &end1Crd = nodePointers[0]->getCrds();
00758     const Vector &end2Crd = nodePointers[1]->getCrds(); 
00759     const Vector &end3Crd = nodePointers[2]->getCrds(); 
00760     const Vector &end4Crd = nodePointers[3]->getCrds(); 
00761 
00762     const Vector &end5Crd = nodePointers[4]->getCrds();
00763     const Vector &end6Crd = nodePointers[5]->getCrds(); 
00764     const Vector &end7Crd = nodePointers[6]->getCrds(); 
00765     const Vector &end8Crd = nodePointers[7]->getCrds(); 
00766 
00767     const Vector &end1Disp = nodePointers[0]->getDisp();
00768     const Vector &end2Disp = nodePointers[1]->getDisp();
00769     const Vector &end3Disp = nodePointers[2]->getDisp();
00770     const Vector &end4Disp = nodePointers[3]->getDisp();
00771 
00772     const Vector &end5Disp = nodePointers[4]->getDisp();
00773     const Vector &end6Disp = nodePointers[5]->getDisp();
00774     const Vector &end7Disp = nodePointers[6]->getDisp();
00775     const Vector &end8Disp = nodePointers[7]->getDisp();
00776 
00777     static Matrix coords(4,3);
00778     static Vector values(4);
00779     static Vector P(24) ;
00780 
00781     values(0) = 1 ;
00782     values(1) = 1 ;
00783     values(2) = 1 ;
00784     values(3) = 1 ;
00785 
00786     if (displayMode < 3 && displayMode > 0)
00787       P = this->getResistingForce();
00788 
00789     int error = 0;
00790 
00791     int i;
00792     for (i = 0; i < 3; i++) {
00793       coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
00794       coords(1,i) = end2Crd(i) + end2Disp(i)*fact;    
00795       coords(2,i) = end3Crd(i) + end3Disp(i)*fact;    
00796       coords(3,i) = end4Crd(i) + end4Disp(i)*fact;
00797     }
00798     error += theViewer.drawPolygon (coords, values);
00799 
00800     for (i = 0; i < 3; i++) {
00801       coords(0,i) = end5Crd(i) + end5Disp(i)*fact;
00802       coords(1,i) = end6Crd(i) + end6Disp(i)*fact;
00803       coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
00804       coords(3,i) = end8Crd(i) + end8Disp(i)*fact;
00805     }
00806     error += theViewer.drawPolygon (coords, values);
00807 
00808     for (i = 0; i < 3; i++) {
00809       coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
00810       coords(1,i) = end4Crd(i) + end4Disp(i)*fact;
00811       coords(2,i) = end8Crd(i) + end8Disp(i)*fact;
00812       coords(3,i) = end5Crd(i) + end5Disp(i)*fact;
00813     }
00814     error += theViewer.drawPolygon (coords, values);
00815 
00816     for (i = 0; i < 3; i++) {
00817       coords(0,i) = end2Crd(i) + end2Disp(i)*fact;
00818       coords(1,i) = end3Crd(i) + end3Disp(i)*fact;
00819       coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
00820       coords(3,i) = end6Crd(i) + end6Disp(i)*fact;
00821     }
00822     error += theViewer.drawPolygon (coords, values);
00823 
00824     for (i = 0; i < 3; i++) {
00825       coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
00826       coords(1,i) = end2Crd(i) + end2Disp(i)*fact;
00827       coords(2,i) = end6Crd(i) + end6Disp(i)*fact;
00828       coords(3,i) = end5Crd(i) + end5Disp(i)*fact;
00829     }
00830     error += theViewer.drawPolygon (coords, values);
00831 
00832     for (i = 0; i < 3; i++) {
00833       coords(0,i) = end4Crd(i) + end4Disp(i)*fact;
00834       coords(1,i) = end3Crd(i) + end3Disp(i)*fact;
00835       coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
00836       coords(3,i) = end8Crd(i) + end8Disp(i)*fact;
00837     }
00838     error += theViewer.drawPolygon (coords, values);
00839 
00840     return error;
00841 }
Copyright Contact Us