BbarBrick.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.18 $
00022 // $Date: 2006/08/04 21:25:56 $
00023 // $Source: /usr/local/cvs/OpenSees/SRC/element/brick/BbarBrick.cpp,v $
00024 
00025 // Ed "C++" Love
00026 //
00027 // Eight node BbarBrick element
00028 //
00029 
00030 #include <stdio.h> 
00031 #include <stdlib.h> 
00032 #include <math.h> 
00033 
00034 #include <ID.h> 
00035 #include <Vector.h>
00036 #include <Matrix.h>
00037 #include <Element.h>
00038 #include <Node.h>
00039 #include <Domain.h>
00040 #include <ErrorHandler.h>
00041 #include <BbarBrick.h>
00042 #include <shp3d.h>
00043 #include <Renderer.h>
00044 #include <ElementResponse.h>
00045 
00046 
00047 #include <Channel.h>
00048 #include <FEM_ObjectBroker.h>
00049 
00050 //static data
00051 double  BbarBrick::xl[3][8] ;
00052 
00053 Matrix  BbarBrick::stiff(24,24) ;
00054 Vector  BbarBrick::resid(24) ;
00055 Matrix  BbarBrick::mass(24,24) ;
00056 
00057     
00058 //quadrature data
00059 const double  BbarBrick::root3 = sqrt(3.0) ;
00060 const double  BbarBrick::one_over_root3 = 1.0 / root3 ;
00061 
00062 const double  BbarBrick::sg[] = { -one_over_root3,  
00063                                one_over_root3  } ;
00064 
00065 const double  BbarBrick::wg[] = { 1.0, 1.0, 1.0, 1.0, 
00066                               1.0, 1.0, 1.0, 1.0  } ;
00067 
00068   
00069 
00070 //null constructor
00071 BbarBrick::BbarBrick( ) :
00072 Element( 0, ELE_TAG_BbarBrick ),
00073 connectedExternalNodes(8), load(0), Ki(0)
00074 { 
00075   for (int i=0; i<8; i++ ) {
00076     materialPointers[i] = 0;
00077     nodePointers[i] = 0;
00078   }
00079   b[0] = 0.0;
00080   b[1] = 0.0; 
00081   b[2] = 0.0;
00082 }
00083 
00084 
00085 //*********************************************************************
00086 //full constructor
00087 BbarBrick::BbarBrick(  int tag, 
00088                          int node1,
00089                          int node2,
00090                          int node3,
00091                          int node4,
00092                          int node5,
00093                          int node6,
00094                          int node7,
00095                          int node8,
00096                          NDMaterial &theMaterial,
00097                          double b1, double b2, double b3) :
00098 Element( tag, ELE_TAG_BbarBrick ),
00099 connectedExternalNodes(8), load(0), Ki(0)
00100 {
00101   connectedExternalNodes(0) = node1 ;
00102   connectedExternalNodes(1) = node2 ;
00103   connectedExternalNodes(2) = node3 ;
00104   connectedExternalNodes(3) = node4 ;
00105 
00106   connectedExternalNodes(4) = node5 ;
00107   connectedExternalNodes(5) = node6 ;
00108   connectedExternalNodes(6) = node7 ;
00109   connectedExternalNodes(7) = node8 ;
00110 
00111   int i ;
00112   for ( i=0; i<8; i++ ) {
00113 
00114       materialPointers[i] = theMaterial.getCopy("ThreeDimensional") ;
00115 
00116       if (materialPointers[i] == 0) {
00117           opserr <<"BbarBrick::constructor - failed to get a material of type: ThreeDimensional\n";
00118           exit(-1);
00119       } //end if
00120       
00121   } //end for i 
00122   
00123   // Body forces
00124   b[0] = b1;
00125   b[1] = b2;
00126   b[2] = b3;
00127 }
00128 //******************************************************************
00129 
00130 
00131 //destructor 
00132 BbarBrick::~BbarBrick( )
00133 {
00134   int i ;
00135   for ( i=0 ; i<8; i++ ) {
00136 
00137     delete materialPointers[i] ;
00138     materialPointers[i] = 0 ; 
00139 
00140     nodePointers[i] = 0 ;
00141 
00142   } //end for i
00143   
00144   if (load != 0)
00145     delete load;
00146   
00147   if (Ki != 0)
00148     delete Ki;
00149 }
00150 
00151 
00152 //set domain
00153 void  BbarBrick::setDomain( Domain *theDomain ) 
00154 {  
00155 
00156   int i ;
00157 
00158   //node pointers
00159   for ( i=0; i<8; i++ ) 
00160      nodePointers[i] = theDomain->getNode( connectedExternalNodes(i) ) ;
00161 
00162   this->DomainComponent::setDomain(theDomain);
00163 
00164 }
00165 
00166 
00167 //get the number of external nodes
00168 int  BbarBrick::getNumExternalNodes( ) const
00169 {
00170   return 8 ;
00171 } 
00172  
00173 
00174 //return connected external nodes
00175 const ID&  BbarBrick::getExternalNodes( ) 
00176 {
00177   return connectedExternalNodes ;
00178 } 
00179 
00180 //return connected external node
00181 Node **  
00182 BbarBrick::getNodePtrs(void) 
00183 {
00184   return nodePointers ;
00185 } 
00186 
00187 
00188 //return number of dofs
00189 int  BbarBrick::getNumDOF( ) 
00190 {
00191   return 24 ;
00192 }
00193 
00194 
00195 //commit state
00196 int  BbarBrick::commitState( )
00197 {
00198   int success = 0 ;
00199 
00200   // call element commitState to do any base class stuff
00201   if ((success = this->Element::commitState()) != 0) {
00202     opserr << "Brick::commitState () - failed in base class";
00203   }    
00204 
00205   for (int i=0; i<8; i++ ) 
00206     success += materialPointers[i]->commitState( ) ;
00207   
00208   return success ;
00209 }
00210  
00211 
00212 
00213 //revert to last commit 
00214 int  BbarBrick::revertToLastCommit( ) 
00215 {
00216   int i ;
00217   int success = 0 ;
00218 
00219   for ( i=0; i<8; i++ ) 
00220     success += materialPointers[i]->revertToLastCommit( ) ;
00221   
00222   return success ;
00223 }
00224     
00225 
00226 //revert to start 
00227 int  BbarBrick::revertToStart( ) 
00228 {
00229   int i ;
00230   int success = 0 ;
00231 
00232   for ( i=0; i<8; i++ ) 
00233     success += materialPointers[i]->revertToStart( ) ;
00234   
00235   return success ;
00236 }
00237 
00238 //print out element data
00239 void  BbarBrick::Print( OPS_Stream &s, int flag )
00240 {
00241   s << endln ;
00242   s << "Volume/Pressure Eight Node BbarBrick \n" ;
00243   s << "Element Number: " << this->getTag() << endln ;
00244   s << "Node 1 : " << connectedExternalNodes(0) << endln ;
00245   s << "Node 2 : " << connectedExternalNodes(1) << endln ;
00246   s << "Node 3 : " << connectedExternalNodes(2) << endln ;
00247   s << "Node 4 : " << connectedExternalNodes(3) << endln ;
00248   s << "Node 5 : " << connectedExternalNodes(4) << endln ;
00249   s << "Node 6 : " << connectedExternalNodes(5) << endln ;
00250   s << "Node 7 : " << connectedExternalNodes(6) << endln ;
00251   s << "Node 8 : " << connectedExternalNodes(7) << endln ;
00252 
00253   s << "Material Information : \n " ;
00254   materialPointers[0]->Print( s, flag ) ;
00255 
00256   s << endln ;
00257 }
00258 
00259 //return stiffness matrix 
00260 const Matrix&  BbarBrick::getTangentStiff( ) 
00261 {
00262   int tang_flag = 1 ; //get the tangent 
00263 
00264   //do tangent and residual here
00265   formResidAndTangent( tang_flag ) ;  
00266 
00267   return stiff ;
00268 }    
00269 
00270 
00271 //return stiffness matrix 
00272 const Matrix&  BbarBrick::getInitialStiff( ) 
00273 {
00274   if (Ki != 0)
00275     return *Ki;
00276 
00277   //strains ordered : eps11, eps22, eps33, 2*eps12, 2*eps23, 2*eps31 
00278 
00279   static const int ndm = 3 ;
00280 
00281   static const int ndf = 3 ; 
00282 
00283   static const int nstress = 6 ;
00284  
00285   static const int numberNodes = 8 ;
00286 
00287   static const int numberGauss = 8 ;
00288 
00289   static const int nShape = 4 ;
00290 
00291   int i, j, k, p, q ;
00292   int jj, kk ;
00293 
00294   static double volume ;
00295 
00296   static double xsj ;  // determinant jacaobian matrix 
00297 
00298   static double dvol[numberGauss] ; //volume element
00299 
00300   static double gaussPoint[ndm] ;
00301 
00302   static Vector strain(nstress) ;  //strain
00303 
00304   static double shp[nShape][numberNodes] ;  //shape functions at a gauss point
00305 
00306   static double Shape[nShape][numberNodes][numberGauss] ; //all the shape functions
00307 
00308   static double shpBar[nShape][numberNodes] ;  //mean value of shape functions
00309 
00310   static Matrix stiffJK(ndf,ndf) ; //nodeJK stiffness 
00311 
00312   static Matrix dd(nstress,nstress) ;  //material tangent
00313 
00314 
00315   //---------B-matrices------------------------------------
00316 
00317     static Matrix BJ(nstress,ndf) ;      // B matrix node J
00318 
00319     static Matrix BJtran(ndf,nstress) ;
00320 
00321     static Matrix BK(nstress,ndf) ;      // B matrix node k
00322 
00323     static Matrix BJtranD(ndf,nstress) ;
00324 
00325   //-------------------------------------------------------
00326 
00327   
00328   //zero stiffness and residual 
00329   stiff.Zero( ) ;
00330 
00331 
00332   //compute basis vectors and local nodal coordinates
00333   computeBasis( ) ;
00334 
00335   //zero mean shape functions
00336   for ( p = 0; p < nShape; p++ ) {
00337     for ( q = 0; q < numberNodes; q++ )
00338       shpBar[p][q] = 0.0 ;
00339   } // end for p
00340 
00341   //zero volume
00342   volume = 0.0 ;
00343 
00344 
00345   //gauss loop to compute and save shape functions 
00346   int count = 0 ;
00347 
00348   for ( i = 0; i < 2; i++ ) {
00349     for ( j = 0; j < 2; j++ ) {
00350       for ( k = 0; k < 2; k++ ) {
00351 
00352         gaussPoint[0] = sg[i] ;        
00353         gaussPoint[1] = sg[j] ;        
00354         gaussPoint[2] = sg[k] ;
00355 
00356         //get shape functions    
00357         shp3d( gaussPoint, xsj, shp, xl ) ;
00358 
00359         //save shape functions
00360         for ( p = 0; p < nShape; p++ ) {
00361           for ( q = 0; q < numberNodes; q++ )
00362             Shape[p][q][count] = shp[p][q] ;
00363         } // end for p
00364 
00365         //volume element to also be saved
00366         dvol[count] = wg[count] * xsj ;  
00367 
00368         //add to volume
00369         volume += dvol[count] ;
00370 
00371         //add to mean shape functions
00372         for ( p = 0; p < nShape; p++ ) {
00373           for ( q = 0; q < numberNodes; q++ )
00374             shpBar[p][q] += ( dvol[count] * shp[p][q] ) ;
00375         } // end for p
00376 
00377         count++ ;
00378 
00379       } //end for k
00380     } //end for j
00381   } // end for i 
00382   
00383 
00384   //mean value of shape functions
00385   for ( p = 0; p < nShape; p++ ) {
00386     for ( q = 0; q < numberNodes; q++ )
00387       shpBar[p][q] /= volume ;
00388   } // end for p
00389 
00390 
00391   //gauss loop 
00392   for ( i = 0; i < numberGauss; i++ ) {
00393 
00394     //extract shape functions from saved array
00395     for ( p = 0; p < nShape; p++ ) {
00396        for ( q = 0; q < numberNodes; q++ )
00397           shp[p][q]  = Shape[p][q][i] ;
00398     } // end for p
00399 
00400     dd = materialPointers[i]->getInitialTangent( ) ;
00401     dd *= dvol[i] ;
00402 
00403     //residual and tangent calculations node loops
00404 
00405     jj = 0 ;
00406     for ( j = 0; j < numberNodes; j++ ) {
00407 
00408       BJ = computeBbar( j, shp, shpBar ) ;
00409    
00410       //transpose 
00411       //BJtran = transpose( nstress, ndf, BJ ) ;
00412       for (p=0; p<ndf; p++) {
00413         for (q=0; q<nstress; q++) 
00414           BJtran(p,q) = BJ(q,p) ;
00415       }//end for p
00416 
00417       //BJtranD = BJtran * dd ;
00418       BJtranD.addMatrixProduct(0.0,  BJtran,dd,1.0);
00419       
00420       kk = 0 ;
00421       for ( k = 0; k < numberNodes; k++ ) {
00422         
00423         BK = computeBbar( k, shp, shpBar ) ;
00424         
00425         //stiffJK =  BJtranD * BK  ;
00426         stiffJK.addMatrixProduct(0.0,  BJtranD,BK,1.0) ;
00427 
00428         for ( p = 0; p < ndf; p++ )  {
00429           for ( q = 0; q < ndf; q++ )
00430             stiff( jj+p, kk+q ) += stiffJK( p, q ) ;
00431         } //end for p
00432 
00433         kk += ndf ;
00434       } // end for k loop
00435 
00436       jj += ndf ;
00437 
00438     } // end for j loop
00439   } //end for i gauss loop 
00440 
00441   Ki = new Matrix(stiff);
00442 
00443   return stiff ;
00444 }    
00445 
00446 
00447 //return mass matrix
00448 const Matrix&  BbarBrick::getMass( ) 
00449 {
00450   int tangFlag = 1 ;
00451 
00452   formInertiaTerms( tangFlag ) ;
00453 
00454   return mass ;
00455 } 
00456 
00457 
00458 void  BbarBrick::zeroLoad( )
00459 {
00460   if (load != 0)
00461     load->Zero();
00462 
00463   return ;
00464 }
00465 
00466 int 
00467 BbarBrick::addLoad(ElementalLoad *theLoad, double loadFactor)
00468 {
00469   opserr << "BbarBrick::addLoad - load type unknown for ele with tag: " << this->getTag() << endln;
00470   return -1;
00471 }
00472 
00473 int
00474 BbarBrick::addInertiaLoadToUnbalance(const Vector &accel)
00475 {
00476   static const int numberNodes = 8 ;
00477   static const int numberGauss = 8 ;
00478   static const int ndf = 3 ; 
00479 
00480   int i;
00481 
00482   // check to see if have mass
00483   int haveRho = 0;
00484   for (i = 0; i < numberGauss; i++) {
00485     if (materialPointers[i]->getRho() != 0.0)
00486       haveRho = 1;
00487   }
00488 
00489   if (haveRho == 0)
00490     return 0;
00491 
00492   // Compute mass matrix
00493   int tangFlag = 1 ;
00494   formInertiaTerms( tangFlag ) ;
00495 
00496   // store computed RV fro nodes in resid vector
00497   int count = 0;
00498   for (i=0; i<numberNodes; i++) {
00499     const Vector &Raccel = nodePointers[i]->getRV(accel);
00500     for (int j=0; j<ndf; j++)
00501       resid(count++) = Raccel(j);
00502   }
00503 
00504   // create the load vector if one does not exist
00505   if (load == 0) 
00506     load = new Vector(numberNodes*ndf);
00507 
00508   // add -M * RV(accel) to the load vector
00509   load->addMatrixVector(1.0, mass, resid, -1.0);
00510   
00511   return 0;
00512 }
00513 
00514 
00515 //get residual
00516 const Vector&  BbarBrick::getResistingForce( ) 
00517 {
00518   int tang_flag = 0 ; //don't get the tangent
00519 
00520   formResidAndTangent( tang_flag ) ;
00521 
00522   if (load != 0)
00523     resid -= *load;
00524 
00525   return resid ;   
00526 }
00527 
00528 
00529 //get residual with inertia terms
00530 const Vector&  BbarBrick::getResistingForceIncInertia( )
00531 {
00532   static Vector res(24);
00533 
00534   int tang_flag = 0 ; //don't get the tangent
00535 
00536   //do tangent and residual here 
00537   formResidAndTangent( tang_flag ) ;
00538 
00539   formInertiaTerms( tang_flag ) ;
00540 
00541   res = resid;
00542 
00543   // add the damping forces if rayleigh damping
00544   if (alphaM != 0.0 || betaK != 0.0 || betaK0 != 0.0 || betaKc != 0.0)
00545       res += this->getRayleighDampingForces();
00546 
00547   if (load != 0)
00548     res -= *load;
00549 
00550   return res ;
00551 }
00552 
00553 
00554 //*********************************************************************
00555 //form inertia terms
00556 
00557 void   BbarBrick::formInertiaTerms( int tangFlag ) 
00558 {
00559 
00560   static const int ndm = 3 ;
00561 
00562   static const int ndf = 3 ; 
00563 
00564   static const int numberNodes = 8 ;
00565 
00566   static const int numberGauss = 8 ;
00567 
00568   static const int nShape = 4 ;
00569 
00570   static const int massIndex = nShape - 1 ;
00571 
00572   double xsj ;  // determinant jacaobian matrix 
00573 
00574   double dvol[numberGauss] ; //volume element
00575 
00576   static double shp[nShape][numberNodes] ;  //shape functions at a gauss point
00577 
00578   static double Shape[nShape][numberNodes][numberGauss] ; //all the shape functions
00579 
00580   static double gaussPoint[ndm] ;
00581 
00582   static Vector momentum(ndf) ;
00583 
00584   int i, j, k, p, q ;
00585   int jj, kk ;
00586 
00587   double temp, rho, massJK ;
00588 
00589 
00590   //zero mass 
00591   mass.Zero( ) ;
00592 
00593   //compute basis vectors and local nodal coordinates
00594   computeBasis( ) ;
00595 
00596   //gauss loop to compute and save shape functions 
00597 
00598   int count = 0 ;
00599 
00600   for ( i = 0; i < 2; i++ ) {
00601     for ( j = 0; j < 2; j++ ) {
00602       for ( k = 0; k < 2; k++ ) {
00603 
00604         gaussPoint[0] = sg[i] ;        
00605         gaussPoint[1] = sg[j] ;        
00606         gaussPoint[2] = sg[k] ;
00607 
00608         //get shape functions    
00609         shp3d( gaussPoint, xsj, shp, xl ) ;
00610 
00611         //save shape functions
00612         for ( p = 0; p < nShape; p++ ) {
00613           for ( q = 0; q < numberNodes; q++ )
00614             Shape[p][q][count] = shp[p][q] ;
00615         } // end for p
00616 
00617         //volume element to also be saved
00618         dvol[count] = wg[count] * xsj ;  
00619 
00620         count++ ;
00621 
00622       } //end for k
00623     } //end for j
00624   } // end for i 
00625   
00626 
00627 
00628   //gauss loop 
00629   for ( i = 0; i < numberGauss; i++ ) {
00630 
00631     //extract shape functions from saved array
00632     for ( p = 0; p < nShape; p++ ) {
00633        for ( q = 0; q < numberNodes; q++ )
00634           shp[p][q]  = Shape[p][q][i] ;
00635     } // end for p
00636 
00637 
00638     //node loop to compute acceleration
00639     momentum.Zero( ) ;
00640     for ( j = 0; j < numberNodes; j++ ) 
00641       //momentum += shp[massIndex][j] * ( nodePointers[j]->getTrialAccel()  ) ; 
00642       momentum.addVector( 1.0,
00643                           nodePointers[j]->getTrialAccel(),
00644                           shp[massIndex][j] ) ;
00645 
00646 
00647     //density
00648     rho = materialPointers[i]->getRho() ;
00649 
00650     //multiply acceleration by density to form momentum
00651     momentum *= rho ;
00652 
00653 
00654     //residual and tangent calculations node loops
00655     jj = 0 ;
00656     for ( j = 0; j < numberNodes; j++ ) {
00657 
00658       temp = shp[massIndex][j] * dvol[i] ;
00659 
00660       for ( p = 0; p < ndf; p++ )
00661         resid( jj+p ) += ( temp * momentum(p) )  ;
00662 
00663       
00664       if ( tangFlag == 1 ) {
00665 
00666          //multiply by density
00667          temp *= rho ;
00668 
00669          //node-node mass
00670          kk = 0 ;
00671          for ( k = 0; k < numberNodes; k++ ) {
00672 
00673             massJK = temp * shp[massIndex][k] ;
00674 
00675             for ( p = 0; p < ndf; p++ )  
00676               mass( jj+p, kk+p ) += massJK ;
00677             
00678             kk += ndf ;
00679           } // end for k loop
00680 
00681       } // end if tang_flag 
00682 
00683       jj += ndf ;
00684     } // end for j loop
00685 
00686 
00687   } //end for i gauss loop 
00688 
00689 }
00690 
00691 //*********************************************************************
00692 //form residual and tangent
00693 void  BbarBrick::formResidAndTangent( int tang_flag ) 
00694 {
00695 
00696   //strains ordered : eps11, eps22, eps33, 2*eps12, 2*eps23, 2*eps31 
00697 
00698   static const int ndm = 3 ;
00699 
00700   static const int ndf = 3 ; 
00701 
00702   static const int nstress = 6 ;
00703  
00704   static const int numberNodes = 8 ;
00705 
00706   static const int numberGauss = 8 ;
00707 
00708   static const int nShape = 4 ;
00709 
00710   int i, j, k, p, q ;
00711   int jj, kk ;
00712 
00713   int success ;
00714   
00715   static double volume ;
00716 
00717   static double xsj ;  // determinant jacaobian matrix 
00718 
00719   static double dvol[numberGauss] ; //volume element
00720 
00721   static double gaussPoint[ndm] ;
00722 
00723   static Vector strain(nstress) ;  //strain
00724 
00725   static double shp[nShape][numberNodes] ;  //shape functions at a gauss point
00726 
00727   static double Shape[nShape][numberNodes][numberGauss] ; //all the shape functions
00728 
00729   static double shpBar[nShape][numberNodes] ;  //mean value of shape functions
00730 
00731   static Vector residJ(ndf) ; //nodeJ residual 
00732 
00733   static Matrix stiffJK(ndf,ndf) ; //nodeJK stiffness 
00734 
00735   static Vector stress(nstress) ;  //stress
00736 
00737   static Matrix dd(nstress,nstress) ;  //material tangent
00738 
00739 
00740   //---------B-matrices------------------------------------
00741 
00742     static Matrix BJ(nstress,ndf) ;      // B matrix node J
00743 
00744     static Matrix BJtran(ndf,nstress) ;
00745 
00746     static Matrix BK(nstress,ndf) ;      // B matrix node k
00747 
00748     static Matrix BJtranD(ndf,nstress) ;
00749 
00750   //-------------------------------------------------------
00751 
00752   
00753   //zero stiffness and residual 
00754   stiff.Zero( ) ;
00755   resid.Zero( ) ;
00756 
00757   //compute basis vectors and local nodal coordinates
00758   computeBasis( ) ;
00759 
00760 
00761   //zero mean shape functions
00762   for ( p = 0; p < nShape; p++ ) {
00763     for ( q = 0; q < numberNodes; q++ )
00764       shpBar[p][q] = 0.0 ;
00765   } // end for p
00766 
00767   //zero volume
00768   volume = 0.0 ;
00769 
00770 
00771   //gauss loop to compute and save shape functions 
00772   int count = 0 ;
00773 
00774   for ( i = 0; i < 2; i++ ) {
00775     for ( j = 0; j < 2; j++ ) {
00776       for ( k = 0; k < 2; k++ ) {
00777 
00778         gaussPoint[0] = sg[i] ;        
00779         gaussPoint[1] = sg[j] ;        
00780         gaussPoint[2] = sg[k] ;
00781 
00782         //get shape functions    
00783         shp3d( gaussPoint, xsj, shp, xl ) ;
00784 
00785         //save shape functions
00786         for ( p = 0; p < nShape; p++ ) {
00787           for ( q = 0; q < numberNodes; q++ )
00788             Shape[p][q][count] = shp[p][q] ;
00789         } // end for p
00790 
00791         //volume element to also be saved
00792         dvol[count] = wg[count] * xsj ;  
00793 
00794         //add to volume
00795         volume += dvol[count] ;
00796 
00797         //add to mean shape functions
00798         for ( p = 0; p < nShape; p++ ) {
00799           for ( q = 0; q < numberNodes; q++ )
00800             shpBar[p][q] += ( dvol[count] * shp[p][q] ) ;
00801         } // end for p
00802 
00803         count++ ;
00804 
00805       } //end for k
00806     } //end for j
00807   } // end for i 
00808   
00809 
00810   //mean value of shape functions
00811   for ( p = 0; p < nShape; p++ ) {
00812     for ( q = 0; q < numberNodes; q++ )
00813       shpBar[p][q] /= volume ;
00814   } // end for p
00815 
00816 
00817   //gauss loop 
00818   for ( i = 0; i < numberGauss; i++ ) {
00819 
00820     //extract shape functions from saved array
00821     for ( p = 0; p < nShape; p++ ) {
00822        for ( q = 0; q < numberNodes; q++ )
00823           shp[p][q]  = Shape[p][q][i] ;
00824     } // end for p
00825 
00826 
00827     //zero the strains
00828     strain.Zero( ) ;
00829 
00830 
00831     // j-node loop to compute strain 
00832     for ( j = 0; j < numberNodes; j++ )  {
00833 
00834       //compute B matrix 
00835 
00836       BJ = computeBbar( j, shp, shpBar ) ;
00837       
00838       //nodal displacements 
00839       const Vector &ul = nodePointers[j]->getTrialDisp( ) ;
00840 
00841       //compute the strain
00842       //strain += (BJ*ul) ; 
00843       strain.addMatrixVector(1.0,  BJ,ul,1.0 ) ;
00844 
00845     } // end for j
00846   
00847 
00848 
00849     //send the strain to the material 
00850     success = materialPointers[i]->setTrialStrain( strain ) ;
00851 
00852     //compute the stress
00853     stress = materialPointers[i]->getStress( ) ;
00854 
00855 
00856     //multiply by volume element
00857     stress  *= dvol[i] ;
00858 
00859     if ( tang_flag == 1 ) {
00860       dd = materialPointers[i]->getTangent( ) ;
00861       dd *= dvol[i] ;
00862     } //end if tang_flag
00863 
00864 
00865     //residual and tangent calculations node loops
00866 
00867     jj = 0 ;
00868     for ( j = 0; j < numberNodes; j++ ) {
00869 
00870       BJ = computeBbar( j, shp, shpBar ) ;
00871    
00872       //transpose 
00873       //BJtran = transpose( nstress, ndf, BJ ) ;
00874       for (p=0; p<ndf; p++) {
00875         for (q=0; q<nstress; q++) 
00876           BJtran(p,q) = BJ(q,p) ;
00877       }//end for p
00878 
00879 
00880       //residual
00881       //residJ = BJtran * stress ;
00882       residJ.addMatrixVector(0.0,  BJtran,stress,1.0);
00883 
00884       //residual 
00885       for ( p = 0; p < ndf; p++ ) {
00886         resid( jj + p ) += residJ(p)  ;
00887             resid( jj + p ) -= dvol[i]*b[p]*shp[3][j];
00888       }
00889 
00890       if ( tang_flag == 1 ) {
00891 
00892         //BJtranD = BJtran * dd ;
00893         BJtranD.addMatrixProduct(0.0,  BJtran,dd,1.0);
00894 
00895          kk = 0 ;
00896          for ( k = 0; k < numberNodes; k++ ) {
00897 
00898             BK = computeBbar( k, shp, shpBar ) ;
00899   
00900  
00901             //stiffJK =  BJtranD * BK  ;
00902             stiffJK.addMatrixProduct(0.0,  BJtranD,BK,1.0) ;
00903 
00904             for ( p = 0; p < ndf; p++ )  {
00905                for ( q = 0; q < ndf; q++ )
00906                   stiff( jj+p, kk+q ) += stiffJK( p, q ) ;
00907             } //end for p
00908 
00909             kk += ndf ;
00910           } // end for k loop
00911 
00912       } // end if tang_flag 
00913 
00914       jj += ndf ;
00915     } // end for j loop
00916 
00917 
00918   } //end for i gauss loop 
00919 
00920   
00921   return ;
00922 }
00923 
00924 
00925 //************************************************************************
00926 //compute local coordinates and basis
00927 
00928 void   BbarBrick::computeBasis( ) 
00929 {
00930 
00931   //nodal coordinates 
00932 
00933   int i ;
00934   for ( i = 0; i < 8; i++ ) {
00935 
00936        const Vector &coorI = nodePointers[i]->getCrds( ) ;
00937 
00938        xl[0][i] = coorI(0) ;
00939        xl[1][i] = coorI(1) ;
00940        xl[2][i] = coorI(2) ;
00941 
00942   }  //end for i 
00943 
00944 }
00945 
00946 //*************************************************************************
00947 //compute B
00948 
00949 const Matrix&   
00950 BbarBrick::computeBbar( int node, 
00951                                  const double shp[4][8], 
00952                                  const double shpBar[4][8] )
00953 {
00954 
00955   static Matrix Bbar(6,3) ;
00956 
00957   //static Matrix Bdev(3,3) ;
00958   static double Bdev[3][3] ;
00959 
00960   //static Matrix BbarVol(3,3) ;
00961   static double BbarVol[3][3] ;
00962 
00963   static const double one3 = 1.0/3.0 ;
00964 
00965 
00966 //---B Matrices in standard {1,2,3} mechanics notation---------
00967 //
00968 //                -                        -
00969 //               |  2N,1    -N,2     -N,3   | 
00970 // Bdev =  (1/3) |  -N,1    2N,2     -N,3   |  (3x3)
00971 //               |  -N,1    -N,2     2N,3   |   
00972 //                -                        -
00973 //
00974 //                -                       -
00975 //               |  N,1      N,2     N,3   | 
00976 // Bvol =  (1/3) |  N,1      N,2     N.3   |  (3x3)
00977 //               |  N,1      N,2     N,3   |   
00978 //                -                       -
00979 //
00980 //                -                   -
00981 //               |                     |
00982 //               |    Bdev + Bvol      |
00983 //   B       =   |                     | 
00984 //               |---------------------|   (6x3)
00985 //               | N,2     N,1     0   |
00986 //               |   0     N,3    N,2  |
00987 //               | N,3      0     N,1  |
00988 //                -                   -       
00989 //
00990 //---------------------------------------------------------------
00991 
00992 
00993   Bbar.Zero( ) ;
00994 
00995   //deviatoric
00996   Bdev[0][0] = 2.0*shp[0][node] ;
00997   Bdev[0][1] =    -shp[1][node] ;
00998   Bdev[0][2] =    -shp[2][node] ;
00999 
01000   Bdev[1][0] =    -shp[0][node] ;
01001   Bdev[1][1] = 2.0*shp[1][node] ;
01002   Bdev[1][2] =    -shp[2][node] ;
01003 
01004   Bdev[2][0] =    -shp[0][node] ;
01005   Bdev[2][1] =    -shp[1][node] ;
01006   Bdev[2][2] = 2.0*shp[2][node] ;
01007 
01008   //volumetric 
01009   BbarVol[0][0] = shpBar[0][node] ;
01010   BbarVol[0][1] = shpBar[1][node] ;
01011   BbarVol[0][2] = shpBar[2][node] ;
01012 
01013   BbarVol[1][0] = shpBar[0][node] ;
01014   BbarVol[1][1] = shpBar[1][node] ;
01015   BbarVol[1][2] = shpBar[2][node] ;
01016 
01017   BbarVol[2][0] = shpBar[0][node] ;
01018   BbarVol[2][1] = shpBar[1][node] ;
01019   BbarVol[2][2] = shpBar[2][node] ;
01020 
01021 
01022 
01023   //extensional terms
01024   for ( int i=0; i<3; i++ ){
01025     for ( int j=0; j<3; j++ ) 
01026       Bbar(i,j) = one3*( Bdev[i][j] + BbarVol[i][j] ) ;
01027   }//end for i
01028 
01029 
01030   //shear terms
01031   Bbar(3,0) = shp[1][node] ;
01032   Bbar(3,1) = shp[0][node] ;
01033 
01034   Bbar(4,1) = shp[2][node] ;
01035   Bbar(4,2) = shp[1][node] ;
01036 
01037   Bbar(5,0) = shp[2][node] ;
01038   Bbar(5,2) = shp[0][node] ;
01039 
01040   return Bbar ;
01041 
01042 }
01043 
01044 //***********************************************************************
01045 
01046 Matrix  BbarBrick::transpose( int dim1, 
01047                                        int dim2, 
01048                                        const Matrix &M ) 
01049 {
01050   int i ;
01051   int j ;
01052 
01053   Matrix Mtran( dim2, dim1 ) ;
01054 
01055   for ( i = 0; i < dim1; i++ ) {
01056      for ( j = 0; j < dim2; j++ ) 
01057          Mtran(j,i) = M(i,j) ;
01058   } // end for i
01059 
01060   return Mtran ;
01061 }
01062 
01063 //**********************************************************************
01064 
01065 
01066 
01067 
01068     
01069 int  BbarBrick::sendSelf (int commitTag, Channel &theChannel)
01070 {
01071 
01072   int res = 0;
01073   
01074   // note: we don't check for dataTag == 0 for Element
01075   // objects as that is taken care of in a commit by the Domain
01076   // object - don't want to have to do the check if sending data
01077   int dataTag = this->getDbTag();
01078   
01079   // Quad packs its data into a Vector and sends this to theChannel
01080   // along with its dbTag and the commitTag passed in the arguments
01081 
01082   // Now quad sends the ids of its materials
01083   int matDbTag;
01084   
01085   static ID idData(25);
01086 
01087   idData(24) = this->getTag();
01088   
01089   int i;
01090   for (i = 0; i < 8; i++) {
01091     idData(i) = materialPointers[i]->getClassTag();
01092     matDbTag = materialPointers[i]->getDbTag();
01093     // NOTE: we do have to ensure that the material has a database
01094     // tag if we are sending to a database channel.
01095     if (matDbTag == 0) {
01096       matDbTag = theChannel.getDbTag();
01097                         if (matDbTag != 0)
01098                           materialPointers[i]->setDbTag(matDbTag);
01099     }
01100     idData(i+8) = matDbTag;
01101   }
01102   
01103   idData(16) = connectedExternalNodes(0);
01104   idData(17) = connectedExternalNodes(1);
01105   idData(18) = connectedExternalNodes(2);
01106   idData(19) = connectedExternalNodes(3);
01107   idData(20) = connectedExternalNodes(4);
01108   idData(21) = connectedExternalNodes(5);
01109   idData(22) = connectedExternalNodes(6);
01110   idData(23) = connectedExternalNodes(7);
01111 
01112   res += theChannel.sendID(dataTag, commitTag, idData);
01113   if (res < 0) {
01114     opserr << "WARNING BbarBrick::sendSelf() - " << this->getTag() << "failed to send ID\n";
01115     return res;
01116   }
01117 
01118 
01119   // Finally, quad asks its material objects to send themselves
01120   for (i = 0; i < 8; i++) {
01121     res += materialPointers[i]->sendSelf(commitTag, theChannel);
01122     if (res < 0) {
01123       opserr << "WARNING BbarBrick::sendSelf() - " << this->getTag() << " failed to send its Material\n";
01124       return res;
01125     }
01126   }
01127   
01128   return res;
01129 }
01130 
01131 int  BbarBrick::recvSelf (int commitTag, 
01132                        Channel &theChannel, 
01133                        FEM_ObjectBroker &theBroker)
01134 {
01135   int res = 0;
01136   
01137   int dataTag = this->getDbTag();
01138 
01139   static ID idData(25);
01140   // Quad now receives the tags of its four external nodes
01141   res += theChannel.recvID(dataTag, commitTag, idData);
01142   if (res < 0) {
01143     opserr << "WARNING BbarBrick::recvSelf() - " << this->getTag() << " failed to receive ID\n";
01144     return res;
01145   }
01146 
01147   this->setTag(idData(24));
01148 
01149   connectedExternalNodes(0) = idData(16);
01150   connectedExternalNodes(1) = idData(17);
01151   connectedExternalNodes(2) = idData(18);
01152   connectedExternalNodes(3) = idData(19);
01153   connectedExternalNodes(4) = idData(20);
01154   connectedExternalNodes(5) = idData(21);
01155   connectedExternalNodes(6) = idData(22);
01156   connectedExternalNodes(7) = idData(23);
01157   
01158 
01159   if (materialPointers[0] == 0) {
01160     for (int i = 0; i < 8; i++) {
01161       int matClassTag = idData(i);
01162       int matDbTag = idData(i+8);
01163       // Allocate new material with the sent class tag
01164       materialPointers[i] = theBroker.getNewNDMaterial(matClassTag);
01165       if (materialPointers[i] == 0) {
01166           opserr << "BbarBrick::recvSelf() - Broker could not create NDMaterial of class type" <<
01167             matClassTag << endln;
01168           exit(-1);
01169       }
01170       // Now receive materials into the newly allocated space
01171       materialPointers[i]->setDbTag(matDbTag);
01172       res += materialPointers[i]->recvSelf(commitTag, theChannel, theBroker);
01173       if (res < 0) {
01174         opserr << "NLBeamColumn3d::recvSelf() - material " << 
01175           i << "failed to recv itself\n";
01176         return res;
01177       }
01178     }
01179   }
01180   // materials exist , ensure materials of correct type and recvSelf on them
01181   else {
01182     for (int i = 0; i < 8; i++) {
01183       int matClassTag = idData(i);
01184       int matDbTag = idData(i+8);
01185       // Check that material is of the right type; if not,
01186       // delete it and create a new one of the right type
01187       if (materialPointers[i]->getClassTag() != matClassTag) {
01188         delete materialPointers[i];
01189         materialPointers[i] = theBroker.getNewNDMaterial(matClassTag);
01190         if (materialPointers[i] == 0) {
01191           opserr << "BbarBrick::recvSelf() - Broker could not create NDMaterial of class type" <<
01192             matClassTag << endln;
01193           exit(-1);
01194         }
01195       materialPointers[i]->setDbTag(matDbTag);
01196       }
01197       // Receive the material
01198 
01199       res += materialPointers[i]->recvSelf(commitTag, theChannel, theBroker);
01200       if (res < 0) {
01201         opserr << "NLBeamColumn3d::recvSelf() - material " << 
01202           i << "failed to recv itself\n";
01203         return res;
01204       }
01205     }
01206   }
01207 
01208   return res;
01209 }
01210 //**************************************************************************
01211 
01212 int
01213 BbarBrick::displaySelf(Renderer &theViewer, int displayMode, float fact)
01214 {
01215 
01216     const Vector &end1Crd = nodePointers[0]->getCrds();
01217     const Vector &end2Crd = nodePointers[1]->getCrds(); 
01218     const Vector &end3Crd = nodePointers[2]->getCrds(); 
01219     const Vector &end4Crd = nodePointers[3]->getCrds(); 
01220 
01221     const Vector &end5Crd = nodePointers[4]->getCrds();
01222     const Vector &end6Crd = nodePointers[5]->getCrds(); 
01223     const Vector &end7Crd = nodePointers[6]->getCrds(); 
01224     const Vector &end8Crd = nodePointers[7]->getCrds(); 
01225 
01226     static Matrix coords(4,3);
01227     static Vector values(4);
01228     static Vector P(24) ;
01229 
01230     values(0) = 1 ;
01231     values(1) = 1 ;
01232     values(2) = 1 ;
01233     values(3) = 1 ;
01234 
01235     int error = 0;
01236     int i;
01237 
01238     if (displayMode >= 0) {
01239       const Vector &end1Disp = nodePointers[0]->getDisp();
01240       const Vector &end2Disp = nodePointers[1]->getDisp();
01241       const Vector &end3Disp = nodePointers[2]->getDisp();
01242       const Vector &end4Disp = nodePointers[3]->getDisp();
01243       
01244       const Vector &end5Disp = nodePointers[4]->getDisp();
01245       const Vector &end6Disp = nodePointers[5]->getDisp();
01246       const Vector &end7Disp = nodePointers[6]->getDisp();
01247       const Vector &end8Disp = nodePointers[7]->getDisp();
01248       
01249       if (displayMode < 3 && displayMode > 0)
01250         P = this->getResistingForce();
01251       
01252       for (i = 0; i < 3; i++) {
01253         coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
01254         coords(1,i) = end2Crd(i) + end2Disp(i)*fact;
01255         coords(2,i) = end3Crd(i) + end3Disp(i)*fact;
01256         coords(3,i) = end4Crd(i) + end4Disp(i)*fact;
01257       }
01258       error += theViewer.drawPolygon (coords, values);
01259       
01260       for (i = 0; i < 3; i++) {
01261         coords(0,i) = end5Crd(i) + end5Disp(i)*fact;
01262         coords(1,i) = end6Crd(i) + end6Disp(i)*fact;
01263         coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
01264         coords(3,i) = end8Crd(i) + end8Disp(i)*fact;
01265       }
01266       error += theViewer.drawPolygon (coords, values);
01267       
01268       for (i = 0; i < 3; i++) {
01269         coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
01270         coords(1,i) = end4Crd(i) + end4Disp(i)*fact;
01271         coords(2,i) = end8Crd(i) + end8Disp(i)*fact;
01272         coords(3,i) = end5Crd(i) + end5Disp(i)*fact;
01273       }
01274       error += theViewer.drawPolygon (coords, values);
01275       
01276       for (i = 0; i < 3; i++) {
01277         coords(0,i) = end2Crd(i) + end2Disp(i)*fact;
01278         coords(1,i) = end3Crd(i) + end3Disp(i)*fact;
01279         coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
01280         coords(3,i) = end6Crd(i) + end6Disp(i)*fact;
01281       }
01282       error += theViewer.drawPolygon (coords, values);
01283       
01284       
01285       for (i = 0; i < 3; i++) {
01286         coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
01287         coords(1,i) = end2Crd(i) + end2Disp(i)*fact;
01288         coords(2,i) = end6Crd(i) + end6Disp(i)*fact;
01289         coords(3,i) = end5Crd(i) + end5Disp(i)*fact;
01290       }
01291       error += theViewer.drawPolygon (coords, values);
01292       
01293       for (i = 0; i < 3; i++) {
01294         coords(0,i) = end4Crd(i) + end4Disp(i)*fact;
01295         coords(1,i) = end3Crd(i) + end3Disp(i)*fact;
01296         coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
01297         coords(3,i) = end8Crd(i) + end8Disp(i)*fact;
01298       }
01299       error += theViewer.drawPolygon (coords, values);
01300 
01301     } else {
01302       int mode = displayMode  *  -1;
01303 
01304       const Matrix &eigen1 = nodePointers[0]->getEigenvectors();
01305       const Matrix &eigen2 = nodePointers[1]->getEigenvectors();
01306       const Matrix &eigen3 = nodePointers[2]->getEigenvectors();
01307       const Matrix &eigen4 = nodePointers[3]->getEigenvectors();
01308       const Matrix &eigen5 = nodePointers[4]->getEigenvectors();
01309       const Matrix &eigen6 = nodePointers[5]->getEigenvectors();
01310       const Matrix &eigen7 = nodePointers[6]->getEigenvectors();
01311       const Matrix &eigen8 = nodePointers[7]->getEigenvectors();
01312       
01313       if (eigen1.noCols() >= mode) {
01314 
01315         for (i = 0; i < 3; i++) {
01316           coords(0,i) = end1Crd(i) + eigen1(i,mode-1)*fact;
01317           coords(1,i) = end2Crd(i) + eigen2(i,mode-1)*fact;    
01318           coords(2,i) = end3Crd(i) + eigen3(i,mode-1)*fact;    
01319           coords(3,i) = end4Crd(i) + eigen4(i,mode-1)*fact;
01320         }
01321         
01322         error += theViewer.drawPolygon (coords, values);
01323         
01324         for (i = 0; i < 3; i++) {
01325           coords(0,i) = end5Crd(i) + eigen5(i,mode-1)*fact;
01326           coords(1,i) = end6Crd(i) + eigen6(i,mode-1)*fact;
01327           coords(2,i) = end7Crd(i) + eigen7(i,mode-1)*fact;
01328           coords(3,i) = end8Crd(i) + eigen8(i,mode-1)*fact;
01329         }
01330         
01331         error += theViewer.drawPolygon (coords, values);
01332         
01333         for (i = 0; i < 3; i++) {
01334           coords(0,i) = end1Crd(i) + eigen1(i,mode-1)*fact;
01335           coords(1,i) = end4Crd(i) + eigen4(i,mode-1)*fact;
01336           coords(2,i) = end8Crd(i) + eigen8(i,mode-1)*fact;
01337           coords(3,i) = end5Crd(i) + eigen5(i,mode-1)*fact;
01338         }
01339         
01340         error += theViewer.drawPolygon (coords, values);
01341         
01342         for (i = 0; i < 3; i++) {
01343           coords(0,i) = end2Crd(i) + eigen2(i,mode-1)*fact;
01344           coords(1,i) = end3Crd(i) + eigen3(i,mode-1)*fact;
01345           coords(2,i) = end7Crd(i) + eigen7(i,mode-1)*fact;
01346           coords(3,i) = end6Crd(i) + eigen6(i,mode-1)*fact;
01347         }
01348         
01349         error += theViewer.drawPolygon (coords, values);
01350         
01351         
01352         for (i = 0; i < 3; i++) {
01353           coords(0,i) = end1Crd(i) + eigen1(i,mode-1)*fact;
01354           coords(1,i) = end2Crd(i) + eigen2(i,mode-1)*fact;
01355           coords(2,i) = end6Crd(i) + eigen6(i,mode-1)*fact;
01356           coords(3,i) = end5Crd(i) + eigen5(i,mode-1)*fact;
01357         }
01358         
01359         error += theViewer.drawPolygon (coords, values);
01360         
01361         for (i = 0; i < 3; i++) {
01362           coords(0,i) = end4Crd(i) + eigen4(i,mode-1)*fact;
01363           coords(1,i) = end3Crd(i) + eigen3(i,mode-1)*fact;
01364           coords(2,i) = end7Crd(i) + eigen7(i,mode-1)*fact;
01365           coords(3,i) = end8Crd(i) + eigen8(i,mode-1)*fact;
01366         }
01367         
01368         error += theViewer.drawPolygon (coords, values);
01369       } else {
01370         values.Zero();
01371         for (i = 0; i < 3; i++) {
01372           coords(0,i) = end1Crd(i); 
01373           coords(1,i) = end2Crd(i); 
01374           coords(2,i) = end3Crd(i); 
01375           coords(3,i) = end4Crd(i); 
01376         }
01377         
01378         error += theViewer.drawPolygon (coords, values);
01379         
01380         for (i = 0; i < 3; i++) {
01381           coords(0,i) = end5Crd(i); 
01382           coords(1,i) = end6Crd(i); 
01383           coords(2,i) = end7Crd(i); 
01384           coords(3,i) = end8Crd(i); 
01385         }
01386         
01387         error += theViewer.drawPolygon (coords, values);
01388         
01389         for (i = 0; i < 3; i++) {
01390           coords(0,i) = end1Crd(i); 
01391           coords(1,i) = end4Crd(i); 
01392           coords(2,i) = end8Crd(i); 
01393           coords(3,i) = end5Crd(i); 
01394         }
01395         
01396         error += theViewer.drawPolygon (coords, values);
01397         
01398         for (i = 0; i < 3; i++) {
01399           coords(0,i) = end2Crd(i); 
01400           coords(1,i) = end3Crd(i); 
01401           coords(2,i) = end7Crd(i); 
01402           coords(3,i) = end6Crd(i); 
01403         }
01404         
01405         error += theViewer.drawPolygon (coords, values);
01406         
01407         
01408         for (i = 0; i < 3; i++) {
01409           coords(0,i) = end1Crd(i); 
01410           coords(1,i) = end2Crd(i); 
01411           coords(2,i) = end6Crd(i); 
01412           coords(3,i) = end5Crd(i); 
01413         }
01414         
01415         error += theViewer.drawPolygon (coords, values);
01416         
01417         for (i = 0; i < 3; i++) {
01418           coords(0,i) = end4Crd(i); 
01419           coords(1,i) = end3Crd(i); 
01420           coords(2,i) = end7Crd(i); 
01421           coords(3,i) = end8Crd(i); 
01422         }
01423         
01424         error += theViewer.drawPolygon (coords, values);
01425       }      
01426     }
01427 
01428     return error;
01429 }
01430 
01431 Response*
01432 BbarBrick::setResponse(const char **argv, int argc, Information &eleInfo, OPS_Stream &output)
01433 {
01434   Response *theResponse = 0;
01435 
01436   char outputData[32];
01437 
01438   output.tag("ElementOutput");
01439   output.attr("eleType","BbarBrick");
01440   output.attr("eleTag",this->getTag());
01441   for (int i=1; i<=8; i++) {
01442     sprintf(outputData,"node%d",i);
01443     output.attr(outputData,nodePointers[i-1]->getTag());
01444   }
01445 
01446   if (strcmp(argv[0],"force") == 0 || strcmp(argv[0],"forces") == 0) {
01447 
01448     char outputData[10];
01449     for (int i=1; i<=8; i++) {
01450       sprintf(outputData,"P1_%d",i);
01451       output.tag("ResponseType",outputData);
01452       sprintf(outputData,"P2_%d",i);
01453       output.tag("ResponseType",outputData);
01454       sprintf(outputData,"P3_%d",i);
01455       output.tag("ResponseType",outputData);
01456     }
01457 
01458     theResponse = new ElementResponse(this, 1, resid);
01459   
01460   }   else if (strcmp(argv[0],"material") == 0 || strcmp(argv[0],"integrPoint") == 0) {
01461 
01462     int pointNum = atoi(argv[1]);
01463     if (pointNum > 0 && pointNum <= 8) {
01464 
01465       output.tag("GaussPoint");
01466       output.attr("number",pointNum);
01467 
01468       theResponse =  materialPointers[pointNum-1]->setResponse(&argv[2], argc-2, eleInfo, output);
01469 
01470       output.endTag(); // GaussPoint
01471     }
01472 
01473 
01474   } else if (strcmp(argv[0],"stresses") ==0) {
01475 
01476     for (int i=0; i<8; i++) {
01477       output.tag("GaussPoint");
01478       output.attr("number",i+1);
01479       output.tag("NdMaterialOutput");
01480       output.attr("classType", materialPointers[i]->getClassTag());
01481       output.attr("tag", materialPointers[i]->getTag());
01482 
01483       output.tag("ResponseType","sigma11");
01484       output.tag("ResponseType","sigma22");
01485       output.tag("ResponseType","sigma33");
01486       output.tag("ResponseType","sigma12");
01487       output.tag("ResponseType","sigma13");
01488       output.tag("ResponseType","sigma23");      
01489 
01490       output.endTag(); // NdMaterialOutput
01491       output.endTag(); // GaussPoint
01492     }
01493     theResponse =  new ElementResponse(this, 3, Vector(48));
01494   }
01495   
01496   output.endTag(); // ElementOutput
01497   return theResponse;
01498 }
01499 
01500 int 
01501 BbarBrick::getResponse(int responseID, Information &eleInfo)
01502 {
01503   static Vector stresses(48);
01504 
01505   if (responseID == 1)
01506     return eleInfo.setVector(this->getResistingForce());
01507 
01508   else if (responseID == 2)
01509     return eleInfo.setMatrix(this->getTangentStiff());
01510     
01511   else if (responseID == 3) {
01512     
01513     // Loop over the integration points
01514     int cnt = 0;
01515     for (int i = 0; i < 8; i++) {
01516       
01517       // Get material stress response
01518       const Vector &sigma = materialPointers[i]->getStress();
01519       stresses(cnt++) = sigma(0);
01520       stresses(cnt++) = sigma(1);
01521       stresses(cnt++) = sigma(2);
01522       stresses(cnt++) = sigma(3);
01523       stresses(cnt++) = sigma(4);
01524       stresses(cnt++) = sigma(5);
01525     }
01526     return eleInfo.setVector(stresses);
01527     
01528   }
01529   else
01530     
01531     return -1;
01532 }

Generated on Mon Oct 23 15:05:05 2006 for OpenSees by doxygen 1.5.0