xc
ShellMITC9.h
1 /* ****************************************************************** **
2 ** OpenSees - Open System for Earthquake Engineering Simulation **
3 ** Pacific Earthquake Engineering Research Center **
4 ** **
5 ** **
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 // Written: Leopoldo Tesser, Diego Talledo
22 // 9-node lagrangian shell element with membrane and drill
23 //
24 #ifndef ShellMITC9_h
25 #define ShellMITC9_h
26 
27 #include "domain/mesh/element/plane/QuadBase9N.h"
28 #include "domain/mesh/element/utils/physical_properties/SectionFDPhysicalProperties.h"
29 #include "domain/mesh/element/utils/coordTransformation/ShellLinearCrdTransf3d.h"
30 #include "domain/mesh/element/utils/fvectors/FVectorShell.h"
31 
32 namespace XC {
33 
35 //
37 class ShellMITC9: public QuadBase9N<SectionFDPhysicalProperties>
38  {
39  private :
40  double Ktt;
41  ShellLinearCrdTransf3d theCoordTransf;
42  mutable Matrix Ki;
43 
44  double xl[2][9];
45 
46  FVectorShell p0;
47 
48  //static data
49  static Matrix stiff;
50  static Vector resid;
51  static Matrix mass;
52  static Matrix damping;
53 
54 
55  void computeBasis(void);
56  void formInertiaTerms(int tangFlag) const;
57  void formResidAndTangent(int tang_flag) const;
58 
59  //compute Jacobian matrix and inverse at point {L1,L2}
60  //void computeJacobian( double L1, double L2,const double x[2][9],
61  // Matrix &JJ,Matrix &JJinv );
62 
63  double* computeBdrill(int node,const double shp[3][9]) const;
64  const Matrix& assembleB(const Matrix &Bmembrane, const Matrix &Bbend, const Matrix &Bshear) const;
65  const Matrix& computeBmembrane(int node, const double shp[3][9]) const;
66  const Matrix& computeBbend(int node, const double shp[3][9]) const;
67  const Matrix& computeBshear(int node, const double shp[3][9]) const;
68 
69  //Matrix transpose
70  Matrix transpose( int dim1, int dim2, const Matrix &M);
71 
72  //shape function routine for four node quads
73  static void shape2d( double ss, double tt, const double x[2][9], double shp[3][9], double &xsj );
74  protected:
75  DbTagData &getDbTagData(void) const;
76  int sendData(Communicator &);
77  int recvData(const Communicator &);
78 
79  public:
80  //null constructor
81  ShellMITC9(void);
82  //full constructor
83  ShellMITC9(int tag,const SectionForceDeformation *theMaterial);
84  Element *getCopy(void) const;
85 
86  int getNumDOF(void) const;
87 
88  void setDomain(Domain *theDomain);
89 
90  //print out element data
91  void Print(std::ostream &, int flag) const;
92 
93  //return stiffness matrix
94  const Matrix &getTangentStiff(void) const;
95  const Matrix &getInitialStiff() const;
96  const Matrix &getMass() const;
97 
98  const GaussModel &getGaussModel(void) const;
99 
100  void alive(void);
101  void zeroLoad(void);
102  // methods for applying loads
103  int addLoad(ElementalLoad *theLoad, double loadFactor);
104  int addInertiaLoadToUnbalance(const Vector &accel);
105 
106  //get residual
107  const Vector &getResistingForce(void) const;
108  const Vector &getResistingForceIncInertia(void) const;
109 
110  virtual ShellCrdTransf3dBase *getCoordTransf(void);
111  virtual const ShellCrdTransf3dBase *getCoordTransf(void) const;
112  double getArea(bool initialGeometry= true) const;
113 
114  int sendSelf(Communicator &);
115  int recvSelf(const Communicator &);
116 
117  };
118 } // end of XC namespace
119 
120 #endif
121 
122 
123 
Base class for force deformation section models.
Definition: SectionForceDeformation.h:88
void Print(std::ostream &, int flag) const
print out element data
Definition: ShellMITC9.cpp:1032
void alive(void)
Reactivates the element.
Definition: ShellMITC9.cpp:297
Float vector abstraction.
Definition: Vector.h:94
Lagrangian shell element with membrane and drill.
Definition: ShellMITC9.h:37
Communication parameters between processes.
Definition: Communicator.h:66
virtual ShellCrdTransf3dBase * getCoordTransf(void)
Returns a pointer a la coordinate transformation.
Definition: ShellMITC9.cpp:1024
const Matrix & getInitialStiff() const
return secant matrix
Definition: ShellMITC9.cpp:124
Vector that stores the dbTags of the class members.
Definition: DbTagData.h:44
double getArea(bool initialGeometry=true) const
Returns element area.
Definition: ShellMITC9.cpp:120
const GaussModel & getGaussModel(void) const
Return the Gauss points of the element.
Definition: ShellMITC9.cpp:66
int sendData(Communicator &)
Send members through the communicator argument.
Definition: ShellMITC9.cpp:998
Base class for the finite elements.
Definition: Element.h:112
Base class for 3D coordinate transformations.
Definition: ShellCrdTransf3dBase.h:49
Base class for nine node quads.
Definition: QuadBase9N.h:45
int recvData(const Communicator &)
Receives members through the communicator argument.
Definition: ShellMITC9.cpp:1011
Ingernal forces for a shell element.
Definition: FVectorShell.h:41
void setDomain(Domain *theDomain)
set domain
Definition: ShellMITC9.cpp:78
int sendSelf(Communicator &)
Send the object.
Definition: ShellMITC9.cpp:1085
const Vector & getResistingForceIncInertia(void) const
get residual with inertia terms
Definition: ShellMITC9.cpp:375
ShellMITC9(void)
null constructor
Definition: ShellMITC9.cpp:62
int recvSelf(const Communicator &)
Receive the object.
Definition: ShellMITC9.cpp:1098
const Matrix & getTangentStiff(void) const
return stiffness matrix
Definition: ShellMITC9.cpp:110
Element * getCopy(void) const
Virtual constructor.
Definition: ShellMITC9.cpp:74
Base class for Gauss integration models.
Definition: GaussModel.h:41
void zeroLoad(void)
Zeroes the element load vector.
Definition: ShellMITC9.cpp:309
Base class for small displacement 3D coordinate transformations.
Definition: ShellLinearCrdTransf3d.h:42
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
int getNumDOF(void) const
return number of dofs
Definition: ShellMITC9.cpp:106
Matrix of floats.
Definition: Matrix.h:111
const Matrix & getMass() const
return mass matrix
Definition: ShellMITC9.cpp:287
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
DbTagData & getDbTagData(void) const
Returns a vector to store the dbTags of the class members.
Definition: ShellMITC9.cpp:991
const Vector & getResistingForce(void) const
get residual
Definition: ShellMITC9.cpp:361