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

EnhancedQuad.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.1 $
00022 // $Date: 2001/07/11 22:59:25 $
00023 // $Source: /usr/local/cvs/OpenSees/SRC/element/fourNodeQuad/EnhancedQuad.cpp,v $
00024 
00025 #include <iostream.h>
00026 #include <stdio.h> 
00027 #include <stdlib.h> 
00028 #include <math.h> 
00029 
00030 #include <ID.h> 
00031 #include <Vector.h>
00032 #include <Matrix.h>
00033 #include <Element.h>
00034 #include <Node.h>
00035 #include <NDMaterial.h>
00036 #include <Domain.h>
00037 #include <ErrorHandler.h>
00038 #include <EnhancedQuad.h>
00039 #include <Renderer.h>
00040 
00041 
00042 //static data
00043 double  EnhancedQuad::xl[2][4] ;
00044 Matrix  EnhancedQuad::stiff(8,8) ;
00045 Vector  EnhancedQuad::resid(8) ;
00046 Matrix  EnhancedQuad::mass(8,8) ;
00047 Matrix  EnhancedQuad::damping(8,8) ;
00048 
00049 double  EnhancedQuad::stressData[3][4] ;
00050 double  EnhancedQuad::tangentData[3][3][4] ;
00051 
00052     
00053 //quadrature data
00054 const double  EnhancedQuad::root3 = sqrt(3.0) ;
00055 const double  EnhancedQuad::one_over_root3 = 1.0 / root3 ;
00056 
00057 const double  EnhancedQuad::sg[] = { -one_over_root3,  
00058          one_over_root3, 
00059          one_over_root3, 
00060          -one_over_root3 } ;
00061 
00062 const double  EnhancedQuad::tg[] = { -one_over_root3, 
00063          -one_over_root3, 
00064          one_over_root3,  
00065          one_over_root3 } ;
00066 
00067 const double  EnhancedQuad::wg[] = { 1.0, 1.0, 1.0, 1.0 } ;
00068 
00069   
00070 
00071 //null constructor
00072 EnhancedQuad::EnhancedQuad( ) :
00073 Element( 0, ELE_TAG_EnhancedQuad ),
00074 connectedExternalNodes(4),
00075 alpha(4) 
00076 { 
00077 
00078 }
00079 
00080 
00081 //full constructor
00082 EnhancedQuad::EnhancedQuad(  int tag, 
00083                          int node1,
00084                          int node2,
00085                      int node3,
00086                          int node4,
00087                   NDMaterial &theMaterial,
00088     const char *type ) 
00089   :
00090 Element( tag, ELE_TAG_EnhancedQuad ),
00091 connectedExternalNodes(4),
00092 alpha(4)
00093 {
00094 
00095   connectedExternalNodes(0) = node1 ;
00096   connectedExternalNodes(1) = node2 ;
00097   connectedExternalNodes(2) = node3 ;
00098   connectedExternalNodes(3) = node4 ;
00099 
00100   int i ;
00101   for ( i = 0 ;  i < 4; i++ ) {
00102 
00103       materialPointers[i] = theMaterial.getCopy(type) ;
00104 
00105       if (materialPointers[i] == 0) {
00106 
00107    g3ErrorHandler->fatal("EnhancedQuad::constructor %s",
00108   "- failed to get a material of type: ShellSection");
00109       } //end if
00110       
00111   } //end for i 
00112 
00113 }
00114 
00115 
00116 //destructor 
00117 EnhancedQuad::~EnhancedQuad( )
00118 {
00119   int i ;
00120   for ( i = 0 ;  i < 4; i++ ) {
00121 
00122     delete materialPointers[i] ;
00123     materialPointers[i] = 0 ; 
00124 
00125     nodePointers[i] = 0 ;
00126 
00127   } //end for i
00128 }
00129 
00130 
00131 //set domain
00132 void  EnhancedQuad::setDomain( Domain *theDomain ) 
00133 {  
00134 
00135   int i ;
00136 
00137   //zero enhanced parameters
00138   alpha.Zero( ) ;
00139 
00140   //node pointers
00141   for ( i = 0; i < 4; i++ ) 
00142      nodePointers[i] = theDomain->getNode( connectedExternalNodes(i) ) ;
00143 
00144   this->DomainComponent::setDomain(theDomain) ;
00145 }
00146 
00147 
00148 //get the number of external nodes
00149 int  EnhancedQuad::getNumExternalNodes( ) const
00150 {
00151   return 4 ;
00152 } 
00153  
00154 
00155 //return connected external nodes
00156 const ID&  EnhancedQuad::getExternalNodes( ) 
00157 {
00158   return connectedExternalNodes ;
00159 } 
00160 
00161 
00162 //return number of dofs
00163 int  EnhancedQuad::getNumDOF( ) 
00164 {
00165   return 8 ;
00166 }
00167 
00168 
00169 //commit state
00170 int  EnhancedQuad::commitState( )
00171 {
00172   int i ;
00173   int success = 0 ;
00174 
00175   for ( i = 0; i < 4; i++ ) 
00176     success += materialPointers[i]->commitState( ) ;
00177   
00178   return success ;
00179 }
00180  
00181 
00182 
00183 //revert to last commit 
00184 int  EnhancedQuad::revertToLastCommit( ) 
00185 {
00186   int i ;
00187   int success = 0 ;
00188 
00189   for ( i = 0; i < 4; i++ ) 
00190     success += materialPointers[i]->revertToLastCommit( ) ;
00191   
00192   return success ;
00193 }
00194     
00195 
00196 //revert to start 
00197 int  EnhancedQuad::revertToStart( ) 
00198 {
00199   int i ;
00200   int success = 0 ;
00201 
00202   //zero enhanced parameters
00203   this->alpha.Zero( ) ;
00204 
00205   for ( i = 0; i < 4; i++ ) 
00206     success += materialPointers[i]->revertToStart( ) ;
00207   
00208   return success ;
00209 }
00210 
00211 //print out element data
00212 void  EnhancedQuad::Print( ostream &s, int flag )
00213 {
00214   s << '\n' ;
00215   s << "Enhanced Strain Four Node Quad \n" ;
00216   s << "Node 1 : " << connectedExternalNodes(0) << '\n' ;
00217   s << "Node 2 : " << connectedExternalNodes(1) << '\n' ;
00218   s << "Node 3 : " << connectedExternalNodes(2) << '\n' ;
00219   s << "Node 4 : " << connectedExternalNodes(3) << '\n' ;
00220 
00221   s << "Material Information : \n " ;
00222   materialPointers[0]->Print( s, flag ) ;
00223 
00224   s << endl ;
00225 }
00226 
00227 //return stiffness matrix 
00228 const Matrix&  EnhancedQuad::getTangentStiff( ) 
00229 {
00230   int tang_flag = 1 ; //get the tangent 
00231 
00232   //do tangent and residual here
00233   formResidAndTangent( tang_flag ) ;  
00234 
00235   return stiff ;
00236 }    
00237 
00238 
00239 //return secant matrix 
00240 const Matrix&  EnhancedQuad::getSecantStiff( ) 
00241 {
00242    int tang_flag = 1 ; //get the tangent
00243 
00244   //do tangent and residual here
00245   formResidAndTangent( tang_flag ) ;  
00246     
00247   return stiff ;
00248 }
00249     
00250 
00251 //return damping matrix because frank is a dumb ass 
00252 const Matrix&  EnhancedQuad::getDamp( ) 
00253 {
00254   //not supported
00255   return damping ;
00256 }    
00257 
00258 
00259 //return mass matrix
00260 const Matrix&  EnhancedQuad::getMass( ) 
00261 {
00262 
00263   int tangFlag = 1 ;
00264 
00265   formInertiaTerms( tangFlag ) ;
00266 
00267   return mass ;
00268 
00269 } 
00270 
00271 
00272 //zero the load -- what load?
00273 void  EnhancedQuad::zeroLoad( )
00274 {
00275   return ;
00276 }
00277 
00278 //add load -- what load?
00279 int  EnhancedQuad::addLoad( const Vector &addP )
00280 {
00281   return -1 ;
00282 }
00283 
00284 
00285 //get residual
00286 const Vector&  EnhancedQuad::getResistingForce( ) 
00287 {
00288   int tang_flag = 0 ; //don't get the tangent
00289 
00290   formResidAndTangent( tang_flag ) ;
00291 
00292   return resid ;   
00293 }
00294 
00295 
00296 //get residual with inertia terms
00297 const Vector&  EnhancedQuad::getResistingForceIncInertia( )
00298 {
00299   int tang_flag = 0 ; //don't get the tangent
00300 
00301   //do tangent and residual here 
00302   formResidAndTangent( tang_flag ) ;
00303 
00304   //inertia terms
00305   formInertiaTerms( tang_flag ) ;
00306 
00307   return resid ;
00308 }
00309 
00310 
00311 //*********************************************************************
00312 //form inertia terms
00313 
00314 void   EnhancedQuad::formInertiaTerms( int tangFlag ) 
00315 {
00316 
00317   static const int ndm = 2 ;
00318 
00319   static const int ndf = 2 ; 
00320 
00321   static const int numberNodes = 4 ;
00322 
00323   static const int numberGauss = 4 ;
00324 
00325   static const int nShape = 3 ;
00326 
00327   static const int massIndex = nShape - 1 ;
00328 
00329   double xsj ;  // determinant jacaobian matrix 
00330 
00331   double dvol ; //volume element
00332 
00333   static double shp[nShape][numberNodes] ;  //shape functions at a gauss point
00334 
00335   static Vector momentum(ndf) ;
00336 
00337   int i, j, k, p ;
00338   int jj, kk ;
00339 
00340   double temp, rho, massJK ;
00341 
00342 
00343   //zero mass 
00344   mass.Zero( ) ;
00345 
00346   //compute basis vectors and local nodal coordinates
00347   computeBasis( ) ;
00348 
00349 
00350   //gauss loop 
00351   for ( i = 0; i < numberGauss; i++ ) {
00352 
00353     //get shape functions    
00354     shape2d( sg[i], tg[i], xl, shp, xsj ) ;
00355     
00356     //volume element
00357     dvol = wg[i] * xsj ;
00358 
00359     //node loop to compute acceleration
00360     momentum.Zero( ) ;
00361     for ( j = 0; j < numberNodes; j++ ) 
00362       momentum += shp[massIndex][j] * ( nodePointers[j]->getTrialAccel()  ) ; 
00363 
00364     //density
00365     rho = materialPointers[i]->getRho() ;
00366 
00367     //multiply acceleration by density to form momentum
00368     momentum *= rho ;
00369 
00370 
00371     //residual and tangent calculations node loops
00372     jj = 0 ;
00373     for ( j = 0; j < numberNodes; j++ ) {
00374 
00375       temp = shp[massIndex][j] * dvol ;
00376 
00377       for ( p = 0; p < ndf; p++ )
00378         resid( jj+p ) += ( temp * momentum(p) )  ;
00379 
00380       
00381       if ( tangFlag == 1 ) {
00382 
00383   //multiply by density
00384   temp *= rho ;
00385 
00386   //node-node mass
00387          kk = 0 ;
00388          for ( k = 0; k < numberNodes; k++ ) {
00389 
00390      massJK = temp * shp[massIndex][k] ;
00391 
00392             for ( p = 0; p < ndf; p++ )  
00393        mass( jj+p, kk+p ) += massJK ;
00394             
00395             kk += ndf ;
00396           } // end for k loop
00397 
00398       } // end if tang_flag 
00399 
00400       jj += ndf ;
00401     } // end for j loop
00402 
00403 
00404   } //end for i gauss loop 
00405 
00406 
00407 
00408 }
00409 
00410 //*********************************************************************
00411 //form residual and tangent
00412 void  EnhancedQuad::formResidAndTangent( int tang_flag ) 
00413 {
00414 
00415   static const double tolerance = 1.0e-08 ;
00416 
00417   static const int nIterations = 10 ;
00418 
00419   static const int ndm = 2 ;
00420 
00421   static const int ndf = 2 ; 
00422 
00423   static const int nstress = 3 ;
00424  
00425   static const int numberNodes = 4 ;
00426 
00427   static const int numberGauss = 4 ;
00428 
00429   static const int nShape = 3 ;
00430 
00431   static const int nEnhanced = 4 ;
00432   
00433   static const int nModes = 2 ;
00434 
00435   static const int numberDOF = 8 ;
00436 
00437 
00438   int i, j, k, p, q ;
00439   int jj, kk ;
00440 
00441   int success ;
00442   
00443 
00444   static double xsj[numberGauss] ;  // determinant jacaobian matrix 
00445 
00446   static double dvol[numberGauss] ; //volume element
00447 
00448   static Vector strain(nstress) ;  //strain
00449 
00450   static double shp[nShape][numberNodes] ;  //shape functions at a gauss point
00451 
00452   static double Shape[nShape][numberNodes][numberGauss] ; //all the shape functions
00453 
00454   static Vector residJ(ndf) ; //nodeJ residual 
00455 
00456   static Matrix stiffJK(ndf,ndf) ; //nodeJK stiffness 
00457 
00458   static Matrix stiffKJ(ndf,ndf) ; //nodeKJ stiffness
00459 
00460   static Vector stress(nstress) ;  //stress
00461 
00462   static Matrix dd(nstress,nstress) ;  //material tangent
00463 
00464   static Matrix J0(ndm,ndm) ; //Jacobian matrix at center of element
00465 
00466   static Matrix J0inv(ndm,ndm) ; //inverse of above
00467 
00468 
00469   static Matrix Kee(nEnhanced,nEnhanced) ;
00470 
00471   static Vector residE(nEnhanced) ;
00472 
00473   static Vector Umode(ndf) ;
00474 
00475   static Vector dalpha(nEnhanced) ;
00476 
00477   static Matrix Kue(numberDOF,nEnhanced) ;
00478 
00479   static Matrix Keu(nEnhanced,numberDOF) ;
00480 
00481   static Matrix KeeInvKeu(nEnhanced,numberDOF) ;
00482 
00483 
00484   //---------B-matrices------------------------------------
00485 
00486     static Matrix BJ(nstress,ndf) ;      // B matrix node J
00487 
00488     static Matrix BJtran(ndf,nstress) ;
00489 
00490     static Matrix BK(nstress,ndf) ;      // B matrix node k
00491 
00492     static Matrix BKtran(ndf,nstress) ;
00493 
00494     static Matrix BJtranD(ndf,nstress) ;
00495 
00496   //-------------------------------------------------------
00497 
00498   
00499   //zero stiffness and residual 
00500   stiff.Zero( ) ;
00501   resid.Zero( ) ;
00502 
00503   Kee.Zero( ) ;
00504   residE.Zero( ) ;
00505 
00506   Kue.Zero( ) ;
00507   Keu.Zero( ) ;
00508 
00509 
00510   //compute basis vectors and local nodal coordinates
00511   computeBasis( ) ;
00512 
00513   //compute Jacobian and inverse at center
00514   double L1 = 0.0 ;
00515   double L2 = 0.0 ;
00516   computeJacobian( L1, L2, xl, J0, J0inv ) ; 
00517 
00518 
00519   //gauss loop to compute and save shape functions 
00520   double det ;
00521   for ( i = 0; i < numberGauss; i++ ) {
00522 
00523     //get shape functions    
00524     shape2d( sg[i], tg[i], xl, shp, det ) ;
00525 
00526     //save shape functions
00527     for ( p = 0; p < nShape; p++ ) {
00528        for ( q = 0; q < numberNodes; q++ )
00529    Shape[p][q][i] = shp[p][q] ;
00530     } // end for p
00531 
00532     //save jacobian determinant
00533     xsj[i] = det ;
00534 
00535     //volume element to also be saved
00536     dvol[i] = wg[i] * det ;  
00537 
00538   } // end for i gauss loop
00539 
00540 
00541   //-------------------------------------------------------------------
00542   //newton loop to solve for enhanced strain parameters
00543 
00544   int count = 0 ;
00545   do {
00546 
00547     residE.Zero( ) ;
00548     Kee.Zero( ) ;
00549 
00550     //gauss loop
00551     for ( i = 0; i < numberGauss; i++ ) {
00552 
00553       //extract shape functions from saved array
00554       for ( p = 0; p < nShape; p++ ) {
00555  for ( q = 0; q < numberNodes; q++ )
00556    shp[p][q]  = Shape[p][q][i] ;
00557       } // end for p
00558 
00559 
00560       //zero the strain
00561       strain.Zero( ) ;
00562 
00563       // j-node loop to compute nodal strain contributions
00564       for ( j = 0; j < numberNodes; j++ )  {
00565 
00566         //compute B matrix 
00567         BJ = computeB( j, shp ) ;
00568       
00569         //nodal displacements 
00570         const Vector &ul = nodePointers[j]->getTrialDisp( ) ;
00571 
00572         //compute the strain
00573         strain += (BJ*ul) ; 
00574 
00575       } // end for j
00576 
00577 
00578       // j-node loop to compute enhanced strain contributions
00579       for ( j = 0; j < nModes; j++ )  {
00580 
00581         //compute B matrix 
00582         BJ = computeBenhanced( j, sg[i], tg[i], xsj[i], J0inv ) ; 
00583       
00584         //enhanced "displacements" 
00585         Umode(0) = this->alpha( 2*j     ) ;
00586         Umode(1) = this->alpha( 2*j + 1 ) ;
00587 
00588         //compute the strain
00589         strain += (BJ*Umode) ; 
00590 
00591       } // end for j
00592 
00593 
00594       //send the strain to the material 
00595       success = materialPointers[i]->setTrialStrain( strain ) ;
00596 
00597       //compute the stress
00598       stress = materialPointers[i]->getStress( ) ;
00599 
00600       //multiply by volume element
00601       stress  *= dvol[i] ;
00602 
00603       //tangent 
00604       dd = materialPointers[i]->getTangent( ) ;
00605 
00606       //multiply by volume element
00607       dd *= dvol[i] ;
00608 
00609 
00610       //save stress and tangent (already multiplied by volume element)
00611       saveData( i, stress, dd ) ; 
00612 
00613 
00614       //enhanced residual and tangent calculations loops
00615 
00616       jj = 0 ;
00617       for ( j = 0; j < nModes; j++ ) {
00618 
00619         //compute B matrix 
00620         BJ = computeBenhanced( j, sg[i], tg[i], xsj[i], J0inv ) ; 
00621 
00622  //transpose 
00623  BJtran = transpose( nstress, ndf, BJ ) ;
00624 
00625 
00626  //residual
00627  residJ = (BJtran * stress) ;
00628  residJ *= (-1.0) ;
00629 
00630  //residual 
00631  for ( p = 0; p < ndf; p++ )
00632    residE( jj+p ) += residJ(p)  ;
00633 
00634 
00635         BJtranD = BJtran * dd ;
00636  
00637  kk = 0 ;
00638  for ( k = 0; k < nModes; k++ ) {
00639 
00640           BK = computeBenhanced( k, sg[i], tg[i], xsj[i], J0inv ) ;
00641   
00642           stiffJK =  BJtranD * BK  ;
00643 
00644           for ( p = 0; p < ndf; p++ )  {
00645              for ( q = 0; q < ndf; q++ )
00646                 Kee( jj+p, kk+q ) += stiffJK( p, q ) ;
00647           } //end for p
00648 
00649             kk += ndf ;
00650           } // end for k loop
00651 
00652  jj += ndf ;
00653       } // end for j loop
00654 
00655 
00656     } // end for i gauss loop 
00657 
00658     //cerr << "residE = " << residE << endl ;
00659 
00660     //solve for enhanced strain parameters
00661 
00662     dalpha.Zero( ) ;
00663 
00664     Kee.Solve( residE, dalpha ) ;
00665 
00666     this->alpha += dalpha ;
00667 
00668     count++ ;
00669     if ( count > nIterations ) {
00670       cerr << "Exceeded " << nIterations
00671     << " solving for enhanced strain parameters " << endl ;
00672       break ;
00673     } //end if 
00674 
00675 
00676     //do at least 2 iterations so saved data is good
00677   } while ( residE.Norm() > tolerance  ||  count < 2 ) ;
00678 
00679   //cerr << "alpha = " << alpha << endl ;
00680 
00681   //end enhanced strain parameters newton loop
00682   //-------------------------------------------------------------------
00683 
00684 
00685 
00686   //gauss loop 
00687   for ( i = 0; i < numberGauss; i++ ) {
00688 
00689     //extract shape functions from saved array
00690     for ( p = 0; p < nShape; p++ ) {
00691        for ( q = 0; q < numberNodes; q++ )
00692    shp[p][q]  = Shape[p][q][i] ;
00693     } // end for p
00694 
00695 
00696     //recover stress and tangent from saved data
00697     getData( i, stress, dd ) ;
00698 
00699 
00700     //residual and tangent calculations node loops
00701 
00702     jj = 0 ;
00703     for ( j = 0; j < numberNodes; j++ ) {
00704 
00705       BJ = computeB( j, shp ) ;
00706    
00707       //transpose 
00708       BJtran = transpose( nstress, ndf, BJ ) ;
00709 
00710       //residual
00711       residJ = BJtran * stress ;
00712 
00713       for ( p = 0; p < ndf; p++ )
00714         resid( jj+p ) += residJ(p)  ;
00715 
00716 
00717       if ( tang_flag == 1 ) {
00718 
00719          BJtranD = BJtran * dd ;
00720 
00721   //node-node stiffness
00722          kk = 0 ;
00723          for ( k = 0; k < numberNodes; k++ ) {
00724 
00725             BK = computeB( k, shp ) ;
00726   
00727             stiffJK =  BJtranD * BK  ;
00728 
00729             for ( p = 0; p < ndf; p++ )  {
00730                for ( q = 0; q < ndf; q++ )
00731                   stiff( jj+p, kk+q ) += stiffJK( p, q ) ;
00732             } //end for p
00733 
00734             kk += ndf ;
00735           } // end for k loop
00736 
00737 
00738   //node-enhanced stiffness Kue 
00739          kk = 0 ;
00740          for ( k = 0; k < nModes; k++ ) {
00741 
00742             BK = computeBenhanced( k, sg[i], tg[i], xsj[i], J0inv ) ;
00743   
00744             stiffJK =  BJtranD * BK  ;
00745 
00746             for ( p = 0; p < ndf; p++ )  {
00747                for ( q = 0; q < ndf; q++ )
00748                   Kue( jj+p, kk+q ) += stiffJK( p, q ) ;
00749             } //end for p 
00750 
00751             kk += ndf ;
00752           } // end for k loop
00753 
00754 
00755   //enhanced-node stiffness Keu 
00756          kk = 0 ;
00757          for ( k = 0; k < nModes; k++ ) {
00758 
00759             BK = computeBenhanced( k, sg[i], tg[i], xsj[i], J0inv ) ;
00760 
00761      //transpose 
00762      BKtran = transpose( nstress, ndf, BK ) ;
00763   
00764             stiffKJ = (BKtran*dd)*BJ ;
00765 
00766             for ( p = 0; p < ndf; p++ )  {
00767                for ( q = 0; q < ndf; q++ )
00768                   Keu( kk+p, jj+q ) += stiffKJ( p, q ) ;
00769             } //end for p  
00770 
00771             kk += ndf ;
00772           } // end for k loop
00773 
00774 
00775 
00776       } // end if tang_flag 
00777 
00778       jj += ndf ;
00779     } // end for j loop
00780 
00781 
00782   } //end for i gauss loop 
00783 
00784 
00785   //static condensation of enhanced parameters
00786 
00787   if ( tang_flag == 1 ) {
00788     
00789      Kee.Solve( Keu, KeeInvKeu ) ;
00790 
00791      stiff -= ( Kue * KeeInvKeu ) ;
00792 
00793   } //end if 
00794 
00795   return ;
00796 }
00797   
00798 
00799 //************************************************************************
00800 void   EnhancedQuad::saveData( int gp, 
00801           const Vector &stress,
00802           const Matrix &tangent ) 
00803 {
00804 
00805   int i, j ;
00806 
00807   //save stress
00808   for ( i=0; i<3; i++ )
00809     stressData[i][gp] = stress(i) ;
00810 
00811   //save tangent
00812   for ( i=0; i<3; i++ ) {
00813     for (j=0; j<3; j++ ) 
00814       tangentData[i][j][gp] = tangent(i,j) ;
00815   } //end for i
00816 
00817   return ;
00818 }
00819 
00820 
00821 //************************************************************************
00822 //recover stress and tangent data
00823 
00824 void  EnhancedQuad::getData( int gp,
00825         Vector &stress,
00826         Matrix &tangent ) 
00827 {
00828 
00829   int i, j ;
00830 
00831   //get stress
00832   for ( i=0; i<3; i++ )
00833     stress(i) = stressData[i][gp] ;
00834 
00835   //get tangent
00836   for ( i=0; i<3; i++ ) {
00837     for ( j=0; j<3; j++ ) 
00838       tangent(i,j) = tangentData[i][j][gp] ;
00839   } //end for i
00840 
00841   return ;
00842 
00843 }
00844 
00845 //************************************************************************
00846 //compute local coordinates and basis
00847 
00848 void   EnhancedQuad::computeBasis( ) 
00849 {
00850   //nodal coordinates 
00851 
00852   int i ;
00853   for ( i = 0; i < 4; i++ ) {
00854 
00855        const Vector &coorI = nodePointers[i]->getCrds( ) ;
00856 
00857        xl[0][i] = coorI(0) ;
00858        xl[1][i] = coorI(1) ;
00859 
00860   }  //end for i 
00861 
00862 }
00863 
00864 //*************************************************************************
00865 //compute B matrix
00866 
00867 Matrix   
00868 EnhancedQuad::computeB( int node, const double shp[3][4] ) 
00869 {
00870 
00871   static Matrix B(3,2) ;
00872 
00873 //---B Matrix in standard {1,2,3} mechanics notation---------------
00874 //
00875 //                -             -
00876 //               | +N,1      0   | 
00877 // B    =        |   0     +N,2  |    (3x2)
00878 //               | +N,2    +N,1  |
00879 //                -             -  
00880 //
00881 //-------------------------------------------------------------------
00882 
00883   B.Zero( ) ;
00884 
00885   B(0,0) = shp[0][node] ;
00886   B(1,1) = shp[1][node] ;
00887   B(2,0) = shp[1][node] ;
00888   B(2,1) = shp[0][node] ;
00889 
00890   return B ;
00891 
00892 }
00893 
00894 //***********************************************************************
00895 //compute enhanced strain B-matrices
00896 
00897 Matrix   EnhancedQuad::computeBenhanced( int node, 
00898       double L1,
00899       double L2,
00900       double j, 
00901       const Matrix &Jinv )
00902 {
00903   static Matrix B(3,2) ;
00904 
00905   static Matrix JinvTran(2,2) ;
00906 
00907   static double shape[2] ;
00908 
00909   static double parameter ;
00910 
00911 
00912   //compute JinvTran
00913   JinvTran(0,0) = Jinv(0,0) ;
00914   JinvTran(1,1) = Jinv(1,1) ;
00915   JinvTran(0,1) = Jinv(1,0) ;
00916   JinvTran(1,0) = Jinv(0,1) ;      //residual 
00917 
00918 
00919 
00920   if ( node == 0 ) {
00921 
00922     //first column of JinvTran 
00923     shape[0] = JinvTran(0,0) ;
00924     shape[1] = JinvTran(1,0) ;
00925 
00926     parameter = L1 / j ;
00927 
00928   }
00929   else if ( node == 1 ) {
00930 
00931     //second column of JinvTran 
00932     shape[0] = JinvTran(0,1) ;
00933     shape[1] = JinvTran(1,1) ;
00934 
00935     parameter = L2 / j ;
00936 
00937   } //end if 
00938 
00939   shape[0] *= parameter ;
00940   shape[1] *= parameter ; 
00941 
00942 
00943   B.Zero( ) ;
00944 
00945   B(0,0) = shape[0] ;
00946   B(1,1) = shape[1] ;
00947   B(2,0) = shape[1] ;
00948   B(2,1) = shape[0] ;
00949   
00950   return B ;
00951 
00952 }
00953 
00954 //***********************************************************************
00955 //compute Jacobian matrix and inverse at point {L1,L2} 
00956 
00957 void   EnhancedQuad::computeJacobian( double L1, double L2,
00958           const double x[2][4], 
00959                                       Matrix &JJ, 
00960                                       Matrix &JJinv )
00961 {
00962   int i, j, k ;
00963      
00964   static const double s[] = { -0.5,  0.5, 0.5, -0.5 } ;
00965   static const double t[] = { -0.5, -0.5, 0.5,  0.5 } ;
00966 
00967   static double shp[2][4] ;
00968 
00969   double ss = L1 ;
00970   double tt = L2 ;
00971 
00972   for ( i = 0; i < 4; i++ ) {
00973       shp[0][i] = s[i] * ( 0.5 + t[i]*tt ) ;
00974       shp[1][i] = t[i] * ( 0.5 + s[i]*ss ) ;
00975   } // end for i
00976 
00977   
00978   // Construct jacobian and its inverse
00979   
00980   JJ.Zero( ) ;
00981   for ( i = 0; i < 2; i++ ) {
00982     for ( j = 0; j < 2; j++ ) {
00983 
00984       for ( k = 0; k < 4; k++ )
00985    JJ(i,j) +=  x[i][k] * shp[j][k] ;
00986 
00987     } //end for j
00988   }  // end for i 
00989 
00990   double xsj = JJ(0,0)*JJ(1,1) - JJ(0,1)*JJ(1,0) ;
00991 
00992   //inverse jacobian
00993   JJinv(0,0) =  JJ(1,1) / xsj ;
00994   JJinv(1,1) =  JJ(0,0) / xsj ;
00995   JJinv(0,1) = -JJ(0,1) / xsj ; 
00996   JJinv(1,0) = -JJ(1,0) / xsj ; 
00997 
00998   return ;
00999 
01000 }
01001 
01002 //************************************************************************
01003 //shape function routine for four node quads
01004 
01005 void  EnhancedQuad::shape2d( double ss, double tt, 
01006              const double x[2][4], 
01007              double shp[3][4], 
01008              double &xsj            )
01009 
01010 { 
01011 
01012   int i, j, k ;
01013 
01014   double temp ;
01015      
01016   static const double s[] = { -0.5,  0.5, 0.5, -0.5 } ;
01017   static const double t[] = { -0.5, -0.5, 0.5,  0.5 } ;
01018 
01019   static Matrix xs(2,2) ;
01020   static Matrix sx(2,2) ;
01021 
01022   for ( i = 0; i < 4; i++ ) {
01023       shp[2][i] = ( 0.5 + s[i]*ss )*( 0.5 + t[i]*tt ) ;
01024       shp[0][i] = s[i] * ( 0.5 + t[i]*tt ) ;
01025       shp[1][i] = t[i] * ( 0.5 + s[i]*ss ) ;
01026   } // end for i
01027 
01028   
01029   // Construct jacobian and its inverse
01030   
01031   xs.Zero( ) ;
01032   for ( i = 0; i < 2; i++ ) {
01033     for ( j = 0; j < 2; j++ ) {
01034 
01035       for ( k = 0; k < 4; k++ )
01036    xs(i,j) +=  x[i][k] * shp[j][k] ;
01037 
01038     } //end for j
01039   }  // end for i 
01040 
01041   xsj = xs(0,0)*xs(1,1) - xs(0,1)*xs(1,0) ;
01042 
01043   //inverse jacobian
01044   sx(0,0) =  xs(1,1) / xsj ;
01045   sx(1,1) =  xs(0,0) / xsj ;
01046   sx(0,1) = -xs(0,1) / xsj ; 
01047   sx(1,0) = -xs(1,0) / xsj ; 
01048 
01049 
01050   //form global derivatives 
01051   
01052   for ( i = 0; i < 4; i++ ) {
01053     temp      = shp[0][i]*sx(0,0) + shp[1][i]*sx(1,0) ;
01054     shp[1][i] = shp[0][i]*sx(0,1) + shp[1][i]*sx(1,1) ;
01055     shp[0][i] = temp ;
01056   } // end for i
01057 
01058   return ;
01059 }
01060     
01061 //**********************************************************************
01062 
01063 Matrix  EnhancedQuad::transpose( int dim1, 
01064                                  int dim2, 
01065                    const Matrix &M ) 
01066 {
01067   int i ;
01068   int j ;
01069 
01070   Matrix Mtran( dim2, dim1 ) ;
01071 
01072   for ( i = 0; i < dim1; i++ ) {
01073      for ( j = 0; j < dim2; j++ ) 
01074          Mtran(j,i) = M(i,j) ;
01075   } // end for i
01076 
01077   return Mtran ;
01078 }
01079 
01080 //**********************************************************************
01081 
01082 int  EnhancedQuad::sendSelf (int commitTag, Channel &theChannel)
01083 {
01084     return -1;
01085 }
01086     
01087 int  EnhancedQuad::recvSelf (int commitTag, 
01088          Channel &theChannel, 
01089          FEM_ObjectBroker &theBroker)
01090 {
01091     return -1;
01092 }
01093 //**************************************************************************
01094 
01095 int
01096 EnhancedQuad::displaySelf(Renderer &theViewer, int displayMode, float fact)
01097 {
01098     // first determine the end points of the quad based on
01099     // the display factor (a measure of the distorted image)
01100     // store this information in 4 3d vectors v1 through v4
01101     const Vector &end1Crd = nodePointers[0]->getCrds();
01102     const Vector &end2Crd = nodePointers[1]->getCrds(); 
01103     const Vector &end3Crd = nodePointers[2]->getCrds(); 
01104     const Vector &end4Crd = nodePointers[3]->getCrds(); 
01105 
01106     const Vector &end1Disp = nodePointers[0]->getDisp();
01107     const Vector &end2Disp = nodePointers[1]->getDisp();
01108     const Vector &end3Disp = nodePointers[2]->getDisp();
01109     const Vector &end4Disp = nodePointers[3]->getDisp();
01110 
01111     static Matrix coords(4,3) ;
01112     static Vector values(4) ;
01113     static Vector P(8) ;
01114 
01115     coords.Zero( ) ;
01116 
01117     values(0) = 1 ;
01118     values(1) = 1 ;
01119     values(2) = 1 ;
01120     values(3) = 1 ;
01121 
01122     if (displayMode < 3 && displayMode > 0)
01123       P = this->getResistingForce();
01124 
01125     for (int i = 0; i < 2; i++) {
01126       coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
01127       coords(1,i) = end2Crd(i) + end2Disp(i)*fact;    
01128       coords(2,i) = end3Crd(i) + end3Disp(i)*fact;    
01129       coords(3,i) = end4Crd(i) + end4Disp(i)*fact;    
01130       /*      if (displayMode < 3 && displayMode > 0)
01131  values(i) = P(displayMode*2+i);
01132       else
01133       values(i) = 1;  */
01134     }
01135 
01136     //cerr << coords;
01137     int error = 0;
01138 
01139     error += theViewer.drawPolygon (coords, values);
01140 
01141     return error;
01142 }
01143 
01144 //***************************************************************************
01145 /* 
01146 
01147     //zero the strains
01148     strain.Zero( ) ;
01149 
01150 
01151     // j-node loop to compute strain 
01152     for ( j = 0; j < numberNodes; j++ )  {
01153 
01154       //compute B matrix 
01155       BJ = computeB( j, shp ) ;
01156       
01157       //nodal displacements 
01158       const Vector &ul = nodePointers[j]->getTrialDisp( ) ;
01159 
01160       //compute the strain
01161       strain += (BJ*ul) ; 
01162 
01163     } // end for j
01164 
01165 
01166     // j-node loop to compute enhanced strain contributions
01167     for ( j = 0; j < nModes; j++ )  {
01168 
01169       //compute B matrix 
01170       BJ = computeBenhanced( j, sg[i], tg[i], xsj[i], J0inv ) ; 
01171       
01172       //enhanced "displacements" 
01173       Umode(0) = this->alpha( 2*j     ) ;
01174       Umode(1) = this->alpha( 2*j + 1 ) ;
01175 
01176       //compute the strain
01177       strain += (BJ*Umode) ; 
01178 
01179     } // end for j
01180 
01181 
01182     //send the strain to the material 
01183     success = materialPointers[i]->setTrialStrain( strain ) ;
01184 
01185     //compute the stress
01186     stress = materialPointers[i]->getStress( ) ;
01187 
01188     //multiply by volume element
01189     stress  *= dvol[i] ;
01190 
01191     if ( tang_flag == 1 ) {
01192       dd = materialPointers[i]->getTangent( ) ;
01193       dd *= dvol[i] ;
01194     } //end if tang_flag
01195     
01196 */
Copyright Contact Us