J2Plasticity.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 ** ****************************************************************** */
00015                                                                         
00016 // $Revision: 1.9 $
00017 // $Date: 2005/03/25 00:32:11 $
00018 // $Source: /usr/local/cvs/OpenSees/SRC/material/nD/J2Plasticity.cpp,v $
00019 
00020 // Written: Ed "C++" Love
00021 //
00022 // J2 isotropic hardening material class
00023 // 
00024 //  Elastic Model
00025 //  sigma = K*trace(epsilion_elastic) + (2*G)*dev(epsilon_elastic)
00026 //
00027 //  Yield Function
00028 //  phi(sigma,q) = || dev(sigma) ||  - sqrt(2/3)*q(xi) 
00029 //
00030 //  Saturation Isotropic Hardening with linear term
00031 //  q(xi) = simga_infty + (sigma_0 - sigma_infty)*exp(-delta*xi) + H*xi 
00032 //
00033 //  Flow Rules
00034 //  \dot{epsilon_p} =  gamma * d_phi/d_sigma
00035 //  \dot{xi}        = -gamma * d_phi/d_q 
00036 //
00037 //  Linear Viscosity 
00038 //  gamma = phi / eta  ( if phi > 0 ) 
00039 //
00040 //  Backward Euler Integration Routine 
00041 //  Yield condition enforced at time n+1 
00042 //
00043 //  set eta := 0 for rate independent case
00044 //
00045 
00046 #include <J2Plasticity.h>
00047 #include <J2PlaneStress.h>
00048 #include <J2PlaneStrain.h>
00049 #include <J2AxiSymm.h>
00050 #include <J2PlateFiber.h>
00051 #include <J2ThreeDimensional.h> 
00052 
00053 #include <Channel.h>
00054 #include <FEM_ObjectBroker.h>
00055 
00056 //this is mike's problem
00057 Tensor J2Plasticity :: rank2(2, def_dim_2, 0.0 ) ;
00058 Tensor J2Plasticity :: rank4(2, def_dim_2, 0.0 ) ;
00059 
00060 //parameters
00061 const double J2Plasticity :: one3   = 1.0 / 3.0 ;
00062 const double J2Plasticity :: two3   = 2.0 / 3.0 ;
00063 const double J2Plasticity :: four3  = 4.0 / 3.0 ;
00064 const double J2Plasticity :: root23 = sqrt( 2.0 / 3.0 ) ;
00065 
00066 double J2Plasticity::initialTangent[3][3][3][3] ;   //material tangent
00067 double J2Plasticity::IIdev[3][3][3][3] ; //rank 4 deviatoric 
00068 double J2Plasticity::IbunI[3][3][3][3] ; //rank 4 I bun I 
00069 
00070 //zero internal variables
00071 void J2Plasticity :: zero ( ) 
00072 {
00073   xi_n = 0.0 ;
00074   xi_nplus1 = 0.0 ;
00075   
00076   epsilon_p_n.Zero( ) ;
00077   epsilon_p_nplus1.Zero( ) ;
00078 
00079   stress.Zero();
00080   strain.Zero();
00081 }
00082 
00083 
00084 //null constructor
00085 J2Plasticity ::  J2Plasticity( ) : 
00086 NDMaterial( ),
00087 epsilon_p_n(3,3),
00088 epsilon_p_nplus1(3,3),
00089 stress(3,3),
00090 strain(3,3)
00091 { 
00092   bulk        = 0.0 ;
00093   shear       = 0.0 ;
00094   sigma_0     = 0.0 ;
00095   sigma_infty = 0.0 ;
00096   delta       = 0.0 ;
00097   Hard        = 0.0 ;
00098   eta         = 0.0 ;
00099 
00100   this->zero( ) ;     // or (*this).zero( ) 
00101 
00102   int i, j, k, l ;
00103 
00104   //zero rank4 IIdev and IbunI 
00105   for ( i = 0; i < 3; i++ ) {
00106     for ( j = 0; j < 3; j++ )  {
00107       for ( k = 0; k < 3; k++ ) {
00108         for ( l = 0; l < 3; l++)  { 
00109 
00110           IbunI[i][j][k][l] = 0.0 ;
00111 
00112           IIdev[i][j][k][l] = 0.0 ;
00113 
00114         } // end for l
00115       } // end for k
00116     } // end for j
00117   } // end for i
00118 
00119 
00120   //form rank4 IbunI 
00121 
00122   IbunI [0][0] [0][0] = 1.0 ;
00123   IbunI [0][0] [1][1] = 1.0 ;
00124   IbunI [0][0] [2][2] = 1.0 ;
00125   IbunI [1][1] [0][0] = 1.0 ;
00126   IbunI [1][1] [1][1] = 1.0 ;
00127   IbunI [1][1] [2][2] = 1.0 ;
00128   IbunI [2][2] [0][0] = 1.0 ;
00129   IbunI [2][2] [1][1] = 1.0 ;
00130   IbunI [2][2] [2][2] = 1.0 ;
00131 
00132   //form rank4 IIdev
00133 
00134   IIdev [0][0] [0][0] =  two3 ; // 0.666667 
00135   IIdev [0][0] [1][1] = -one3 ; //-0.333333 
00136   IIdev [0][0] [2][2] = -one3 ; //-0.333333 
00137   IIdev [0][1] [0][1] = 0.5 ;
00138   IIdev [0][1] [1][0] = 0.5 ;
00139   IIdev [0][2] [0][2] = 0.5 ;
00140   IIdev [0][2] [2][0] = 0.5 ;
00141   IIdev [1][0] [0][1] = 0.5 ;
00142   IIdev [1][0] [1][0] = 0.5 ;
00143   IIdev [1][1] [0][0] = -one3 ; //-0.333333 
00144   IIdev [1][1] [1][1] =  two3 ; // 0.666667 
00145   IIdev [1][1] [2][2] = -one3 ; //-0.333333 
00146   IIdev [1][2] [1][2] = 0.5 ;
00147   IIdev [1][2] [2][1] = 0.5 ;
00148   IIdev [2][0] [0][2] = 0.5 ;
00149   IIdev [2][0] [2][0] = 0.5 ;
00150   IIdev [2][1] [1][2] = 0.5 ;
00151   IIdev [2][1] [2][1] = 0.5 ;
00152   IIdev [2][2] [0][0] = -one3 ; //-0.333333 
00153   IIdev [2][2] [1][1] = -one3 ; //-0.333333 
00154   IIdev [2][2] [2][2] =  two3 ; // 0.666667 
00155 }
00156 
00157 
00158 //full constructor
00159 J2Plasticity :: J2Plasticity(int    tag,
00160                              int classTag,
00161                              double K,
00162                              double G,
00163                              double yield0,
00164                              double yield_infty,
00165                              double d,
00166                              double H,
00167                              double viscosity) 
00168 : 
00169   NDMaterial(tag, classTag),
00170   epsilon_p_n(3,3),
00171   epsilon_p_nplus1(3,3),
00172   stress(3,3),
00173   strain(3,3)
00174 {
00175   bulk        = K ;
00176   shear       = G ;
00177   sigma_0     = yield0 ;
00178   sigma_infty = yield_infty ;
00179   delta       = d ;
00180   Hard        = H ;
00181   eta         = viscosity ;
00182 
00183   this->zero( ) ;
00184 
00185   int i, j, k, l ;
00186 
00187   //zero rank4 IIdev and IbunI 
00188   for ( i = 0; i < 3; i++ ) {
00189     for ( j = 0; j < 3; j++ )  {
00190       for ( k = 0; k < 3; k++ ) {
00191         for ( l = 0; l < 3; l++)  { 
00192 
00193           IbunI[i][j][k][l] = 0.0 ;
00194 
00195           IIdev[i][j][k][l] = 0.0 ;
00196 
00197         } // end for l
00198       } // end for k
00199     } // end for j
00200   } // end for i
00201 
00202 
00203   //form rank4 IbunI 
00204 
00205   IbunI [0][0] [0][0] = 1.0 ;
00206   IbunI [0][0] [1][1] = 1.0 ;
00207   IbunI [0][0] [2][2] = 1.0 ;
00208   IbunI [1][1] [0][0] = 1.0 ;
00209   IbunI [1][1] [1][1] = 1.0 ;
00210   IbunI [1][1] [2][2] = 1.0 ;
00211   IbunI [2][2] [0][0] = 1.0 ;
00212   IbunI [2][2] [1][1] = 1.0 ;
00213   IbunI [2][2] [2][2] = 1.0 ;
00214 
00215   //form rank4 IIdev
00216 
00217   IIdev [0][0] [0][0] =  two3 ; // 0.666667 
00218   IIdev [0][0] [1][1] = -one3 ; //-0.333333 
00219   IIdev [0][0] [2][2] = -one3 ; //-0.333333 
00220   IIdev [0][1] [0][1] = 0.5 ;
00221   IIdev [0][1] [1][0] = 0.5 ;
00222   IIdev [0][2] [0][2] = 0.5 ;
00223   IIdev [0][2] [2][0] = 0.5 ;
00224   IIdev [1][0] [0][1] = 0.5 ;
00225   IIdev [1][0] [1][0] = 0.5 ;
00226   IIdev [1][1] [0][0] = -one3 ; //-0.333333 
00227   IIdev [1][1] [1][1] =  two3 ; // 0.666667 
00228   IIdev [1][1] [2][2] = -one3 ; //-0.333333 
00229   IIdev [1][2] [1][2] = 0.5 ;
00230   IIdev [1][2] [2][1] = 0.5 ;
00231   IIdev [2][0] [0][2] = 0.5 ;
00232   IIdev [2][0] [2][0] = 0.5 ;
00233   IIdev [2][1] [1][2] = 0.5 ;
00234   IIdev [2][1] [2][1] = 0.5 ;
00235   IIdev [2][2] [0][0] = -one3 ; //-0.333333 
00236   IIdev [2][2] [1][1] = -one3 ; //-0.333333 
00237   IIdev [2][2] [2][2] =  two3 ; // 0.666667 
00238 }
00239 
00240 
00241 //elastic constructor
00242 J2Plasticity :: 
00243 J2Plasticity(   int    tag, 
00244                 int  classTag,
00245                 double K, 
00246                 double G ) :
00247 NDMaterial(tag, classTag),
00248 epsilon_p_n(3,3),
00249 epsilon_p_nplus1(3,3),
00250 stress(3,3),
00251 strain(3,3)
00252 {
00253   bulk        = K ;
00254   shear       = G ; 
00255   sigma_0     = 1.0e16*shear ;
00256   sigma_infty = sigma_0 ;
00257   delta       = 0.0 ;
00258   Hard        = 0.0 ;
00259   eta         = 0.0 ;
00260 
00261   this->zero( ) ;
00262 
00263   int i, j, k, l ;
00264 
00265   //zero rank4 IIdev and IbunI 
00266   for ( i = 0; i < 3; i++ ) {
00267     for ( j = 0; j < 3; j++ )  {
00268       for ( k = 0; k < 3; k++ ) {
00269         for ( l = 0; l < 3; l++)  { 
00270 
00271           IbunI[i][j][k][l] = 0.0 ;
00272 
00273           IIdev[i][j][k][l] = 0.0 ;
00274 
00275         } // end for l
00276       } // end for k
00277     } // end for j
00278   } // end for i
00279 
00280 
00281   //form rank4 IbunI 
00282 
00283   IbunI [0][0] [0][0] = 1.0 ;
00284   IbunI [0][0] [1][1] = 1.0 ;
00285   IbunI [0][0] [2][2] = 1.0 ;
00286   IbunI [1][1] [0][0] = 1.0 ;
00287   IbunI [1][1] [1][1] = 1.0 ;
00288   IbunI [1][1] [2][2] = 1.0 ;
00289   IbunI [2][2] [0][0] = 1.0 ;
00290   IbunI [2][2] [1][1] = 1.0 ;
00291   IbunI [2][2] [2][2] = 1.0 ;
00292 
00293   //form rank4 IIdev
00294 
00295   IIdev [0][0] [0][0] =  two3 ; // 0.666667 
00296   IIdev [0][0] [1][1] = -one3 ; //-0.333333 
00297   IIdev [0][0] [2][2] = -one3 ; //-0.333333 
00298   IIdev [0][1] [0][1] = 0.5 ;
00299   IIdev [0][1] [1][0] = 0.5 ;
00300   IIdev [0][2] [0][2] = 0.5 ;
00301   IIdev [0][2] [2][0] = 0.5 ;
00302   IIdev [1][0] [0][1] = 0.5 ;
00303   IIdev [1][0] [1][0] = 0.5 ;
00304   IIdev [1][1] [0][0] = -one3 ; //-0.333333 
00305   IIdev [1][1] [1][1] =  two3 ; // 0.666667 
00306   IIdev [1][1] [2][2] = -one3 ; //-0.333333 
00307   IIdev [1][2] [1][2] = 0.5 ;
00308   IIdev [1][2] [2][1] = 0.5 ;
00309   IIdev [2][0] [0][2] = 0.5 ;
00310   IIdev [2][0] [2][0] = 0.5 ;
00311   IIdev [2][1] [1][2] = 0.5 ;
00312   IIdev [2][1] [2][1] = 0.5 ;
00313   IIdev [2][2] [0][0] = -one3 ; //-0.333333 
00314   IIdev [2][2] [1][1] = -one3 ; //-0.333333 
00315   IIdev [2][2] [2][2] =  two3 ; // 0.666667 
00316 }
00317 
00318 
00319 //destructor
00320 J2Plasticity :: ~J2Plasticity( ) 
00321 {  } 
00322 
00323 
00324 
00325 NDMaterial*
00326 J2Plasticity :: getCopy (const char *type)
00327 {
00328     if (strcmp(type,"PlaneStress2D") == 0 || strcmp(type,"PlaneStress") == 0)
00329     {
00330         J2PlaneStress  *clone ;
00331         clone = new J2PlaneStress(this->getTag(), bulk, shear, sigma_0,
00332                                   sigma_infty, delta, Hard, eta) ;
00333         return clone ;
00334     }
00335     else if (strcmp(type,"PlaneStrain2D") == 0 || strcmp(type,"PlaneStrain") == 0)
00336     {
00337         J2PlaneStrain  *clone ;
00338         clone = new J2PlaneStrain(this->getTag(), bulk, shear, sigma_0,
00339                                   sigma_infty, delta, Hard, eta) ;
00340         return clone ;
00341     }
00342     else if (strcmp(type,"AxiSymmetric2D") == 0 || strcmp(type,"AxiSymmetric") == 0)
00343     {
00344         J2AxiSymm  *clone ;
00345         clone = new J2AxiSymm(this->getTag(), bulk, shear, sigma_0,
00346                               sigma_infty, delta, Hard, eta) ;
00347         return clone ;  
00348     }
00349     else if ((strcmp(type,"ThreeDimensional") == 0) ||
00350              (strcmp(type,"3D") == 0))
00351     {
00352         J2ThreeDimensional  *clone ;
00353         clone = new J2ThreeDimensional(this->getTag(), bulk, shear, sigma_0,
00354                               sigma_infty, delta, Hard, eta) ;
00355         return clone ;  
00356     }
00357     else if ( (strcmp(type,"PlateFiber") == 0) )
00358     {
00359         J2PlateFiber  *clone ;
00360         clone = new J2PlateFiber(this->getTag(), bulk, shear, sigma_0,
00361                               sigma_infty, delta, Hard, eta) ;
00362         return clone ;  
00363     }
00364 
00365  
00366     // Handle other cases
00367     else
00368     {
00369       opserr << "J2Plasticity::getModel failed to get model: " << type << endln;
00370       return 0;
00371     }
00372 }
00373 
00374 //print out material data
00375 void J2Plasticity :: Print( OPS_Stream &s, int flag )
00376 {
00377   s << endln ;
00378   s << "J2-Plasticity : " ; 
00379   s << this->getType( ) << endln ;
00380   s << "Bulk Modulus =   " << bulk        << endln ;
00381   s << "Shear Modulus =  " << shear       << endln ;
00382   s << "Sigma_0 =        " << sigma_0     << endln ;
00383   s << "Sigma_infty =    " << sigma_infty << endln ;
00384   s << "Delta =          " << delta       << endln ;
00385   s << "H =              " << Hard        << endln ;
00386   s << "Eta =            " << eta         << endln ;
00387   s << endln ;
00388 }
00389 
00390 
00391 //--------------------Plasticity-------------------------------------
00392 
00393 //plasticity integration routine
00394 void J2Plasticity :: plastic_integrator( )
00395 {
00396   const double tolerance = (1.0e-8)*sigma_0 ;
00397 
00398   const double dt = ops_Dt ; //time step
00399 
00400   static Matrix dev_strain(3,3) ; //deviatoric strain
00401 
00402   static Matrix dev_stress(3,3) ; //deviatoric stress
00403  
00404   static Matrix normal(3,3) ;     //normal to yield surface
00405 
00406   double NbunN ; //normal bun normal 
00407 
00408   double norm_tau = 0.0 ;   //norm of deviatoric stress 
00409   double inv_norm_tau = 0.0 ;
00410 
00411   double phi = 0.0 ; //trial value of yield function
00412 
00413   double trace = 0.0 ; //trace of strain
00414 
00415   double gamma = 0.0 ; //consistency parameter
00416 
00417   double resid = 1.0 ; 
00418   double tang  = 0.0 ;
00419   
00420   double theta = 0.0 ; 
00421   double theta_inv = 0.0 ;
00422 
00423   double c1 = 0.0 ; 
00424   double c2 = 0.0 ;
00425   double c3 = 0.0 ;
00426 
00427   int i,j,k,l;
00428   int ii, jj ; 
00429 
00430   int iteration_counter ;
00431   const int max_iterations = 25 ;
00432 
00433   //compute the deviatoric strains
00434 
00435   trace = strain(0,0) + strain(1,1) + strain(2,2) ;
00436  
00437   dev_strain = strain ;
00438   for ( i = 0; i < 3; i++ )
00439     dev_strain(i,i) -= ( one3*trace ) ;
00440    
00441   //compute the trial deviatoric stresses
00442 
00443   //   dev_stress = (2.0*shear) * ( dev_strain - epsilon_p_n ) ;
00444   dev_stress = dev_strain;
00445   dev_stress -= epsilon_p_n;
00446   dev_stress *= 2.0 * shear;
00447 
00448   //compute norm of deviatoric stress
00449 
00450   norm_tau = 0.0 ;
00451   for ( i = 0; i < 3; i++ ){
00452     for ( j = 0; j < 3; j++ ) 
00453       norm_tau += dev_stress(i,j)*dev_stress(i,j) ;
00454   } //end for i 
00455 
00456   norm_tau = sqrt( norm_tau ) ;
00457 
00458   if ( norm_tau > tolerance ) {
00459     inv_norm_tau = 1.0 / norm_tau ;
00460     normal =  inv_norm_tau * dev_stress ;
00461   }
00462   else {
00463     normal.Zero( ) ;
00464     inv_norm_tau = 0.0 ;
00465   } //end if 
00466 
00467   //compute trial value of yield function
00468 
00469   phi = norm_tau -  root23 * q(xi_n) ;
00470 
00471   // check if phi > 0 
00472   
00473   if ( phi > 0.0 ) { //plastic
00474 
00475      //solve for gamma 
00476      gamma = 0.0 ;
00477      resid = 1.0 ;
00478      iteration_counter = 0 ;
00479      while ( fabs(resid) > tolerance ) {
00480 
00481         resid = norm_tau 
00482               - (2.0*shear) * gamma 
00483               - root23 * q( xi_n + root23*gamma ) 
00484               - (eta/dt) * gamma ;
00485 
00486         tang =  - (2.0*shear)  
00487                 - two3 * qprime( xi_n + root23*gamma )
00488                 - (eta/dt) ;
00489 
00490         gamma -= ( resid / tang ) ;
00491 
00492         iteration_counter++ ;
00493 
00494         if ( iteration_counter > max_iterations ) {
00495             opserr << "More than " << max_iterations ;
00496             opserr << " iterations in constituive subroutine J2-plasticity \n" ;
00497             break ;
00498         } //end if 
00499         
00500      } //end while resid
00501 
00502      gamma *= (1.0 - 1e-08) ;
00503 
00504      //update plastic internal variables
00505 
00506      epsilon_p_nplus1 = epsilon_p_n + gamma*normal ;
00507 
00508      xi_nplus1 = xi_n + root23*gamma ;
00509 
00510      //recompute deviatoric stresses 
00511 
00512      dev_stress = (2.0*shear) * ( dev_strain - epsilon_p_nplus1 ) ;
00513 
00514      //compute the terms for plastic part of tangent
00515 
00516      theta =  (2.0*shear)  
00517            +  two3 * qprime( xi_nplus1 )
00518            +  (eta/dt) ;
00519 
00520      theta_inv = 1.0/theta ;
00521 
00522   }
00523   else { //elastic 
00524 
00525     //update history variables -- they remain unchanged
00526 
00527     epsilon_p_nplus1 = epsilon_p_n ;
00528 
00529     xi_nplus1 = xi_n ;
00530 
00531     //no extra tangent terms to compute 
00532     
00533     gamma = 0.0 ; 
00534     theta = 0.0 ;
00535     theta_inv = 0.0 ;
00536 
00537   } //end if phi > 0
00538 
00539 
00540   //add on bulk part of stress
00541 
00542   stress = dev_stress ;
00543   for ( i = 0; i < 3; i++ )
00544      stress(i,i) += bulk*trace ;
00545 
00546   //compute the tangent
00547 
00548   c1 = -4.0 * shear * shear ;
00549   c2 = c1 * theta_inv ;
00550   c3 = c1 * gamma * inv_norm_tau ;
00551 
00552   for ( ii = 0; ii < 6; ii++ ) {
00553     for ( jj = 0; jj < 6; jj++ )  {
00554 
00555           index_map( ii, i, j ) ;
00556           index_map( jj, k, l ) ;
00557 
00558           NbunN  = normal(i,j)*normal(k,l) ; 
00559 
00560           //elastic terms
00561           tangent[i][j][k][l]  = bulk * IbunI[i][j][k][l] ;
00562 
00563           tangent[i][j][k][l] += (2.0*shear) * IIdev[i][j][k][l] ;
00564 
00565           //plastic terms 
00566           tangent[i][j][k][l] += c2 * NbunN ;
00567 
00568           tangent[i][j][k][l] += c3 * (  IIdev[i][j][k][l] - NbunN ) ;
00569 
00570           //minor symmetries 
00571           tangent [j][i][k][l] = tangent[i][j][k][l] ;
00572           tangent [i][j][l][k] = tangent[i][j][k][l] ;
00573           tangent [j][i][l][k] = tangent[i][j][k][l] ;
00574 
00575     } // end for jj
00576   } // end for ii
00577 
00578   return ;
00579 } 
00580 
00581 
00582 
00583 
00584 
00585 // set up for initial elastic
00586 void J2Plasticity :: doInitialTangent( )
00587 {
00588   int ii,jj,i,j,k,l;
00589 
00590   //compute the deviatoric strains
00591   for ( ii = 0; ii < 6; ii++ ) {
00592     for ( jj = 0; jj < 6; jj++ )  {
00593 
00594           index_map( ii, i, j ) ;
00595           index_map( jj, k, l ) ;
00596 
00597           //elastic terms
00598           initialTangent[i][j][k][l]  = bulk * IbunI[i][j][k][l] ;
00599           initialTangent[i][j][k][l] += (2.0*shear) * IIdev[i][j][k][l] ;
00600 
00601           //minor symmetries 
00602           //minor symmetries 
00603           initialTangent [j][i][k][l] = initialTangent[i][j][k][l] ;
00604           initialTangent [i][j][l][k] = initialTangent[i][j][k][l] ;
00605           initialTangent [j][i][l][k] = initialTangent[i][j][k][l] ;
00606 
00607     } // end for jj
00608   } // end for ii
00609 
00610   return ;
00611 } 
00612 
00613 
00614 
00615 //hardening function
00616 double J2Plasticity :: q( double xi ) 
00617 {
00618 //  q(xi) = simga_infty + (sigma_0 - sigma_infty)*exp(-delta*xi) + H*xi 
00619 
00620  return    sigma_infty
00621          + (sigma_0 - sigma_infty)*exp(-delta*xi)
00622          + Hard*xi ;
00623 }
00624 
00625 
00626 //hardening function derivative
00627 double J2Plasticity :: qprime( double xi )
00628 {
00629   return  (sigma_0 - sigma_infty) * (-delta) * exp(-delta*xi)
00630          + Hard ;
00631 }
00632 
00633 
00634 //matrix_index ---> tensor indices i,j
00635 void J2Plasticity :: index_map( int matrix_index, int &i, int &j )
00636 {
00637   switch ( matrix_index+1 ) { //add 1 for standard tensor indices
00638 
00639     case 1 :
00640       i = 1 ; 
00641       j = 1 ;
00642       break ;
00643  
00644     case 2 :
00645       i = 2 ;
00646       j = 2 ; 
00647       break ;
00648 
00649     case 3 :
00650       i = 3 ;
00651       j = 3 ;
00652       break ;
00653 
00654     case 4 :
00655       i = 1 ;
00656       j = 2 ;
00657       break ;
00658 
00659     case 5 :
00660       i = 2 ;
00661       j = 3 ;
00662       break ;
00663 
00664     case 6 :
00665       i = 3 ;
00666       j = 1 ;
00667       break ;
00668 
00669 
00670     default :
00671       i = 1 ;
00672       j = 1 ;
00673       break ;
00674 
00675   } //end switch
00676 
00677 i-- ; //subtract 1 for C-indexing
00678 j-- ;
00679 
00680 return ; 
00681 }
00682 
00683 
00684 NDMaterial*
00685 J2Plasticity::getCopy (void)
00686 {
00687   opserr << "J2Plasticity::getCopy -- subclass responsibility\n"; 
00688   exit(-1);
00689   return 0;
00690 }
00691 
00692 const char*
00693 J2Plasticity::getType (void) const
00694 {
00695     opserr << "J2Plasticity::getType -- subclass responsibility\n";
00696     exit(-1);
00697     return 0;
00698 }
00699 
00700 int
00701 J2Plasticity::getOrder (void) const
00702 {
00703     opserr << "J2Plasticity::getOrder -- subclass responsibility\n";
00704     exit(-1);
00705     return 0;
00706 }
00707 
00708 
00709 int 
00710 J2Plasticity::commitState( ) 
00711 {
00712   epsilon_p_n = epsilon_p_nplus1 ;
00713   xi_n        = xi_nplus1 ;
00714 
00715   return 0;
00716 }
00717 
00718 int 
00719 J2Plasticity::revertToLastCommit( ) 
00720 {
00721   return 0;
00722 }
00723 
00724 
00725 int 
00726 J2Plasticity::revertToStart( ) {
00727 
00728   this->zero( ) ;
00729 
00730   return 0;
00731 }
00732 
00733 int
00734 J2Plasticity::sendSelf(int commitTag, Channel &theChannel)
00735 {
00736   // we place all the data needed to define material and it's state
00737   // int a vector object
00738   static Vector data(9);
00739   int cnt = 0;
00740   data(cnt++) = this->getTag();
00741   data(cnt++) = bulk;
00742   data(cnt++) = shear;
00743   data(cnt++) = sigma_0;
00744   data(cnt++) = sigma_infty;
00745   data(cnt++) = delta;
00746   data(cnt++) = Hard;
00747   data(cnt++) = eta;
00748   data(cnt++) = xi_n;
00749 
00750   // send the vector object to the channel
00751   if (theChannel.sendVector(this->getDbTag(), commitTag, data) < 0) {
00752     opserr << "J2Plasticity::sendSelf - failed to send vector to channel\n";
00753     return -1;
00754   }
00755 
00756   return 0;
00757 }
00758 
00759 int
00760 J2Plasticity::recvSelf (int commitTag, Channel &theChannel, 
00761                          FEM_ObjectBroker &theBroker)
00762 {
00763 
00764   // recv the vector object from the channel which defines material param and state
00765   static Vector data(9);
00766   if (theChannel.recvVector(this->getDbTag(), commitTag, data) < 0) {
00767     opserr << "J2Plasticity::recvSelf - failed to recv vector from channel\n";
00768     return -1;
00769   }
00770 
00771   // set the material parameters and state variables
00772   int cnt = 0;
00773   this->setTag(data(cnt++));
00774   bulk = data(cnt++);
00775   shear = data(cnt++);
00776   sigma_0 = data(cnt++);
00777   sigma_infty = data(cnt++);
00778   delta = data(cnt++);
00779   Hard = data(cnt++);
00780   eta = data(cnt++);
00781   xi_n = data(cnt++);
00782 
00783   epsilon_p_nplus1 = epsilon_p_n;
00784   xi_nplus1        = xi_n;
00785 
00786   return 0;
00787 }

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