00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033 #include "beam2d02.h"
00034 #include <Domain.h>
00035 #include <Channel.h>
00036 #include <FEM_ObjectBroker.h>
00037 #include <Renderer.h>
00038
00039 #include <math.h>
00040 #include <stdlib.h>
00041 #include <CrdTransf2d.h>
00042
00043 beam2d02::beam2d02()
00044 :Element(0,ELE_TAG_beam2d02),
00045 A(0.0), E(0.0), I(0.0), M(0.0), L(0),
00046 connectedExternalNodes(2),
00047 Kd(3,3), m(6,6),
00048 q(3), rForce(6), load(6),
00049 theCoordTrans(0)
00050 {
00051
00052 }
00053
00054
00055
00056
00057 beam2d02::beam2d02(int tag, double a, double e, double i, int Nd1, int Nd2,
00058 CrdTransf2d &theTrans, double rho)
00059 :Element(tag,ELE_TAG_beam2d02),
00060 A(a), E(e), I(i), M(rho), L(0),
00061 connectedExternalNodes(2),
00062 Kd(3,3), m(6,6),
00063 q(3), rForce(6), load(6),
00064 theCoordTrans(0)
00065 {
00066 connectedExternalNodes(0) = Nd1;
00067 connectedExternalNodes(1) = Nd2;
00068
00069 theNodes[0] = 0;
00070 theNodes[1] = 0;
00071
00072 theCoordTrans = theTrans.getCopy();
00073 }
00074
00075
00076
00077
00078
00079
00080 beam2d02::~beam2d02()
00081 {
00082 if (theCoordTrans != 0)
00083 delete theCoordTrans;
00084 }
00085
00086
00087 int
00088 beam2d02::getNumExternalNodes(void) const
00089 {
00090 return 2;
00091 }
00092
00093 const ID &
00094 beam2d02::getExternalNodes(void)
00095 {
00096 return connectedExternalNodes;
00097 }
00098
00099 Node **
00100 beam2d02::getNodePtrs(void)
00101 {
00102 return theNodes;
00103 }
00104
00105 int
00106 beam2d02::getNumDOF(void) {
00107 return 6;
00108 }
00109
00110
00111 void
00112 beam2d02::setDomain(Domain *theDomain)
00113 {
00114
00115 int Nd1 = connectedExternalNodes(0);
00116 int Nd2 = connectedExternalNodes(1);
00117 theNodes[0] = theDomain->getNode(Nd1);
00118 theNodes[1] = theDomain->getNode(Nd2);
00119 if (theNodes[0] == 0) {
00120 opserr << "WARNING beam2d02::setDomain(): Nd1: ";
00121 opserr << Nd1 << "does not exist in model for beam \n" << *this;
00122 return;
00123 }
00124 if (theNodes[1] == 0) {
00125 opserr << "WARNING beam2d02::setDomain(): Nd2: ";
00126 opserr << Nd2 << "does not exist in model for beam\n" << *this;
00127 return;
00128 }
00129
00130
00131 int dofNd1 = theNodes[0]->getNumberDOF();
00132 int dofNd2 = theNodes[1]->getNumberDOF();
00133 if (dofNd1 != 3 && dofNd2 != 3) {
00134 opserr << "WARNING beam2d02::setDomain(): node " << Nd1;
00135 opserr << " and/or node " << Nd2 << " have/has incorrect number ";
00136 opserr << "of dof's at end for beam\n " << *this;
00137 return;
00138 }
00139
00140
00141 this->DomainComponent::setDomain(theDomain);
00142
00143
00144 double dx,dy;
00145 const Vector &end1Crd = theNodes[0]->getCrds();
00146 const Vector &end2Crd = theNodes[1]->getCrds();
00147
00148 dx = end2Crd(0)-end1Crd(0);
00149 dy = end2Crd(1)-end1Crd(1);
00150
00151 L = sqrt(dx*dx + dy*dy);
00152 if (L == 0.0) {
00153 opserr << "WARNING beam2d02::setDomain(): beam " << this->getTag();
00154 opserr << " has zero length for beam\n" << *this;
00155 return;
00156 }
00157
00158 Kd(0,0) = E*A/L;
00159 Kd(0,1) = 0.0;
00160 Kd(0,2) = 0.0;
00161
00162 Kd(1,0) = 0.0;
00163 Kd(1,1) = 4.0*E*I/L;
00164 Kd(1,2) = 2.0*E*I/L;
00165
00166 Kd(2,0) = 0.0;
00167 Kd(2,1) = 2.0*E*I/L;
00168 Kd(2,2) = 4.0*E*I/L;
00169
00170 cs = dx/L;
00171 sn = dy/L;
00172
00173
00174 M = 0.5*M*A*L;
00175 }
00176
00177
00178 int
00179 beam2d02::commitState()
00180 {
00181 int retVal = 0;
00182
00183
00184 if ((retVal = this->Element::commitState()) != 0) {
00185 opserr << "beam2d02::commitState () - failed in base class";
00186 }
00187
00188 retVal = theCoordTrans->commitState();
00189 return retVal;
00190 }
00191
00192 int
00193 beam2d02::revertToLastCommit()
00194 {
00195 return theCoordTrans->revertToLastCommit();
00196 }
00197
00198 int
00199 beam2d02::revertToStart()
00200 {
00201 return theCoordTrans->revertToStart();
00202 }
00203
00204 const Matrix &
00205 beam2d02::getTangentStiff(void)
00206 {
00207 return this->getStiff();
00208 }
00209
00210 const Matrix &
00211 beam2d02::getInitialStiff(void)
00212 {
00213 return this->getStiff();
00214 }
00215
00216
00217
00218
00219
00220 const Matrix &
00221 beam2d02::getStiff(void)
00222 {
00223 const Vector &v = theCoordTrans->getBasicTrialDisp();
00224 q.addMatrixVector(0.0,Kd,v,1.0);
00225
00226 return theCoordTrans->getGlobalStiffMatrix(Kd, q);
00227
00228 }
00229
00230 const Matrix &
00231 beam2d02::getMass(void)
00232 {
00233
00234 m(0,0) = M;
00235 m(1,1) = M;
00236 m(3,3) = M;
00237 m(4,4) = M;
00238 return m;
00239 }
00240
00241
00242
00243 void
00244 beam2d02::zeroLoad(void)
00245 {
00246 load.Zero();
00247 }
00248
00249
00250 int
00251 beam2d02::addLoad(ElementalLoad *theLoad, double loadFactor)
00252 {
00253 opserr << "ElasticBeam2d::addLoad() - beam " << this->getTag() << ", does not handle ele loads\n";
00254 return -1;
00255 }
00256
00257
00258 int
00259 beam2d02::addInertiaLoadToUnbalance(const Vector &accel)
00260 {
00261 if (M == 0.0)
00262 return 0;
00263
00264
00265 const Vector &Raccel1 = theNodes[0]->getRV(accel);
00266 const Vector &Raccel2 = theNodes[1]->getRV(accel);
00267
00268
00269
00270 load(0) -= M * Raccel1(0);
00271 load(1) -= M * Raccel1(1);
00272
00273 load(3) -= M * Raccel2(0);
00274 load(4) -= M * Raccel2(1);
00275
00276 return 0;
00277 }
00278
00279
00280 const Vector &
00281 beam2d02::getResistingForce()
00282 {
00283
00284 const Vector &v = theCoordTrans->getBasicTrialDisp();
00285 q.addMatrixVector(0.0,Kd,v,1.0);
00286
00287 static Vector uniLoad(2);
00288
00289 rForce = theCoordTrans->getGlobalResistingForce(q, uniLoad);
00290
00291
00292 rForce -= load;
00293 return rForce;
00294 }
00295
00296
00297 const Vector &
00298 beam2d02::getResistingForceIncInertia()
00299 {
00300 this->getResistingForce();
00301
00302 if (M != 0.0) {
00303 const Vector &end1Accel = theNodes[0]->getTrialAccel();
00304 const Vector &end2Accel = theNodes[1]->getTrialAccel();
00305 Vector inertia(6);
00306 rForce(0) += M*end1Accel(0);
00307 rForce(1) += M*end1Accel(1);
00308 rForce(3) += M*end2Accel(0);
00309 rForce(4) += M*end2Accel(1);
00310
00311
00312 if (alphaM != 0.0 || betaK != 0.0 || betaK0 != 0.0 || betaKc != 0.0)
00313 rForce += this->getRayleighDampingForces();
00314
00315 } else {
00316
00317
00318 if (betaK != 0.0 || betaK0 != 0.0 || betaKc != 0.0)
00319 rForce += this->getRayleighDampingForces();
00320
00321 }
00322
00323 return rForce;
00324 }
00325
00326
00327
00328 int
00329 beam2d02::sendSelf(int commitTag, Channel &theChannel)
00330 {
00331 int dataTag = this->getDbTag();
00332
00333 Vector data(4);
00334 data(0) = A; data(1) = E; data(2) = I; data(3) = this->getTag();
00335
00336 int result = 0;
00337 result = theChannel.sendVector(dataTag, commitTag, data);
00338 if (result < 0) {
00339 opserr << "beam2d02::sendSelf - failed to send data\n";
00340 return -1;
00341 }
00342
00343 result = theChannel.sendID(dataTag, commitTag, connectedExternalNodes);
00344 if (result < 0) {
00345 opserr << "beam2d02::sendSelf - failed to send data\n";
00346 return -1;
00347 }
00348
00349 return 0;
00350 }
00351
00352 int
00353 beam2d02::recvSelf(int commitTag, Channel &theChannel, FEM_ObjectBroker &theBroker)
00354 {
00355 Vector data(4);
00356 int result = 0;
00357 int dataTag = this->getDbTag();
00358
00359 result = theChannel.recvVector(dataTag, commitTag, data);
00360 if (result < 0) {
00361 opserr << "beam2d02::recvSelf - failed to recv data\n";
00362 return -1;
00363 }
00364
00365 A = data(0); E = data(1); I=data(2);
00366 this->setTag((int)data(3));
00367
00368 result = theChannel.recvID(dataTag, commitTag, connectedExternalNodes);
00369 if (result < 0) {
00370 opserr << "beam2d02::recvSelf - failed to recv data\n";
00371 return -1;
00372 }
00373
00374 return 0;
00375 }
00376
00377 int
00378 beam2d02::displaySelf(Renderer &theViewer, int displayMode, float fact)
00379 {
00380
00381
00382
00383 const Vector &end1Crd = theNodes[0]->getCrds();
00384 const Vector &end2Crd = theNodes[1]->getCrds();
00385 const Vector &end1Disp = theNodes[0]->getDisp();
00386 const Vector &end2Disp = theNodes[1]->getDisp();
00387
00388 if (displayMode == 1 || displayMode == 2) {
00389 Vector v1(3);
00390 Vector v2(3);
00391 for (int i=0; i<2; i++) {
00392 v1(i) = end1Crd(i)+end1Disp(i)*fact;
00393 v2(i) = end2Crd(i)+end2Disp(i)*fact;
00394 }
00395
00396 return theViewer.drawLine(v1, v2, 1.0,1.0);
00397 } else
00398 return 0;
00399 }
00400
00401 void
00402 beam2d02::Print(OPS_Stream &s, int flag)
00403 {
00404
00405 this->getResistingForce();
00406
00407 s << "Element: " << this->getTag();
00408 s << " type: beam2d02 iNode: " << connectedExternalNodes(0);
00409 s << " jNode: " << connectedExternalNodes(1);
00410 s << " Area: " << A << " E: " << E << " I: " << I << endln;
00411 s << " resisting Force: " << rForce;
00412 }
00413
00414
00415
00416
00417
00418
00419