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

DispBeamColumn2d.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/13 21:11:24 $
00023 // $Source: /usr/local/cvs/OpenSees/SRC/element/dispBeamColumn/DispBeamColumn2d.cpp,v $
00024 
00025 // Written: MHS
00026 // Created: Feb 2001
00027 //
00028 // Description: This file contains the class definition for DispBeamColumn2d.
00029 
00030 #include <DispBeamColumn2d.h>
00031 #include <Node.h>
00032 #include <SectionForceDeformation.h>
00033 #include <GaussQuadRule1d01.h>
00034 #include <CrdTransf2d.h>
00035 #include <Matrix.h>
00036 #include <Vector.h>
00037 #include <ID.h>
00038 #include <Renderer.h>
00039 #include <Domain.h>
00040 #include <string.h>
00041 #include <Information.h>
00042 #include <Channel.h>
00043 #include <FEM_ObjectBroker.h>
00044 #include <ElementResponse.h>
00045 
00046 #include <G3Globals.h>
00047 
00048 Matrix DispBeamColumn2d::K(6,6);
00049 Vector DispBeamColumn2d::P(6);
00050 double DispBeamColumn2d::workArea[100];
00051 
00052 DispBeamColumn2d::DispBeamColumn2d(int tag, int nd1, int nd2,
00053   int numSec, SectionForceDeformation &s,
00054   CrdTransf2d &coordTransf, double r)
00055 :Element (tag, ELE_TAG_DispBeamColumn2d), rho(r),
00056  Q(6), q(3), connectedExternalNodes(2), L(0.0),
00057  numSections(numSec), theSections(0), crdTransf(0)
00058 {
00059     // Allocate arrays of pointers to SectionForceDeformations
00060     theSections = new SectionForceDeformation *[numSections];
00061     
00062  if (theSections == 0)
00063      g3ErrorHandler->fatal("%s - failed to allocate section model pointer",
00064    "DispBeamColumn2d::DispBeamColumn2d");
00065 
00066     for (int i = 0; i < numSections; i++) {
00067 
00068   // Get copies of the material model for each integration point
00069   theSections[i] = s.getCopy();
00070    
00071   // Check allocation
00072   if (theSections[i] == 0)
00073    g3ErrorHandler->fatal("%s -- failed to get a copy of section model",
00074     "DispBeamColumn2d::DispBeamColumn2d");
00075  }
00076 
00077  crdTransf = coordTransf.getCopy();
00078 
00079  if (crdTransf == 0)
00080      g3ErrorHandler->fatal("%s - failed to copy coordinate transformation",
00081    "DispBeamColumn2d::DispBeamColumn2d");
00082 
00083  // Set connected external node IDs
00084     connectedExternalNodes(0) = nd1;
00085     connectedExternalNodes(1) = nd2;
00086 }
00087 
00088 DispBeamColumn2d::DispBeamColumn2d()
00089 :Element (0, ELE_TAG_DispBeamColumn2d), rho(0.0),
00090  Q(6), q(3), connectedExternalNodes(2), L(0.0),
00091  numSections(0), theSections(0), crdTransf(0)
00092 {
00093 
00094 }
00095 
00096 DispBeamColumn2d::~DispBeamColumn2d()
00097 {    
00098     for (int i = 0; i < numSections; i++) {
00099   if (theSections[i])
00100    delete theSections[i];
00101  }
00102 
00103     // Delete the array of pointers to SectionForceDeformation pointer arrays
00104     if (theSections)
00105   delete [] theSections;
00106 }
00107 
00108 int
00109 DispBeamColumn2d::getNumExternalNodes() const
00110 {
00111     return 2;
00112 }
00113 
00114 const ID&
00115 DispBeamColumn2d::getExternalNodes()
00116 {
00117     return connectedExternalNodes;
00118 }
00119 
00120 int
00121 DispBeamColumn2d::getNumDOF()
00122 {
00123     return 6;
00124 }
00125 
00126 void
00127 DispBeamColumn2d::setDomain(Domain *theDomain)
00128 {
00129  // Check Domain is not null - invoked when object removed from a domain
00130     if (theDomain == 0) {
00131  nd1Ptr = 0;
00132  nd2Ptr = 0;
00133  return;
00134     }
00135 
00136     int Nd1 = connectedExternalNodes(0);
00137     int Nd2 = connectedExternalNodes(1);
00138 
00139     nd1Ptr = theDomain->getNode(Nd1);
00140     nd2Ptr = theDomain->getNode(Nd2);
00141 
00142     if (nd1Ptr == 0 || nd2Ptr == 0) {
00143  //g3ErrorHandler->fatal("FATAL ERROR DispBeamColumn2d (tag: %d), node not found in domain",
00144  // this->getTag());
00145  
00146  return;
00147     }
00148 
00149     int dofNd1 = nd1Ptr->getNumberDOF();
00150     int dofNd2 = nd2Ptr->getNumberDOF();
00151     
00152     if (dofNd1 != 3 || dofNd2 != 3) {
00153  //g3ErrorHandler->fatal("FATAL ERROR DispBeamColumn2d (tag: %d), has differing number of DOFs at its nodes",
00154  // this->getTag());
00155  
00156  return;
00157     }
00158 
00159  if (crdTransf->initialize(nd1Ptr, nd2Ptr)) {
00160   // Add some error check
00161  }
00162 
00163  L = crdTransf->getInitialLength();
00164 
00165  if (L == 0.0) {
00166   // Add some error check
00167  }
00168 
00169     this->DomainComponent::setDomain(theDomain);
00170 
00171     this->update();
00172 }
00173 
00174 int
00175 DispBeamColumn2d::commitState()
00176 {
00177     int retVal = 0;
00178 
00179     // Loop over the integration points and commit the material states
00180     for (int i = 0; i < numSections; i++)
00181   retVal += theSections[i]->commitState();
00182 
00183     retVal += crdTransf->commitState();
00184 
00185     return retVal;
00186 }
00187 
00188 int
00189 DispBeamColumn2d::revertToLastCommit()
00190 {
00191     int retVal = 0;
00192 
00193     // Loop over the integration points and revert to last committed state
00194     for (int i = 0; i < numSections; i++)
00195   retVal += theSections[i]->revertToLastCommit();
00196 
00197     retVal += crdTransf->revertToLastCommit();
00198 
00199     return retVal;
00200 }
00201 
00202 int
00203 DispBeamColumn2d::revertToStart()
00204 {
00205     int retVal = 0;
00206 
00207     // Loop over the integration points and revert states to start
00208     for (int i = 0; i < numSections; i++)
00209   retVal += theSections[i]->revertToStart();
00210 
00211     retVal += crdTransf->revertToStart();
00212 
00213     return retVal;
00214 }
00215 
00216 int
00217 DispBeamColumn2d::update(void)
00218 {
00219  // Update the transformation
00220  crdTransf->update();
00221 
00222     // Get basic deformations
00223     static Vector v(3);
00224  v = crdTransf->getBasicTrialDisp();
00225 
00226  double oneOverL = 1.0/L;
00227  GaussQuadRule1d01 quadrat(numSections);
00228  const Matrix &pts = quadrat.getIntegrPointCoords();
00229  const Vector &wts = quadrat.getIntegrPointWeights();
00230 
00231  // Assuming member is prismatic ... have to move inside
00232  // the loop if it is not prismatic
00233  int order = theSections[0]->getOrder();
00234  const ID &code = theSections[0]->getType();
00235  Vector e(workArea, order);
00236 
00237  // Loop over the integration points
00238  for (int i = 0; i < numSections; i++) {
00239 
00240   double xi6 = 6.0*pts(i,0);
00241 
00242   int j;
00243   for (j = 0; j < order; j++) {
00244    switch(code(j)) {
00245    case SECTION_RESPONSE_P:
00246     e(j) = oneOverL*v(0); break;
00247    case SECTION_RESPONSE_MZ:
00248     e(j) = oneOverL*((xi6-4.0)*v(1) + (xi6-2.0)*v(2)); break;
00249    default:
00250     e(j) = 0.0; break;
00251    }
00252   }
00253 
00254   // Set the section deformations
00255   theSections[i]->setTrialSectionDeformation(e);
00256  }
00257 
00258  return 0;
00259 }
00260 
00261 const Matrix&
00262 DispBeamColumn2d::getTangentStiff()
00263 {
00264  static Matrix kb(3,3);
00265 
00266  // Zero for integral
00267  kb.Zero();
00268  q.Zero();
00269 
00270  GaussQuadRule1d01 quadrat(numSections);
00271  const Matrix &pts = quadrat.getIntegrPointCoords();
00272  const Vector &wts = quadrat.getIntegrPointWeights();
00273 
00274  // Assuming member is prismatic ... have to move inside
00275  // the loop if it is not prismatic
00276  int order = theSections[0]->getOrder();
00277  const ID &code = theSections[0]->getType();
00278 
00279  double oneOverL = 1.0/L;
00280  Matrix ka(workArea, order, 3);
00281 
00282  // Loop over the integration points
00283  for (int i = 0; i < numSections; i++) {
00284 
00285   // Get the section tangent stiffness and stress resultant
00286   const Matrix &ks = theSections[i]->getSectionTangent();
00287   const Vector &s = theSections[i]->getStressResultant();
00288 
00289   double xi6 = 6.0*pts(i,0);
00290   ka.Zero();
00291 
00292   // Perform numerical integration
00293   //kb.addMatrixTripleProduct(1.0, *B, ks, wts(i)/L);
00294   double wti = wts(i)*oneOverL;
00295   double tmp;
00296   int j, k;
00297   for (j = 0; j < order; j++) {
00298    switch(code(j)) {
00299    case SECTION_RESPONSE_P:
00300     for (k = 0; k < order; k++)
00301      ka(k,0) += ks(k,j)*wti;
00302     break;
00303    case SECTION_RESPONSE_MZ:
00304     for (k = 0; k < order; k++) {
00305      tmp = ks(k,j)*wti;
00306      ka(k,1) += (xi6-4.0)*tmp;
00307      ka(k,2) += (xi6-2.0)*tmp;
00308     }
00309     break;
00310    default:
00311     break;
00312    }
00313   }
00314   for (j = 0; j < order; j++) {
00315    switch (code(j)) {
00316    case SECTION_RESPONSE_P:
00317     for (k = 0; k < 3; k++)
00318      kb(0,k) += ka(j,k);
00319     break;
00320    case SECTION_RESPONSE_MZ:
00321     for (k = 0; k < 3; k++) {
00322      tmp = ka(j,k);
00323      kb(1,k) += (xi6-4.0)*tmp;
00324      kb(2,k) += (xi6-2.0)*tmp;
00325     }
00326     break;
00327    default:
00328     break;
00329    }
00330   }
00331 
00332   //q.addMatrixTransposeVector(1.0, *B, s, wts(i));
00333   double si;
00334   for (j = 0; j < order; j++) {
00335    si = s(j)*wts(i);
00336    switch(code(j)) {
00337    case SECTION_RESPONSE_P:
00338     q(0) += si; break;
00339    case SECTION_RESPONSE_MZ:
00340     q(1) += (xi6-4.0)*si; q(2) += (xi6-2.0)*si; break;
00341    default:
00342     break;
00343    }
00344   }
00345 
00346  }
00347 
00348  // Transform to global stiffness
00349  K = crdTransf->getGlobalStiffMatrix(kb, q);
00350 
00351  return K;
00352 }
00353 
00354 const Matrix&
00355 DispBeamColumn2d::getSecantStiff()
00356 {
00357  return this->getTangentStiff();
00358 }
00359 
00360 const Matrix&
00361 DispBeamColumn2d::getDamp()
00362 {
00363  K.Zero();
00364  
00365  return K;
00366 }
00367 
00368 const Matrix&
00369 DispBeamColumn2d::getMass()
00370 {
00371  K.Zero();
00372 
00373  if (rho == 0.0)
00374   return K;
00375 
00376  double m = 0.5*rho*L;
00377 
00378  K(0,0) = K(1,1) = K(3,3) = K(4,4) = m;
00379 
00380  return K;
00381 }
00382 
00383 void
00384 DispBeamColumn2d::zeroLoad(void)
00385 {
00386  Q(0) = 0.0;
00387  Q(1) = 0.0;
00388  Q(2) = 0.0;
00389  Q(3) = 0.0;
00390  Q(4) = 0.0;
00391  Q(5) = 0.0;
00392 
00393  return;
00394 }
00395 
00396 int 
00397 DispBeamColumn2d::addLoad(const Vector &addLoad)
00398 {
00399  if (addLoad.Size() != 6) {
00400   g3ErrorHandler->warning("DispBeamColumn2d::addLoad %s\n",
00401     "Vector not of correct size");
00402   return -1;
00403  }
00404 
00405  // Add to the external nodal loads
00406  //Q.addVector(1.0, addLoad, 1.0);
00407  Q(0) += addLoad(0);
00408  Q(1) += addLoad(1);
00409  Q(2) += addLoad(2);
00410  Q(3) += addLoad(3);
00411  Q(4) += addLoad(4);
00412  Q(5) += addLoad(5);
00413 
00414  return 0;
00415 }
00416 
00417 int 
00418 DispBeamColumn2d::addInertiaLoadToUnbalance(const Vector &accel)
00419 {
00420  // Check for a quick return
00421  if (rho == 0.0) 
00422   return 0;
00423 
00424  // Get R * accel from the nodes
00425  const Vector &Raccel1 = nd1Ptr->getRV(accel);
00426  const Vector &Raccel2 = nd2Ptr->getRV(accel);
00427 
00428     if (3 != Raccel1.Size() || 3 != Raccel2.Size()) {
00429   g3ErrorHandler->warning("DispBeamColumn2d::addInertiaLoadToUnbalance %s\n",
00430     "matrix and vector sizes are incompatable");
00431   return -1;
00432     }
00433 
00434  double m = 0.5*rho*L;
00435 
00436     // Want to add ( - fact * M R * accel ) to unbalance
00437  // Take advantage of lumped mass matrix
00438  Q(0) -= m*Raccel1(0);
00439  Q(1) -= m*Raccel1(1);
00440  Q(3) -= m*Raccel2(0);
00441  Q(4) -= m*Raccel2(1);
00442 
00443     return 0;
00444 }
00445 
00446 const Vector&
00447 DispBeamColumn2d::getResistingForce()
00448 {
00449  GaussQuadRule1d01 quadrat(numSections);
00450  const Matrix &pts = quadrat.getIntegrPointCoords();
00451  const Vector &wts = quadrat.getIntegrPointWeights();
00452 
00453  // Assuming member is prismatic ... have to move inside
00454  // the loop if it is not prismatic
00455  int order = theSections[0]->getOrder();
00456  const ID &code = theSections[0]->getType();
00457 
00458  // Zero for integration
00459  q.Zero();
00460 
00461  // Loop over the integration points
00462  for (int i = 0; i < numSections; i++) {
00463 
00464   double xi6 = 6.0*pts(i,0);
00465 
00466   // Get section stress resultant
00467   const Vector &s = theSections[i]->getStressResultant();
00468 
00469   // Perform numerical integration on internal force
00470   //q.addMatrixTransposeVector(1.0, *B, s, wts(i));
00471 
00472   double si;
00473   for (int j = 0; j < order; j++) {
00474    si = s(j)*wts(i);
00475    switch(code(j)) {
00476    case SECTION_RESPONSE_P:
00477     q(0) += si; break;
00478    case SECTION_RESPONSE_MZ:
00479     q(1) += (xi6-4.0)*si; q(2) += (xi6-2.0)*si; break;
00480    default:
00481     break;
00482    }
00483   }
00484 
00485  }
00486 
00487  // Transform forces
00488  static Vector dummy(2);  // No distributed loads
00489  P = crdTransf->getGlobalResistingForce(q,dummy);
00490 
00491  // Subtract other external nodal loads ... P_res = P_int - P_ext
00492  //P.addVector(1.0, Q, -1.0);
00493  P(0) -= Q(0);
00494  P(1) -= Q(1);
00495  P(2) -= Q(2);
00496  P(3) -= Q(3);
00497  P(4) -= Q(4);
00498  P(5) -= Q(5);
00499 
00500  return P;
00501 }
00502 
00503 const Vector&
00504 DispBeamColumn2d::getResistingForceIncInertia()
00505 {
00506  // Check for a quick return
00507  if (rho == 0.0)
00508   return this->getResistingForce();
00509 
00510  const Vector &accel1 = nd1Ptr->getTrialAccel();
00511  const Vector &accel2 = nd2Ptr->getTrialAccel();
00512  
00513  // Compute the current resisting force
00514  this->getResistingForce();
00515 
00516  double m = 0.5*rho*L;
00517 
00518  P(0) += m*accel1(0);
00519  P(1) += m*accel1(1);
00520  P(3) += m*accel2(0);
00521  P(4) += m*accel2(1);
00522 
00523  return P;
00524 }
00525 
00526 int
00527 DispBeamColumn2d::sendSelf(int commitTag, Channel &theChannel)
00528 {
00529  return -1;
00530 }
00531 
00532 int
00533 DispBeamColumn2d::recvSelf(int commitTag, Channel &theChannel,
00534       FEM_ObjectBroker &theBroker)
00535 {
00536  return -1;
00537 }
00538 
00539 void
00540 DispBeamColumn2d::Print(ostream &s, int flag)
00541 {
00542  s << "\nDispBeamColumn2d, element id:  " << this->getTag() << endl;
00543  s << "\tConnected external nodes:  " << connectedExternalNodes;
00544  s << "\tmass density:  " << rho << endl;
00545  theSections[0]->Print(s,flag);
00546 }
00547 
00548 
00549 int
00550 DispBeamColumn2d::displaySelf(Renderer &theViewer, int displayMode, float fact)
00551 {
00552     // first determine the end points of the quad based on
00553     // the display factor (a measure of the distorted image)
00554     const Vector &end1Crd = nd1Ptr->getCrds();
00555     const Vector &end2Crd = nd2Ptr->getCrds(); 
00556 
00557     const Vector &end1Disp = nd1Ptr->getDisp();
00558     const Vector &end2Disp = nd2Ptr->getDisp();
00559 
00560  static Vector v1(3);
00561  static Vector v2(3);
00562 
00563  for (int i = 0; i < 2; i++) {
00564   v1(i) = end1Crd(i) + end1Disp(i)*fact;
00565   v2(i) = end2Crd(i) + end2Disp(i)*fact;    
00566  }
00567  
00568  return theViewer.drawLine (v1, v2, 1.0, 1.0);
00569 }
00570 
00571 Response*
00572 DispBeamColumn2d::setResponse(char **argv, int argc, Information &eleInfo)
00573 {
00574     // global force - 
00575     if (strcmp(argv[0],"forces") == 0 || strcmp(argv[0],"force") == 0
00576   || strcmp(argv[0],"globalForce") == 0 || strcmp(argv[0],"globalForces") == 0)
00577   return new ElementResponse(this, 1, P);
00578 
00579     // local force -
00580     else if (strcmp(argv[0],"localForce") == 0 || strcmp(argv[0],"localForces") == 0)
00581   return new ElementResponse(this, 2, P);
00582     
00583     // section response -
00584     else if (strcmp(argv[0],"section") ==0) {
00585   if (argc <= 2)
00586    return 0;
00587  
00588   int sectionNum = atoi(argv[1]);
00589   if (sectionNum > 0 && sectionNum <= numSections)
00590    return theSections[sectionNum-1]->setResponse(&argv[2], argc-2, eleInfo);
00591   else
00592    return 0;
00593  }
00594     
00595  else
00596   return 0;
00597 }
00598 
00599 int 
00600 DispBeamColumn2d::getResponse(int responseID, Information &eleInfo)
00601 {
00602   double V;
00603 
00604   switch (responseID) {
00605     case 1:  // global forces
00606       return eleInfo.setVector(P);
00607 
00608     case 2:
00609       P(3) = q(0);
00610       P(0) = -q(0);
00611       P(2) = q(1);
00612       P(5) = q(2);
00613       V = (q(1)+q(2))/L;
00614       P(1) = V;
00615       P(4) = -V;
00616       return eleInfo.setVector(P);
00617       
00618     default: 
00619    return -1;
00620   }
00621 }
Copyright Contact Us