00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030 #include <stdio.h>
00031 #include <stdlib.h>
00032 #include <math.h>
00033
00034 #include <ID.h>
00035 #include <Vector.h>
00036 #include <Matrix.h>
00037 #include <Element.h>
00038 #include <Node.h>
00039 #include <Domain.h>
00040 #include <ErrorHandler.h>
00041 #include <Brick.h>
00042 #include <shp3d.h>
00043 #include <Renderer.h>
00044 #include <ElementResponse.h>
00045
00046 #include <Channel.h>
00047 #include <FEM_ObjectBroker.h>
00048
00049
00050 double Brick::xl[3][8] ;
00051
00052 Matrix Brick::stiff(24,24) ;
00053 Vector Brick::resid(24) ;
00054 Matrix Brick::mass(24,24) ;
00055
00056
00057
00058 const double Brick::root3 = sqrt(3.0) ;
00059 const double Brick::one_over_root3 = 1.0 / root3 ;
00060
00061 const double Brick::sg[] = { -one_over_root3,
00062 one_over_root3 } ;
00063
00064 const double Brick::wg[] = { 1.0, 1.0, 1.0, 1.0,
00065 1.0, 1.0, 1.0, 1.0 } ;
00066
00067
00068 static Matrix B(6,3) ;
00069
00070
00071 Brick::Brick( )
00072 :Element( 0, ELE_TAG_Brick ),
00073 connectedExternalNodes(8), load(0), Ki(0)
00074 {
00075 B.Zero();
00076
00077 for (int i=0; i<8; i++ ) {
00078 materialPointers[i] = 0;
00079 nodePointers[i] = 0;
00080 }
00081
00082 b[0] = 0.0;
00083 b[1] = 0.0;
00084 b[2] = 0.0;
00085 }
00086
00087
00088
00089
00090 Brick::Brick(int tag,
00091 int node1,
00092 int node2,
00093 int node3,
00094 int node4,
00095 int node5,
00096 int node6,
00097 int node7,
00098 int node8,
00099 NDMaterial &theMaterial,
00100 double b1, double b2, double b3)
00101 :Element(tag, ELE_TAG_Brick),
00102 connectedExternalNodes(8) , load(0), Ki(0)
00103 {
00104 B.Zero();
00105
00106 connectedExternalNodes(0) = node1 ;
00107 connectedExternalNodes(1) = node2 ;
00108 connectedExternalNodes(2) = node3 ;
00109 connectedExternalNodes(3) = node4 ;
00110
00111 connectedExternalNodes(4) = node5 ;
00112 connectedExternalNodes(5) = node6 ;
00113 connectedExternalNodes(6) = node7 ;
00114 connectedExternalNodes(7) = node8 ;
00115
00116 for (int i=0; i<8; i++ ) {
00117 materialPointers[i] = theMaterial.getCopy("ThreeDimensional") ;
00118 if (materialPointers[i] == 0) {
00119 opserr << "Brick::constructor - failed to get a material of type: ThreeDimensional\n";
00120 exit(-1);
00121 }
00122 nodePointers[i] = 0;
00123 }
00124
00125
00126 b[0] = b1;
00127 b[1] = b2;
00128 b[2] = b3;
00129 }
00130
00131
00132
00133
00134 Brick::~Brick( )
00135 {
00136
00137 for (int i=0 ; i<8; i++ ) {
00138 delete materialPointers[i] ;
00139 }
00140
00141 if (load != 0)
00142 delete load;
00143
00144 if (Ki != 0)
00145 delete Ki;
00146
00147 }
00148
00149
00150
00151 void Brick::setDomain( Domain *theDomain )
00152 {
00153
00154 int i ;
00155
00156
00157 for ( i=0; i<8; i++ )
00158 nodePointers[i] = theDomain->getNode( connectedExternalNodes(i) ) ;
00159
00160 this->DomainComponent::setDomain(theDomain);
00161
00162 }
00163
00164
00165
00166 int Brick::getNumExternalNodes( ) const
00167 {
00168 return 8 ;
00169 }
00170
00171
00172
00173 const ID& Brick::getExternalNodes( )
00174 {
00175 return connectedExternalNodes ;
00176 }
00177
00178 Node **
00179 Brick::getNodePtrs(void)
00180 {
00181 return nodePointers ;
00182 }
00183
00184
00185 int Brick::getNumDOF( )
00186 {
00187 return 24 ;
00188 }
00189
00190
00191
00192 int Brick::commitState( )
00193 {
00194 int success = 0 ;
00195
00196
00197 if ((success = this->Element::commitState()) != 0) {
00198 opserr << "Brick::commitState () - failed in base class";
00199 }
00200
00201
00202 for (int i=0; i<8; i++ )
00203 success += materialPointers[i]->commitState( ) ;
00204
00205 return success ;
00206 }
00207
00208
00209
00210
00211 int Brick::revertToLastCommit( )
00212 {
00213 int i ;
00214 int success = 0 ;
00215
00216 for ( i=0; i<8; i++ )
00217 success += materialPointers[i]->revertToLastCommit( ) ;
00218
00219 return success ;
00220 }
00221
00222
00223
00224 int Brick::revertToStart( )
00225 {
00226 int i ;
00227 int success = 0 ;
00228
00229 for ( i=0; i<8; i++ )
00230 success += materialPointers[i]->revertToStart( ) ;
00231
00232 return success ;
00233 }
00234
00235
00236 void Brick::Print( OPS_Stream &s, int flag )
00237 {
00238
00239 if (flag == 2) {
00240
00241 s << "#Brick\n";
00242
00243 int i;
00244 const int numNodes = 8;
00245 const int nstress = 6 ;
00246
00247 for (i=0; i<numNodes; i++) {
00248 const Vector &nodeCrd = nodePointers[i]->getCrds();
00249 const Vector &nodeDisp = nodePointers[i]->getDisp();
00250 s << "#NODE " << nodeCrd(0) << " " << nodeCrd(1) << " " << nodeCrd(2)
00251 << " " << nodeDisp(0) << " " << nodeDisp(1) << " " << nodeDisp(2) << endln;
00252 }
00253
00254
00255 const int numMaterials = 8;
00256
00257 static Vector avgStress(nstress);
00258 static Vector avgStrain(nstress);
00259 avgStress.Zero();
00260 avgStrain.Zero();
00261 for (i=0; i<numMaterials; i++) {
00262 avgStress += materialPointers[i]->getCommittedStress();
00263 avgStrain += materialPointers[i]->getCommittedStrain();
00264 }
00265 avgStress /= numMaterials;
00266 avgStrain /= numMaterials;
00267
00268 s << "#AVERAGE_STRESS ";
00269 for (i=0; i<nstress; i++)
00270 s << avgStress(i) << " " ;
00271 s << endln;
00272
00273 s << "#AVERAGE_STRAIN ";
00274 for (i=0; i<nstress; i++)
00275 s << avgStrain(i) << " " ;
00276 s << endln;
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286 } else {
00287
00288 s << "Standard Eight Node Brick \n" ;
00289 s << "Element Number: " << this->getTag() << endln ;
00290 s << "Nodes: " << connectedExternalNodes;
00291
00292 s << "Material Information : \n " ;
00293 materialPointers[0]->Print( s, flag ) ;
00294
00295 s << endln ;
00296 s << this->getTag() << " " <<connectedExternalNodes(0)
00297 << " " <<connectedExternalNodes(1)
00298 << " " <<connectedExternalNodes(2)
00299 << " " <<connectedExternalNodes(3)
00300 << " " <<connectedExternalNodes(4)
00301 << " " <<connectedExternalNodes(5)
00302 << " " <<connectedExternalNodes(6)
00303 << " " <<connectedExternalNodes(7)
00304 << endln ;
00305
00306 s << "Resisting Force (no inertia): " << this->getResistingForce();
00307 }
00308 }
00309
00310
00311
00312 const Matrix& Brick::getTangentStiff( )
00313 {
00314 int tang_flag = 1 ;
00315
00316
00317 formResidAndTangent( tang_flag ) ;
00318
00319 return stiff ;
00320 }
00321
00322
00323
00324
00325
00326 const Matrix& Brick::getInitialStiff( )
00327
00328 {
00329 if (Ki != 0)
00330 return *Ki;
00331
00332
00333 static const int ndm = 3 ;
00334 static const int ndf = 3 ;
00335 static const int nstress = 6 ;
00336 static const int numberNodes = 8 ;
00337 static const int numberGauss = 8 ;
00338 static const int nShape = 4 ;
00339
00340 int i, j, k, p, q ;
00341 int jj, kk ;
00342
00343
00344 static double volume ;
00345 static double xsj ;
00346 static double dvol[numberGauss] ;
00347 static double gaussPoint[ndm] ;
00348 static Vector strain(nstress) ;
00349 static double shp[nShape][numberNodes] ;
00350 static double Shape[nShape][numberNodes][numberGauss] ;
00351 static Matrix stiffJK(ndf,ndf) ;
00352 static Matrix dd(nstress,nstress) ;
00353
00354
00355
00356
00357 static Matrix BJ(nstress,ndf) ;
00358
00359 static Matrix BJtran(ndf,nstress) ;
00360
00361 static Matrix BK(nstress,ndf) ;
00362
00363 static Matrix BJtranD(ndf,nstress) ;
00364
00365
00366
00367
00368
00369 stiff.Zero( ) ;
00370
00371
00372 computeBasis( ) ;
00373
00374
00375
00376 int count = 0 ;
00377 volume = 0.0 ;
00378
00379 for ( i = 0; i < 2; i++ ) {
00380 for ( j = 0; j < 2; j++ ) {
00381 for ( k = 0; k < 2; k++ ) {
00382
00383 gaussPoint[0] = sg[i] ;
00384 gaussPoint[1] = sg[j] ;
00385 gaussPoint[2] = sg[k] ;
00386
00387
00388 shp3d( gaussPoint, xsj, shp, xl ) ;
00389
00390
00391 for ( p = 0; p < nShape; p++ ) {
00392 for ( q = 0; q < numberNodes; q++ )
00393 Shape[p][q][count] = shp[p][q] ;
00394 }
00395
00396
00397
00398 dvol[count] = wg[count] * xsj ;
00399
00400
00401
00402 count++ ;
00403
00404 }
00405 }
00406 }
00407
00408
00409
00410 for ( i = 0; i < numberGauss; i++ ) {
00411
00412
00413 for ( p = 0; p < nShape; p++ ) {
00414 for ( q = 0; q < numberNodes; q++ )
00415 shp[p][q] = Shape[p][q][i] ;
00416 }
00417
00418
00419 dd = materialPointers[i]->getInitialTangent( ) ;
00420 dd *= dvol[i] ;
00421
00422 jj = 0;
00423 for ( j = 0; j < numberNodes; j++ ) {
00424
00425 BJ = computeB( j, shp ) ;
00426
00427
00428
00429 for (p=0; p<ndf; p++) {
00430 for (q=0; q<nstress; q++)
00431 BJtran(p,q) = BJ(q,p) ;
00432 }
00433
00434
00435 BJtranD.addMatrixProduct(0.0, BJtran, dd, 1.0) ;
00436
00437 kk = 0 ;
00438 for ( k = 0; k < numberNodes; k++ ) {
00439
00440 BK = computeB( k, shp ) ;
00441
00442
00443
00444 stiffJK.addMatrixProduct(0.0, BJtranD, BK, 1.0) ;
00445
00446 for ( p = 0; p < ndf; p++ ) {
00447 for ( q = 0; q < ndf; q++ )
00448 stiff( jj+p, kk+q ) += stiffJK( p, q ) ;
00449 }
00450
00451 kk += ndf ;
00452
00453 }
00454
00455 jj += ndf ;
00456
00457 }
00458 }
00459
00460 Ki = new Matrix(stiff);
00461
00462 return stiff ;
00463 }
00464
00465
00466
00467 const Matrix& Brick::getMass( )
00468 {
00469 int tangFlag = 1 ;
00470
00471 formInertiaTerms( tangFlag ) ;
00472
00473 return mass ;
00474 }
00475
00476
00477
00478 void Brick::zeroLoad( )
00479 {
00480 if (load != 0)
00481 load->Zero();
00482
00483 return ;
00484 }
00485
00486
00487 int
00488 Brick::addLoad(ElementalLoad *theLoad, double loadFactor)
00489 {
00490 opserr << "Brick::addLoad - load type unknown for truss with tag: " << this->getTag() << endln;
00491 return -1;
00492 }
00493
00494 int
00495 Brick::addInertiaLoadToUnbalance(const Vector &accel)
00496 {
00497 static const int numberNodes = 8 ;
00498 static const int numberGauss = 8 ;
00499 static const int ndf = 3 ;
00500
00501 int i;
00502
00503
00504 int haveRho = 0;
00505 for (i = 0; i < numberGauss; i++) {
00506 if (materialPointers[i]->getRho() != 0.0)
00507 haveRho = 1;
00508 }
00509
00510 if (haveRho == 0)
00511 return 0;
00512
00513
00514 int tangFlag = 1 ;
00515 formInertiaTerms( tangFlag ) ;
00516
00517
00518 int count = 0;
00519 for (i=0; i<numberNodes; i++) {
00520 const Vector &Raccel = nodePointers[i]->getRV(accel);
00521 for (int j=0; j<ndf; j++)
00522 resid(count++) = Raccel(j);
00523 }
00524
00525
00526 if (load == 0)
00527 load = new Vector(numberNodes*ndf);
00528
00529
00530 load->addMatrixVector(1.0, mass, resid, -1.0);
00531
00532 return 0;
00533 }
00534
00535
00536
00537 const Vector& Brick::getResistingForce( )
00538 {
00539 int tang_flag = 0 ;
00540
00541 formResidAndTangent( tang_flag ) ;
00542
00543 if (load != 0)
00544 resid -= *load;
00545
00546 return resid ;
00547 }
00548
00549
00550
00551 const Vector& Brick::getResistingForceIncInertia( )
00552 {
00553 static Vector res(24);
00554
00555 int tang_flag = 0 ;
00556
00557
00558 formResidAndTangent( tang_flag ) ;
00559
00560 formInertiaTerms( tang_flag ) ;
00561
00562 res = resid;
00563
00564
00565 if (alphaM != 0.0 || betaK != 0.0 || betaK0 != 0.0 || betaKc != 0.0)
00566 res += this->getRayleighDampingForces();
00567
00568 if (load != 0)
00569 res -= *load;
00570
00571 return res;
00572 }
00573
00574
00575
00576
00577
00578 void Brick::formInertiaTerms( int tangFlag )
00579 {
00580
00581 static const int ndm = 3 ;
00582
00583 static const int ndf = 3 ;
00584
00585 static const int numberNodes = 8 ;
00586
00587 static const int numberGauss = 8 ;
00588
00589 static const int nShape = 4 ;
00590
00591 static const int massIndex = nShape - 1 ;
00592
00593 double xsj ;
00594
00595 double dvol[numberGauss] ;
00596
00597 static double shp[nShape][numberNodes] ;
00598
00599 static double Shape[nShape][numberNodes][numberGauss] ;
00600
00601 static double gaussPoint[ndm] ;
00602
00603 static Vector momentum(ndf) ;
00604
00605 int i, j, k, p, q ;
00606 int jj, kk ;
00607
00608 double temp, rho, massJK ;
00609
00610
00611
00612 mass.Zero( ) ;
00613
00614
00615 computeBasis( ) ;
00616
00617
00618
00619 int count = 0 ;
00620
00621 for ( i = 0; i < 2; i++ ) {
00622 for ( j = 0; j < 2; j++ ) {
00623 for ( k = 0; k < 2; k++ ) {
00624
00625 gaussPoint[0] = sg[i] ;
00626 gaussPoint[1] = sg[j] ;
00627 gaussPoint[2] = sg[k] ;
00628
00629
00630 shp3d( gaussPoint, xsj, shp, xl ) ;
00631
00632
00633 for ( p = 0; p < nShape; p++ ) {
00634 for ( q = 0; q < numberNodes; q++ )
00635 Shape[p][q][count] = shp[p][q] ;
00636 }
00637
00638
00639
00640 dvol[count] = wg[count] * xsj ;
00641
00642 count++ ;
00643
00644 }
00645 }
00646 }
00647
00648
00649
00650
00651 for ( i = 0; i < numberGauss; i++ ) {
00652
00653
00654 for ( p = 0; p < nShape; p++ ) {
00655 for ( q = 0; q < numberNodes; q++ )
00656 shp[p][q] = Shape[p][q][i] ;
00657 }
00658
00659
00660
00661 momentum.Zero( ) ;
00662 for ( j = 0; j < numberNodes; j++ )
00663
00664 momentum.addVector( 1.0,
00665 nodePointers[j]->getTrialAccel(),
00666 shp[massIndex][j] ) ;
00667
00668
00669
00670 rho = materialPointers[i]->getRho() ;
00671
00672
00673
00674 momentum *= rho ;
00675
00676
00677
00678 jj = 0 ;
00679 for ( j = 0; j < numberNodes; j++ ) {
00680
00681 temp = shp[massIndex][j] * dvol[i] ;
00682
00683 for ( p = 0; p < ndf; p++ )
00684 resid( jj+p ) += ( temp * momentum(p) ) ;
00685
00686
00687 if ( tangFlag == 1 ) {
00688
00689
00690 temp *= rho ;
00691
00692
00693 kk = 0 ;
00694 for ( k = 0; k < numberNodes; k++ ) {
00695
00696 massJK = temp * shp[massIndex][k] ;
00697
00698 for ( p = 0; p < ndf; p++ )
00699 mass( jj+p, kk+p ) += massJK ;
00700
00701 kk += ndf ;
00702 }
00703
00704 }
00705
00706 jj += ndf ;
00707 }
00708
00709
00710 }
00711
00712 }
00713
00714
00715
00716 int
00717 Brick::update(void)
00718 {
00719
00720
00721
00722 static const int ndm = 3 ;
00723
00724 static const int ndf = 3 ;
00725
00726 static const int nstress = 6 ;
00727
00728 static const int numberNodes = 8 ;
00729
00730 static const int numberGauss = 8 ;
00731
00732 static const int nShape = 4 ;
00733
00734 int i, j, k, p, q ;
00735 int success ;
00736
00737 static double volume ;
00738
00739 static double xsj ;
00740
00741 static double dvol[numberGauss] ;
00742
00743 static double gaussPoint[ndm] ;
00744
00745 static Vector strain(nstress) ;
00746
00747 static double shp[nShape][numberNodes] ;
00748
00749 static double Shape[nShape][numberNodes][numberGauss] ;
00750
00751
00752
00753 static Matrix BJ(nstress,ndf) ;
00754
00755 static Matrix BJtran(ndf,nstress) ;
00756
00757 static Matrix BK(nstress,ndf) ;
00758
00759 static Matrix BJtranD(ndf,nstress) ;
00760
00761
00762
00763
00764
00765 computeBasis( ) ;
00766
00767
00768
00769 int count = 0 ;
00770 volume = 0.0 ;
00771
00772 for ( i = 0; i < 2; i++ ) {
00773 for ( j = 0; j < 2; j++ ) {
00774 for ( k = 0; k < 2; k++ ) {
00775
00776 gaussPoint[0] = sg[i] ;
00777 gaussPoint[1] = sg[j] ;
00778 gaussPoint[2] = sg[k] ;
00779
00780
00781 shp3d( gaussPoint, xsj, shp, xl ) ;
00782
00783
00784 for ( p = 0; p < nShape; p++ ) {
00785 for ( q = 0; q < numberNodes; q++ )
00786 Shape[p][q][count] = shp[p][q] ;
00787 }
00788
00789
00790
00791 dvol[count] = wg[count] * xsj ;
00792
00793
00794
00795 count++ ;
00796
00797 }
00798 }
00799 }
00800
00801
00802
00803 for ( i = 0; i < numberGauss; i++ ) {
00804
00805
00806 for ( p = 0; p < nShape; p++ ) {
00807 for ( q = 0; q < numberNodes; q++ )
00808 shp[p][q] = Shape[p][q][i] ;
00809 }
00810
00811
00812
00813 strain.Zero( ) ;
00814
00815
00816 for ( j = 0; j < numberNodes; j++ ) {
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832
00833
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848 double b00 = shp[0][j];
00849 double b11 = shp[1][j];
00850 double b22 = shp[2][j];
00851 double b30 = shp[1][j];
00852 double b31 = shp[0][j];
00853 double b41 = shp[2][j];
00854 double b42 = shp[1][j];
00855 double b50 = shp[2][j];
00856 double b52 = shp[0][j];
00857
00858 const Vector &ul = nodePointers[j]->getTrialDisp();
00859
00860 double ul0 = ul(0);
00861 double ul1 = ul(1);
00862 double ul2 = ul(2);
00863
00864 strain(0) += b00 * ul0;
00865 strain(1) += b11 * ul1;
00866 strain(2) += b22 * ul2;
00867 strain(3) += b30 * ul0 + b31 * ul1;
00868 strain(4) += b41 * ul1 + b42 * ul2;
00869 strain(5) += b50 * ul0 + b52 * ul2;
00870
00871 }
00872
00873
00874 success = materialPointers[i]->setTrialStrain( strain ) ;
00875
00876 }
00877
00878 return 0;
00879 }
00880
00881
00882
00883
00884 void Brick::formResidAndTangent( int tang_flag )
00885 {
00886
00887
00888
00889 static const int ndm = 3 ;
00890
00891 static const int ndf = 3 ;
00892
00893 static const int nstress = 6 ;
00894
00895 static const int numberNodes = 8 ;
00896
00897 static const int numberGauss = 8 ;
00898
00899 static const int nShape = 4 ;
00900
00901 int i, j, k, p, q ;
00902
00903 int success ;
00904
00905 static double volume ;
00906
00907 static double xsj ;
00908
00909 static double dvol[numberGauss] ;
00910
00911 static double gaussPoint[ndm] ;
00912
00913 static double shp[nShape][numberNodes] ;
00914
00915 static double Shape[nShape][numberNodes][numberGauss] ;
00916
00917 static Vector residJ(ndf) ;
00918
00919 static Matrix stiffJK(ndf,ndf) ;
00920
00921 static Vector stress(nstress) ;
00922
00923 static Matrix dd(nstress,nstress) ;
00924
00925
00926
00927
00928 static Matrix BJ(nstress,ndf) ;
00929
00930 static Matrix BJtran(ndf,nstress) ;
00931
00932 static Matrix BK(nstress,ndf) ;
00933
00934 static Matrix BJtranD(ndf,nstress) ;
00935
00936
00937
00938
00939
00940 stiff.Zero( ) ;
00941 resid.Zero( ) ;
00942
00943
00944 computeBasis( ) ;
00945
00946
00947
00948 int count = 0 ;
00949 volume = 0.0 ;
00950
00951 for ( i = 0; i < 2; i++ ) {
00952 for ( j = 0; j < 2; j++ ) {
00953 for ( k = 0; k < 2; k++ ) {
00954
00955 gaussPoint[0] = sg[i] ;
00956 gaussPoint[1] = sg[j] ;
00957 gaussPoint[2] = sg[k] ;
00958
00959
00960 shp3d( gaussPoint, xsj, shp, xl ) ;
00961
00962
00963 for ( p = 0; p < nShape; p++ ) {
00964 for ( q = 0; q < numberNodes; q++ )
00965 Shape[p][q][count] = shp[p][q] ;
00966 }
00967
00968
00969
00970 dvol[count] = wg[count] * xsj ;
00971
00972
00973
00974 count++ ;
00975
00976 }
00977 }
00978 }
00979
00980
00981
00982 for ( i = 0; i < numberGauss; i++ ) {
00983
00984
00985 for ( p = 0; p < nShape; p++ ) {
00986 for ( q = 0; q < numberNodes; q++ )
00987 shp[p][q] = Shape[p][q][i] ;
00988 }
00989
00990
00991
00992 stress = materialPointers[i]->getStress( ) ;
00993
00994
00995
00996 stress *= dvol[i] ;
00997
00998 if ( tang_flag == 1 ) {
00999 dd = materialPointers[i]->getTangent( ) ;
01000 dd *= dvol[i] ;
01001 }
01002
01003
01004 double stress0 = stress(0);
01005 double stress1 = stress(1);
01006 double stress2 = stress(2);
01007 double stress3 = stress(3);
01008 double stress4 = stress(4);
01009 double stress5 = stress(5);
01010
01011
01012
01013 int jj = 0 ;
01014 for ( j = 0; j < numberNodes; j++ ) {
01015
01016
01017
01018
01019
01020
01021
01022
01023
01024
01025
01026
01027
01028
01029
01030
01031
01032
01033
01034
01035
01036 double b00 = shp[0][j];
01037 double b11 = shp[1][j];
01038 double b22 = shp[2][j];
01039 double b30 = shp[1][j];
01040 double b31 = shp[0][j];
01041 double b41 = shp[2][j];
01042 double b42 = shp[1][j];
01043 double b50 = shp[2][j];
01044 double b52 = shp[0][j];
01045
01046 residJ(0) = b00 * stress0 + b30 * stress3 + b50 * stress5;
01047 residJ(1) = b11 * stress1 + b31 * stress3 + b41 * stress4;
01048 residJ(2) = b22 * stress2 + b42 * stress4 + b52 * stress5;
01049
01050 BJ = computeB( j, shp ) ;
01051
01052
01053
01054 for (p=0; p<ndf; p++) {
01055 for (q=0; q<nstress; q++)
01056 BJtran(p,q) = BJ(q,p) ;
01057 }
01058
01059
01060
01061 for ( p = 0; p < ndf; p++ ) {
01062 resid( jj + p ) += residJ(p) ;
01063 resid( jj + p ) -= dvol[i]*b[p]*shp[3][j];
01064 }
01065
01066 if ( tang_flag == 1 ) {
01067
01068
01069 BJtranD.addMatrixProduct(0.0, BJtran,dd,1.0) ;
01070
01071 int kk = 0 ;
01072 for ( k = 0; k < numberNodes; k++ ) {
01073
01074 BK = computeB( k, shp ) ;
01075
01076
01077
01078 stiffJK.addMatrixProduct(0.0, BJtranD,BK,1.0) ;
01079
01080 for ( p = 0; p < ndf; p++ ) {
01081 for ( q = 0; q < ndf; q++ )
01082 stiff( jj+p, kk+q ) += stiffJK( p, q ) ;
01083 }
01084
01085 kk += ndf ;
01086 }
01087
01088 }
01089
01090 jj += ndf ;
01091 }
01092
01093
01094 }
01095
01096
01097 return ;
01098 }
01099
01100
01101
01102
01103
01104 void Brick::computeBasis( )
01105 {
01106
01107
01108
01109 int i ;
01110 for ( i = 0; i < 8; i++ ) {
01111
01112 const Vector &coorI = nodePointers[i]->getCrds( ) ;
01113
01114 xl[0][i] = coorI(0) ;
01115 xl[1][i] = coorI(1) ;
01116 xl[2][i] = coorI(2) ;
01117
01118 }
01119
01120 }
01121
01122
01123
01124
01125 const Matrix&
01126 Brick::computeB( int node, const double shp[4][8] )
01127 {
01128
01129
01130
01131
01132
01133
01134
01135
01136
01137
01138
01139
01140
01141
01142 B(0,0) = shp[0][node] ;
01143 B(1,1) = shp[1][node] ;
01144 B(2,2) = shp[2][node] ;
01145
01146 B(3,0) = shp[1][node] ;
01147 B(3,1) = shp[0][node] ;
01148
01149 B(4,1) = shp[2][node] ;
01150 B(4,2) = shp[1][node] ;
01151
01152 B(5,0) = shp[2][node] ;
01153 B(5,2) = shp[0][node] ;
01154
01155 return B ;
01156
01157 }
01158
01159
01160
01161 Matrix Brick::transpose( int dim1,
01162 int dim2,
01163 const Matrix &M )
01164 {
01165 int i ;
01166 int j ;
01167
01168 Matrix Mtran( dim2, dim1 ) ;
01169
01170 for ( i = 0; i < dim1; i++ ) {
01171 for ( j = 0; j < dim2; j++ )
01172 Mtran(j,i) = M(i,j) ;
01173 }
01174
01175 return Mtran ;
01176 }
01177
01178
01179
01180 int Brick::sendSelf (int commitTag, Channel &theChannel)
01181 {
01182 int res = 0;
01183
01184
01185
01186
01187 int dataTag = this->getDbTag();
01188
01189
01190
01191
01192
01193 int matDbTag;
01194
01195 static ID idData(25);
01196
01197 idData(24) = this->getTag();
01198
01199 int i;
01200 for (i = 0; i < 8; i++) {
01201 idData(i) = materialPointers[i]->getClassTag();
01202 matDbTag = materialPointers[i]->getDbTag();
01203
01204
01205 if (matDbTag == 0) {
01206 matDbTag = theChannel.getDbTag();
01207 if (matDbTag != 0)
01208 materialPointers[i]->setDbTag(matDbTag);
01209 }
01210 idData(i+8) = matDbTag;
01211 }
01212
01213 idData(16) = connectedExternalNodes(0);
01214 idData(17) = connectedExternalNodes(1);
01215 idData(18) = connectedExternalNodes(2);
01216 idData(19) = connectedExternalNodes(3);
01217 idData(20) = connectedExternalNodes(4);
01218 idData(21) = connectedExternalNodes(5);
01219 idData(22) = connectedExternalNodes(6);
01220 idData(23) = connectedExternalNodes(7);
01221
01222 res += theChannel.sendID(dataTag, commitTag, idData);
01223 if (res < 0) {
01224 opserr << "WARNING Brick::sendSelf() - " << this->getTag() << " failed to send ID\n";
01225 return res;
01226 }
01227
01228
01229
01230 for (i = 0; i < 8; i++) {
01231 res += materialPointers[i]->sendSelf(commitTag, theChannel);
01232 if (res < 0) {
01233 opserr << "WARNING Brick::sendSelf() - " << this->getTag() << " failed to send its Material\n";
01234 return res;
01235 }
01236 }
01237
01238 return res;
01239
01240 }
01241
01242 int Brick::recvSelf (int commitTag,
01243 Channel &theChannel,
01244 FEM_ObjectBroker &theBroker)
01245 {
01246 int res = 0;
01247
01248 int dataTag = this->getDbTag();
01249
01250 static ID idData(25);
01251 res += theChannel.recvID(dataTag, commitTag, idData);
01252 if (res < 0) {
01253 opserr << "WARNING Brick::recvSelf() - " << this->getTag() << " failed to receive ID\n";
01254 return res;
01255 }
01256
01257 this->setTag(idData(24));
01258
01259 connectedExternalNodes(0) = idData(16);
01260 connectedExternalNodes(1) = idData(17);
01261 connectedExternalNodes(2) = idData(18);
01262 connectedExternalNodes(3) = idData(19);
01263 connectedExternalNodes(4) = idData(20);
01264 connectedExternalNodes(5) = idData(21);
01265 connectedExternalNodes(6) = idData(22);
01266 connectedExternalNodes(7) = idData(23);
01267
01268
01269 if (materialPointers[0] == 0) {
01270 for (int i = 0; i < 8; i++) {
01271 int matClassTag = idData(i);
01272 int matDbTag = idData(i+8);
01273
01274 materialPointers[i] = theBroker.getNewNDMaterial(matClassTag);
01275 if (materialPointers[i] == 0) {
01276 opserr << "Brick::recvSelf() - Broker could not create NDMaterial of class type " << matClassTag << endln;
01277 return -1;
01278 }
01279
01280 materialPointers[i]->setDbTag(matDbTag);
01281 res += materialPointers[i]->recvSelf(commitTag, theChannel, theBroker);
01282 if (res < 0) {
01283 opserr << "NLBeamColumn3d::recvSelf() - material " << i << "failed to recv itself\n";
01284 return res;
01285 }
01286 }
01287 }
01288
01289 else {
01290 for (int i = 0; i < 8; i++) {
01291 int matClassTag = idData(i);
01292 int matDbTag = idData(i+8);
01293
01294
01295 if (materialPointers[i]->getClassTag() != matClassTag) {
01296 delete materialPointers[i];
01297 materialPointers[i] = theBroker.getNewNDMaterial(matClassTag);
01298 if (materialPointers[i] == 0) {
01299 opserr << "Brick::recvSelf() - Broker could not create NDMaterial of class type " <<
01300 matClassTag << endln;
01301 exit(-1);
01302 }
01303 materialPointers[i]->setDbTag(matDbTag);
01304 }
01305
01306
01307 res += materialPointers[i]->recvSelf(commitTag, theChannel, theBroker);
01308 if (res < 0) {
01309 opserr << "Brick::recvSelf() - material " << i << "failed to recv itself\n";
01310 return res;
01311 }
01312 }
01313 }
01314
01315 return res;
01316 }
01317
01318
01319 int
01320 Brick::displaySelf(Renderer &theViewer, int displayMode, float fact)
01321 {
01322
01323 const Vector &end1Crd = nodePointers[0]->getCrds();
01324 const Vector &end2Crd = nodePointers[1]->getCrds();
01325 const Vector &end3Crd = nodePointers[2]->getCrds();
01326 const Vector &end4Crd = nodePointers[3]->getCrds();
01327
01328 const Vector &end5Crd = nodePointers[4]->getCrds();
01329 const Vector &end6Crd = nodePointers[5]->getCrds();
01330 const Vector &end7Crd = nodePointers[6]->getCrds();
01331 const Vector &end8Crd = nodePointers[7]->getCrds();
01332
01333 static Matrix coords(4,3);
01334 static Vector values(4);
01335 static Vector P(24) ;
01336
01337 values(0) = 1 ;
01338 values(1) = 1 ;
01339 values(2) = 1 ;
01340 values(3) = 1 ;
01341
01342 int error = 0;
01343 int i;
01344
01345
01346 if (displayMode >= 0) {
01347
01348 const Vector &end1Disp = nodePointers[0]->getDisp();
01349 const Vector &end2Disp = nodePointers[1]->getDisp();
01350 const Vector &end3Disp = nodePointers[2]->getDisp();
01351 const Vector &end4Disp = nodePointers[3]->getDisp();
01352 const Vector &end5Disp = nodePointers[4]->getDisp();
01353 const Vector &end6Disp = nodePointers[5]->getDisp();
01354 const Vector &end7Disp = nodePointers[6]->getDisp();
01355 const Vector &end8Disp = nodePointers[7]->getDisp();
01356
01357 const Vector &stress1 = materialPointers[0]->getStress();
01358 const Vector &stress2 = materialPointers[1]->getStress();
01359 const Vector &stress3 = materialPointers[2]->getStress();
01360 const Vector &stress4 = materialPointers[3]->getStress();
01361 const Vector &stress5 = materialPointers[4]->getStress();
01362 const Vector &stress6 = materialPointers[5]->getStress();
01363 const Vector &stress7 = materialPointers[6]->getStress();
01364 const Vector &stress8 = materialPointers[7]->getStress();
01365
01366
01367
01368
01369
01370
01371
01372 for (i = 0; i < 3; i++) {
01373 coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
01374 coords(1,i) = end2Crd(i) + end2Disp(i)*fact;
01375 coords(2,i) = end3Crd(i) + end3Disp(i)*fact;
01376 coords(3,i) = end4Crd(i) + end4Disp(i)*fact;
01377 }
01378
01379 if (displayMode < 3 && displayMode > 0) {
01380 int index = displayMode - 1;
01381 values(0) = stress1(index);
01382 values(1) = stress2(index);
01383 values(2) = stress3(index);
01384 values(3) = stress4(index);
01385 }
01386
01387 error += theViewer.drawPolygon (coords, values);
01388
01389 for (i = 0; i < 3; i++) {
01390 coords(0,i) = end5Crd(i) + end5Disp(i)*fact;
01391 coords(1,i) = end6Crd(i) + end6Disp(i)*fact;
01392 coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
01393 coords(3,i) = end8Crd(i) + end8Disp(i)*fact;
01394 }
01395
01396 if (displayMode < 3 && displayMode > 0) {
01397 int index = displayMode - 1;
01398 values(0) = stress5(index);
01399 values(1) = stress6(index);
01400 values(2) = stress7(index);
01401 values(3) = stress8(index);
01402 }
01403
01404 error += theViewer.drawPolygon (coords, values);
01405
01406 for (i = 0; i < 3; i++) {
01407 coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
01408 coords(1,i) = end4Crd(i) + end4Disp(i)*fact;
01409 coords(2,i) = end8Crd(i) + end8Disp(i)*fact;
01410 coords(3,i) = end5Crd(i) + end5Disp(i)*fact;
01411 }
01412
01413 if (displayMode < 3 && displayMode > 0) {
01414 int index = displayMode - 1;
01415 values(0) = stress1(index);
01416 values(1) = stress4(index);
01417 values(2) = stress8(index);
01418 values(3) = stress5(index);
01419 }
01420
01421 error += theViewer.drawPolygon (coords, values);
01422
01423 for (i = 0; i < 3; i++) {
01424 coords(0,i) = end2Crd(i) + end2Disp(i)*fact;
01425 coords(1,i) = end3Crd(i) + end3Disp(i)*fact;
01426 coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
01427 coords(3,i) = end6Crd(i) + end6Disp(i)*fact;
01428 }
01429 if (displayMode < 3 && displayMode > 0) {
01430 int index = displayMode - 1;
01431 values(0) = stress2(index);
01432 values(1) = stress3(index);
01433 values(2) = stress7(index);
01434 values(3) = stress6(index);
01435 }
01436
01437 error += theViewer.drawPolygon (coords, values);
01438
01439
01440 for (i = 0; i < 3; i++) {
01441 coords(0,i) = end1Crd(i) + end1Disp(i)*fact;
01442 coords(1,i) = end2Crd(i) + end2Disp(i)*fact;
01443 coords(2,i) = end6Crd(i) + end6Disp(i)*fact;
01444 coords(3,i) = end5Crd(i) + end5Disp(i)*fact;
01445 }
01446
01447 if (displayMode < 3 && displayMode > 0) {
01448 int index = displayMode - 1;
01449 values(0) = stress1(index);
01450 values(1) = stress2(index);
01451 values(2) = stress6(index);
01452 values(3) = stress5(index);
01453 }
01454
01455 error += theViewer.drawPolygon (coords, values);
01456
01457 for (i = 0; i < 3; i++) {
01458 coords(0,i) = end4Crd(i) + end4Disp(i)*fact;
01459 coords(1,i) = end3Crd(i) + end3Disp(i)*fact;
01460 coords(2,i) = end7Crd(i) + end7Disp(i)*fact;
01461 coords(3,i) = end8Crd(i) + end8Disp(i)*fact;
01462 }
01463
01464 if (displayMode < 3 && displayMode > 0) {
01465 int index = displayMode - 1;
01466 values(0) = stress3(index);
01467 values(1) = stress4(index);
01468 values(2) = stress7(index);
01469 values(3) = stress8(index);
01470 }
01471
01472 error += theViewer.drawPolygon (coords, values);
01473
01474 } else {
01475
01476 int mode = displayMode * -1;
01477
01478 const Matrix &eigen1 = nodePointers[0]->getEigenvectors();
01479 const Matrix &eigen2 = nodePointers[1]->getEigenvectors();
01480 const Matrix &eigen3 = nodePointers[2]->getEigenvectors();
01481 const Matrix &eigen4 = nodePointers[3]->getEigenvectors();
01482 const Matrix &eigen5 = nodePointers[4]->getEigenvectors();
01483 const Matrix &eigen6 = nodePointers[5]->getEigenvectors();
01484 const Matrix &eigen7 = nodePointers[6]->getEigenvectors();
01485 const Matrix &eigen8 = nodePointers[7]->getEigenvectors();
01486
01487 if (eigen1.noCols() >= mode) {
01488
01489 for (i = 0; i < 3; i++) {
01490 coords(0,i) = end1Crd(i) + eigen1(i,mode-1)*fact;
01491 coords(1,i) = end2Crd(i) + eigen2(i,mode-1)*fact;
01492 coords(2,i) = end3Crd(i) + eigen3(i,mode-1)*fact;
01493 coords(3,i) = end4Crd(i) + eigen4(i,mode-1)*fact;
01494 }
01495
01496 error += theViewer.drawPolygon (coords, values);
01497
01498 for (i = 0; i < 3; i++) {
01499 coords(0,i) = end5Crd(i) + eigen5(i,mode-1)*fact;
01500 coords(1,i) = end6Crd(i) + eigen6(i,mode-1)*fact;
01501 coords(2,i) = end7Crd(i) + eigen7(i,mode-1)*fact;
01502 coords(3,i) = end8Crd(i) + eigen8(i,mode-1)*fact;
01503 }
01504
01505 error += theViewer.drawPolygon (coords, values);
01506
01507 for (i = 0; i < 3; i++) {
01508 coords(0,i) = end1Crd(i) + eigen1(i,mode-1)*fact;
01509 coords(1,i) = end4Crd(i) + eigen4(i,mode-1)*fact;
01510 coords(2,i) = end8Crd(i) + eigen8(i,mode-1)*fact;
01511 coords(3,i) = end5Crd(i) + eigen5(i,mode-1)*fact;
01512 }
01513
01514 error += theViewer.drawPolygon (coords, values);
01515
01516 for (i = 0; i < 3; i++) {
01517 coords(0,i) = end2Crd(i) + eigen2(i,mode-1)*fact;
01518 coords(1,i) = end3Crd(i) + eigen3(i,mode-1)*fact;
01519 coords(2,i) = end7Crd(i) + eigen7(i,mode-1)*fact;
01520 coords(3,i) = end6Crd(i) + eigen6(i,mode-1)*fact;
01521 }
01522
01523 error += theViewer.drawPolygon (coords, values);
01524
01525
01526 for (i = 0; i < 3; i++) {
01527 coords(0,i) = end1Crd(i) + eigen1(i,mode-1)*fact;
01528 coords(1,i) = end2Crd(i) + eigen2(i,mode-1)*fact;
01529 coords(2,i) = end6Crd(i) + eigen6(i,mode-1)*fact;
01530 coords(3,i) = end5Crd(i) + eigen5(i,mode-1)*fact;
01531 }
01532
01533 error += theViewer.drawPolygon (coords, values);
01534
01535 for (i = 0; i < 3; i++) {
01536 coords(0,i) = end4Crd(i) + eigen4(i,mode-1)*fact;
01537 coords(1,i) = end3Crd(i) + eigen3(i,mode-1)*fact;
01538 coords(2,i) = end7Crd(i) + eigen7(i,mode-1)*fact;
01539 coords(3,i) = end8Crd(i) + eigen8(i,mode-1)*fact;
01540 }
01541
01542 error += theViewer.drawPolygon (coords, values);
01543 } else {
01544 values.Zero();
01545 for (i = 0; i < 3; i++) {
01546 coords(0,i) = end1Crd(i);
01547 coords(1,i) = end2Crd(i);
01548 coords(2,i) = end3Crd(i);
01549 coords(3,i) = end4Crd(i);
01550 }
01551
01552 error += theViewer.drawPolygon (coords, values);
01553
01554 for (i = 0; i < 3; i++) {
01555 coords(0,i) = end5Crd(i);
01556 coords(1,i) = end6Crd(i);
01557 coords(2,i) = end7Crd(i);
01558 coords(3,i) = end8Crd(i);
01559 }
01560
01561 error += theViewer.drawPolygon (coords, values);
01562
01563 for (i = 0; i < 3; i++) {
01564 coords(0,i) = end1Crd(i);
01565 coords(1,i) = end4Crd(i);
01566 coords(2,i) = end8Crd(i);
01567 coords(3,i) = end5Crd(i);
01568 }
01569
01570 error += theViewer.drawPolygon (coords, values);
01571
01572 for (i = 0; i < 3; i++) {
01573 coords(0,i) = end2Crd(i);
01574 coords(1,i) = end3Crd(i);
01575 coords(2,i) = end7Crd(i);
01576 coords(3,i) = end6Crd(i);
01577 }
01578
01579 error += theViewer.drawPolygon (coords, values);
01580
01581
01582 for (i = 0; i < 3; i++) {
01583 coords(0,i) = end1Crd(i);
01584 coords(1,i) = end2Crd(i);
01585 coords(2,i) = end6Crd(i);
01586 coords(3,i) = end5Crd(i);
01587 }
01588
01589 error += theViewer.drawPolygon (coords, values);
01590
01591 for (i = 0; i < 3; i++) {
01592 coords(0,i) = end4Crd(i);
01593 coords(1,i) = end3Crd(i);
01594 coords(2,i) = end7Crd(i);
01595 coords(3,i) = end8Crd(i);
01596 }
01597
01598 error += theViewer.drawPolygon (coords, values);
01599 }
01600 }
01601
01602
01603 return error;
01604 }
01605
01606 Response*
01607 Brick::setResponse(const char **argv, int argc, Information &eleInfo, OPS_Stream &output)
01608 {
01609 Response *theResponse = 0;
01610
01611 char outputData[32];
01612
01613 output.tag("ElementOutput");
01614 output.attr("eleType","Brick");
01615 output.attr("eleTag",this->getTag());
01616 for (int i=1; i<=8; i++) {
01617 sprintf(outputData,"node%d",i);
01618 output.attr(outputData,nodePointers[i-1]->getTag());
01619 }
01620
01621
01622 if (strcmp(argv[0],"force") == 0 || strcmp(argv[0],"forces") == 0) {
01623
01624 for (int i=1; i<=8; i++) {
01625 sprintf(outputData,"P1_%d",i);
01626 output.tag("ResponseType",outputData);
01627 sprintf(outputData,"P2_%d",i);
01628 output.tag("ResponseType",outputData);
01629 sprintf(outputData,"P3_%d",i);
01630 output.tag("ResponseType",outputData);
01631 }
01632
01633 theResponse = new ElementResponse(this, 1, resid);
01634
01635 } else if (strcmp(argv[0],"material") == 0 || strcmp(argv[0],"integrPoint") == 0) {
01636
01637 int pointNum = atoi(argv[1]);
01638 if (pointNum > 0 && pointNum <= 8) {
01639
01640 output.tag("GaussPoint");
01641 output.attr("number",pointNum);
01642
01643 theResponse = materialPointers[pointNum-1]->setResponse(&argv[2], argc-2, eleInfo, output);
01644
01645 output.endTag();
01646 }
01647
01648
01649 } else if (strcmp(argv[0],"stresses") ==0) {
01650
01651 for (int i=0; i<8; i++) {
01652 output.tag("GaussPoint");
01653 output.attr("number",i+1);
01654 output.tag("NdMaterialOutput");
01655 output.attr("classType", materialPointers[i]->getClassTag());
01656 output.attr("tag", materialPointers[i]->getTag());
01657
01658 output.tag("ResponseType","sigma11");
01659 output.tag("ResponseType","sigma22");
01660 output.tag("ResponseType","sigma33");
01661 output.tag("ResponseType","sigma12");
01662 output.tag("ResponseType","sigma13");
01663 output.tag("ResponseType","sigma23");
01664
01665 output.endTag();
01666 output.endTag();
01667 }
01668 theResponse = new ElementResponse(this, 3, Vector(48));
01669 }
01670
01671 output.endTag();
01672 return theResponse;
01673 }
01674
01675 int
01676 Brick::getResponse(int responseID, Information &eleInfo)
01677 {
01678 static Vector stresses(48);
01679
01680 if (responseID == 1)
01681 return eleInfo.setVector(this->getResistingForce());
01682
01683 else if (responseID == 2)
01684 return eleInfo.setMatrix(this->getTangentStiff());
01685
01686 else if (responseID == 3) {
01687
01688
01689 int cnt = 0;
01690 for (int i = 0; i < 8; i++) {
01691
01692
01693 const Vector &sigma = materialPointers[i]->getStress();
01694 stresses(cnt++) = sigma(0);
01695 stresses(cnt++) = sigma(1);
01696 stresses(cnt++) = sigma(2);
01697 stresses(cnt++) = sigma(3);
01698 stresses(cnt++) = sigma(4);
01699 stresses(cnt++) = sigma(5);
01700 }
01701 return eleInfo.setVector(stresses);
01702
01703 }
01704 else
01705
01706 return -1;
01707 }