Rev 4864 | Details | Compare with Previous | Last modification | View Log | RSS feed
| Rev | Author | Line No. | Line |
|---|---|---|---|
| 4617 | fmk | 1 | /* ****************************************************************** ** |
| 2 | ** OpenSees - Open System for Earthquake Engineering Simulation ** |
||
| 3 | ** Pacific Earthquake Engineering Research Center ** |
||
| 4 | ** ** |
||
| 5141 | fmk | 5 | ** 1 ** |
| 4617 | fmk | 6 | ** (C) Copyright 1999, The Regents of the University of California ** |
| 7 | ** All Rights Reserved. ** |
||
| 8 | ** ** |
||
| 9 | ** Commercial use of this program without express permission of the ** |
||
| 10 | ** University of California, Berkeley, is strictly prohibited. See ** |
||
| 11 | ** file 'COPYRIGHT' in main directory for information on usage and ** |
||
| 12 | ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. ** |
||
| 13 | ** ** |
||
| 14 | ** Developed by: ** |
||
| 15 | ** Frank McKenna (fmckenna@ce.berkeley.edu) ** |
||
| 16 | ** Gregory L. Fenves (fenves@ce.berkeley.edu) ** |
||
| 17 | ** Filip C. Filippou (filippou@ce.berkeley.edu) ** |
||
| 18 | ** ** |
||
| 19 | ** ****************************************************************** */ |
||
| 20 | |||
| 21 | // $Revision: 1.10 $ |
||
| 22 | // $Date: 2011/03/10 22:51:21 $ |
||
| 23 | // $Source: /usr/local/cvs/OpenSees/SRC/element/shell/ShellMITC4.cpp,v $ |
||
| 24 | |||
| 25 | // Written: Leopoldo Tesser, Diego A. Talledo, Véronique Le Corvec |
||
| 26 | // |
||
| 27 | // Bathe MITC 4 four node shell element with membrane and drill |
||
| 28 | // Ref: Dvorkin,Bathe, A continuum mechanics based four node shell |
||
| 29 | // element for general nonlinear analysis, |
||
| 30 | // Eng.Comput.,1,77-88,1984 |
||
| 31 | |||
| 32 | #include <stdio.h> |
||
| 33 | #include <stdlib.h> |
||
| 34 | #include <math.h> |
||
| 35 | |||
| 36 | #include <ID.h> |
||
| 37 | #include <Vector.h> |
||
| 38 | #include <Matrix.h> |
||
| 39 | #include <Element.h> |
||
| 40 | #include <Node.h> |
||
| 41 | #include <SectionForceDeformation.h> |
||
| 42 | #include <Domain.h> |
||
| 43 | #include <ErrorHandler.h> |
||
| 44 | #include <ShellMITC4.h> |
||
| 45 | #include <R3vectors.h> |
||
| 46 | #include <Renderer.h> |
||
| 47 | #include <ElementResponse.h> |
||
| 48 | |||
| 49 | #include <Channel.h> |
||
| 50 | #include <FEM_ObjectBroker.h> |
||
| 51 | #include <elementAPI.h> |
||
| 52 | |||
| 53 | #define min(a,b) ( (a)<(b) ? (a):(b) ) |
||
| 54 | |||
| 55 | static int numShellMITC4 = 0; |
||
| 56 | |||
| 57 | void * |
||
| 58 | OPS_NewShellMITC4(void) |
||
| 59 | { |
||
| 60 | if (numShellMITC4 == 0) { |
||
| 5141 | fmk | 61 | opserr << "Using ShellMITC4 - Developed by: Leopoldo Tesser, Diego A. Talledo, Veronique Le Corvec\n"; |
| 4617 | fmk | 62 | numShellMITC4++; |
| 63 | } |
||
| 64 | |||
| 65 | Element *theElement = 0; |
||
| 66 | |||
| 67 | int numArgs = OPS_GetNumRemainingInputArgs(); |
||
| 68 | |||
| 69 | if (numArgs < 6) { |
||
| 70 | opserr << "Want: element ShellMITC4 $tag $iNode $jNoe $kNode $lNode $secTag"; |
||
| 71 | return 0; |
||
| 72 | } |
||
| 73 | |||
| 74 | int iData[6]; |
||
| 75 | int numData = 6; |
||
| 76 | if (OPS_GetInt(&numData, iData) != 0) { |
||
| 77 | opserr << "WARNING invalid integer tag: element ShellMITC4 \n"; |
||
| 78 | return 0; |
||
| 79 | } |
||
| 80 | |||
| 81 | SectionForceDeformation *theSection = OPS_GetSectionForceDeformation(iData[5]); |
||
| 82 | |||
| 83 | if (theSection == 0) { |
||
| 4864 | fmk | 84 | opserr << "ERROR: element ShellMITC4 " << iData[0] << "section " << iData[5] << " not found\n"; |
| 4617 | fmk | 85 | return 0; |
| 86 | } |
||
| 87 | |||
| 88 | theElement = new ShellMITC4(iData[0], iData[1], iData[2], iData[3], |
||
| 89 | iData[4], *theSection); |
||
| 90 | |||
| 91 | return theElement; |
||
| 92 | } |
||
| 93 | |||
| 94 | |||
| 95 | //static data |
||
| 96 | Matrix ShellMITC4::stiff(24,24) ; |
||
| 97 | Vector ShellMITC4::resid(24) ; |
||
| 98 | Matrix ShellMITC4::mass(24,24) ; |
||
| 99 | |||
| 100 | //quadrature data |
||
| 101 | const double ShellMITC4::root3 = sqrt(3.0) ; |
||
| 102 | const double ShellMITC4::one_over_root3 = 1.0 / root3 ; |
||
| 103 | |||
| 104 | double ShellMITC4::sg[4] ; |
||
| 105 | double ShellMITC4::tg[4] ; |
||
| 106 | double ShellMITC4::wg[4] ; |
||
| 107 | |||
| 108 | |||
| 109 | |||
| 110 | //null constructor |
||
| 111 | ShellMITC4::ShellMITC4( ) : |
||
| 112 | Element( 0, ELE_TAG_ShellMITC4 ), |
||
| 113 | connectedExternalNodes(4), load(0), Ki(0) |
||
| 114 | { |
||
| 115 | for (int i = 0 ; i < 4; i++ ) |
||
| 116 | materialPointers[i] = 0; |
||
| 117 | |||
| 118 | sg[0] = -one_over_root3; |
||
| 119 | sg[1] = one_over_root3; |
||
| 120 | sg[2] = one_over_root3; |
||
| 121 | sg[3] = -one_over_root3; |
||
| 122 | |||
| 123 | tg[0] = -one_over_root3; |
||
| 124 | tg[1] = -one_over_root3; |
||
| 125 | tg[2] = one_over_root3; |
||
| 126 | tg[3] = one_over_root3; |
||
| 127 | |||
| 128 | wg[0] = 1.0; |
||
| 129 | wg[1] = 1.0; |
||
| 130 | wg[2] = 1.0; |
||
| 131 | wg[3] = 1.0; |
||
| 132 | } |
||
| 133 | |||
| 134 | |||
| 135 | //********************************************************************* |
||
| 136 | //full constructor |
||
| 137 | ShellMITC4::ShellMITC4( int tag, |
||
| 138 | int node1, |
||
| 139 | int node2, |
||
| 140 | int node3, |
||
| 141 | int node4, |
||
| 142 | SectionForceDeformation &theMaterial ) : |
||
| 143 | Element( tag, ELE_TAG_ShellMITC4 ), |
||
| 144 | connectedExternalNodes(4), load(0), Ki(0) |
||
| 145 | { |
||
| 146 | int i; |
||
| 147 | |||
| 148 | connectedExternalNodes(0) = node1 ; |
||
| 149 | connectedExternalNodes(1) = node2 ; |
||
| 150 | connectedExternalNodes(2) = node3 ; |
||
| 151 | connectedExternalNodes(3) = node4 ; |
||
| 152 | |||
| 153 | for ( i = 0 ; i < 4; i++ ) { |
||
| 154 | |||
| 155 | materialPointers[i] = theMaterial.getCopy( ) ; |
||
| 156 | |||
| 157 | if (materialPointers[i] == 0) { |
||
| 158 | opserr << "ShellMITC4::constructor - failed to get a material of type: ShellSection\n"; |
||
| 159 | } //end if |
||
| 160 | |||
| 161 | } //end for i |
||
| 162 | |||
| 163 | sg[0] = -one_over_root3; |
||
| 164 | sg[1] = one_over_root3; |
||
| 165 | sg[2] = one_over_root3; |
||
| 166 | sg[3] = -one_over_root3; |
||
| 167 | |||
| 168 | tg[0] = -one_over_root3; |
||
| 169 | tg[1] = -one_over_root3; |
||
| 170 | tg[2] = one_over_root3; |
||
| 171 | tg[3] = one_over_root3; |
||
| 172 | |||
| 173 | wg[0] = 1.0; |
||
| 174 | wg[1] = 1.0; |
||
| 175 | wg[2] = 1.0; |
||
| 176 | wg[3] = 1.0; |
||
| 177 | |||
| 178 | } |
||
| 179 | //****************************************************************** |
||
| 180 | |||
| 181 | //destructor |
||
| 182 | ShellMITC4::~ShellMITC4( ) |
||
| 183 | { |
||
| 184 | int i ; |
||
| 185 | for ( i = 0 ; i < 4; i++ ) { |
||
| 186 | |||
| 187 | delete materialPointers[i] ; |
||
| 188 | materialPointers[i] = 0 ; |
||
| 189 | |||
| 190 | nodePointers[i] = 0 ; |
||
| 191 | |||
| 192 | } //end for i |
||
| 193 | |||
| 194 | if (load != 0) |
||
| 195 | delete load; |
||
| 196 | |||
| 197 | if (Ki != 0) |
||
| 198 | delete Ki; |
||
| 199 | } |
||
| 200 | //************************************************************************** |
||
| 201 | |||
| 202 | |||
| 203 | //set domain |
||
| 204 | void ShellMITC4::setDomain( Domain *theDomain ) |
||
| 205 | { |
||
| 206 | int i, j ; |
||
| 207 | static Vector eig(3) ; |
||
| 208 | static Matrix ddMembrane(3,3) ; |
||
| 209 | |||
| 210 | //node pointers |
||
| 211 | for ( i = 0; i < 4; i++ ) { |
||
| 212 | nodePointers[i] = theDomain->getNode( connectedExternalNodes(i) ) ; |
||
| 213 | if (nodePointers[i] == 0) { |
||
| 214 | opserr << "ShellMITC4::setDomain - no node " << connectedExternalNodes(i); |
||
| 215 | opserr << " exists in the model\n"; |
||
| 216 | } |
||
| 4787 | fmk | 217 | const Vector &nodeDisp=nodePointers[i]->getTrialDisp(); |
| 218 | if (nodeDisp.Size() != 6) { |
||
| 219 | opserr << "ShellMITC4::setDomain - node " << connectedExternalNodes(i); |
||
| 220 | opserr << " NEEDS 6 dof - GARBAGE RESULTS or SEGMENTATION FAULT WILL FOLLOW\n"; |
||
| 221 | } |
||
| 4617 | fmk | 222 | } |
| 223 | |||
| 224 | //compute drilling stiffness penalty parameter |
||
| 225 | const Matrix &dd = materialPointers[0]->getInitialTangent( ) ; |
||
| 226 | |||
| 227 | //assemble ddMembrane ; |
||
| 228 | for ( i = 0; i < 3; i++ ) { |
||
| 229 | for ( j = 0; j < 3; j++ ) |
||
| 230 | ddMembrane(i,j) = dd(i,j) ; |
||
| 231 | } //end for i |
||
| 232 | |||
| 233 | //eigenvalues of ddMembrane |
||
| 234 | eig = LovelyEig( ddMembrane ) ; |
||
| 235 | |||
| 236 | //set ktt |
||
| 237 | //Ktt = dd(2,2) ; //shear modulus |
||
| 238 | Ktt = min( eig(2), min( eig(0), eig(1) ) ) ; |
||
| 239 | //Ktt = dd(2,2); |
||
| 240 | |||
| 241 | //basis vectors and local coordinates |
||
| 242 | computeBasis( ) ; |
||
| 243 | |||
| 244 | this->DomainComponent::setDomain(theDomain); |
||
| 5141 | fmk | 245 | |
| 4617 | fmk | 246 | } |
| 247 | |||
| 248 | |||
| 249 | //get the number of external nodes |
||
| 250 | int ShellMITC4::getNumExternalNodes( ) const |
||
| 251 | { |
||
| 252 | return 4 ; |
||
| 253 | } |
||
| 254 | |||
| 255 | |||
| 256 | //return connected external nodes |
||
| 257 | const ID& ShellMITC4::getExternalNodes( ) |
||
| 258 | { |
||
| 259 | return connectedExternalNodes ; |
||
| 260 | } |
||
| 261 | |||
| 262 | |||
| 263 | Node ** |
||
| 264 | ShellMITC4::getNodePtrs(void) |
||
| 265 | { |
||
| 266 | return nodePointers; |
||
| 267 | } |
||
| 268 | |||
| 269 | //return number of dofs |
||
| 270 | int ShellMITC4::getNumDOF( ) |
||
| 271 | { |
||
| 272 | return 24 ; |
||
| 273 | } |
||
| 274 | |||
| 275 | |||
| 276 | //commit state |
||
| 277 | int ShellMITC4::commitState( ) |
||
| 278 | { |
||
| 279 | int success = 0 ; |
||
| 280 | |||
| 281 | // call element commitState to do any base class stuff |
||
| 282 | if ((success = this->Element::commitState()) != 0) { |
||
| 283 | opserr << "ShellMITC4::commitState () - failed in base class"; |
||
| 284 | } |
||
| 285 | |||
| 286 | for (int i = 0; i < 4; i++ ) |
||
| 287 | success += materialPointers[i]->commitState( ) ; |
||
| 288 | |||
| 289 | return success ; |
||
| 290 | } |
||
| 291 | |||
| 292 | |||
| 293 | |||
| 294 | //revert to last commit |
||
| 295 | int ShellMITC4::revertToLastCommit( ) |
||
| 296 | { |
||
| 297 | int i ; |
||
| 298 | int success = 0 ; |
||
| 299 | |||
| 300 | for ( i = 0; i < 4; i++ ) |
||
| 301 | success += materialPointers[i]->revertToLastCommit( ) ; |
||
| 302 | |||
| 303 | return success ; |
||
| 304 | } |
||
| 305 | |||
| 306 | |||
| 307 | //revert to start |
||
| 308 | int ShellMITC4::revertToStart( ) |
||
| 309 | { |
||
| 310 | int i ; |
||
| 311 | int success = 0 ; |
||
| 312 | |||
| 313 | for ( i = 0; i < 4; i++ ) |
||
| 314 | success += materialPointers[i]->revertToStart( ) ; |
||
| 315 | |||
| 316 | return success ; |
||
| 317 | } |
||
| 318 | |||
| 319 | //print out element data |
||
| 320 | void ShellMITC4::Print( OPS_Stream &s, int flag ) |
||
| 321 | { |
||
| 322 | if (flag == -1) { |
||
| 323 | int eleTag = this->getTag(); |
||
| 324 | s << "EL_ShellMITC4\t" << eleTag << "\t"; |
||
| 325 | s << eleTag << "\t" << 1; |
||
| 326 | s << "\t" << connectedExternalNodes(0) << "\t" << connectedExternalNodes(1); |
||
| 327 | s << "\t" << connectedExternalNodes(2) << "\t" << connectedExternalNodes(3) << "\t0.00"; |
||
| 328 | s << endln; |
||
| 329 | s << "PROP_3D\t" << eleTag << "\t"; |
||
| 330 | s << eleTag << "\t" << 1; |
||
| 331 | s << "\t" << -1 << "\tSHELL\t1.0\0.0"; |
||
| 332 | s << endln; |
||
| 333 | } else if (flag < -1) { |
||
| 334 | |||
| 335 | int counter = (flag + 1) * -1; |
||
| 336 | int eleTag = this->getTag(); |
||
| 337 | int i,j; |
||
| 338 | for ( i = 0; i < 4; i++ ) { |
||
| 339 | const Vector &stress = materialPointers[i]->getStressResultant(); |
||
| 340 | |||
| 341 | s << "STRESS\t" << eleTag << "\t" << counter << "\t" << i << "\tTOP"; |
||
| 342 | for (j=0; j<6; j++) |
||
| 343 | s << "\t" << stress(j); |
||
| 344 | s << endln; |
||
| 345 | } |
||
| 346 | |||
| 347 | } else { |
||
| 348 | s << endln ; |
||
| 349 | s << "MITC4 Non-Locking Four Node Shell \n" ; |
||
| 350 | s << "Element Number: " << this->getTag() << endln ; |
||
| 351 | s << "Node 1 : " << connectedExternalNodes(0) << endln ; |
||
| 352 | s << "Node 2 : " << connectedExternalNodes(1) << endln ; |
||
| 353 | s << "Node 3 : " << connectedExternalNodes(2) << endln ; |
||
| 354 | s << "Node 4 : " << connectedExternalNodes(3) << endln ; |
||
| 355 | |||
| 356 | s << "Material Information : \n " ; |
||
| 357 | materialPointers[0]->Print( s, flag ) ; |
||
| 358 | |||
| 359 | s << endln ; |
||
| 360 | } |
||
| 361 | } |
||
| 362 | |||
| 363 | Response* |
||
| 364 | ShellMITC4::setResponse(const char **argv, int argc, OPS_Stream &output) |
||
| 365 | { |
||
| 366 | Response *theResponse = 0; |
||
| 367 | |||
| 368 | output.tag("ElementOutput"); |
||
| 369 | output.attr("eleType", "ShellMITC4"); |
||
| 370 | output.attr("eleTag",this->getTag()); |
||
| 371 | int numNodes = this->getNumExternalNodes(); |
||
| 372 | const ID &nodes = this->getExternalNodes(); |
||
| 373 | static char nodeData[32]; |
||
| 374 | |||
| 375 | for (int i=0; i<numNodes; i++) { |
||
| 376 | sprintf(nodeData,"node%d",i+1); |
||
| 377 | output.attr(nodeData,nodes(i)); |
||
| 378 | } |
||
| 379 | |||
| 380 | if (strcmp(argv[0],"force") == 0 || strcmp(argv[0],"forces") == 0 || |
||
| 381 | strcmp(argv[0],"globalForce") == 0 || strcmp(argv[0],"globalForces") == 0) { |
||
| 382 | const Vector &force = this->getResistingForce(); |
||
| 383 | int size = force.Size(); |
||
| 384 | for (int i=0; i<size; i++) { |
||
| 385 | sprintf(nodeData,"P%d",i+1); |
||
| 386 | output.tag("ResponseType",nodeData); |
||
| 387 | } |
||
| 388 | theResponse = new ElementResponse(this, 1, this->getResistingForce()); |
||
| 389 | } |
||
| 390 | |||
| 391 | else if (strcmp(argv[0],"material") == 0 || strcmp(argv[0],"Material") == 0) { |
||
| 392 | if (argc < 2) { |
||
| 393 | opserr << "ShellMITC4::setResponse() - need to specify more data\n"; |
||
| 394 | return 0; |
||
| 395 | } |
||
| 396 | int pointNum = atoi(argv[1]); |
||
| 397 | if (pointNum > 0 && pointNum <= 4) { |
||
| 398 | |||
| 399 | output.tag("GaussPoint"); |
||
| 400 | output.attr("number",pointNum); |
||
| 401 | output.attr("eta",sg[pointNum-1]); |
||
| 402 | output.attr("neta",tg[pointNum-1]); |
||
| 403 | |||
| 404 | theResponse = materialPointers[pointNum-1]->setResponse(&argv[2], argc-2, output); |
||
| 405 | |||
| 406 | output.endTag(); |
||
| 407 | } |
||
| 408 | |||
| 409 | } else if (strcmp(argv[0],"stresses") ==0) { |
||
| 410 | |||
| 411 | for (int i=0; i<4; i++) { |
||
| 412 | output.tag("GaussPoint"); |
||
| 413 | output.attr("number",i+1); |
||
| 414 | output.attr("eta",sg[i]); |
||
| 415 | output.attr("neta",tg[i]); |
||
| 416 | |||
| 417 | output.tag("SectionForceDeformation"); |
||
| 418 | output.attr("classType", materialPointers[i]->getClassTag()); |
||
| 419 | output.attr("tag", materialPointers[i]->getTag()); |
||
| 420 | |||
| 421 | output.tag("ResponseType","p11"); |
||
| 422 | output.tag("ResponseType","p22"); |
||
| 423 | output.tag("ResponseType","p1212"); |
||
| 424 | output.tag("ResponseType","m11"); |
||
| 425 | output.tag("ResponseType","m22"); |
||
| 426 | output.tag("ResponseType","m12"); |
||
| 427 | output.tag("ResponseType","q1"); |
||
| 428 | output.tag("ResponseType","q2"); |
||
| 429 | |||
| 430 | output.endTag(); // GaussPoint |
||
| 431 | output.endTag(); // NdMaterialOutput |
||
| 432 | } |
||
| 433 | |||
| 434 | theResponse = new ElementResponse(this, 2, Vector(32)); |
||
| 435 | } |
||
| 436 | |||
| 437 | else if (strcmp(argv[0],"strains") ==0) { |
||
| 438 | |||
| 439 | for (int i=0; i<4; i++) { |
||
| 440 | output.tag("GaussPoint"); |
||
| 441 | output.attr("number",i+1); |
||
| 442 | output.attr("eta",sg[i]); |
||
| 443 | output.attr("neta",tg[i]); |
||
| 444 | |||
| 445 | output.tag("SectionForceDeformation"); |
||
| 446 | output.attr("classType", materialPointers[i]->getClassTag()); |
||
| 447 | output.attr("tag", materialPointers[i]->getTag()); |
||
| 448 | |||
| 449 | output.tag("ResponseType","eps11"); |
||
| 450 | output.tag("ResponseType","eps22"); |
||
| 451 | output.tag("ResponseType","gamma12"); |
||
| 452 | output.tag("ResponseType","theta11"); |
||
| 453 | output.tag("ResponseType","theta22"); |
||
| 454 | output.tag("ResponseType","theta33"); |
||
| 455 | output.tag("ResponseType","gamma13"); |
||
| 456 | output.tag("ResponseType","gamma23"); |
||
| 457 | |||
| 458 | output.endTag(); // GaussPoint |
||
| 459 | output.endTag(); // NdMaterialOutput |
||
| 460 | } |
||
| 461 | |||
| 462 | theResponse = new ElementResponse(this, 3, Vector(32)); |
||
| 463 | } |
||
| 464 | |||
| 465 | output.endTag(); |
||
| 466 | return theResponse; |
||
| 467 | } |
||
| 468 | |||
| 469 | int |
||
| 470 | ShellMITC4::getResponse(int responseID, Information &eleInfo) |
||
| 471 | { |
||
| 472 | int cnt = 0; |
||
| 473 | static Vector stresses(32); |
||
| 474 | static Vector strains(32); |
||
| 475 | |||
| 476 | switch (responseID) { |
||
| 477 | case 1: // global forces |
||
| 478 | return eleInfo.setVector(this->getResistingForce()); |
||
| 479 | break; |
||
| 480 | |||
| 481 | case 2: // stresses |
||
| 482 | for (int i = 0; i < 4; i++) { |
||
| 483 | |||
| 484 | // Get material stress response |
||
| 485 | const Vector &sigma = materialPointers[i]->getStressResultant(); |
||
| 486 | stresses(cnt) = sigma(0); |
||
| 487 | stresses(cnt+1) = sigma(1); |
||
| 488 | stresses(cnt+2) = sigma(2); |
||
| 489 | stresses(cnt+3) = sigma(3); |
||
| 490 | stresses(cnt+4) = sigma(4); |
||
| 491 | stresses(cnt+5) = sigma(5); |
||
| 492 | stresses(cnt+6) = sigma(6); |
||
| 493 | stresses(cnt+7) = sigma(7); |
||
| 494 | cnt += 8; |
||
| 495 | } |
||
| 496 | return eleInfo.setVector(stresses); |
||
| 497 | break; |
||
| 498 | case 3: //strain |
||
| 499 | for (int i = 0; i < 4; i++) { |
||
| 500 | |||
| 501 | // Get section deformation |
||
| 502 | const Vector &deformation = materialPointers[i]->getSectionDeformation(); |
||
| 503 | strains(cnt) = deformation(0); |
||
| 504 | strains(cnt+1) = deformation(1); |
||
| 505 | strains(cnt+2) = deformation(2); |
||
| 506 | strains(cnt+3) = deformation(3); |
||
| 507 | strains(cnt+4) = deformation(4); |
||
| 508 | strains(cnt+5) = deformation(5); |
||
| 509 | strains(cnt+6) = deformation(6); |
||
| 510 | strains(cnt+7) = deformation(7); |
||
| 511 | cnt += 8; |
||
| 512 | } |
||
| 513 | return eleInfo.setVector(strains); |
||
| 514 | break; |
||
| 515 | default: |
||
| 516 | return -1; |
||
| 517 | } |
||
| 518 | cnt=0; |
||
| 519 | } |
||
| 520 | |||
| 521 | |||
| 522 | //return stiffness matrix |
||
| 523 | const Matrix& ShellMITC4::getTangentStiff( ) |
||
| 524 | { |
||
| 525 | int tang_flag = 1 ; //get the tangent |
||
| 526 | |||
| 527 | //do tangent and residual here |
||
| 528 | formResidAndTangent( tang_flag ) ; |
||
| 529 | |||
| 530 | return stiff ; |
||
| 531 | } |
||
| 532 | |||
| 533 | //return secant matrix |
||
| 534 | const Matrix& ShellMITC4::getInitialStiff( ) |
||
| 535 | { |
||
| 536 | if (Ki != 0) |
||
| 537 | return *Ki; |
||
| 538 | |||
| 539 | static const int ndf = 6 ; //two membrane plus three bending plus one drill |
||
| 540 | |||
| 541 | static const int nstress = 8 ; //three membrane, three moment, two shear |
||
| 542 | |||
| 543 | static const int ngauss = 4 ; |
||
| 544 | |||
| 545 | static const int numnodes = 4 ; |
||
| 546 | |||
| 547 | int i, j, k, p, q ; |
||
| 548 | int jj, kk ; |
||
| 549 | |||
| 550 | double volume = 0.0 ; |
||
| 551 | |||
| 552 | static double xsj ; // determinant jacaobian matrix |
||
| 553 | |||
| 554 | static double dvol[ngauss] ; //volume element |
||
| 555 | |||
| 556 | static double shp[3][numnodes] ; //shape functions at a gauss point |
||
| 557 | |||
| 558 | // static double Shape[3][numnodes][ngauss] ; //all the shape functions |
||
| 559 | |||
| 560 | static Matrix stiffJK(ndf,ndf) ; //nodeJK stiffness |
||
| 561 | |||
| 562 | static Matrix dd(nstress,nstress) ; //material tangent |
||
| 563 | |||
| 564 | static Matrix J0(2,2) ; //Jacobian at center |
||
| 565 | |||
| 566 | static Matrix J0inv(2,2) ; //inverse of Jacobian at center |
||
| 567 | |||
| 568 | //---------B-matrices------------------------------------ |
||
| 569 | |||
| 570 | static Matrix BJ(nstress,ndf) ; // B matrix node J |
||
| 571 | |||
| 572 | static Matrix BJtran(ndf,nstress) ; |
||
| 573 | |||
| 574 | static Matrix BK(nstress,ndf) ; // B matrix node k |
||
| 575 | |||
| 576 | static Matrix BJtranD(ndf,nstress) ; |
||
| 577 | |||
| 578 | |||
| 579 | static Matrix Bbend(3,3) ; // bending B matrix |
||
| 580 | |||
| 581 | static Matrix Bshear(2,3) ; // shear B matrix |
||
| 582 | |||
| 583 | static Matrix Bmembrane(3,2) ; // membrane B matrix |
||
| 584 | |||
| 585 | |||
| 586 | static double BdrillJ[ndf] ; //drill B matrix |
||
| 587 | |||
| 588 | static double BdrillK[ndf] ; |
||
| 589 | |||
| 590 | double *drillPointer ; |
||
| 591 | |||
| 592 | static double saveB[nstress][ndf][numnodes] ; |
||
| 593 | |||
| 594 | //------------------------------------------------------- |
||
| 595 | |||
| 596 | stiff.Zero( ) ; |
||
| 597 | |||
| 598 | |||
| 599 | double dx34 = xl[0][2]-xl[0][3]; |
||
| 600 | double dy34 = xl[1][2]-xl[1][3]; |
||
| 601 | |||
| 602 | double dx21 = xl[0][1]-xl[0][0]; |
||
| 603 | double dy21 = xl[1][1]-xl[1][0]; |
||
| 604 | |||
| 605 | double dx32 = xl[0][2]-xl[0][1]; |
||
| 606 | double dy32 = xl[1][2]-xl[1][1]; |
||
| 607 | |||
| 608 | double dx41 = xl[0][3]-xl[0][0]; |
||
| 609 | double dy41 = xl[1][3]-xl[1][0]; |
||
| 610 | |||
| 611 | Matrix G(4,12); |
||
| 612 | G.Zero(); |
||
| 613 | double one_over_four = 0.25; |
||
| 614 | G(0,0)=-0.5; |
||
| 615 | G(0,1)=-dy41*one_over_four; |
||
| 616 | G(0,2)=dx41*one_over_four; |
||
| 617 | G(0,9)=0.5; |
||
| 618 | G(0,10)=-dy41*one_over_four; |
||
| 619 | G(0,11)=dx41*one_over_four; |
||
| 620 | G(1,0)=-0.5; |
||
| 621 | G(1,1)=-dy21*one_over_four; |
||
| 622 | G(1,2)=dx21*one_over_four; |
||
| 623 | G(1,3)=0.5; |
||
| 624 | G(1,4)=-dy21*one_over_four; |
||
| 625 | G(1,5)=dx21*one_over_four; |
||
| 626 | G(2,3)=-0.5; |
||
| 627 | G(2,4)=-dy32*one_over_four; |
||
| 628 | G(2,5)=dx32*one_over_four; |
||
| 629 | G(2,6)=0.5; |
||
| 630 | G(2,7)=-dy32*one_over_four; |
||
| 631 | G(2,8)=dx32*one_over_four; |
||
| 632 | G(3,6)=0.5; |
||
| 633 | G(3,7)=-dy34*one_over_four; |
||
| 634 | G(3,8)=dx34*one_over_four; |
||
| 635 | G(3,9)=-0.5; |
||
| 636 | G(3,10)=-dy34*one_over_four; |
||
| 637 | G(3,11)=dx34*one_over_four; |
||
| 638 | |||
| 639 | Matrix Ms(2,4); |
||
| 640 | Ms.Zero(); |
||
| 641 | Matrix Bsv(2,12); |
||
| 642 | Bsv.Zero(); |
||
| 643 | |||
| 644 | double Ax = -xl[0][0]+xl[0][1]+xl[0][2]-xl[0][3]; |
||
| 645 | double Bx = xl[0][0]-xl[0][1]+xl[0][2]-xl[0][3]; |
||
| 646 | double Cx = -xl[0][0]-xl[0][1]+xl[0][2]+xl[0][3]; |
||
| 647 | |||
| 648 | double Ay = -xl[1][0]+xl[1][1]+xl[1][2]-xl[1][3]; |
||
| 649 | double By = xl[1][0]-xl[1][1]+xl[1][2]-xl[1][3]; |
||
| 650 | double Cy = -xl[1][0]-xl[1][1]+xl[1][2]+xl[1][3]; |
||
| 651 | |||
| 652 | double alph = atan(Ay/Ax); |
||
| 653 | double beta = 3.141592653589793/2-atan(Cx/Cy); |
||
| 654 | Matrix Rot(2,2); |
||
| 655 | Rot.Zero(); |
||
| 656 | Rot(0,0)=sin(beta); |
||
| 657 | Rot(0,1)=-sin(alph); |
||
| 658 | Rot(1,0)=-cos(beta); |
||
| 659 | Rot(1,1)=cos(alph); |
||
| 660 | Matrix Bs(2,12); |
||
| 661 | |||
| 662 | double r1 = 0; |
||
| 663 | double r2 = 0; |
||
| 664 | double r3 = 0; |
||
| 665 | |||
| 666 | //gauss loop |
||
| 667 | for ( i = 0; i < ngauss; i++ ) { |
||
| 668 | |||
| 669 | r1 = Cx + sg[i]*Bx; |
||
| 670 | r3 = Cy + sg[i]*By; |
||
| 671 | r1 = r1*r1 + r3*r3; |
||
| 672 | r1 = sqrt (r1); |
||
| 673 | r2 = Ax + tg[i]*Bx; |
||
| 674 | r3 = Ay + tg[i]*By; |
||
| 675 | r2 = r2*r2 + r3*r3; |
||
| 676 | r2 = sqrt (r2); |
||
| 677 | |||
| 678 | //get shape functions |
||
| 679 | shape2d( sg[i], tg[i], xl, shp, xsj ) ; |
||
| 680 | //volume element to also be saved |
||
| 681 | dvol[i] = wg[i] * xsj ; |
||
| 682 | volume += dvol[i] ; |
||
| 683 | |||
| 684 | Ms(1,0)=1-sg[i]; |
||
| 685 | Ms(0,1)=1-tg[i]; |
||
| 686 | Ms(1,2)=1+sg[i]; |
||
| 687 | Ms(0,3)=1+tg[i]; |
||
| 688 | Bsv = Ms*G; |
||
| 689 | |||
| 690 | for ( j = 0; j < 12; j++ ) { |
||
| 691 | Bsv(0,j)=Bsv(0,j)*r1/(8*xsj); |
||
| 692 | Bsv(1,j)=Bsv(1,j)*r2/(8*xsj); |
||
| 693 | } |
||
| 694 | Bs=Rot*Bsv; |
||
| 695 | |||
| 696 | // j-node loop to compute strain |
||
| 697 | for ( j = 0; j < numnodes; j++ ) { |
||
| 698 | |||
| 699 | //compute B matrix |
||
| 700 | |||
| 701 | Bmembrane = computeBmembrane( j, shp ) ; |
||
| 702 | |||
| 703 | Bbend = computeBbend( j, shp ) ; |
||
| 704 | |||
| 705 | for ( p = 0; p < 3; p++) { |
||
| 706 | Bshear(0,p) = Bs(0,j*3+p); |
||
| 707 | Bshear(1,p) = Bs(1,j*3+p); |
||
| 708 | }//end for p |
||
| 709 | |||
| 710 | BJ = assembleB( Bmembrane, Bbend, Bshear ) ; |
||
| 711 | |||
| 712 | //save the B-matrix |
||
| 713 | for (p=0; p<nstress; p++) { |
||
| 714 | for (q=0; q<ndf; q++ ) |
||
| 715 | saveB[p][q][j] = BJ(p,q) ; |
||
| 716 | }//end for p |
||
| 717 | |||
| 718 | //drilling B matrix |
||
| 719 | drillPointer = computeBdrill( j, shp ) ; |
||
| 720 | for (p=0; p<ndf; p++ ) { |
||
| 721 | //BdrillJ[p] = *drillPointer++ ; |
||
| 722 | BdrillJ[p] = *drillPointer ; //set p-th component |
||
| 723 | drillPointer++ ; //pointer arithmetic |
||
| 724 | }//end for p |
||
| 725 | } // end for j |
||
| 726 | |||
| 727 | |||
| 728 | dd = materialPointers[i]->getInitialTangent( ) ; |
||
| 729 | dd *= dvol[i] ; |
||
| 730 | |||
| 731 | //residual and tangent calculations node loops |
||
| 732 | |||
| 733 | jj = 0 ; |
||
| 734 | for ( j = 0; j < numnodes; j++ ) { |
||
| 735 | |||
| 736 | //extract BJ |
||
| 737 | for (p=0; p<nstress; p++) { |
||
| 738 | for (q=0; q<ndf; q++ ) |
||
| 739 | BJ(p,q) = saveB[p][q][j] ; |
||
| 740 | }//end for p |
||
| 741 | |||
| 742 | //multiply bending terms by (-1.0) for correct statement |
||
| 743 | // of equilibrium |
||
| 744 | for ( p = 3; p < 6; p++ ) { |
||
| 745 | for ( q = 3; q < 6; q++ ) |
||
| 746 | BJ(p,q) *= (-1.0) ; |
||
| 747 | } //end for p |
||
| 748 | |||
| 749 | |||
| 750 | //transpose |
||
| 751 | //BJtran = transpose( 8, ndf, BJ ) ; |
||
| 752 | for (p=0; p<ndf; p++) { |
||
| 753 | for (q=0; q<nstress; q++) |
||
| 754 | BJtran(p,q) = BJ(q,p) ; |
||
| 755 | }//end for p |
||
| 756 | |||
| 757 | //drilling B matrix |
||
| 758 | drillPointer = computeBdrill( j, shp ) ; |
||
| 759 | for (p=0; p<ndf; p++ ) { |
||
| 760 | BdrillJ[p] = *drillPointer ; |
||
| 761 | drillPointer++ ; |
||
| 762 | }//end for p |
||
| 763 | |||
| 764 | //BJtranD = BJtran * dd ; |
||
| 765 | BJtranD.addMatrixProduct(0.0, BJtran,dd,1.0 ) ; |
||
| 766 | |||
| 767 | for (p=0; p<ndf; p++) |
||
| 768 | BdrillJ[p] *= ( Ktt*dvol[i] ) ; |
||
| 769 | |||
| 770 | kk = 0 ; |
||
| 771 | for ( k = 0; k < numnodes; k++ ) { |
||
| 772 | |||
| 773 | //extract BK |
||
| 774 | for (p=0; p<nstress; p++) { |
||
| 775 | for (q=0; q<ndf; q++ ) |
||
| 776 | BK(p,q) = saveB[p][q][k] ; |
||
| 777 | }//end for p |
||
| 778 | |||
| 779 | |||
| 780 | //drilling B matrix |
||
| 781 | drillPointer = computeBdrill( k, shp ) ; |
||
| 782 | for (p=0; p<ndf; p++ ) { |
||
| 783 | BdrillK[p] = *drillPointer ; |
||
| 784 | drillPointer++ ; |
||
| 785 | }//end for p |
||
| 786 | |||
| 787 | //stiffJK = BJtranD * BK ; |
||
| 788 | // + transpose( 1,ndf,BdrillJ ) * BdrillK ; |
||
| 789 | stiffJK.addMatrixProduct(0.0, BJtranD,BK,1.0 ) ; |
||
| 790 | |||
| 791 | for ( p = 0; p < ndf; p++ ) { |
||
| 792 | for ( q = 0; q < ndf; q++ ) { |
||
| 793 | stiff( jj+p, kk+q ) += stiffJK(p,q) |
||
| 794 | + ( BdrillJ[p]*BdrillK[q] ) ; |
||
| 795 | }//end for q |
||
| 796 | }//end for p |
||
| 797 | |||
| 798 | kk += ndf ; |
||
| 799 | } // end for k loop |
||
| 800 | |||
| 801 | jj += ndf ; |
||
| 802 | } // end for j loop |
||
| 803 | |||
| 804 | } //end for i gauss loop |
||
| 805 | |||
| 806 | Ki = new Matrix(stiff); |
||
| 807 | |||
| 808 | return stiff ; |
||
| 809 | } |
||
| 810 | |||
| 811 | |||
| 812 | //return mass matrix |
||
| 813 | const Matrix& ShellMITC4::getMass( ) |
||
| 814 | { |
||
| 815 | int tangFlag = 1 ; |
||
| 816 | |||
| 817 | formInertiaTerms( tangFlag ) ; |
||
| 818 | |||
| 819 | return mass ; |
||
| 820 | } |
||
| 821 | |||
| 822 | |||
| 823 | void ShellMITC4::zeroLoad( ) |
||
| 824 | { |
||
| 825 | if (load != 0) |
||
| 826 | load->Zero(); |
||
| 827 | |||
| 828 | return ; |
||
| 829 | } |
||
| 830 | |||
| 831 | |||
| 832 | int |
||
| 833 | ShellMITC4::addLoad(ElementalLoad *theLoad, double loadFactor) |
||
| 834 | { |
||
| 835 | opserr << "ShellMITC4::addLoad - load type unknown for ele with tag: " << this->getTag() << endln; |
||
| 836 | return -1; |
||
| 837 | } |
||
| 838 | |||
| 839 | |||
| 840 | |||
| 841 | int |
||
| 842 | ShellMITC4::addInertiaLoadToUnbalance(const Vector &accel) |
||
| 843 | { |
||
| 844 | int tangFlag = 1 ; |
||
| 845 | |||
| 846 | int i; |
||
| 847 | |||
| 848 | int allRhoZero = 0; |
||
| 849 | for (i=0; i<4; i++) { |
||
| 850 | if (materialPointers[i]->getRho() != 0.0) |
||
| 851 | allRhoZero = 1; |
||
| 852 | } |
||
| 853 | |||
| 854 | if (allRhoZero == 0) |
||
| 855 | return 0; |
||
| 856 | |||
| 857 | int count = 0; |
||
| 858 | for (i=0; i<4; i++) { |
||
| 859 | const Vector &Raccel = nodePointers[i]->getRV(accel); |
||
| 860 | for (int j=0; j<6; j++) |
||
| 861 | resid(count++) = Raccel(i); |
||
| 862 | } |
||
| 863 | |||
| 864 | formInertiaTerms( tangFlag ) ; |
||
| 865 | if (load == 0) |
||
| 866 | load = new Vector(24); |
||
| 867 | load->addMatrixVector(1.0, mass, resid, -1.0); |
||
| 868 | |||
| 869 | return 0; |
||
| 870 | } |
||
| 871 | |||
| 872 | |||
| 873 | |||
| 874 | //get residual |
||
| 875 | const Vector& ShellMITC4::getResistingForce( ) |
||
| 876 | { |
||
| 877 | int tang_flag = 0 ; //don't get the tangent |
||
| 878 | |||
| 879 | formResidAndTangent( tang_flag ) ; |
||
| 880 | |||
| 881 | // subtract external loads |
||
| 882 | if (load != 0) |
||
| 883 | resid -= *load; |
||
| 884 | |||
| 885 | return resid ; |
||
| 886 | } |
||
| 887 | |||
| 888 | |||
| 889 | //get residual with inertia terms |
||
| 890 | const Vector& ShellMITC4::getResistingForceIncInertia( ) |
||
| 891 | { |
||
| 892 | static Vector res(24); |
||
| 893 | int tang_flag = 0 ; //don't get the tangent |
||
| 894 | |||
| 895 | //do tangent and residual here |
||
| 896 | formResidAndTangent( tang_flag ) ; |
||
| 897 | |||
| 898 | formInertiaTerms( tang_flag ) ; |
||
| 899 | |||
| 900 | res = resid; |
||
| 901 | // add the damping forces if rayleigh damping |
||
| 902 | if (alphaM != 0.0 || betaK != 0.0 || betaK0 != 0.0 || betaKc != 0.0) |
||
| 903 | res += this->getRayleighDampingForces(); |
||
| 904 | |||
| 905 | // subtract external loads |
||
| 906 | if (load != 0) |
||
| 907 | res -= *load; |
||
| 908 | |||
| 909 | return res; |
||
| 910 | } |
||
| 911 | |||
| 912 | //********************************************************************* |
||
| 913 | //form inertia terms |
||
| 914 | |||
| 915 | void |
||
| 916 | ShellMITC4::formInertiaTerms( int tangFlag ) |
||
| 917 | { |
||
| 918 | |||
| 919 | //translational mass only |
||
| 920 | //rotational inertia terms are neglected |
||
| 921 | |||
| 922 | |||
| 923 | static const int ndf = 6 ; |
||
| 924 | |||
| 925 | static const int numberNodes = 4 ; |
||
| 926 | |||
| 927 | static const int numberGauss = 4 ; |
||
| 928 | |||
| 929 | static const int nShape = 3 ; |
||
| 930 | |||
| 931 | static const int massIndex = nShape - 1 ; |
||
| 932 | |||
| 933 | double xsj ; // determinant jacaobian matrix |
||
| 934 | |||
| 935 | double dvol ; //volume element |
||
| 936 | |||
| 937 | static double shp[nShape][numberNodes] ; //shape functions at a gauss point |
||
| 938 | |||
| 939 | static Vector momentum(ndf) ; |
||
| 940 | |||
| 941 | |||
| 942 | int i, j, k, p; |
||
| 943 | int jj, kk ; |
||
| 944 | |||
| 945 | double temp, rhoH, massJK ; |
||
| 946 | |||
| 947 | |||
| 948 | //zero mass |
||
| 949 | mass.Zero( ) ; |
||
| 950 | |||
| 951 | //gauss loop |
||
| 952 | for ( i = 0; i < numberGauss; i++ ) { |
||
| 953 | |||
| 954 | //get shape functions |
||
| 955 | shape2d( sg[i], tg[i], xl, shp, xsj ) ; |
||
| 956 | |||
| 957 | //volume element to also be saved |
||
| 958 | dvol = wg[i] * xsj ; |
||
| 959 | |||
| 960 | |||
| 961 | //node loop to compute accelerations |
||
| 962 | momentum.Zero( ) ; |
||
| 963 | for ( j = 0; j < numberNodes; j++ ) |
||
| 964 | //momentum += ( shp[massIndex][j] * nodePointers[j]->getTrialAccel() ) ; |
||
| 965 | momentum.addVector(1.0, |
||
| 966 | nodePointers[j]->getTrialAccel(), |
||
| 967 | shp[massIndex][j] ) ; |
||
| 968 | |||
| 969 | |||
| 970 | //density |
||
| 971 | rhoH = materialPointers[i]->getRho() ; |
||
| 972 | |||
| 973 | //multiply acceleration by density to form momentum |
||
| 974 | momentum *= rhoH ; |
||
| 975 | |||
| 976 | |||
| 977 | //residual and tangent calculations node loops |
||
| 978 | for ( j=0, jj=0; j<numberNodes; j++, jj+=ndf ) { |
||
| 979 | |||
| 980 | temp = shp[massIndex][j] * dvol ; |
||
| 981 | |||
| 982 | for ( p = 0; p < 3; p++ ) |
||
| 983 | resid( jj+p ) += ( temp * momentum(p) ) ; |
||
| 984 | |||
| 985 | |||
| 986 | if ( tangFlag == 1 && rhoH != 0.0) { |
||
| 987 | |||
| 988 | //multiply by density |
||
| 989 | temp *= rhoH ; |
||
| 990 | |||
| 991 | //node-node translational mass |
||
| 992 | for ( k=0, kk=0; k<numberNodes; k++, kk+=ndf ) { |
||
| 993 | |||
| 994 | massJK = temp * shp[massIndex][k] ; |
||
| 995 | |||
| 996 | for ( p = 0; p < 3; p++ ) |
||
| 997 | mass( jj+p, kk+p ) += massJK ; |
||
| 998 | |||
| 999 | } // end for k loop |
||
| 1000 | |||
| 1001 | } // end if tang_flag |
||
| 1002 | |||
| 1003 | } // end for j loop |
||
| 1004 | |||
| 1005 | } //end for i gauss loop |
||
| 1006 | |||
| 1007 | } |
||
| 1008 | |||
| 1009 | //********************************************************************* |
||
| 1010 | |||
| 1011 | //form residual and tangent |
||
| 1012 | void |
||
| 1013 | ShellMITC4::formResidAndTangent( int tang_flag ) |
||
| 1014 | { |
||
| 1015 | // |
||
| 1016 | // six(6) nodal dof's ordered : |
||
| 1017 | // |
||
| 1018 | // - - |
||
| 1019 | // | u1 | <---plate membrane |
||
| 1020 | // | u2 | |
||
| 1021 | // |----------| |
||
| 1022 | // | w = u3 | <---plate bending |
||
| 1023 | // | theta1 | |
||
| 1024 | // | theta2 | |
||
| 1025 | // |----------| |
||
| 1026 | // | theta3 | <---drill |
||
| 1027 | // - - |
||
| 1028 | // |
||
| 1029 | // membrane strains ordered : |
||
| 1030 | // |
||
| 1031 | // strain(0) = eps00 i.e. (11)-strain |
||
| 1032 | // strain(1) = eps11 i.e. (22)-strain |
||
| 1033 | // strain(2) = gamma01 i.e. (12)-shear |
||
| 1034 | // |
||
| 1035 | // curvatures and shear strains ordered : |
||
| 1036 | // |
||
| 1037 | // strain(3) = kappa00 i.e. (11)-curvature |
||
| 1038 | // strain(4) = kappa11 i.e. (22)-curvature |
||
| 1039 | // strain(5) = 2*kappa01 i.e. 2*(12)-curvature |
||
| 1040 | // |
||
| 1041 | // strain(6) = gamma02 i.e. (13)-shear |
||
| 1042 | // strain(7) = gamma12 i.e. (23)-shear |
||
| 1043 | // |
||
| 1044 | // same ordering for moments/shears but no 2 |
||
| 1045 | // |
||
| 1046 | // Then, |
||
| 1047 | // epsilon00 = -z * kappa00 + eps00_membrane |
||
| 1048 | // epsilon11 = -z * kappa11 + eps11_membrane |
||
| 1049 | // gamma01 = 2*epsilon01 = -z * (2*kappa01) + gamma01_membrane |
||
| 1050 | // |
||
| 1051 | // Shear strains gamma02, gamma12 constant through cross section |
||
| 1052 | // |
||
| 1053 | |||
| 1054 | static const int ndf = 6 ; //two membrane plus three bending plus one drill |
||
| 1055 | |||
| 1056 | static const int nstress = 8 ; //three membrane, three moment, two shear |
||
| 1057 | |||
| 1058 | static const int ngauss = 4 ; |
||
| 1059 | |||
| 1060 | static const int numnodes = 4 ; |
||
| 1061 | |||
| 1062 | int i, j, k, p, q ; |
||
| 1063 | int jj, kk ; |
||
| 1064 | |||
| 1065 | int success ; |
||
| 1066 | |||
| 1067 | double volume = 0.0 ; |
||
| 1068 | |||
| 1069 | static double xsj ; // determinant jacaobian matrix |
||
| 1070 | |||
| 1071 | static double dvol[ngauss] ; //volume element |
||
| 1072 | |||
| 1073 | static Vector strain(nstress) ; //strain |
||
| 1074 | |||
| 1075 | static double shp[3][numnodes] ; //shape functions at a gauss point |
||
| 1076 | |||
| 1077 | // static double Shape[3][numnodes][ngauss] ; //all the shape functions |
||
| 1078 | |||
| 1079 | static Vector residJ(ndf) ; //nodeJ residual |
||
| 1080 | |||
| 1081 | static Matrix stiffJK(ndf,ndf) ; //nodeJK stiffness |
||
| 1082 | |||
| 1083 | static Vector stress(nstress) ; //stress resultants |
||
| 1084 | |||
| 1085 | static Matrix dd(nstress,nstress) ; //material tangent |
||
| 1086 | |||
| 1087 | static Matrix J0(2,2) ; //Jacobian at center |
||
| 1088 | |||
| 1089 | static Matrix J0inv(2,2) ; //inverse of Jacobian at center |
||
| 1090 | |||
| 1091 | double epsDrill = 0.0 ; //drilling "strain" |
||
| 1092 | |||
| 1093 | double tauDrill = 0.0 ; //drilling "stress" |
||
| 1094 | |||
| 1095 | //---------B-matrices------------------------------------ |
||
| 1096 | |||
| 1097 | static Matrix BJ(nstress,ndf) ; // B matrix node J |
||
| 1098 | |||
| 1099 | static Matrix BJtran(ndf,nstress) ; |
||
| 1100 | |||
| 1101 | static Matrix BK(nstress,ndf) ; // B matrix node k |
||
| 1102 | |||
| 1103 | static Matrix BJtranD(ndf,nstress) ; |
||
| 1104 | |||
| 1105 | |||
| 1106 | static Matrix Bbend(3,3) ; // bending B matrix |
||
| 1107 | |||
| 1108 | static Matrix Bshear(2,3) ; // shear B matrix |
||
| 1109 | |||
| 1110 | static Matrix Bmembrane(3,2) ; // membrane B matrix |
||
| 1111 | |||
| 1112 | |||
| 1113 | static double BdrillJ[ndf] ; //drill B matrix |
||
| 1114 | |||
| 1115 | static double BdrillK[ndf] ; |
||
| 1116 | |||
| 1117 | double *drillPointer ; |
||
| 1118 | |||
| 1119 | static double saveB[nstress][ndf][numnodes] ; |
||
| 1120 | |||
| 1121 | //------------------------------------------------------- |
||
| 1122 | |||
| 1123 | //zero stiffness and residual |
||
| 1124 | stiff.Zero( ) ; |
||
| 1125 | resid.Zero( ) ; |
||
| 1126 | |||
| 1127 | double dx34 = xl[0][2]-xl[0][3]; |
||
| 1128 | double dy34 = xl[1][2]-xl[1][3]; |
||
| 1129 | |||
| 1130 | double dx21 = xl[0][1]-xl[0][0]; |
||
| 1131 | double dy21 = xl[1][1]-xl[1][0]; |
||
| 1132 | |||
| 1133 | double dx32 = xl[0][2]-xl[0][1]; |
||
| 1134 | double dy32 = xl[1][2]-xl[1][1]; |
||
| 1135 | |||
| 1136 | double dx41 = xl[0][3]-xl[0][0]; |
||
| 1137 | double dy41 = xl[1][3]-xl[1][0]; |
||
| 1138 | |||
| 1139 | Matrix G(4,12); |
||
| 1140 | G.Zero(); |
||
| 1141 | double one_over_four = 0.25; |
||
| 1142 | G(0,0)=-0.5; |
||
| 1143 | G(0,1)=-dy41*one_over_four; |
||
| 1144 | G(0,2)=dx41*one_over_four; |
||
| 1145 | G(0,9)=0.5; |
||
| 1146 | G(0,10)=-dy41*one_over_four; |
||
| 1147 | G(0,11)=dx41*one_over_four; |
||
| 1148 | G(1,0)=-0.5; |
||
| 1149 | G(1,1)=-dy21*one_over_four; |
||
| 1150 | G(1,2)=dx21*one_over_four; |
||
| 1151 | G(1,3)=0.5; |
||
| 1152 | G(1,4)=-dy21*one_over_four; |
||
| 1153 | G(1,5)=dx21*one_over_four; |
||
| 1154 | G(2,3)=-0.5; |
||
| 1155 | G(2,4)=-dy32*one_over_four; |
||
| 1156 | G(2,5)=dx32*one_over_four; |
||
| 1157 | G(2,6)=0.5; |
||
| 1158 | G(2,7)=-dy32*one_over_four; |
||
| 1159 | G(2,8)=dx32*one_over_four; |
||
| 1160 | G(3,6)=0.5; |
||
| 1161 | G(3,7)=-dy34*one_over_four; |
||
| 1162 | G(3,8)=dx34*one_over_four; |
||
| 1163 | G(3,9)=-0.5; |
||
| 1164 | G(3,10)=-dy34*one_over_four; |
||
| 1165 | G(3,11)=dx34*one_over_four; |
||
| 1166 | |||
| 1167 | Matrix Ms(2,4); |
||
| 1168 | Ms.Zero(); |
||
| 1169 | Matrix Bsv(2,12); |
||
| 1170 | Bsv.Zero(); |
||
| 1171 | |||
| 1172 | double Ax = -xl[0][0]+xl[0][1]+xl[0][2]-xl[0][3]; |
||
| 1173 | double Bx = xl[0][0]-xl[0][1]+xl[0][2]-xl[0][3]; |
||
| 1174 | double Cx = -xl[0][0]-xl[0][1]+xl[0][2]+xl[0][3]; |
||
| 1175 | |||
| 1176 | double Ay = -xl[1][0]+xl[1][1]+xl[1][2]-xl[1][3]; |
||
| 1177 | double By = xl[1][0]-xl[1][1]+xl[1][2]-xl[1][3]; |
||
| 1178 | double Cy = -xl[1][0]-xl[1][1]+xl[1][2]+xl[1][3]; |
||
| 1179 | |||
| 1180 | double alph = atan(Ay/Ax); |
||
| 1181 | double beta = 3.141592653589793/2-atan(Cx/Cy); |
||
| 1182 | Matrix Rot(2,2); |
||
| 1183 | Rot.Zero(); |
||
| 1184 | Rot(0,0)=sin(beta); |
||
| 1185 | Rot(0,1)=-sin(alph); |
||
| 1186 | Rot(1,0)=-cos(beta); |
||
| 1187 | Rot(1,1)=cos(alph); |
||
| 1188 | Matrix Bs(2,12); |
||
| 1189 | |||
| 1190 | double r1 = 0; |
||
| 1191 | double r2 = 0; |
||
| 1192 | double r3 = 0; |
||
| 1193 | |||
| 1194 | //gauss loop |
||
| 1195 | for ( i = 0; i < ngauss; i++ ) { |
||
| 1196 | |||
| 1197 | r1 = Cx + sg[i]*Bx; |
||
| 1198 | r3 = Cy + sg[i]*By; |
||
| 1199 | r1 = r1*r1 + r3*r3; |
||
| 1200 | r1 = sqrt (r1); |
||
| 1201 | r2 = Ax + tg[i]*Bx; |
||
| 1202 | r3 = Ay + tg[i]*By; |
||
| 1203 | r2 = r2*r2 + r3*r3; |
||
| 1204 | r2 = sqrt (r2); |
||
| 1205 | |||
| 1206 | //get shape functions |
||
| 1207 | shape2d( sg[i], tg[i], xl, shp, xsj ) ; |
||
| 1208 | //volume element to also be saved |
||
| 1209 | dvol[i] = wg[i] * xsj ; |
||
| 1210 | volume += dvol[i] ; |
||
| 1211 | |||
| 1212 | Ms(1,0)=1-sg[i]; |
||
| 1213 | Ms(0,1)=1-tg[i]; |
||
| 1214 | Ms(1,2)=1+sg[i]; |
||
| 1215 | Ms(0,3)=1+tg[i]; |
||
| 1216 | Bsv = Ms*G; |
||
| 1217 | |||
| 1218 | for ( j = 0; j < 12; j++ ) { |
||
| 1219 | Bsv(0,j)=Bsv(0,j)*r1/(8*xsj); |
||
| 1220 | Bsv(1,j)=Bsv(1,j)*r2/(8*xsj); |
||
| 1221 | } |
||
| 1222 | Bs=Rot*Bsv; |
||
| 1223 | |||
| 1224 | //zero the strains |
||
| 1225 | strain.Zero( ) ; |
||
| 1226 | epsDrill = 0.0 ; |
||
| 1227 | |||
| 1228 | |||
| 1229 | // j-node loop to compute strain |
||
| 1230 | for ( j = 0; j < numnodes; j++ ) { |
||
| 1231 | |||
| 1232 | //compute B matrix |
||
| 1233 | |||
| 1234 | Bmembrane = computeBmembrane( j, shp ) ; |
||
| 1235 | |||
| 1236 | Bbend = computeBbend( j, shp ) ; |
||
| 1237 | |||
| 1238 | for ( p = 0; p < 3; p++) { |
||
| 1239 | Bshear(0,p) = Bs(0,j*3+p); |
||
| 1240 | Bshear(1,p) = Bs(1,j*3+p); |
||
| 1241 | }//end for p |
||
| 1242 | |||
| 1243 | BJ = assembleB( Bmembrane, Bbend, Bshear ) ; |
||
| 1244 | |||
| 1245 | //save the B-matrix |
||
| 1246 | for (p=0; p<nstress; p++) { |
||
| 1247 | for (q=0; q<ndf; q++ ) |
||
| 1248 | saveB[p][q][j] = BJ(p,q) ; |
||
| 1249 | }//end for p |
||
| 1250 | |||
| 1251 | |||
| 1252 | //nodal "displacements" |
||
| 1253 | const Vector &ul = nodePointers[j]->getTrialDisp( ) ; |
||
| 1254 | |||
| 1255 | //compute the strain |
||
| 1256 | //strain += (BJ*ul) ; |
||
| 1257 | strain.addMatrixVector(1.0, BJ,ul,1.0 ) ; |
||
| 1258 | |||
| 1259 | //drilling B matrix |
||
| 1260 | drillPointer = computeBdrill( j, shp ) ; |
||
| 1261 | for (p=0; p<ndf; p++ ) { |
||
| 1262 | BdrillJ[p] = *drillPointer ; //set p-th component |
||
| 1263 | drillPointer++ ; //pointer arithmetic |
||
| 1264 | }//end for p |
||
| 1265 | |||
| 1266 | //drilling "strain" |
||
| 1267 | for ( p = 0; p < ndf; p++ ) |
||
| 1268 | epsDrill += BdrillJ[p]*ul(p) ; |
||
| 1269 | } // end for j |
||
| 1270 | |||
| 1271 | |||
| 1272 | //send the strain to the material |
||
| 1273 | success = materialPointers[i]->setTrialSectionDeformation( strain ) ; |
||
| 1274 | |||
| 1275 | //compute the stress |
||
| 1276 | stress = materialPointers[i]->getStressResultant( ) ; |
||
| 1277 | |||
| 1278 | //drilling "stress" |
||
| 1279 | tauDrill = Ktt * epsDrill ; |
||
| 1280 | |||
| 1281 | //multiply by volume element |
||
| 1282 | stress *= dvol[i] ; |
||
| 1283 | tauDrill *= dvol[i] ; |
||
| 1284 | |||
| 1285 | if ( tang_flag == 1 ) { |
||
| 1286 | dd = materialPointers[i]->getSectionTangent( ) ; |
||
| 5141 | fmk | 1287 | |
| 4617 | fmk | 1288 | dd *= dvol[i] ; |
| 1289 | } //end if tang_flag |
||
| 1290 | |||
| 1291 | |||
| 1292 | //residual and tangent calculations node loops |
||
| 1293 | |||
| 1294 | jj = 0 ; |
||
| 1295 | for ( j = 0; j < numnodes; j++ ) { |
||
| 1296 | |||
| 1297 | //extract BJ |
||
| 1298 | for (p=0; p<nstress; p++) { |
||
| 1299 | for (q=0; q<ndf; q++ ) |
||
| 1300 | BJ(p,q) = saveB[p][q][j] ; |
||
| 1301 | }//end for p |
||
| 1302 | |||
| 1303 | //multiply bending terms by (-1.0) for correct statement |
||
| 1304 | // of equilibrium |
||
| 1305 | for ( p = 3; p < 6; p++ ) { |
||
| 1306 | for ( q = 3; q < 6; q++ ) |
||
| 1307 | BJ(p,q) *= (-1.0) ; |
||
| 1308 | } //end for p |
||
| 1309 | |||
| 1310 | //transpose |
||
| 1311 | for (p=0; p<ndf; p++) { |
||
| 1312 | for (q=0; q<nstress; q++) |
||
| 1313 | BJtran(p,q) = BJ(q,p) ; |
||
| 1314 | }//end for p |
||
| 1315 | |||
| 1316 | residJ.addMatrixVector(0.0, BJtran,stress,1.0 ) ; |
||
| 1317 | |||
| 1318 | //drilling B matrix |
||
| 1319 | drillPointer = computeBdrill( j, shp ) ; |
||
| 1320 | for (p=0; p<ndf; p++ ) { |
||
| 1321 | BdrillJ[p] = *drillPointer ; |
||
| 1322 | drillPointer++ ; |
||
| 1323 | }//end for p |
||
| 1324 | |||
| 1325 | //residual including drill |
||
| 1326 | for ( p = 0; p < ndf; p++ ) |
||
| 1327 | resid( jj + p ) += ( residJ(p) + BdrillJ[p]*tauDrill ) ; |
||
| 1328 | |||
| 1329 | if ( tang_flag == 1 ) { |
||
| 1330 | |||
| 1331 | BJtranD.addMatrixProduct(0.0, BJtran,dd,1.0 ) ; |
||
| 1332 | |||
| 1333 | for (p=0; p<ndf; p++) |
||
| 1334 | BdrillJ[p] *= ( Ktt*dvol[i] ) ; |
||
| 1335 | |||
| 1336 | kk = 0 ; |
||
| 1337 | for ( k = 0; k < numnodes; k++ ) { |
||
| 1338 | |||
| 1339 | //extract BK |
||
| 1340 | for (p=0; p<nstress; p++) { |
||
| 1341 | for (q=0; q<ndf; q++ ) |
||
| 1342 | BK(p,q) = saveB[p][q][k] ; |
||
| 1343 | }//end for p |
||
| 1344 | |||
| 1345 | //drilling B matrix |
||
| 1346 | drillPointer = computeBdrill( k, shp ) ; |
||
| 1347 | for (p=0; p<ndf; p++ ) { |
||
| 1348 | BdrillK[p] = *drillPointer ; |
||
| 1349 | drillPointer++ ; |
||
| 1350 | }//end for p |
||
| 1351 | |||
| 1352 | //stiffJK = BJtranD * BK ; |
||
| 1353 | // + transpose( 1,ndf,BdrillJ ) * BdrillK ; |
||
| 1354 | stiffJK.addMatrixProduct(0.0, BJtranD,BK,1.0 ) ; |
||
| 1355 | |||
| 1356 | for ( p = 0; p < ndf; p++ ) { |
||
| 1357 | for ( q = 0; q < ndf; q++ ) { |
||
| 1358 | stiff( jj+p, kk+q ) += stiffJK(p,q) |
||
| 1359 | + ( BdrillJ[p]*BdrillK[q] ) ; |
||
| 1360 | }//end for q |
||
| 1361 | }//end for p |
||
| 1362 | |||
| 1363 | kk += ndf ; |
||
| 1364 | } // end for k loop |
||
| 1365 | |||
| 1366 | } // end if tang_flag |
||
| 1367 | |||
| 1368 | jj += ndf ; |
||
| 1369 | } // end for j loop |
||
| 1370 | |||
| 1371 | } //end for i gauss loop |
||
| 1372 | |||
| 1373 | return ; |
||
| 1374 | } |
||
| 1375 | |||
| 1376 | |||
| 1377 | //************************************************************************ |
||
| 1378 | //compute local coordinates and basis |
||
| 1379 | |||
| 1380 | void |
||
| 1381 | ShellMITC4::computeBasis( ) |
||
| 1382 | { |
||
| 1383 | //could compute derivatives \frac{ \partial {\bf x} }{ \partial L_1 } |
||
| 1384 | // and \frac{ \partial {\bf x} }{ \partial L_2 } |
||
| 1385 | //and use those as basis vectors but this is easier |
||
| 1386 | //and the shell is flat anyway. |
||
| 1387 | |||
| 1388 | static Vector temp(3) ; |
||
| 1389 | |||
| 1390 | static Vector v1(3) ; |
||
| 1391 | static Vector v2(3) ; |
||
| 1392 | static Vector v3(3) ; |
||
| 1393 | |||
| 1394 | //get two vectors (v1, v2) in plane of shell by |
||
| 1395 | // nodal coordinate differences |
||
| 1396 | |||
| 1397 | const Vector &coor0 = nodePointers[0]->getCrds( ) ; |
||
| 1398 | |||
| 1399 | const Vector &coor1 = nodePointers[1]->getCrds( ) ; |
||
| 1400 | |||
| 1401 | const Vector &coor2 = nodePointers[2]->getCrds( ) ; |
||
| 1402 | |||
| 1403 | const Vector &coor3 = nodePointers[3]->getCrds( ) ; |
||
| 1404 | |||
| 1405 | v1.Zero( ) ; |
||
| 1406 | //v1 = 0.5 * ( coor2 + coor1 - coor3 - coor0 ) ; |
||
| 1407 | v1 = coor2 ; |
||
| 1408 | v1 += coor1 ; |
||
| 1409 | v1 -= coor3 ; |
||
| 1410 | v1 -= coor0 ; |
||
| 1411 | v1 *= 0.50 ; |
||
| 1412 | |||
| 1413 | v2.Zero( ) ; |
||
| 1414 | //v2 = 0.5 * ( coor3 + coor2 - coor1 - coor0 ) ; |
||
| 1415 | v2 = coor3 ; |
||
| 1416 | v2 += coor2 ; |
||
| 1417 | v2 -= coor1 ; |
||
| 1418 | v2 -= coor0 ; |
||
| 1419 | v2 *= 0.50 ; |
||
| 1420 | |||
| 1421 | //normalize v1 |
||
| 1422 | //double length = LovelyNorm( v1 ) ; |
||
| 1423 | double length = v1.Norm( ) ; |
||
| 1424 | v1 /= length ; |
||
| 1425 | |||
| 1426 | //Gram-Schmidt process for v2 |
||
| 1427 | |||
| 1428 | //double alpha = LovelyInnerProduct( v2, v1 ) ; |
||
| 1429 | double alpha = v2^v1 ; |
||
| 1430 | |||
| 1431 | //v2 -= alpha*v1 ; |
||
| 1432 | temp = v1 ; |
||
| 1433 | temp *= alpha ; |
||
| 1434 | v2 -= temp ; |
||
| 1435 | |||
| 1436 | //normalize v2 |
||
| 1437 | //length = LovelyNorm( v2 ) ; |
||
| 1438 | length = v2.Norm( ) ; |
||
| 1439 | v2 /= length ; |
||
| 1440 | |||
| 1441 | //cross product for v3 |
||
| 1442 | v3 = LovelyCrossProduct( v1, v2 ) ; |
||
| 1443 | |||
| 1444 | //local nodal coordinates in plane of shell |
||
| 1445 | |||
| 1446 | int i ; |
||
| 1447 | for ( i = 0; i < 4; i++ ) { |
||
| 1448 | |||
| 1449 | const Vector &coorI = nodePointers[i]->getCrds( ) ; |
||
| 1450 | xl[0][i] = coorI^v1 ; |
||
| 1451 | xl[1][i] = coorI^v2 ; |
||
| 1452 | |||
| 1453 | } //end for i |
||
| 1454 | |||
| 1455 | //basis vectors stored as array of doubles |
||
| 1456 | for ( i = 0; i < 3; i++ ) { |
||
| 1457 | g1[i] = v1(i) ; |
||
| 1458 | g2[i] = v2(i) ; |
||
| 1459 | g3[i] = v3(i) ; |
||
| 1460 | } //end for i |
||
| 1461 | |||
| 1462 | } |
||
| 1463 | |||
| 1464 | //************************************************************************* |
||
| 1465 | //compute Bdrill |
||
| 1466 | |||
| 1467 | double* |
||
| 1468 | ShellMITC4::computeBdrill( int node, const double shp[3][4] ) |
||
| 1469 | { |
||
| 1470 | |||
| 1471 | //static Matrix Bdrill(1,6) ; |
||
| 1472 | static double Bdrill[6] ; |
||
| 1473 | static double B1 ; |
||
| 1474 | static double B2 ; |
||
| 1475 | static double B6 ; |
||
| 1476 | |||
| 1477 | |||
| 1478 | //---Bdrill Matrix in standard {1,2,3} mechanics notation--------- |
||
| 1479 | // |
||
| 1480 | // - - |
||
| 1481 | // Bdrill = | -0.5*N,2 +0.5*N,1 0 0 0 -N | (1x6) |
||
| 1482 | // - - |
||
| 1483 | // |
||
| 1484 | //---------------------------------------------------------------- |
||
| 1485 | |||
| 1486 | // Bdrill.Zero( ) ; |
||
| 1487 | |||
| 1488 | //Bdrill(0,0) = -0.5*shp[1][node] ; |
||
| 1489 | |||
| 1490 | //Bdrill(0,1) = +0.5*shp[0][node] ; |
||
| 1491 | |||
| 1492 | //Bdrill(0,5) = -shp[2][node] ; |
||
| 1493 | |||
| 1494 | |||
| 1495 | B1 = -0.5*shp[1][node] ; |
||
| 1496 | |||
| 1497 | B2 = +0.5*shp[0][node] ; |
||
| 1498 | |||
| 1499 | B6 = -shp[2][node] ; |
||
| 1500 | |||
| 1501 | Bdrill[0] = B1*g1[0] + B2*g2[0] ; |
||
| 1502 | Bdrill[1] = B1*g1[1] + B2*g2[1] ; |
||
| 1503 | Bdrill[2] = B1*g1[2] + B2*g2[2] ; |
||
| 1504 | |||
| 1505 | Bdrill[3] = B6*g3[0] ; |
||
| 1506 | Bdrill[4] = B6*g3[1] ; |
||
| 1507 | Bdrill[5] = B6*g3[2] ; |
||
| 1508 | |||
| 1509 | return Bdrill ; |
||
| 1510 | |||
| 1511 | } |
||
| 1512 | |||
| 1513 | |||
| 1514 | //******************************************************************** |
||
| 1515 | //assemble a B matrix |
||
| 1516 | |||
| 1517 | const Matrix& |
||
| 1518 | ShellMITC4::assembleB( const Matrix &Bmembrane, |
||
| 1519 | const Matrix &Bbend, |
||
| 1520 | const Matrix &Bshear ) |
||
| 1521 | { |
||
| 1522 | |||
| 1523 | //Matrix Bbend(3,3) ; // plate bending B matrix |
||
| 1524 | |||
| 1525 | //Matrix Bshear(2,3) ; // plate shear B matrix |
||
| 1526 | |||
| 1527 | //Matrix Bmembrane(3,2) ; // plate membrane B matrix |
||
| 1528 | |||
| 1529 | |||
| 1530 | static Matrix B(8,6) ; |
||
| 1531 | |||
| 1532 | static Matrix BmembraneShell(3,3) ; |
||
| 1533 | |||
| 1534 | static Matrix BbendShell(3,3) ; |
||
| 1535 | |||
| 1536 | static Matrix BshearShell(2,6) ; |
||
| 1537 | |||
| 1538 | static Matrix Gmem(2,3) ; |
||
| 1539 | |||
| 1540 | static Matrix Gshear(3,6) ; |
||
| 1541 | |||
| 1542 | int p, q ; |
||
| 1543 | int pp ; |
||
| 1544 | |||
| 1545 | // |
||
| 1546 | // For Shell : |
||
| 1547 | // |
||
| 1548 | //---B Matrices in standard {1,2,3} mechanics notation--------- |
||
| 1549 | // |
||
| 1550 | // - _ |
||
| 1551 | // | Bmembrane | 0 | |
||
| 1552 | // | --------------------- | |
||
| 1553 | // B = | 0 | Bbend | (8x6) |
||
| 1554 | // | --------------------- | |
||
| 1555 | // | Bshear | |
||
| 1556 | // - - - |
||
| 1557 | // |
||
| 1558 | //------------------------------------------------------------- |
||
| 1559 | // |
||
| 1560 | // |
||
| 1561 | |||
| 1562 | |||
| 1563 | //shell modified membrane terms |
||
| 1564 | |||
| 1565 | Gmem(0,0) = g1[0] ; |
||
| 1566 | Gmem(0,1) = g1[1] ; |
||
| 1567 | Gmem(0,2) = g1[2] ; |
||
| 1568 | |||
| 1569 | Gmem(1,0) = g2[0] ; |
||
| 1570 | Gmem(1,1) = g2[1] ; |
||
| 1571 | Gmem(1,2) = g2[2] ; |
||
| 1572 | |||
| 1573 | //BmembraneShell = Bmembrane * Gmem ; |
||
| 1574 | BmembraneShell.addMatrixProduct(0.0, Bmembrane,Gmem,1.0 ) ; |
||
| 1575 | |||
| 1576 | |||
| 1577 | //shell modified bending terms |
||
| 1578 | |||
| 1579 | Matrix &Gbend = Gmem ; |
||
| 1580 | |||
| 1581 | //BbendShell = Bbend * Gbend ; |
||
| 1582 | BbendShell.addMatrixProduct(0.0, Bbend,Gbend,1.0 ) ; |
||
| 1583 | |||
| 1584 | |||
| 1585 | //shell modified shear terms |
||
| 1586 | |||
| 1587 | Gshear.Zero( ) ; |
||
| 1588 | |||
| 1589 | Gshear(0,0) = g3[0] ; |
||
| 1590 | Gshear(0,1) = g3[1] ; |
||
| 1591 | Gshear(0,2) = g3[2] ; |
||
| 1592 | |||
| 1593 | Gshear(1,3) = g1[0] ; |
||
| 1594 | Gshear(1,4) = g1[1] ; |
||
| 1595 | Gshear(1,5) = g1[2] ; |
||
| 1596 | |||
| 1597 | Gshear(2,3) = g2[0] ; |
||
| 1598 | Gshear(2,4) = g2[1] ; |
||
| 1599 | Gshear(2,5) = g2[2] ; |
||
| 1600 | |||
| 1601 | //BshearShell = Bshear * Gshear ; |
||
| 1602 | BshearShell.addMatrixProduct(0.0, Bshear,Gshear,1.0 ) ; |
||
| 1603 | |||
| 1604 | |||
| 1605 | B.Zero( ) ; |
||
| 1606 | |||
| 1607 | //assemble B from sub-matrices |
||
| 1608 | |||
| 1609 | //membrane terms |
||
| 1610 | for ( p = 0; p < 3; p++ ) { |
||
| 1611 | |||
| 1612 | for ( q = 0; q < 3; q++ ) |
||
| 1613 | B(p,q) = BmembraneShell(p,q) ; |
||
| 1614 | |||
| 1615 | } //end for p |
||
| 1616 | |||
| 1617 | |||
| 1618 | //bending terms |
||
| 1619 | for ( p = 3; p < 6; p++ ) { |
||
| 1620 | pp = p - 3 ; |
||
| 1621 | for ( q = 3; q < 6; q++ ) |
||
| 1622 | B(p,q) = BbendShell(pp,q-3) ; |
||
| 1623 | } // end for p |
||
| 1624 | |||
| 1625 | |||
| 1626 | //shear terms |
||
| 1627 | for ( p = 0; p < 2; p++ ) { |
||
| 1628 | pp = p + 6 ; |
||
| 1629 | |||
| 1630 | for ( q = 0; q < 6; q++ ) { |
||
| 1631 | B(pp,q) = BshearShell(p,q) ; |
||
| 1632 | } // end for q |
||
| 1633 | |||
| 1634 | } //end for p |
||
| 1635 | |||
| 1636 | return B ; |
||
| 1637 | |||
| 1638 | } |
||
| 1639 | |||
| 1640 | //*********************************************************************** |
||
| 1641 | //compute Bmembrane matrix |
||
| 1642 | |||
| 1643 | const Matrix& |
||
| 1644 | ShellMITC4::computeBmembrane( int node, const double shp[3][4] ) |
||
| 1645 | { |
||
| 1646 | |||
| 1647 | static Matrix Bmembrane(3,2) ; |
||
| 1648 | |||
| 1649 | //---Bmembrane Matrix in standard {1,2,3} mechanics notation--------- |
||
| 1650 | // |
||
| 1651 | // - - |
||
| 1652 | // | +N,1 0 | |
||
| 1653 | // Bmembrane = | 0 +N,2 | (3x2) |
||
| 1654 | // | +N,2 +N,1 | |
||
| 1655 | // - - |
||
| 1656 | // |
||
| 1657 | // three(3) strains and two(2) displacements (for plate) |
||
| 1658 | //------------------------------------------------------------------- |
||
| 1659 | |||
| 1660 | Bmembrane.Zero( ) ; |
||
| 1661 | |||
| 1662 | Bmembrane(0,0) = shp[0][node] ; |
||
| 1663 | Bmembrane(1,1) = shp[1][node] ; |
||
| 1664 | Bmembrane(2,0) = shp[1][node] ; |
||
| 1665 | Bmembrane(2,1) = shp[0][node] ; |
||
| 1666 | |||
| 1667 | return Bmembrane ; |
||
| 1668 | |||
| 1669 | } |
||
| 1670 | |||
| 1671 | //*********************************************************************** |
||
| 1672 | //compute Bbend matrix |
||
| 1673 | |||
| 1674 | const Matrix& |
||
| 1675 | ShellMITC4::computeBbend( int node, const double shp[3][4] ) |
||
| 1676 | { |
||
| 1677 | |||
| 1678 | static Matrix Bbend(3,2) ; |
||
| 1679 | |||
| 1680 | //---Bbend Matrix in standard {1,2,3} mechanics notation--------- |
||
| 1681 | // |
||
| 1682 | // - - |
||
| 1683 | // Bbend = | 0 -N,1 | |
||
| 1684 | // | +N,2 0 | (3x2) |
||
| 1685 | // | +N,1 -N,2 | |
||
| 1686 | // - - |
||
| 1687 | // |
||
| 1688 | // three(3) curvatures and two(2) rotations (for plate) |
||
| 1689 | //---------------------------------------------------------------- |
||
| 1690 | |||
| 1691 | Bbend.Zero( ) ; |
||
| 1692 | |||
| 1693 | Bbend(0,1) = -shp[0][node] ; |
||
| 1694 | Bbend(1,0) = shp[1][node] ; |
||
| 1695 | Bbend(2,0) = shp[0][node] ; |
||
| 1696 | Bbend(2,1) = -shp[1][node] ; |
||
| 1697 | |||
| 1698 | return Bbend ; |
||
| 1699 | } |
||
| 1700 | |||
| 1701 | |||
| 1702 | |||
| 1703 | //************************************************************************ |
||
| 1704 | //shape function routine for four node quads |
||
| 1705 | |||
| 1706 | void |
||
| 1707 | ShellMITC4::shape2d( double ss, double tt, |
||
| 1708 | const double x[2][4], |
||
| 1709 | double shp[3][4], |
||
| 1710 | double &xsj ) |
||
| 1711 | |||
| 1712 | { |
||
| 1713 | |||
| 1714 | int i, j, k ; |
||
| 1715 | |||
| 1716 | double temp ; |
||
| 1717 | |||
| 1718 | static const double s[] = { -0.5, 0.5, 0.5, -0.5 } ; |
||
| 1719 | static const double t[] = { -0.5, -0.5, 0.5, 0.5 } ; |
||
| 1720 | |||
| 1721 | static double xs[2][2] ; |
||
| 1722 | static double sx[2][2] ; |
||
| 1723 | |||
| 1724 | for ( i = 0; i < 4; i++ ) { |
||
| 1725 | shp[2][i] = ( 0.5 + s[i]*ss )*( 0.5 + t[i]*tt ) ; |
||
| 1726 | shp[0][i] = s[i] * ( 0.5 + t[i]*tt ) ; |
||
| 1727 | shp[1][i] = t[i] * ( 0.5 + s[i]*ss ) ; |
||
| 1728 | } // end for i |
||
| 1729 | |||
| 1730 | |||
| 1731 | // Construct jacobian and its inverse |
||
| 1732 | |||
| 1733 | for ( i = 0; i < 2; i++ ) { |
||
| 1734 | for ( j = 0; j < 2; j++ ) { |
||
| 1735 | |||
| 1736 | xs[i][j] = 0.0 ; |
||
| 1737 | |||
| 1738 | for ( k = 0; k < 4; k++ ) |
||
| 1739 | xs[i][j] += x[i][k] * shp[j][k] ; |
||
| 1740 | |||
| 1741 | } //end for j |
||
| 1742 | } // end for i |
||
| 1743 | |||
| 1744 | xsj = xs[0][0]*xs[1][1] - xs[0][1]*xs[1][0] ; |
||
| 1745 | |||
| 1746 | //inverse jacobian |
||
| 1747 | double jinv = 1.0 / xsj ; |
||
| 1748 | sx[0][0] = xs[1][1] * jinv ; |
||
| 1749 | sx[1][1] = xs[0][0] * jinv ; |
||
| 1750 | sx[0][1] = -xs[0][1] * jinv ; |
||
| 1751 | sx[1][0] = -xs[1][0] * jinv ; |
||
| 1752 | |||
| 1753 | |||
| 1754 | //form global derivatives |
||
| 1755 | |||
| 1756 | for ( i = 0; i < 4; i++ ) { |
||
| 1757 | temp = shp[0][i]*sx[0][0] + shp[1][i]*sx[1][0] ; |
||
| 1758 | shp[1][i] = shp[0][i]*sx[0][1] + shp[1][i]*sx[1][1] ; |
||
| 1759 | shp[0][i] = temp ; |
||
| 1760 | } // end for i |
||
| 1761 | |||
| 1762 | return ; |
||
| 1763 | } |
||
| 1764 | |||
| 1765 | //********************************************************************** |
||
| 1766 | |||
| 1767 | Matrix |
||
| 1768 | ShellMITC4::transpose( int dim1, |
||
| 1769 | int dim2, |
||
| 1770 | const Matrix &M ) |
||
| 1771 | { |
||
| 1772 | int i ; |
||
| 1773 | int j ; |
||
| 1774 | |||
| 1775 | Matrix Mtran( dim2, dim1 ) ; |
||
| 1776 | |||
| 1777 | for ( i = 0; i < dim1; i++ ) { |
||
| 1778 | for ( j = 0; j < dim2; j++ ) |
||
| 1779 | Mtran(j,i) = M(i,j) ; |
||
| 1780 | } // end for i |
||
| 1781 | |||
| 1782 | return Mtran ; |
||
| 1783 | } |
||
| 1784 | |||
| 1785 | //********************************************************************** |
||
| 1786 | |||
| 1787 | int ShellMITC4::sendSelf (int commitTag, Channel &theChannel) |
||
| 1788 | { |
||
| 1789 | int res = 0; |
||
| 1790 | |||
| 1791 | // note: we don't check for dataTag == 0 for Element |
||
| 1792 | // objects as that is taken care of in a commit by the Domain |
||
| 1793 | // object - don't want to have to do the check if sending data |
||
| 1794 | int dataTag = this->getDbTag(); |
||
| 1795 | |||
| 1796 | |||
| 1797 | // Now quad sends the ids of its materials |
||
| 1798 | int matDbTag; |
||
| 1799 | |||
| 1800 | static ID idData(13); |
||
| 1801 | |||
| 1802 | int i; |
||
| 1803 | for (i = 0; i < 4; i++) { |
||
| 1804 | idData(i) = materialPointers[i]->getClassTag(); |
||
| 1805 | matDbTag = materialPointers[i]->getDbTag(); |
||
| 1806 | // NOTE: we do have to ensure that the material has a database |
||
| 1807 | // tag if we are sending to a database channel. |
||
| 1808 | if (matDbTag == 0) { |
||
| 1809 | matDbTag = theChannel.getDbTag(); |
||
| 1810 | if (matDbTag != 0) |
||
| 1811 | materialPointers[i]->setDbTag(matDbTag); |
||
| 1812 | } |
||
| 1813 | idData(i+4) = matDbTag; |
||
| 1814 | } |
||
| 1815 | |||
| 1816 | idData(8) = this->getTag(); |
||
| 1817 | idData(9) = connectedExternalNodes(0); |
||
| 1818 | idData(10) = connectedExternalNodes(1); |
||
| 1819 | idData(11) = connectedExternalNodes(2); |
||
| 1820 | idData(12) = connectedExternalNodes(3); |
||
| 1821 | |||
| 1822 | res += theChannel.sendID(dataTag, commitTag, idData); |
||
| 1823 | if (res < 0) { |
||
| 1824 | opserr << "WARNING ShellMITC4::sendSelf() - " << this->getTag() << " failed to send ID\n"; |
||
| 1825 | return res; |
||
| 1826 | } |
||
| 1827 | |||
| 1828 | static Vector vectData(5); |
||
| 1829 | vectData(0) = Ktt; |
||
| 1830 | vectData(1) = alphaM; |
||
| 1831 | vectData(2) = betaK; |
||
| 1832 | vectData(3) = betaK0; |
||
| 1833 | vectData(4) = betaKc; |
||
| 1834 | |||
| 1835 | res += theChannel.sendVector(dataTag, commitTag, vectData); |
||
| 1836 | if (res < 0) { |
||
| 1837 | opserr << "WARNING ShellMITC4::sendSelf() - " << this->getTag() << " failed to send ID\n"; |
||
| 1838 | return res; |
||
| 1839 | } |
||
| 1840 | |||
| 1841 | // Finally, quad asks its material objects to send themselves |
||
| 1842 | for (i = 0; i < 4; i++) { |
||
| 1843 | res += materialPointers[i]->sendSelf(commitTag, theChannel); |
||
| 1844 | if (res < 0) { |
||
| 1845 | opserr << "WARNING ShellMITC4::sendSelf() - " << this->getTag() << " failed to send its Material\n"; |
||
| 1846 | return res; |
||
| 1847 | } |
||
| 1848 | } |
||
| 1849 | |||
| 1850 | return res; |
||
| 1851 | } |
||
| 1852 | |||
| 1853 | int ShellMITC4::recvSelf (int commitTag, |
||
| 1854 | Channel &theChannel, |
||
| 1855 | FEM_ObjectBroker &theBroker) |
||
| 1856 | { |
||
| 1857 | int res = 0; |
||
| 1858 | |||
| 1859 | int dataTag = this->getDbTag(); |
||
| 1860 | |||
| 1861 | static ID idData(13); |
||
| 1862 | // Quad now receives the tags of its four external nodes |
||
| 1863 | res += theChannel.recvID(dataTag, commitTag, idData); |
||
| 1864 | if (res < 0) { |
||
| 1865 | opserr << "WARNING ShellMITC4::recvSelf() - " << this->getTag() << " failed to receive ID\n"; |
||
| 1866 | return res; |
||
| 1867 | } |
||
| 1868 | |||
| 1869 | this->setTag(idData(8)); |
||
| 1870 | connectedExternalNodes(0) = idData(9); |
||
| 1871 | connectedExternalNodes(1) = idData(10); |
||
| 1872 | connectedExternalNodes(2) = idData(11); |
||
| 1873 | connectedExternalNodes(3) = idData(12); |
||
| 1874 | |||
| 1875 | static Vector vectData(5); |
||
| 1876 | res += theChannel.recvVector(dataTag, commitTag, vectData); |
||
| 1877 | if (res < 0) { |
||
| 1878 | opserr << "WARNING ShellMITC4::sendSelf() - " << this->getTag() << " failed to send ID\n"; |
||
| 1879 | return res; |
||
| 1880 | } |
||
| 1881 | |||
| 1882 | Ktt = vectData(0); |
||
| 1883 | alphaM = vectData(1); |
||
| 1884 | betaK = vectData(2); |
||
| 1885 | betaK0 = vectData(3); |
||
| 1886 | betaKc = vectData(4); |
||
| 1887 | |||
| 1888 | int i; |
||
| 1889 | |||
| 1890 | if (materialPointers[0] == 0) { |
||
| 1891 | for (i = 0; i < 4; i++) { |
||
| 1892 | int matClassTag = idData(i); |
||
| 1893 | int matDbTag = idData(i+4); |
||
| 1894 | // Allocate new material with the sent class tag |
||
| 1895 | materialPointers[i] = theBroker.getNewSection(matClassTag); |
||
| 1896 | if (materialPointers[i] == 0) { |
||
| 1897 | opserr << "ShellMITC4::recvSelf() - Broker could not create NDMaterial of class type" << matClassTag << endln;; |
||
| 1898 | return -1; |
||
| 1899 | } |
||
| 1900 | // Now receive materials into the newly allocated space |
||
| 1901 | materialPointers[i]->setDbTag(matDbTag); |
||
| 1902 | res += materialPointers[i]->recvSelf(commitTag, theChannel, theBroker); |
||
| 1903 | if (res < 0) { |
||
| 1904 | opserr << "ShellMITC4::recvSelf() - material " << i << "failed to recv itself\n"; |
||
| 1905 | return res; |
||
| 1906 | } |
||
| 1907 | } |
||
| 1908 | } |
||
| 1909 | // Number of materials is the same, receive materials into current space |
||
| 1910 | else { |
||
| 1911 | for (i = 0; i < 4; i++) { |
||
| 1912 | int matClassTag = idData(i); |
||
| 1913 | int matDbTag = idData(i+4); |
||
| 1914 | // Check that material is of the right type; if not, |
||
| 1915 | // delete it and create a new one of the right type |
||
| 1916 | if (materialPointers[i]->getClassTag() != matClassTag) { |
||
| 1917 | delete materialPointers[i]; |
||
| 1918 | materialPointers[i] = theBroker.getNewSection(matClassTag); |
||
| 1919 | if (materialPointers[i] == 0) { |
||
| 1920 | opserr << "ShellMITC4::recvSelf() - Broker could not create NDMaterial of class type" << matClassTag << endln; |
||
| 1921 | exit(-1); |
||
| 1922 | } |
||
| 1923 | } |
||
| 1924 | // Receive the material |
||
| 1925 | materialPointers[i]->setDbTag(matDbTag); |
||
| 1926 | res += materialPointers[i]->recvSelf(commitTag, theChannel, theBroker); |
||
| 1927 | if (res < 0) { |
||
| 1928 | opserr << "ShellMITC4::recvSelf() - material " << i << "failed to recv itself\n"; |
||
| 1929 | return res; |
||
| 1930 | } |
||
| 1931 | } |
||
| 1932 | } |
||
| 1933 | |||
| 1934 | return res; |
||
| 1935 | } |
||
| 1936 | //************************************************************************** |
||
| 1937 | |||
| 1938 | int |
||
| 1939 | ShellMITC4::displaySelf(Renderer &theViewer, int displayMode, float fact) |
||
| 1940 | { |
||
| 1941 | // first determine the end points of the quad based on |
||
| 1942 | // the display factor (a measure of the distorted image) |
||
| 1943 | // store this information in 4 3d vectors v1 through v4 |
||
| 1944 | const Vector &end1Crd = nodePointers[0]->getCrds(); |
||
| 1945 | const Vector &end2Crd = nodePointers[1]->getCrds(); |
||
| 1946 | const Vector &end3Crd = nodePointers[2]->getCrds(); |
||
| 1947 | const Vector &end4Crd = nodePointers[3]->getCrds(); |
||
| 1948 | |||
| 1949 | static Matrix coords(4,3); |
||
| 1950 | static Vector values(4); |
||
| 1951 | static Vector P(24) ; |
||
| 1952 | |||
| 1953 | for (int j=0; j<4; j++) |
||
| 1954 | values(j) = 0.0; |
||
| 1955 | |||
| 1956 | if (displayMode >= 0) { |
||
| 1957 | // Display mode is positive: |
||
| 1958 | // display mode = 0 -> plot no contour |
||
| 1959 | // display mode = 1-8 -> plot 1-8 stress resultant |
||
| 1960 | |||
| 1961 | // Get nodal displacements |
||
| 1962 | const Vector &end1Disp = nodePointers[0]->getDisp(); |
||
| 1963 | const Vector &end2Disp = nodePointers[1]->getDisp(); |
||
| 1964 | const Vector &end3Disp = nodePointers[2]->getDisp(); |
||
| 1965 | const Vector &end4Disp = nodePointers[3]->getDisp(); |
||
| 1966 | |||
| 1967 | // Get stress resultants |
||
| 1968 | if (displayMode <= 8 && displayMode > 0) { |
||
| 1969 | for (int i=0; i<4; i++) { |
||
| 1970 | const Vector &stress = materialPointers[i]->getStressResultant(); |
||
| 1971 | values(i) = stress(displayMode-1); |
||
| 1972 | } |
||
| 1973 | } |
||
| 1974 | |||
| 1975 | // Get nodal absolute position = OriginalPosition + (Displacement*factor) |
||
| 1976 | for (int i = 0; i < 3; i++) { |
||
| 1977 | coords(0,i) = end1Crd(i) + end1Disp(i)*fact; |
||
| 1978 | coords(1,i) = end2Crd(i) + end2Disp(i)*fact; |
||
| 1979 | coords(2,i) = end3Crd(i) + end3Disp(i)*fact; |
||
| 1980 | coords(3,i) = end4Crd(i) + end4Disp(i)*fact; |
||
| 1981 | } |
||
| 1982 | } else { |
||
| 1983 | // Display mode is negative. |
||
| 1984 | // Plot eigenvectors |
||
| 1985 | int mode = displayMode * -1; |
||
| 1986 | const Matrix &eigen1 = nodePointers[0]->getEigenvectors(); |
||
| 1987 | const Matrix &eigen2 = nodePointers[1]->getEigenvectors(); |
||
| 1988 | const Matrix &eigen3 = nodePointers[2]->getEigenvectors(); |
||
| 1989 | const Matrix &eigen4 = nodePointers[3]->getEigenvectors(); |
||
| 1990 | if (eigen1.noCols() >= mode) { |
||
| 1991 | for (int i = 0; i < 3; i++) { |
||
| 1992 | coords(0,i) = end1Crd(i) + eigen1(i,mode-1)*fact; |
||
| 1993 | coords(1,i) = end2Crd(i) + eigen2(i,mode-1)*fact; |
||
| 1994 | coords(2,i) = end3Crd(i) + eigen3(i,mode-1)*fact; |
||
| 1995 | coords(3,i) = end4Crd(i) + eigen4(i,mode-1)*fact; |
||
| 1996 | } |
||
| 1997 | } else { |
||
| 1998 | for (int i = 0; i < 3; i++) { |
||
| 1999 | coords(0,i) = end1Crd(i); |
||
| 2000 | coords(1,i) = end2Crd(i); |
||
| 2001 | coords(2,i) = end3Crd(i); |
||
| 2002 | coords(3,i) = end4Crd(i); |
||
| 2003 | } |
||
| 2004 | } |
||
| 2005 | } |
||
| 2006 | |||
| 2007 | int error = 0; |
||
| 2008 | |||
| 2009 | // Draw a poligon with coordinates coords and values (colors) corresponding to values vector |
||
| 2010 | error += theViewer.drawPolygon (coords, values); |
||
| 2011 | |||
| 2012 | return error; |
||
| 2013 | } |