TotalLagrangianFD8NodeBrick.cpp

Go to the documentation of this file.
00001 //===============================================================================
00002 
00003 //# COPYRIGHT (C): Woody's license (by BJ):
00004 
00005 //                 ``This    source  code is Copyrighted in
00006 
00007 //                 U.S.,  for  an  indefinite  period,  and anybody
00008 
00009 //                 caught  using it without our permission, will be
00010 
00011 //                 mighty good friends of ourn, cause we don't give
00012 
00013 //                 a  darn.  Hack it. Compile it. Debug it. Run it.
00014 
00015 //                 Yodel  it.  Enjoy it. We wrote it, that's all we
00016 
00017 //                 wanted to do.''
00018 
00019 //
00020 
00021 //# PROJECT:           Object Oriented Finite Element Program
00022 
00023 //# PURPOSE:           Finite Deformation Hyper-Elastic classes
00024 
00025 //# CLASS:
00026 
00027 //#
00028 
00029 //# VERSION:           0.6_(1803398874989) (golden section)
00030 
00031 //# LANGUAGE:          C++
00032 
00033 //# TARGET OS:         all...
00034 
00035 //# DESIGN:            Zhao Cheng, Boris Jeremic (jeremic@ucdavis.edu)
00036 
00037 //# PROGRAMMER(S):     Zhao Cheng, Boris Jeremic
00038 
00039 //#
00040 
00041 //#
00042 
00043 //# DATE:              Sept2005              
00044 
00045 //# UPDATE HISTORY:   
00046 
00047 //#                     
00048 
00049 //#                    
00050 
00051 //#
00052 
00053 //===============================================================================
00054 
00055 #ifndef TOTALLAGRANGIANFD8NODEBRICK_CPP
00056 
00057 #define TOTALLAGRANGIANFD8NODEBRICK_CPP
00058 
00059 
00060 
00061 #include <TotalLagrangianFD8NodeBrick.h>
00062 
00063 
00064 
00065 const int TotalLagrangianFD8NodeBrick::NumIntegrationPts = 2;
00066 
00067 const int TotalLagrangianFD8NodeBrick::NumTotalGaussPts = 8;
00068 
00069 const int TotalLagrangianFD8NodeBrick::NumNodes = 8;
00070 
00071 const int TotalLagrangianFD8NodeBrick::NumDof = 3;
00072 
00073 const int TotalLagrangianFD8NodeBrick::NumElemDof = NumNodes*NumDof;
00074 
00075 
00076 
00077 Matrix TotalLagrangianFD8NodeBrick::K(NumElemDof, NumElemDof);
00078 
00079 Matrix TotalLagrangianFD8NodeBrick::M(NumElemDof, NumElemDof);
00080 
00081 Vector TotalLagrangianFD8NodeBrick::P(NumElemDof);
00082 
00083 const double TotalLagrangianFD8NodeBrick::pts[2] = {-0.577350269189626, +0.577350269189626};
00084 
00085 const double TotalLagrangianFD8NodeBrick::wts[2] = {1.0, 1.0};
00086 
00087 
00088 
00089 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00090 
00091 
00092 
00093 TotalLagrangianFD8NodeBrick::TotalLagrangianFD8NodeBrick(int tag,
00094 
00095 int node_numb_1,  int node_numb_2,  int node_numb_3,  int node_numb_4,
00096 
00097 int node_numb_5,  int node_numb_6,  int node_numb_7,  int node_numb_8,
00098 
00099 NDMaterial &m, double b1, double b2, double b3)
00100 
00101 :Element(tag, ELE_TAG_TotalLagrangianFD8NodeBrick ),
00102 
00103  theMaterial(0), connectedExternalNodes(NumNodes), Q(0), bf(NumDof), Ki(0)
00104 
00105 {
00106 
00107       connectedExternalNodes( 0) = node_numb_1;
00108 
00109       connectedExternalNodes( 1) = node_numb_2;
00110 
00111       connectedExternalNodes( 2) = node_numb_3;
00112 
00113       connectedExternalNodes( 3) = node_numb_4;
00114 
00115       connectedExternalNodes( 4) = node_numb_5;
00116 
00117       connectedExternalNodes( 5) = node_numb_6;
00118 
00119       connectedExternalNodes( 6) = node_numb_7;
00120 
00121       connectedExternalNodes( 7) = node_numb_8;
00122 
00123 
00124 
00125       bf(0) = b1;
00126 
00127       bf(1) = b2;
00128 
00129       bf(2) = b3;
00130 
00131 
00132 
00133       theMaterial = new NDMaterial *[NumTotalGaussPts];
00134 
00135 
00136 
00137       if (theMaterial == 0) {
00138 
00139        opserr<<"FiniteDeformationElastic3D::FiniteDeformationElastic3D -- failed allocate material model pointer\n";
00140 
00141        exit(-1);
00142 
00143       }
00144 
00145       
00146 
00147       int i;
00148 
00149       for (i=0; i<NumTotalGaussPts; i++) {
00150 
00151        theMaterial[i] = m.getCopy();
00152 
00153        if (theMaterial[i] == 0) {
00154 
00155         opserr<<"FiniteDeformationElastic3D::FiniteDeformationElastic3D -- failed allocate material model pointer\n";
00156 
00157         exit(-1);
00158 
00159        }
00160 
00161       }
00162 
00163 
00164 
00165       rho = m.getRho();
00166 
00167 
00168 
00169       for (i=0; i<NumNodes; i++) theNodes[i] = 0;
00170 
00171 
00172 
00173 }
00174 
00175 
00176 
00177 //-------------------------------------------------------------------------------------------
00178 
00179 TotalLagrangianFD8NodeBrick::TotalLagrangianFD8NodeBrick ()
00180 
00181 :Element(0, ELE_TAG_TotalLagrangianFD8NodeBrick ),
00182 
00183  theMaterial(0), connectedExternalNodes(NumNodes), Q(0), bf(NumDof), Ki(0)
00184 
00185 {    
00186 
00187          int i;
00188 
00189      for (i=0; i<NumNodes; i++) {  
00190 
00191        theNodes[i] = 0;
00192 
00193      }
00194 
00195 
00196 
00197      bf(0) = 0.0;
00198 
00199      bf(1) = 0.0;
00200 
00201      bf(2) = 0.0;
00202 
00203 
00204 
00205      rho = 0.0;
00206 
00207 }
00208 
00209 
00210 
00211 //-------------------------------------------------------------------------------------------------
00212 
00213 TotalLagrangianFD8NodeBrick::~TotalLagrangianFD8NodeBrick ()
00214 
00215 {   
00216 
00217         int i;
00218 
00219         for (i=0; i<NumTotalGaussPts; i++) {
00220 
00221       if (theMaterial[i]) 
00222 
00223         delete theMaterial[i];
00224 
00225     }
00226 
00227 
00228 
00229     if(theMaterial) 
00230 
00231       delete [] theMaterial;
00232 
00233 
00234 
00235     if(Ki) 
00236 
00237       delete Ki;
00238 
00239 
00240 
00241 }
00242 
00243 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00244 
00245 
00246 
00247 
00248 
00249 //=============================================================================
00250 
00251 int TotalLagrangianFD8NodeBrick::getNumExternalNodes () const
00252 
00253 {
00254 
00255     return NumNodes;
00256 
00257 }
00258 
00259 
00260 
00261 //=============================================================================
00262 
00263 const ID& TotalLagrangianFD8NodeBrick::getExternalNodes ()
00264 
00265 {
00266 
00267     return connectedExternalNodes;
00268 
00269 }
00270 
00271 
00272 
00273 //=============================================================================
00274 
00275 Node **TotalLagrangianFD8NodeBrick::getNodePtrs(void)
00276 
00277 {
00278 
00279     return theNodes;
00280 
00281 }
00282 
00283 
00284 
00285 //=============================================================================
00286 
00287 int TotalLagrangianFD8NodeBrick::getNumDOF ()
00288 
00289 {
00290 
00291     return NumElemDof;
00292 
00293 }
00294 
00295 
00296 
00297 //=============================================================================
00298 
00299 void TotalLagrangianFD8NodeBrick::setDomain (Domain *theDomain)
00300 
00301 {
00302 
00303     int i;
00304 
00305 
00306 
00307     // Check Domain is not null - invoked when object removed from a domain
00308 
00309     if (theDomain == 0) {
00310 
00311       for (i=0; i<NumNodes; i++) {
00312 
00313              theNodes[i] = 0; 
00314 
00315       }
00316 
00317       return;
00318 
00319     }
00320 
00321     
00322 
00323     for (i=0; i<NumNodes; i++) {
00324 
00325           theNodes[i] = theDomain->getNode(connectedExternalNodes(i));
00326 
00327           if ( theNodes[i]==0 ) {
00328 
00329                   opserr << "FATAL ERROR TotalLagrangianFD8NodeBrick (tag: " << this->getTag() <<
00330 
00331           " ), node not found in domain\n";
00332 
00333           exit(-1);               
00334 
00335           }        
00336 
00337     }
00338 
00339 
00340 
00341     this->DomainComponent::setDomain(theDomain);  // Very Important!!
00342 
00343 
00344 
00345     for (i=0; i<NumNodes; i++) {
00346 
00347           if ( theNodes[i]->getNumberDOF() != NumDof ) {
00348 
00349         opserr << "FATAL ERROR TotalLagrangianFD8NodeBrick (tag: " << this->getTag() <<
00350 
00351         "), has differing number of DOFs at its nodes\n";
00352 
00353         exit(-1);                                   
00354 
00355           }           
00356 
00357     }
00358 
00359 
00360 
00361 }
00362 
00363 
00364 
00365 //=============================================================================
00366 
00367 int TotalLagrangianFD8NodeBrick::commitState ()
00368 
00369 {
00370 
00371   int retVal = 0;
00372 
00373 
00374 
00375   if ((retVal = this->Element::commitState()) != 0)
00376 
00377     opserr << "TotalLagrangianFD8NodeBrick::commitState () - failed in base class";
00378 
00379 
00380 
00381   int i;
00382 
00383   for (i=0; i<NumTotalGaussPts; i++)
00384 
00385     retVal += theMaterial[i]->commitState();
00386 
00387 
00388 
00389   return retVal;
00390 
00391 }
00392 
00393 
00394 
00395 //=============================================================================
00396 
00397 int TotalLagrangianFD8NodeBrick::revertToLastCommit ()
00398 
00399 {
00400 
00401     int retVal = 0;
00402 
00403     
00404 
00405     int i;
00406 
00407     for (i=0; i<NumTotalGaussPts; i++)
00408 
00409        retVal += theMaterial[i]->revertToLastCommit();
00410 
00411 
00412 
00413     return retVal;
00414 
00415 }
00416 
00417 
00418 
00419 //=============================================================================
00420 
00421 int TotalLagrangianFD8NodeBrick::revertToStart ()
00422 
00423 {
00424 
00425     int retVal = 0;
00426 
00427     
00428 
00429     int i;
00430 
00431     for (i=0; i<NumTotalGaussPts; i++)
00432 
00433        retVal += theMaterial[i]->revertToStart();
00434 
00435 
00436 
00437     return retVal;
00438 
00439 
00440 
00441 }
00442 
00443 
00444 
00445 //=============================================================================
00446 
00447 int TotalLagrangianFD8NodeBrick::update ()
00448 
00449 {
00450 
00451     int ret = 0;
00452 
00453     tensor dh;
00454 
00455     tensor dH_dX;
00456 
00457     int where = 0;
00458 
00459     int GP_c_r, GP_c_s, GP_c_t;
00460 
00461     double r = 0.0;
00462 
00463     double s = 0.0;
00464 
00465     double t = 0.0;
00466 
00467     tensor I_ij("I", 2, def_dim_2);
00468 
00469     tensor currentF;
00470 
00471     tensor updatedF;
00472 
00473 
00474 
00475     tensor InitialNodesCrds = this->getNodesCrds();
00476 
00477     tensor CurrentNodesDisp = this->getNodesDisp();
00478 
00479 
00480 
00481     for( GP_c_r = 0 ; GP_c_r < NumIntegrationPts ; GP_c_r++ ) {
00482 
00483       r = pts[GP_c_r ];
00484 
00485       for( GP_c_s = 0 ; GP_c_s < NumIntegrationPts ; GP_c_s++ ) {
00486 
00487         s = pts[GP_c_s ];
00488 
00489         for( GP_c_t = 0 ; GP_c_t < NumIntegrationPts ; GP_c_t++ ) {
00490 
00491           t = pts[GP_c_t ];
00492 
00493           where =(GP_c_r * NumIntegrationPts + GP_c_s) * NumIntegrationPts + GP_c_t;
00494 
00495           //dh = shapeFunctionDerivative(r,s,t);
00496 
00497           dH_dX = this->dh_Global(r,s,t);
00498 
00499           currentF = CurrentNodesDisp("Ia") * dH_dX("Ib");
00500 
00501             currentF.null_indices();
00502 
00503           updatedF = currentF + I_ij;
00504 
00505           ret += theMaterial[where]->setTrialF(updatedF);
00506 
00507         }
00508 
00509       }
00510 
00511     }
00512 
00513     return ret;
00514 
00515 }
00516 
00517 
00518 
00519 //======================================================================
00520 
00521 tensor TotalLagrangianFD8NodeBrick::Jacobian_3D(double x, double y, double z)
00522 
00523   {
00524 
00525      tensor N_C = this->getNodesCrds();
00526 
00527      tensor dh = this->shapeFunctionDerivative(x, y, z);
00528 
00529      
00530 
00531      tensor J3D = N_C("ki") * dh("kj");
00532 
00533        J3D.null_indices();
00534 
00535      return J3D;
00536 
00537   }
00538 
00539 
00540 
00541 //======================================================================
00542 
00543 tensor TotalLagrangianFD8NodeBrick::Jacobian_3Dinv(double x, double y, double z)
00544 
00545   {
00546 
00547      tensor J = this->Jacobian_3D(x,y,z);
00548 
00549      return J.inverse();
00550 
00551   }
00552 
00553 
00554 
00555 //======================================================================
00556 
00557 tensor TotalLagrangianFD8NodeBrick::dh_Global(double x, double y, double z)
00558 
00559   {
00560 
00561      tensor JacobianINV0 = this->Jacobian_3Dinv(x, y, z);
00562 
00563      tensor dh = this->shapeFunctionDerivative(x, y, z);
00564 
00565      tensor  dhGlobal = dh("ik") * JacobianINV0("kj");
00566 
00567        dhGlobal.null_indices();
00568 
00569      return dhGlobal;
00570 
00571   }
00572 
00573 
00574 
00575 //======================================================================
00576 
00577 tensor TotalLagrangianFD8NodeBrick::getStiffnessTensor(void)
00578 
00579   {
00580 
00581     tensor tI2("I", 2, def_dim_2);
00582 
00583           
00584 
00585         int K_dim[] = {NumNodes,NumDof,NumDof,NumNodes};
00586 
00587     tensor Kk(4,K_dim,0.0);
00588 
00589 
00590 
00591     double r  = 0.0;
00592 
00593     double rw = 0.0;
00594 
00595     double s  = 0.0;
00596 
00597     double sw = 0.0;
00598 
00599     double t  = 0.0;
00600 
00601     double tw = 0.0;
00602 
00603    
00604 
00605     int where = 0;
00606 
00607     int GP_c_r, GP_c_s, GP_c_t;
00608 
00609     double weight = 0.0;
00610 
00611 
00612 
00613     int dh_dim[] = {NumNodes,NumDof};
00614 
00615     tensor dh(2, dh_dim, 0.0);
00616 
00617     stresstensor PK2Stress;
00618 
00619     tensor L2;
00620 
00621 
00622 
00623     double det_of_Jacobian = 0.0;
00624 
00625 
00626 
00627     tensor Jacobian;
00628 
00629     tensor dhGlobal;
00630 
00631     tensor nodesDisp;
00632 
00633     tensor F;
00634 
00635     //tensor temp01;
00636 
00637     tensor temp02;
00638 
00639     tensor temp03;
00640 
00641     tensor temp04; 
00642 
00643     tensor temp05;
00644 
00645     tensor temp06;
00646 
00647 
00648 
00649     nodesDisp = this->getNodesDisp( );
00650 
00651 
00652 
00653     for( GP_c_r = 0 ; GP_c_r < NumIntegrationPts ; GP_c_r++ ) {
00654 
00655       r = pts[GP_c_r ];
00656 
00657       rw = wts[GP_c_r ];
00658 
00659       for( GP_c_s = 0 ; GP_c_s < NumIntegrationPts ; GP_c_s++ ) {
00660 
00661         s = pts[GP_c_s ];
00662 
00663         sw = wts[GP_c_s ];
00664 
00665         for( GP_c_t = 0 ; GP_c_t < NumIntegrationPts ; GP_c_t++ ) {
00666 
00667           t = pts[GP_c_t ];
00668 
00669           tw = wts[GP_c_t ];
00670 
00671           where =(GP_c_r * NumIntegrationPts + GP_c_s) * NumIntegrationPts + GP_c_t;
00672 
00673           //dh = shapeFunctionDerivative(r,s,t);
00674 
00675           Jacobian = this->Jacobian_3D(r,s,t);
00676 
00677           det_of_Jacobian  = Jacobian.determinant();
00678 
00679           dhGlobal = this->dh_Global(r,s,t);
00680 
00681           weight = rw * sw * tw * det_of_Jacobian;
00682 
00683           PK2Stress = theMaterial[where]->getStressTensor();
00684 
00685           L2 = theMaterial[where]->getTangentTensor();
00686 
00687           F = theMaterial[where]->getF();
00688 
00689                         
00690 
00691           //K1
00692 
00693           temp04 = dhGlobal("Pb") * tI2("mn");
00694 
00695             temp04.null_indices(); 
00696 
00697           temp02 = PK2Stress("bd") * dhGlobal("Qd");   
00698 
00699             temp02.null_indices();
00700 
00701           temp06 = temp04("Pbmn") * temp02("bQ") * weight;
00702 
00703             temp06.null_indices(); 
00704 
00705           Kk += temp06;
00706 
00707                         
00708 
00709           //K2
00710 
00711           temp03 =  dhGlobal("Pb") * F("ma");
00712 
00713             temp03.null_indices(); 
00714 
00715           temp04 = F("nc") * L2("abcd");
00716 
00717             temp04.null_indices(); 
00718 
00719           temp05 = temp04("nabd") * dhGlobal("Qd"); 
00720 
00721             temp05.null_indices(); 
00722 
00723           temp06 = temp03("Pbma") * temp05("nabQ") * weight;
00724 
00725             temp06.null_indices(); 
00726 
00727           Kk += temp06;
00728 
00729         }
00730 
00731       }
00732 
00733     }
00734 
00735 
00736 
00737     return Kk;
00738 
00739   }
00740 
00741 
00742 
00743 //======================================================================
00744 
00745 tensor TotalLagrangianFD8NodeBrick::getRtensor(void)
00746 
00747   {
00748 
00749     int R_dim[] = {NumNodes,NumDof};
00750 
00751     tensor Rr(2,R_dim,0.0);
00752 
00753 
00754 
00755     double r  = 0.0;
00756 
00757     double rw = 0.0;
00758 
00759     double s  = 0.0;
00760 
00761     double sw = 0.0;
00762 
00763     double t  = 0.0;
00764 
00765     double tw = 0.0;
00766 
00767 
00768 
00769     int where = 0;
00770 
00771     int GP_c_r, GP_c_s, GP_c_t;
00772 
00773     double weight = 0.0;
00774 
00775 
00776 
00777     int dh_dim[] = {NumNodes,NumDof};
00778 
00779     tensor dh(2, dh_dim, 0.0);
00780 
00781 
00782 
00783     double det_of_Jacobian = 0.0;
00784 
00785 
00786 
00787     tensor Jacobian;
00788 
00789     tensor JacobianINV;
00790 
00791     tensor dhGlobal;
00792 
00793     tensor currentF;
00794 
00795     tensor nodesDisp;
00796 
00797     //stresstensor PK2Stress;
00798 
00799     tensor temp01;
00800 
00801     tensor temp02;
00802 
00803     tensor F;
00804 
00805 
00806 
00807     nodesDisp = this->getNodesDisp( );
00808 
00809 
00810 
00811     for( GP_c_r = 0 ; GP_c_r < NumIntegrationPts ; GP_c_r++ ) {
00812 
00813           r = pts[GP_c_r ];
00814 
00815       rw = wts[GP_c_r ];
00816 
00817       for( GP_c_s = 0 ; GP_c_s < NumIntegrationPts ; GP_c_s++ ) {
00818 
00819         s = pts[GP_c_s ];
00820 
00821         sw = wts[GP_c_s ];
00822 
00823         for( GP_c_t = 0 ; GP_c_t < NumIntegrationPts ; GP_c_t++ ) {
00824 
00825           t = pts[GP_c_t ];
00826 
00827           tw = wts[GP_c_t ];
00828 
00829           where =(GP_c_r * NumIntegrationPts + GP_c_s) * NumIntegrationPts + GP_c_t;
00830 
00831           //dh = shapeFunctionDerivative(r,s,t);
00832 
00833           Jacobian = this->Jacobian_3D(r,s,t);
00834 
00835           det_of_Jacobian  = Jacobian.determinant();
00836 
00837           dhGlobal = this->dh_Global(r,s,t);
00838 
00839           weight = rw * sw * tw * det_of_Jacobian;
00840 
00841           //PK2Stress = theMaterial[where]->getStressTensor();
00842 
00843           //F = theMaterial[where]->getF();
00844 
00845           //temp01 = PK2Stress("ik") * F("jk");
00846 
00847           //  temp01.null_indices();
00848 
00849           temp01 = theMaterial[where]->getPK1StressTensor();
00850 
00851           temp02 = dhGlobal("PJ") * temp01("iJ") * weight;
00852 
00853             temp02.null_indices(); 
00854 
00855           Rr += temp02;
00856 
00857         }
00858 
00859       }
00860 
00861     }
00862 
00863 
00864 
00865     return Rr;
00866 
00867   }
00868 
00869 
00870 
00871 //======================================================================
00872 
00873 tensor TotalLagrangianFD8NodeBrick::getBodyForce(void)
00874 
00875   {
00876 
00877     int B_dim[] = {NumNodes,NumDof};
00878 
00879     tensor Bb(2,B_dim,0.0);
00880 
00881 
00882 
00883     double r  = 0.0;
00884 
00885     double rw = 0.0;
00886 
00887     double s  = 0.0;
00888 
00889     double sw = 0.0;
00890 
00891     double t  = 0.0;
00892 
00893     double tw = 0.0;
00894 
00895 
00896 
00897     int where = 0;
00898 
00899     int GP_c_r, GP_c_s, GP_c_t;
00900 
00901     double weight = 0.0;
00902 
00903 
00904 
00905     int h_dim[] = {20};
00906 
00907     tensor h(1, h_dim, 0.0);
00908 
00909     int dh_dim[] = {NumNodes,NumDof};
00910 
00911     tensor dh(2, dh_dim, 0.0);
00912 
00913     int bodyforce_dim[] = {3};
00914 
00915     tensor bodyforce(1, bodyforce_dim, 0.0);
00916 
00917 
00918 
00919     double det_of_Jacobian = 0.0;
00920 
00921 
00922 
00923     tensor Jacobian;
00924 
00925     tensor JacobianINV;
00926 
00927 
00928 
00929     bodyforce.val(1) = bf(0);
00930 
00931     bodyforce.val(2) = bf(1);
00932 
00933     bodyforce.val(3) = bf(2);
00934 
00935 
00936 
00937     for( GP_c_r = 0 ; GP_c_r < NumIntegrationPts ; GP_c_r++ ) {
00938 
00939       r = pts[GP_c_r ];
00940 
00941       rw = wts[GP_c_r ];
00942 
00943       for( GP_c_s = 0 ; GP_c_s < NumIntegrationPts ; GP_c_s++ ) {
00944 
00945         s = pts[GP_c_s ];
00946 
00947         sw = wts[GP_c_s ];
00948 
00949         for( GP_c_t = 0 ; GP_c_t < NumIntegrationPts ; GP_c_t++ ) {
00950 
00951           t = pts[GP_c_t ];
00952 
00953           tw = wts[GP_c_t ];
00954 
00955           where =(GP_c_r * NumIntegrationPts + GP_c_s) * NumIntegrationPts + GP_c_t;
00956 
00957           h = shapeFunction(r,s,t);
00958 
00959           dh = shapeFunctionDerivative(r,s,t);
00960 
00961           Jacobian = this->Jacobian_3D(r,s,t);
00962 
00963           det_of_Jacobian  = Jacobian.determinant();
00964 
00965           weight = rw * sw * tw * det_of_Jacobian;
00966 
00967           Bb = Bb +  h("P") * bodyforce("i") * rho *weight;
00968 
00969              Bb.null_indices();
00970 
00971         }
00972 
00973       }
00974 
00975     }
00976 
00977     return Bb;
00978 
00979   }
00980 
00981 
00982 
00983 //======================================================================
00984 
00985 tensor TotalLagrangianFD8NodeBrick::getSurfaceForce(void)
00986 
00987   {
00988 
00989     int S_dim[] = {NumNodes, NumDof};
00990 
00991     tensor Ss(2,S_dim,0.0);
00992 
00993     // Need Work Here!
00994 
00995 
00996 
00997     return Ss;
00998 
00999   }
01000 
01001 
01002 
01003 //============================================================================
01004 
01005 tensor TotalLagrangianFD8NodeBrick::getForces(void)
01006 
01007   {
01008 
01009     int F_dim[] = {NumNodes,NumDof};
01010 
01011     tensor Ff(2,F_dim,0.0);
01012 
01013 
01014 
01015     Ff = this->getBodyForce( ) + this->getSurfaceForce( );
01016 
01017 
01018 
01019     return Ff;
01020 
01021   }
01022 
01023 
01024 
01025 //=============================================================================
01026 
01027 const Matrix &TotalLagrangianFD8NodeBrick::getTangentStiff ()
01028 
01029 {
01030 
01031      K.Zero();
01032 
01033 
01034 
01035      tensor stifftensor = this->getStiffnessTensor();
01036 
01037 
01038 
01039      int kki=0;
01040 
01041      int kkj=0;
01042 
01043      
01044 
01045      int i, j, k, l;
01046 
01047      for (i=1 ; i<=NumNodes ; i++ ) {
01048 
01049         for (j=1 ; j<=NumNodes ; j++ ) {
01050 
01051            for (k=1 ; k<=NumDof ; k++ ) {
01052 
01053               for (l=1 ; l<=NumDof ; l++ ) {
01054 
01055                  kki = k + NumDof*(i-1);
01056 
01057                  kkj = l + NumDof*(j-1);
01058 
01059                  K(kki-1 , kkj-1) = stifftensor.cval(i,k,l,j); 
01060 
01061               }
01062 
01063            }
01064 
01065         }
01066 
01067      }
01068 
01069 
01070 
01071      return K;
01072 
01073 }
01074 
01075 
01076 
01077 //=============================================================================
01078 
01079 const Matrix &TotalLagrangianFD8NodeBrick::getInitialStiff ()
01080 
01081 {
01082 
01083      if (Ki != 0) return *Ki;
01084 
01085 
01086 
01087      K.Zero();
01088 
01089      K = this->getTangentStiff ();
01090 
01091 
01092 
01093      Ki = new Matrix(K);
01094 
01095 
01096 
01097      return K;
01098 
01099 }
01100 
01101 
01102 
01103 //=============================================================================
01104 
01105 const Matrix &TotalLagrangianFD8NodeBrick::getMass ()
01106 
01107 {
01108 
01109     // Need Work Here
01110 
01111     M.Zero();
01112 
01113     return M;
01114 
01115 }
01116 
01117 
01118 
01119 //======================================================================
01120 
01121 tensor TotalLagrangianFD8NodeBrick::getNodesCrds(void) 
01122 
01123 {
01124 
01125     const int dimensions[] = {NumNodes, NumDof};
01126 
01127     tensor N_coord(2, dimensions, 0.0);
01128 
01129 
01130 
01131     int i, j;
01132 
01133     for (i=0; i<NumNodes; i++) {
01134 
01135           const Vector &TNodesCrds = theNodes[i]->getCrds();
01136 
01137       for (j=0; j<NumDof; j++) {
01138 
01139         N_coord.val(i+1, j+1) = TNodesCrds(j);
01140 
01141           }                 
01142 
01143     }
01144 
01145     
01146 
01147     return N_coord;
01148 
01149 }
01150 
01151 
01152 
01153 //=============================================================================================
01154 
01155 tensor TotalLagrangianFD8NodeBrick::getNodesDisp(void)
01156 
01157   {
01158 
01159     int i, j;
01160 
01161     int dimU[] = {NumNodes, NumDof};
01162 
01163     tensor total_disp(2, dimU, 0.0);
01164 
01165 
01166 
01167     for (i=0; i<NumNodes; i++) {
01168 
01169       const Vector &TNodesDisp = theNodes[i]->getTrialDisp();
01170 
01171       for (j=0; j<NumDof; j++) {
01172 
01173         total_disp.val(i+1, j+1) = TNodesDisp(j);
01174 
01175       }
01176 
01177     }
01178 
01179 
01180 
01181     return total_disp;
01182 
01183   }
01184 
01185 
01186 
01187 //=============================================================================
01188 
01189 void TotalLagrangianFD8NodeBrick::zeroLoad(void)
01190 
01191 {
01192 
01193 
01194 
01195         if ( Q != 0 )
01196 
01197           Q->Zero();
01198 
01199     
01200 
01201         return;
01202 
01203 }
01204 
01205 
01206 
01207 
01208 
01209 //=============================================================================
01210 
01211 int
01212 
01213 TotalLagrangianFD8NodeBrick::addLoad(ElementalLoad *theLoad, double loadFactor)
01214 
01215 {
01216 
01217     opserr<<"TotalLagrangianFD8NodeBrick::addLoad - load type unknown for ele with tag: "<<this->getTag();          
01218 
01219     return -1;
01220 
01221 }
01222 
01223 
01224 
01225 //=============================================================================
01226 
01227 int TotalLagrangianFD8NodeBrick::addInertiaLoadToUnbalance(const Vector &accel)
01228 
01229 {
01230 
01231     // Check for a quick return
01232 
01233     if (rho == 0.0) return 0;
01234 
01235 
01236 
01237     static Vector ra(NumElemDof);
01238 
01239     int i, j;
01240 
01241 
01242 
01243     for (i=0; i<NumNodes; i++) {
01244 
01245       const Vector &RA = theNodes[i]->getRV(accel);
01246 
01247       if ( RA.Size() != NumDof ) {
01248 
01249         opserr << "TotalLagrangianFD8NodeBrick::addInertiaLoadToUnbalance(): matrix and vector sizes are incompatable \n";
01250 
01251         return (-1);
01252 
01253       }
01254 
01255       
01256 
01257       for (j=0; j<NumDof; j++) {
01258 
01259             ra(i*NumDof +j) = RA(j);
01260 
01261       }
01262 
01263 
01264 
01265     }
01266 
01267 
01268 
01269     this->getMass();
01270 
01271 
01272 
01273     if (Q == 0)  
01274 
01275       Q = new Vector(NumElemDof);
01276 
01277 
01278 
01279     Q->addMatrixVector(1.0, M, ra, -1.0);
01280 
01281 
01282 
01283     return 0;  
01284 
01285     
01286 
01287 }
01288 
01289 
01290 
01291 
01292 
01293 //=============================================================================
01294 
01295 const Vector &TotalLagrangianFD8NodeBrick::getResistingForce ()
01296 
01297 {   
01298 
01299         int i, j;
01300 
01301     int f_dim[] = {NumNodes, NumDof};
01302 
01303     tensor NodalForces_in(2, f_dim, 0.0);
01304 
01305     NodalForces_in = this->getRtensor() - this->getForces();
01306 
01307     
01308 
01309     for (i=0; i<NumNodes; i++) {
01310 
01311       for (j=0; j<NumDof; j++) {
01312 
01313          P(i*NumDof +j) = NodalForces_in.cval(i+1, j+1);
01314 
01315       }
01316 
01317     }
01318 
01319         
01320 
01321     if ( Q != 0 )
01322 
01323       P.addVector(1.0, *Q, -1.0);
01324 
01325     
01326 
01327     return P;
01328 
01329 }
01330 
01331 
01332 
01333 //=============================================================================
01334 
01335 const Vector &TotalLagrangianFD8NodeBrick::getResistingForceIncInertia ()
01336 
01337 {
01338 
01339     int i, j;
01340 
01341     Vector a(NumElemDof);
01342 
01343     
01344 
01345     this->getResistingForce();
01346 
01347 
01348 
01349     if (rho != 0.0)
01350 
01351     {
01352 
01353       for (i=0; i<NumNodes; i++) {
01354 
01355         const Vector &acc = theNodes[i]->getTrialAccel();
01356 
01357         if ( acc.Size() != NumDof ) {
01358 
01359           opserr << "TotalLagrangianFD8NodeBrick::getResistingForceIncInertia matrix and vector sizes are incompatable \n";
01360 
01361           exit(-1);
01362 
01363         }
01364 
01365       for (j=0; j<NumDof; j++) {
01366 
01367         a(i*NumDof +j) = acc(j);
01368 
01369       }
01370 
01371     }
01372 
01373 
01374 
01375     this->getMass();
01376 
01377     P.addMatrixVector(1.0, M, a, 1.0);
01378 
01379 
01380 
01381   }
01382 
01383 
01384 
01385   return P;
01386 
01387 }
01388 
01389 
01390 
01391 //=============================================================================
01392 
01393 int TotalLagrangianFD8NodeBrick::sendSelf (int commitTag, Channel &theChannel)
01394 
01395 {
01396 
01397      // Not implemtented yet
01398 
01399      return 0;
01400 
01401 }
01402 
01403 
01404 
01405 //=============================================================================
01406 
01407 int TotalLagrangianFD8NodeBrick::recvSelf (int commitTag, Channel &theChannel,
01408 
01409 FEM_ObjectBroker &theBroker)
01410 
01411 {
01412 
01413      // Not implemtented yet
01414 
01415      return 0;
01416 
01417 }
01418 
01419 
01420 
01421 
01422 
01423 //=============================================================================
01424 
01425 int TotalLagrangianFD8NodeBrick::displaySelf (Renderer &theViewer, int displayMode, float fact)
01426 
01427 {
01428 
01429      // Not implemtented yet
01430 
01431      return 0;
01432 
01433 }
01434 
01435 
01436 
01437 //=============================================================================
01438 
01439 void TotalLagrangianFD8NodeBrick::Print(OPS_Stream &s, int flag)
01440 
01441 {
01442 
01443     s << "\nTotalLagrangianFD8NodeBrick, element id:  " << this->getTag() << endln;
01444 
01445     s << "\nConnected external nodes:  " << connectedExternalNodes;
01446 
01447     s << "\nBody forces:  " << bf(0) << " " << bf(1) << " " << bf(2) << endln;
01448 
01449 
01450 
01451     theMaterial[0]->Print(s,flag);
01452 
01453 
01454 
01455     tensor sigma;
01456 
01457     Vector P00(6);
01458 
01459     
01460 
01461     int i;
01462 
01463     for (i=0; i<NumTotalGaussPts; i++)
01464 
01465     {
01466 
01467       sigma = theMaterial[i]->getCauchyStressTensor();
01468 
01469       P00(0) = sigma.val(1,1);
01470 
01471       P00(1) = sigma.val(2,2);
01472 
01473       P00(2) = sigma.val(3,3);
01474 
01475       P00(3) = sigma.val(2,3);
01476 
01477       P00(4) = sigma.val(3,1);
01478 
01479       P00(5) = sigma.val(1,2);
01480 
01481 
01482 
01483       s << "\n where = " << i << endln;
01484 
01485       s << " Stress (Cauchy): xx yy zz yz zx xy) " << P00 << endln;
01486 
01487     }
01488 
01489 
01490 
01491 }
01492 
01493 
01494 
01495 //=============================================================================
01496 
01497 Response * TotalLagrangianFD8NodeBrick::setResponse (const char **argv, int argc, Information &eleInfo, OPS_Stream &output)
01498 
01499 {
01500 
01501   Response *theResponse = 0;
01502 
01503 
01504 
01505   char outputData[32];
01506 
01507 
01508 
01509   output.tag("ElementOutput");
01510 
01511   output.attr("eleType","TotalLagrangianFD8NodeBrick");
01512 
01513   output.attr("eleTag",this->getTag());
01514 
01515   for (int i=1; i<=NumNodes; i++) {
01516 
01517     sprintf(outputData,"node%d",i);
01518 
01519     output.attr(outputData,connectedExternalNodes[i-1]);
01520 
01521   }
01522 
01523 
01524 
01525   if (strcmp(argv[0],"force") == 0 || strcmp(argv[0],"forces") == 0) {
01526 
01527     for (int i=1; i<=NumNodes; i++)
01528 
01529       for (int j=1; j<=NumDof; j++) {
01530 
01531         sprintf(outputData,"P%d_%d",j,i);
01532 
01533         output.tag("ResponseType",outputData);
01534 
01535       }
01536 
01537     
01538 
01539     theResponse = new ElementResponse(this, 1, this->getResistingForce());
01540 
01541     
01542 
01543   } else if (strcmp(argv[0],"CauchyStress") == 0 || strcmp(argv[0],"stress") == 0)
01544 
01545       theResponse = new ElementResponse(this, 3, Vector(NumTotalGaussPts*6));
01546 
01547     
01548 
01549     else if (strcmp(argv[0],"PK2Stress") == 0 || strcmp(argv[0],"PK2stress") == 0)
01550 
01551       theResponse = new ElementResponse(this, 4, Vector(NumTotalGaussPts*6));
01552 
01553 
01554 
01555     else if (strcmp(argv[0],"EulerianStrain") == 0 || strcmp(argv[0],"strain") == 0)
01556 
01557       theResponse = new ElementResponse(this, 5, Vector(NumTotalGaussPts*6));
01558 
01559     
01560 
01561     else if (strcmp(argv[0],"LagrangianStrain") == 0 || strcmp(argv[0],"iniStrain") == 0)
01562 
01563       theResponse = new ElementResponse(this, 6, Vector(NumTotalGaussPts*6));
01564 
01565 
01566 
01567     output.endTag(); // ElementOutput
01568 
01569 
01570 
01571   return theResponse; 
01572 
01573 }
01574 
01575 
01576 
01577 //=============================================================================
01578 
01579 int TotalLagrangianFD8NodeBrick::getResponse (int responseID, Information &eleInfo)
01580 
01581 {    
01582 
01583          int i;
01584 
01585      static Vector P0(NumTotalGaussPts*6);
01586 
01587      
01588 
01589      switch (responseID) {
01590 
01591      
01592 
01593      case 1:
01594 
01595           return eleInfo.setVector(this->getResistingForce() );
01596 
01597 
01598 
01599      case 3: { 
01600 
01601         Vector P0(NumTotalGaussPts*6);
01602 
01603         tensor sigma; 
01604 
01605         for (i=0; i<NumTotalGaussPts; i++) {
01606 
01607           sigma = theMaterial[i]->getCauchyStressTensor();
01608 
01609           P0(i*6 +0 ) = sigma.val(1,1);
01610 
01611           P0(i*6 +1 ) = sigma.val(2,2);
01612 
01613           P0(i*6 +2 ) = sigma.val(3,3);
01614 
01615           P0(i*6 +3 ) = sigma.val(2,3);
01616 
01617           P0(i*6 +4 ) = sigma.val(3,1);
01618 
01619           P0(i*6 +5 ) = sigma.val(1,2);
01620 
01621         }
01622 
01623         return eleInfo.setVector(P0);
01624 
01625      }
01626 
01627 
01628 
01629      case 4: { 
01630 
01631         Vector P0(NumTotalGaussPts*6);
01632 
01633         tensor sigma; 
01634 
01635         for (i=0; i<NumTotalGaussPts; i++) {
01636 
01637           sigma = theMaterial[i]->getStressTensor();
01638 
01639           P0(i*6 +0 ) = sigma.val(1,1);
01640 
01641           P0(i*6 +1 ) = sigma.val(2,2);
01642 
01643           P0(i*6 +2 ) = sigma.val(3,3);
01644 
01645           P0(i*6 +3 ) = sigma.val(2,3);
01646 
01647           P0(i*6 +4 ) = sigma.val(3,1);
01648 
01649           P0(i*6 +5 ) = sigma.val(1,2);
01650 
01651         }
01652 
01653         return eleInfo.setVector(P0);
01654 
01655      }
01656 
01657 
01658 
01659      case 5: { 
01660 
01661         Vector P0(NumTotalGaussPts*6);
01662 
01663         tensor e;
01664 
01665         tensor E;
01666 
01667         tensor F;
01668 
01669         tensor tI2("I", 2, def_dim_2); 
01670 
01671         for (i=0; i<NumTotalGaussPts; i++) {
01672 
01673           E = theMaterial[i]->getStrainTensor();
01674 
01675           F = theMaterial[i]->getF();
01676 
01677           F = F.inverse();
01678 
01679           e = F("ki")*F("kj"); e.null_indices();
01680 
01681           e = (tI2-e) *0.5;
01682 
01683           P0(i*6 +0 ) = e.val(1,1);
01684 
01685           P0(i*6 +1 ) = e.val(2,2);
01686 
01687           P0(i*6 +2 ) = e.val(3,3);
01688 
01689           P0(i*6 +3 ) = e.val(2,3);
01690 
01691           P0(i*6 +4 ) = e.val(3,1);
01692 
01693           P0(i*6 +5 ) = e.val(1,2);
01694 
01695         }
01696 
01697         return eleInfo.setVector(P0);
01698 
01699      }
01700 
01701 
01702 
01703      case 6: { 
01704 
01705         Vector P0(NumTotalGaussPts*6);
01706 
01707         tensor E; 
01708 
01709         for (i=0; i<NumTotalGaussPts; i++) {
01710 
01711           E = theMaterial[i]->getStrainTensor();
01712 
01713           P0(i*6 +0 ) = E.val(1,1);
01714 
01715           P0(i*6 +1 ) = E.val(2,2);
01716 
01717           P0(i*6 +2 ) = E.val(3,3);
01718 
01719           P0(i*6 +3 ) = E.val(2,3);
01720 
01721           P0(i*6 +4 ) = E.val(3,1);
01722 
01723           P0(i*6 +5 ) = E.val(1,2);
01724 
01725         }
01726 
01727         return eleInfo.setVector(P0);
01728 
01729      }
01730 
01731     
01732 
01733      default:
01734 
01735      return -1;
01736 
01737 
01738 
01739    }    
01740 
01741 }
01742 
01743 
01744 
01745 
01746 
01747 //#############################################################################
01748 
01749 //===================================================================
01750 
01751 tensor TotalLagrangianFD8NodeBrick::shapeFunction(double r1, double r2, double r3)
01752 
01753   {
01754 
01755 
01756 
01757     int dimension[] = {NumNodes};
01758 
01759     tensor h(1, dimension, 0.0);
01760 
01761 
01762 
01763       // influence of the node number 8
01764 
01765     h.val(8)=(1.0+r1)*(1.0-r2)*(1.0-r3)*0.125;
01766 
01767       // influence of the node number 7
01768 
01769     h.val(7)=(1.0-r1)*(1.0-r2)*(1.0-r3)*0.125;
01770 
01771       // influence of the node number 6
01772 
01773     h.val(6)=(1.0-r1)*(1.0+r2)*(1.0-r3)*0.125;
01774 
01775       // influence of the node number 5
01776 
01777     h.val(5)=(1.0+r1)*(1.0+r2)*(1.0-r3)*0.125;
01778 
01779 
01780 
01781       // influence of the node number 4
01782 
01783     h.val(4)=(1.0+r1)*(1.0-r2)*(1.0+r3)*0.125;
01784 
01785       // influence of the node number 3
01786 
01787     h.val(3)=(1.0-r1)*(1.0-r2)*(1.0+r3)*0.125;
01788 
01789       // influence of the node number 2
01790 
01791     h.val(2)=(1.0-r1)*(1.0+r2)*(1.0+r3)*0.125;
01792 
01793       // influence of the node number 1
01794 
01795     h.val(1)=(1.0+r1)*(1.0+r2)*(1.0+r3)*0.125;
01796 
01797 
01798 
01799     return h;
01800 
01801   }
01802 
01803 
01804 
01805 
01806 
01807 //==============================================================
01808 
01809 tensor TotalLagrangianFD8NodeBrick::shapeFunctionDerivative(double r1, double r2, double r3)
01810 
01811   {
01812 
01813 
01814 
01815     int dimensions[] = {NumNodes, NumDof};
01816 
01817     tensor dh(2, dimensions, 0.0);
01818 
01819 
01820 
01821       // influence of the node number 8
01822 
01823     dh.val(8,1)= (1.0-r2)*(1.0-r3)*0.125;
01824 
01825     dh.val(8,2)=-(1.0+r1)*(1.0-r3)*0.125;
01826 
01827     dh.val(8,3)=-(1.0+r1)*(1.0-r2)*0.125;
01828 
01829       // influence of the node number 7
01830 
01831     dh.val(7,1)=-(1.0-r2)*(1.0-r3)*0.125;
01832 
01833     dh.val(7,2)=-(1.0-r1)*(1.0-r3)*0.125;
01834 
01835     dh.val(7,3)=-(1.0-r1)*(1.0-r2)*0.125;
01836 
01837       // influence of the node number 6
01838 
01839     dh.val(6,1)=-(1.0+r2)*(1.0-r3)*0.125;
01840 
01841     dh.val(6,2)= (1.0-r1)*(1.0-r3)*0.125;
01842 
01843     dh.val(6,3)=-(1.0-r1)*(1.0+r2)*0.125;
01844 
01845       // influence of the node number 5
01846 
01847     dh.val(5,1)= (1.0+r2)*(1.0-r3)*0.125;
01848 
01849     dh.val(5,2)= (1.0+r1)*(1.0-r3)*0.125;
01850 
01851     dh.val(5,3)=-(1.0+r1)*(1.0+r2)*0.125;
01852 
01853 
01854 
01855       // influence of the node number 4
01856 
01857     dh.val(4,1)= (1.0-r2)*(1.0+r3)*0.125;
01858 
01859     dh.val(4,2)=-(1.0+r1)*(1.0+r3)*0.125;
01860 
01861     dh.val(4,3)= (1.0+r1)*(1.0-r2)*0.125;
01862 
01863       // influence of the node number 3
01864 
01865     dh.val(3,1)=-(1.0-r2)*(1.0+r3)*0.125;
01866 
01867     dh.val(3,2)=-(1.0-r1)*(1.0+r3)*0.125;
01868 
01869     dh.val(3,3)= (1.0-r1)*(1.0-r2)*0.125;
01870 
01871       // influence of the node number 2
01872 
01873     dh.val(2,1)=-(1.0+r2)*(1.0+r3)*0.125;
01874 
01875     dh.val(2,2)= (1.0-r1)*(1.0+r3)*0.125;
01876 
01877     dh.val(2,3)= (1.0-r1)*(1.0+r2)*0.125;
01878 
01879       // influence of the node number 1
01880 
01881     dh.val(1,1)= (1.0+r2)*(1.0+r3)*0.125;
01882 
01883     dh.val(1,2)= (1.0+r1)*(1.0+r3)*0.125;
01884 
01885     dh.val(1,3)= (1.0+r1)*(1.0+r2)*0.125;
01886 
01887 
01888 
01889     return dh;
01890 
01891   }
01892 
01893 
01894 
01895 
01896 
01897 #endif
01898 
01899 
01900 

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