|
OOFEM
2.1
|
00001 /* 00002 * 00003 * ##### ##### ###### ###### ### ### 00004 * ## ## ## ## ## ## ## ### ## 00005 * ## ## ## ## #### #### ## # ## 00006 * ## ## ## ## ## ## ## ## 00007 * ## ## ## ## ## ## ## ## 00008 * ##### ##### ## ###### ## ## 00009 * 00010 * 00011 * OOFEM : Object Oriented Finite Element Code 00012 * 00013 * Copyright (C) 1993 - 2013 Borek Patzak 00014 * 00015 * 00016 * 00017 * Czech Technical University, Faculty of Civil Engineering, 00018 * Department of Structural Mechanics, 166 29 Prague, Czech Republic 00019 * 00020 * This program is free software; you can redistribute it and/or modify 00021 * it under the terms of the GNU General Public License as published by 00022 * the Free Software Foundation; either version 2 of the License, or 00023 * (at your option) any later version. 00024 * 00025 * This program is distributed in the hope that it will be useful, 00026 * but WITHOUT ANY WARRANTY; without even the implied warranty of 00027 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00028 * GNU General Public License for more details. 00029 * 00030 * You should have received a copy of the GNU General Public License 00031 * along with this program; if not, write to the Free Software 00032 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. 00033 */ 00034 00035 #include "libeam3dnl.h" 00036 #include "node.h" 00037 #include "material.h" 00038 #include "crosssection.h" 00039 #include "gausspnt.h" 00040 #include "gaussintegrationrule.h" 00041 #include "structuralms.h" 00042 #include "flotmtrx.h" 00043 #include "intarray.h" 00044 #include "flotarry.h" 00045 #include "mathfem.h" 00046 #include "timestep.h" 00047 00048 #ifdef __OOFEG 00049 #include "oofeggraphiccontext.h" 00050 #endif 00051 00052 namespace oofem { 00053 LIBeam3dNL :: LIBeam3dNL(int n, Domain *aDomain) : NLStructuralElement(n, aDomain), tc(3, 3), tempTc(3, 3) //, kappa (3) 00054 { 00055 numberOfDofMans = 2; 00056 l0 = 0.; 00057 tempTcCounter = 0; 00058 referenceNode = 0; 00059 // init kappa vector at centre 00060 // kappa.zero(); 00061 } 00062 00063 00064 void 00065 LIBeam3dNL :: computeSMtrx(FloatMatrix &answer, FloatArray &vec) 00066 { 00067 if ( vec.giveSize() != 3 ) { 00068 _error("computeSMtrx: vec param size mismatch"); 00069 } 00070 00071 answer.resize(3, 3); 00072 00073 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 0.; 00074 answer.at(1, 2) = -vec.at(3); 00075 answer.at(1, 3) = vec.at(2); 00076 answer.at(2, 1) = vec.at(3); 00077 answer.at(2, 3) = -vec.at(1); 00078 answer.at(3, 1) = -vec.at(2); 00079 answer.at(3, 2) = vec.at(1); 00080 } 00081 00082 00083 void 00084 LIBeam3dNL :: computeRotMtrx(FloatMatrix &answer, FloatArray &psi) 00085 { 00086 FloatMatrix S(3, 3), SS(3, 3); 00087 double psiSize; 00088 00089 if ( psi.giveSize() != 3 ) { 00090 _error("computeSMtrx: psi param size mismatch"); 00091 } 00092 00093 answer.resize(3, 3); 00094 answer.zero(); 00095 00096 psiSize = psi.computeNorm(); 00097 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 1.; 00098 00099 if ( psiSize <= 1.e-40 ) { 00100 return; 00101 } 00102 00103 this->computeSMtrx(S, psi); 00104 SS.beProductOf(S, S); 00105 S.times(sin(psiSize) / psiSize); 00106 SS.times( ( 1. - cos(psiSize) ) / ( psiSize * psiSize ) ); 00107 00108 answer.add(S); 00109 answer.add(SS); 00110 } 00111 00112 00113 void 00114 LIBeam3dNL :: updateTempTriad(TimeStep *tStep) 00115 { 00116 // test if not previously done 00117 if ( tStep->giveSolutionStateCounter() == tempTcCounter ) { 00118 return; 00119 } 00120 00121 FloatArray u, centreSpin(3); 00122 FloatMatrix dR(3, 3); 00123 00124 // ask element's displacement increments 00125 this->computeVectorOf(EID_MomentumBalance, VM_Incremental, tStep, u); 00126 00127 // interpolate spin at the centre 00128 centreSpin.at(1) = 0.5 * ( u.at(4) + u.at(10) ); 00129 centreSpin.at(2) = 0.5 * ( u.at(5) + u.at(11) ); 00130 centreSpin.at(3) = 0.5 * ( u.at(6) + u.at(12) ); 00131 00132 // compute rotation matrix from centreSpin pseudovector 00133 this->computeRotMtrx(dR, centreSpin); 00134 00135 // update triad 00136 tempTc.beProductOf(dR, tc); 00137 // remember timestamp 00138 tempTcCounter = tStep->giveSolutionStateCounter(); 00139 } 00140 00141 00142 void 00143 LIBeam3dNL :: computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep) 00144 { 00145 FloatArray xd(3), eps(3), curv(3); 00146 00147 // update temp triad 00148 this->updateTempTriad(tStep); 00149 00150 this->computeXdVector(xd, tStep); 00151 00152 // compute eps_xl, gamma_l2, gamma_l3 00153 eps.beTProductOf(tempTc, xd); 00154 eps.times(1. / this->l0); 00155 eps.at(1) -= 1.0; 00156 00157 // update curvature at midpoint 00158 this->computeTempCurv(curv, tStep); 00159 00160 answer.resize(6); 00161 answer.at(1) = eps.at(1); // eps_xl 00162 answer.at(2) = eps.at(2); // gamma_l2 00163 answer.at(3) = eps.at(3); // gamma_l3 00164 answer.at(4) = curv.at(1); // kappa_1 00165 answer.at(5) = curv.at(2); // kappa_2 00166 answer.at(6) = curv.at(3); // kappa_3 00167 } 00168 00169 00170 void 00171 LIBeam3dNL :: computeXMtrx(FloatMatrix &answer, TimeStep *tStep) { 00172 int i, j; 00173 FloatArray xd(3); 00174 FloatMatrix s(3, 3); 00175 00176 this->computeXdVector(xd, tStep); 00177 this->computeSMtrx(s, xd); 00178 00179 answer.resize(12, 6); 00180 answer.zero(); 00181 00182 for ( i = 1; i < 4; i++ ) { 00183 answer.at(i, i) = -1.0; 00184 answer.at(i + 6, i) = 1.0; 00185 answer.at(i + 3, i + 3) = -1.0; 00186 answer.at(i + 9, i + 3) = 1.0; 00187 00188 for ( j = 1; j < 4; j++ ) { 00189 answer.at(i + 3, j) = answer.at(i + 9, j) = 0.5 * s.at(j, i); 00190 } 00191 } 00192 } 00193 00194 00195 void 00196 LIBeam3dNL :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord) 00197 { 00198 int i, j; 00199 Material *mat = this->giveMaterial(); 00200 IntegrationRule *iRule = integrationRulesArray [ giveDefaultIntegrationRule() ]; 00201 GaussPoint *gp = iRule->getIntegrationPoint(0); 00202 FloatArray nm(6), TotalStressVector(6); 00203 FloatMatrix x; 00204 double s1, s2; 00205 00206 // update temp triad 00207 this->updateTempTriad(tStep); 00208 00209 if ( useUpdatedGpRecord == 1 ) { 00210 TotalStressVector = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) ) 00211 ->giveStressVector(); 00212 } else { 00213 this->computeStressVector(TotalStressVector, gp, tStep); 00214 } 00215 00216 for ( i = 1; i <= 3; i++ ) { 00217 s1 = s2 = 0.0; 00218 for ( j = 1; j <= 3; j++ ) { 00219 s1 += tempTc.at(i, j) * TotalStressVector.at(j); 00220 s2 += tempTc.at(i, j) * TotalStressVector.at(j + 3); 00221 } 00222 00223 nm.at(i) = s1; 00224 nm.at(i + 3) = s2; 00225 } 00226 00227 this->computeXMtrx(x, tStep); 00228 answer.beProductOf(x, nm); 00229 } 00230 00231 00232 void 00233 LIBeam3dNL :: computeXdVector(FloatArray &answer, TimeStep *tStep) 00234 { 00235 FloatArray u(3); 00236 00237 answer.resize(3); 00238 // ask element's displacements 00239 this->computeVectorOf(EID_MomentumBalance, VM_Total, tStep, u); 00240 00241 answer.at(1) = ( this->giveNode(2)->giveCoordinate(1) + u.at(7) ) - 00242 ( this->giveNode(1)->giveCoordinate(1) + u.at(1) ); 00243 answer.at(2) = ( this->giveNode(2)->giveCoordinate(2) + u.at(8) ) - 00244 ( this->giveNode(1)->giveCoordinate(2) + u.at(2) ); 00245 answer.at(3) = ( this->giveNode(2)->giveCoordinate(3) + u.at(9) ) - 00246 ( this->giveNode(1)->giveCoordinate(3) + u.at(3) ); 00247 } 00248 00249 00250 void 00251 LIBeam3dNL :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) 00252 { 00253 int i, j, k; 00254 double s1, s2; 00255 FloatMatrix d, x, xt(12, 6), dxt, sn, sm, sxd, y; 00256 FloatArray n(3), m(3), xd(3), TotalStressVector; 00257 IntegrationRule *iRule = integrationRulesArray [ giveDefaultIntegrationRule() ]; 00258 GaussPoint *gp = iRule->getIntegrationPoint(0); 00259 00260 answer.resize( this->computeNumberOfDofs(EID_MomentumBalance), this->computeNumberOfDofs(EID_MomentumBalance) ); 00261 answer.zero(); 00262 00263 // linear part 00264 00265 this->updateTempTriad(tStep); // update temp triad 00266 this->computeXMtrx(x, tStep); 00267 xt.zero(); 00268 for ( i = 1; i <= 12; i++ ) { 00269 for ( j = 1; j <= 3; j++ ) { 00270 for ( k = 1; k <= 3; k++ ) { 00271 // compute x*Tbar, taking into account sparsity of Tbar 00272 xt.at(i, j) += x.at(i, k) * tempTc.at(k, j); 00273 xt.at(i, j + 3) += x.at(i, k + 3) * tempTc.at(k, j); 00274 } 00275 } 00276 } 00277 00278 this->computeConstitutiveMatrixAt(d, rMode, gp, tStep); 00279 dxt.beProductTOf(d, xt); 00280 answer.beProductOf(xt, dxt); 00281 answer.times(1. / this->l0); 00282 00283 // geometric stiffness ks = ks1+ks2 00284 // ks1 00285 this->computeStressVector(TotalStressVector, gp, tStep); 00286 00287 for ( i = 1; i <= 3; i++ ) { 00288 s1 = s2 = 0.0; 00289 for ( j = 1; j <= 3; j++ ) { 00290 s1 += tempTc.at(i, j) * TotalStressVector.at(j); 00291 s2 += tempTc.at(i, j) * TotalStressVector.at(j + 3); 00292 } 00293 00294 n.at(i) = s1; 00295 m.at(i) = s2; 00296 } 00297 00298 this->computeSMtrx(sn, n); 00299 this->computeSMtrx(sm, m); 00300 00301 for ( i = 1; i <= 3; i++ ) { 00302 for ( j = 1; j <= 3; j++ ) { 00303 answer.at(i, j + 3) += sn.at(i, j); 00304 answer.at(i, j + 9) += sn.at(i, j); 00305 answer.at(i + 3, j + 3) += sm.at(i, j); 00306 answer.at(i + 3, j + 9) += sm.at(i, j); 00307 00308 answer.at(i + 6, j + 3) -= sn.at(i, j); 00309 answer.at(i + 6, j + 9) -= sn.at(i, j); 00310 answer.at(i + 9, j + 3) -= sm.at(i, j); 00311 answer.at(i + 9, j + 9) -= sm.at(i, j); 00312 } 00313 } 00314 00315 // ks2 00316 this->computeXdVector(xd, tStep); 00317 this->computeSMtrx(sxd, xd); 00318 00319 y.beProductOf(sxd, sn); 00320 y.times(0.5); 00321 00322 for ( i = 1; i <= 3; i++ ) { 00323 for ( j = 1; j <= 3; j++ ) { 00324 answer.at(i + 3, j) -= sn.at(i, j); 00325 answer.at(i + 3, j + 3) += y.at(i, j); 00326 answer.at(i + 3, j + 6) += sn.at(i, j); 00327 answer.at(i + 3, j + 9) += y.at(i, j); 00328 00329 answer.at(i + 9, j) -= sn.at(i, j); 00330 answer.at(i + 9, j + 3) += y.at(i, j); 00331 answer.at(i + 9, j + 6) += sn.at(i, j); 00332 answer.at(i + 9, j + 9) += y.at(i, j); 00333 } 00334 } 00335 } 00336 00337 00338 void 00339 LIBeam3dNL :: computeGaussPoints() 00340 // Sets up the array of Gauss Points of the receiver. 00341 { 00342 if ( !integrationRulesArray ) { 00343 numberOfIntegrationRules = 1; 00344 integrationRulesArray = new IntegrationRule * [ 1 ]; 00345 integrationRulesArray [ 0 ] = new GaussIntegrationRule(1, this, 1, 2); 00346 integrationRulesArray [ 0 ]->setUpIntegrationPoints(_Line, 1, _3dBeam); 00347 } 00348 } 00349 00350 00351 IRResultType 00352 LIBeam3dNL :: initializeFrom(InputRecord *ir) 00353 { 00354 const char *__proc = "initializeFrom"; // Required by IR_GIVE_FIELD macro 00355 IRResultType result; // Required by IR_GIVE_FIELD macro 00356 00357 // first call parent 00358 NLStructuralElement :: initializeFrom(ir); 00359 00360 IR_GIVE_FIELD(ir, referenceNode, IFT_LIBeam3dNL_refnode, "refnode"); // Macro 00361 if ( referenceNode == 0 ) { 00362 _error("instanciateFrom: wrong reference node specified"); 00363 } 00364 00365 /* 00366 * if (this->hasString (initString, "dofstocondense")) { 00367 * dofsToCondense = this->ReadIntArray (initString, "dofstocondense"); 00368 * if (dofsToCondense->giveSize() >= 12) 00369 * _error ("instanciateFrom: wrong input data for condensed dofs"); 00370 * } else { 00371 * dofsToCondense = NULL; 00372 * } 00373 */ 00374 00375 // compute initial triad at centre - requires nodal coordinates 00376 FloatMatrix lcs; 00377 this->giveLocalCoordinateSystem(lcs); 00378 this->tc.beTranspositionOf(lcs); 00379 00380 this->computeGaussPoints(); 00381 00382 return IRRT_OK; 00383 } 00384 00385 00386 double 00387 LIBeam3dNL :: giveLength() 00388 // Returns the original length (l0) of the receiver. 00389 { 00390 double dx, dy, dz; 00391 Node *nodeA, *nodeB; 00392 00393 if ( l0 == 0. ) { 00394 nodeA = this->giveNode(1); 00395 nodeB = this->giveNode(2); 00396 dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1); 00397 dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2); 00398 dz = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3); 00399 l0 = sqrt(dx * dx + dy * dy + dz * dz); 00400 } 00401 00402 return l0; 00403 } 00404 00405 00406 void 00407 LIBeam3dNL :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) 00408 // Returns the lumped mass matrix of the receiver. This expression is 00409 // valid in both local and global axes. 00410 { 00411 Material *mat; 00412 double halfMass; 00413 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0); 00414 00415 mat = this->giveMaterial(); 00416 halfMass = mat->give('d', gp) * this->giveCrossSection()->give(CS_Area) * this->giveLength() / 2.; 00417 answer.resize(12, 12); 00418 answer.zero(); 00419 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = halfMass; 00420 answer.at(7, 7) = answer.at(8, 8) = answer.at(9, 9) = halfMass; 00421 00422 double Ik, Iy, Iz; 00423 Ik = this->giveCrossSection()->give(CS_TorsionMomentX); 00424 Iy = this->giveCrossSection()->give(CS_InertiaMomentY); 00425 Iz = this->giveCrossSection()->give(CS_InertiaMomentZ); 00426 halfMass = mat->give('d', gp) * this->giveLength() / 2.; 00427 answer.at(4, 4) = answer.at(10, 10) = Ik * halfMass; 00428 answer.at(5, 5) = answer.at(11, 11) = Iy * halfMass; 00429 answer.at(6, 6) = answer.at(12, 12) = Iz * halfMass; 00430 } 00431 00432 00433 void 00434 LIBeam3dNL :: computeNmatrixAt(GaussPoint *aGaussPoint, FloatMatrix &answer) 00435 // Returns the displacement interpolation matrix {N} of the receiver, eva- 00436 // luated at aGaussPoint. 00437 { 00438 double ksi, n1, n2; 00439 00440 ksi = aGaussPoint->giveCoordinate(1); 00441 n1 = ( 1. - ksi ) * 0.5; 00442 n2 = ( 1. + ksi ) * 0.5; 00443 00444 answer.resize(6, 12); 00445 answer.zero(); 00446 00447 // u 00448 answer.at(1, 1) = n1; 00449 answer.at(1, 7) = n2; 00450 // v 00451 answer.at(2, 2) = n1; 00452 answer.at(2, 8) = n2; 00453 // w 00454 answer.at(3, 3) = n1; 00455 answer.at(3, 9) = n2; 00456 // fi_x 00457 answer.at(4, 4) = n1; 00458 answer.at(4, 10) = n2; 00459 // fi_y 00460 answer.at(5, 5) = n1; 00461 answer.at(5, 11) = n2; 00462 // fi_z 00463 answer.at(6, 6) = n1; 00464 answer.at(6, 12) = n2; 00465 } 00466 00467 00468 double 00469 LIBeam3dNL :: computeVolumeAround(GaussPoint *aGaussPoint) 00470 // Returns the length of the receiver. This method is valid only if 1 00471 // Gauss point is used. 00472 { 00473 double weight = aGaussPoint->giveWeight(); 00474 return weight * 0.5 * this->giveLength(); 00475 } 00476 00477 00478 void 00479 LIBeam3dNL :: giveDofManDofIDMask(int inode, EquationID, IntArray &answer) const 00480 { 00481 answer.setValues(6, D_u, D_v, D_w, R_u, R_v, R_w); 00482 } 00483 00484 00485 int 00486 LIBeam3dNL :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords) 00487 { 00488 double ksi, n1, n2; 00489 00490 ksi = lcoords.at(1); 00491 n1 = ( 1. - ksi ) * 0.5; 00492 n2 = ( 1. + ksi ) * 0.5; 00493 00494 answer.resize(3); 00495 answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1); 00496 answer.at(2) = n1 * this->giveNode(1)->giveCoordinate(2) + n2 *this->giveNode(2)->giveCoordinate(2); 00497 answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3); 00498 00499 return 1; 00500 } 00501 00502 00503 void 00504 LIBeam3dNL :: computeEgdeNMatrixAt(FloatMatrix &answer, int iedge, GaussPoint *aGaussPoint) 00505 { 00506 /* 00507 * 00508 * computes interpolation matrix for element edge. 00509 * we assemble locally this matrix for only nonzero 00510 * shape functions. 00511 * (for example only two nonzero shape functions for 2 dofs are 00512 * necessary for linear plane stress tringle edge). 00513 * These nonzero shape functions are then mapped to 00514 * global element functions. 00515 * 00516 * Uso of mapping technique will allow to assemble shape functions 00517 * without regarding particular side 00518 */ 00519 00520 this->computeNmatrixAt(aGaussPoint, answer); 00521 } 00522 00523 00524 void 00525 LIBeam3dNL :: giveEdgeDofMapping(IntArray &answer, int iEdge) const 00526 { 00527 /* 00528 * provides dof mapping of local edge dofs (only nonzero are taken into account) 00529 * to global element dofs 00530 */ 00531 int i; 00532 00533 if ( iEdge != 1 ) { 00534 _error("giveEdgeDofMapping: wrong edge number"); 00535 } 00536 00537 answer.resize(12); 00538 for ( i = 1; i <= 12; i++ ) { 00539 answer.at(i) = i; 00540 } 00541 } 00542 00543 00544 double 00545 LIBeam3dNL :: computeEdgeVolumeAround(GaussPoint *aGaussPoint, int iEdge) 00546 { 00547 if ( iEdge != 1 ) { // edge between nodes 1 2 00548 _error("computeEdgeVolumeAround: wrong egde number"); 00549 } 00550 00551 double weight = aGaussPoint->giveWeight(); 00552 return 0.5 * this->giveLength() * weight; 00553 } 00554 00555 00556 int 00557 LIBeam3dNL :: giveLocalCoordinateSystem(FloatMatrix &answer) 00558 // 00559 // returns a unit vectors of local coordinate system at element 00560 // stored rowwise (mainly used by some materials with ortho and anisotrophy) 00561 // 00562 { 00563 FloatArray lx(3), ly(3), lz(3), help(3); 00564 double length = this->giveLength(); 00565 Node *nodeA, *nodeB, *refNode; 00566 int i; 00567 00568 answer.resize(3, 3); 00569 answer.zero(); 00570 nodeA = this->giveNode(1); 00571 nodeB = this->giveNode(2); 00572 refNode = this->giveDomain()->giveNode(this->referenceNode); 00573 00574 for ( i = 1; i <= 3; i++ ) { 00575 lx.at(i) = ( nodeB->giveCoordinate(i) - nodeA->giveCoordinate(i) ) / length; 00576 help.at(i) = ( refNode->giveCoordinate(i) - nodeA->giveCoordinate(i) ); 00577 } 00578 00579 lz.beVectorProductOf(lx, help); 00580 lz.normalize(); 00581 ly.beVectorProductOf(lz, lx); 00582 ly.normalize(); 00583 00584 for ( i = 1; i <= 3; i++ ) { 00585 answer.at(1, i) = lx.at(i); 00586 answer.at(2, i) = ly.at(i); 00587 answer.at(3, i) = lz.at(i); 00588 } 00589 00590 return 1; 00591 } 00592 00593 00594 /* 00595 * int 00596 * LIBeam3dNL :: computeGtoLRotationMatrix (FloatMatrix& answer) 00597 * // Returns the rotation matrix of the receiver. 00598 * // rotation matrix is computed for original configuration only 00599 * { 00600 * FloatMatrix lcs; 00601 * int i,j; 00602 * 00603 * answer.resize(12,12); 00604 * answer.zero(); 00605 * 00606 * this->giveLocalCoordinateSystem (lcs); 00607 * 00608 * for (i=1; i <= 3; i++) 00609 * for (j=1; j <= 3; j++) { 00610 * answer.at(i,j) = lcs.at(i,j); 00611 * answer.at(i+3, j+3) = lcs.at(i,j); 00612 * answer.at(i+6, j+6) = lcs.at(i,j); 00613 * answer.at(i+9, j+9) = lcs.at(i,j); 00614 * } 00615 * 00616 * return 1 ; 00617 * } 00618 */ 00619 00620 int 00621 LIBeam3dNL :: computeLoadGToLRotationMtrx(FloatMatrix &answer) 00622 { 00623 /* 00624 * Returns transformation matrix from global coordinate system to local 00625 * element coordinate system for element load vector components. 00626 * If no transformation is necessary, answer is empty matrix (default); 00627 * 00628 * Does not support follower load 00629 */ 00630 00631 FloatMatrix lcs; 00632 int i, j; 00633 00634 answer.resize(6, 6); 00635 answer.zero(); 00636 00637 this->giveLocalCoordinateSystem(lcs); 00638 00639 for ( i = 1; i <= 3; i++ ) { 00640 for ( j = 1; j <= 3; j++ ) { 00641 answer.at(i, j) = lcs.at(i, j); 00642 answer.at(3 + i, 3 + j) = lcs.at(i, j); 00643 } 00644 } 00645 00646 return 1; 00647 } 00648 00649 00650 int 00651 LIBeam3dNL :: computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp) 00652 { 00653 // returns transformation matrix from 00654 // edge local coordinate system 00655 // to element local c.s 00656 // (same as global c.s in this case) 00657 // 00658 // i.e. f(element local) = T * f(edge local) 00659 // 00660 answer.beEmptyMtrx(); 00661 return 0; 00662 } 00663 00664 00665 void 00666 LIBeam3dNL :: computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode) 00667 { 00668 NLStructuralElement::computeBodyLoadVectorAt(answer, load, tStep, mode); 00669 answer.times(this->giveCrossSection()->give(CS_Area)); 00670 } 00671 00672 00673 00674 void LIBeam3dNL :: updateYourself(TimeStep *tStep) 00675 // Updates the receiver at end of step. 00676 { 00677 NLStructuralElement :: updateYourself(tStep); 00678 00679 // update triad 00680 this->updateTempTriad(tStep); 00681 tc = tempTc; 00682 00683 // update curvature 00684 //FloatArray curv; 00685 //this->computeTempCurv (curv, tStep); 00686 //kappa = curv; 00687 } 00688 00689 void 00690 LIBeam3dNL :: initForNewStep() 00691 // initializes receiver to new time step or can be used 00692 // if current time step must be restarted 00693 // 00694 // call material->initGpForNewStep() for all GPs. 00695 // 00696 { 00697 NLStructuralElement :: initForNewStep(); 00698 tempTc = tc; 00699 } 00700 00701 00702 00703 void 00704 LIBeam3dNL :: computeTempCurv(FloatArray &answer, TimeStep *tStep) 00705 { 00706 Material *mat = this->giveMaterial(); 00707 IntegrationRule *iRule = integrationRulesArray [ giveDefaultIntegrationRule() ]; 00708 GaussPoint *gp = iRule->getIntegrationPoint(0); 00709 ; 00710 FloatArray ui(3), xd(3), curv(3), ac(3), PrevEpsilon; 00711 FloatMatrix sc(3, 3), tmid(3, 3); 00712 00713 answer.resize(3); 00714 00715 // update curvature at midpoint 00716 // first, compute Tmid 00717 // ask increments 00718 this->computeVectorOf(EID_MomentumBalance, VM_Incremental, tStep, ui); 00719 00720 ac.at(1) = 0.5 * ( ui.at(10) - ui.at(4) ); 00721 ac.at(2) = 0.5 * ( ui.at(11) - ui.at(5) ); 00722 ac.at(3) = 0.5 * ( ui.at(12) - ui.at(6) ); 00723 this->computeSMtrx(sc, ac); 00724 sc.times(1. / 2.); 00725 // compute I+sc 00726 sc.at(1, 1) += 1.0; 00727 sc.at(2, 2) += 1.0; 00728 sc.at(3, 3) += 1.0; 00729 tmid.beProductOf(sc, this->tc); 00730 00731 // update curvature at centre 00732 ac.at(1) = ( ui.at(10) - ui.at(4) ); 00733 ac.at(2) = ( ui.at(11) - ui.at(5) ); 00734 ac.at(3) = ( ui.at(12) - ui.at(6) ); 00735 00736 answer.beTProductOf(tmid, ac); 00737 answer.times(1 / this->l0); 00738 00739 // ask for previous kappa 00740 PrevEpsilon = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) )->giveStrainVector(); 00741 if ( PrevEpsilon.giveSize() ) { 00742 answer.at(1) += PrevEpsilon.at(4); 00743 answer.at(2) += PrevEpsilon.at(5); 00744 answer.at(3) += PrevEpsilon.at(6); 00745 } 00746 } 00747 00748 00749 #ifdef __OOFEG 00750 void LIBeam3dNL :: drawRawGeometry(oofegGraphicContext &gc) 00751 { 00752 GraphicObj *go; 00753 00754 if ( !gc.testElementGraphicActivity(this) ) { 00755 return; 00756 } 00757 00758 // if (!go) { // create new one 00759 WCRec p [ 2 ]; /* poin */ 00760 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH); 00761 EASValsSetColor( gc.getElementColor() ); 00762 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER); 00763 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1); 00764 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2); 00765 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3); 00766 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1); 00767 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2); 00768 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3); 00769 go = CreateLine3D(p); 00770 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go); 00771 EGAttachObject(go, ( EObjectP ) this); 00772 EMAddGraphicsToModel(ESIModel(), go); 00773 } 00774 00775 00776 void LIBeam3dNL :: drawDeformedGeometry(oofegGraphicContext &gc, UnknownType type) 00777 { 00778 GraphicObj *go; 00779 00780 if ( !gc.testElementGraphicActivity(this) ) { 00781 return; 00782 } 00783 00784 TimeStep *tStep = domain->giveEngngModel()->giveCurrentStep(); 00785 double defScale = gc.getDefScale(); 00786 // if (!go) { // create new one 00787 WCRec p [ 2 ]; /* poin */ 00788 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH); 00789 EASValsSetColor( gc.getDeformedElementColor() ); 00790 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER); 00791 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, EID_MomentumBalance, defScale); 00792 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, EID_MomentumBalance, defScale); 00793 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, EID_MomentumBalance, defScale); 00794 00795 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, EID_MomentumBalance, defScale); 00796 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, EID_MomentumBalance, defScale); 00797 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, EID_MomentumBalance, defScale); 00798 go = CreateLine3D(p); 00799 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go); 00800 EMAddGraphicsToModel(ESIModel(), go); 00801 } 00802 #endif 00803 } // end namespace oofem