XC Open source finite element analysis program
ShellNL.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 ShellNL_h
25 #define ShellNL_h
26 
27 #include "domain/mesh/element/plane/QuadBase9N.h"
28 #include "domain/mesh/element/utils/physical_properties/SectionFDPhysicalProperties.h"
29 #include "ShellLinearCrdTransf3d.h"
30 #include "domain/mesh/element/utils/fvectors/FVectorShell.h"
31 
32 namespace XC {
33 
35 //
37 class ShellNL : public QuadBase9N<SectionFDPhysicalProperties>
38  {
39  private :
40  double Ktt;
41  ShellLinearCrdTransf3d theCoordTransf;
42  mutable Matrix *Ki;
43  //local nodal coordinates, two coordinates for each of nine nodes
44  double xl[2][9];
45  //shell basis vectors
46 
47  FVectorShell p0; // Reactions in the basic system due to element loads
48 
49  //static data
50  static Matrix stiff;
51  static Vector resid;
52  static Matrix mass;
53  static Matrix damping;
54 
55  //quadrature data
56  static const double root3;
57  static const double root3_over_root5;
58 
59 
60 
61  void computeBasis(void);
62  void formInertiaTerms(int tangFlag) const;
63  void formResidAndTangent(int tang_flag) const;
64 
65  //compute Jacobian matrix and inverse at point {L1,L2}
66  //void computeJacobian( double L1, double L2,const double x[2][9],
67  // Matrix &JJ,Matrix &JJinv );
68 
69  double* computeBdrill(int node,const double shp[3][9]) const;
70  const Matrix& assembleB(const Matrix &Bmembrane, const Matrix &Bbend, const Matrix &Bshear) const;
71  const Matrix& computeBmembrane(int node, const double shp[3][9]) const;
72  const Matrix& computeBbend(int node, const double shp[3][9]) const;
73  const Matrix& computeBshear(int node, const double shp[3][9]) const;
74 
75  //Matrix transpose
76  Matrix transpose( int dim1, int dim2, const Matrix &M);
77 
78  //shape function routine for four node quads
79  static void shape2d( double ss, double tt, const double x[2][9], double shp[3][9], double &xsj );
80  protected:
81  DbTagData &getDbTagData(void) const;
82  int sendData(CommParameters &);
83  int recvData(const CommParameters &);
84 
85  public:
86  //null constructor
87  ShellNL(void);
88  //full constructor
89  ShellNL(int tag,const SectionForceDeformation *theMaterial);
90  Element *getCopy(void) const;
91  virtual ~ShellNL(void);
92 
93  int getNumDOF(void) const;
94 
95  void setDomain(Domain *theDomain);
96 
97  int revertToLastCommit(void);
98  int revertToStart(void);
99 
100  //print out element data
101  void Print(std::ostream &, int flag);
102 
103  //return stiffness matrix
104  const Matrix &getTangentStiff(void) const;
105  const Matrix &getInitialStiff() const;
106  const Matrix &getMass() const;
107 
108  const GaussModel &getGaussModel(void) const;
109 
110  // methods for applying loads
111  void zeroLoad(void);
112  int addLoad(ElementalLoad *theLoad, double loadFactor);
113  int addInertiaLoadToUnbalance(const Vector &accel);
114 
115  //get residual
116  const Vector &getResistingForce(void) const;
117  const Vector &getResistingForceIncInertia(void) const;
118 
119  virtual ShellCrdTransf3dBase *getCoordTransf(void);
120  virtual const ShellCrdTransf3dBase *getCoordTransf(void) const;
121 
122  int sendSelf(CommParameters &);
123  int recvSelf(const CommParameters &);
124 
125  };
126 } // end of XC namespace
127 
128 #endif
129 
130 
131 
DbTagData & getDbTagData(void) const
Returns a vector para almacenar los dbTags de los miembros de la clase.
Definition: ShellNL.cpp:996
const Matrix & getMass() const
return mass matrix
Definition: ShellNL.cpp:302
const Matrix & getInitialStiff() const
return secant matrix
Definition: ShellNL.cpp:139
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
const GaussModel & getGaussModel(void) const
Returns the puntos de Gauss of the element.
Definition: ShellNL.cpp:67
Base class for Gauss integration models.
Definition: GaussModel.h:40
int getNumDOF(void) const
return number of dofs
Definition: ShellNL.cpp:111
Definition: Vector.h:82
Vector que almacena los dbTags de los miembros de la clase.
Definition: DbTagData.h:43
virtual ~ShellNL(void)
destructor
Definition: ShellNL.cpp:79
Lagrangian shell element with membrane and drill.
Definition: ShellNL.h:37
virtual ShellCrdTransf3dBase * getCoordTransf(void)
Returns a pointer a la coordinate transformation.
Definition: ShellNL.cpp:1029
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
const Matrix & getTangentStiff(void) const
return stiffness matrix
Definition: ShellNL.cpp:128
Ingernal forces for a shell element.
Definition: FVectorShell.h:40
Element * getCopy(void) const
Virtual constructor.
Definition: ShellNL.cpp:75
void Print(std::ostream &, int flag)
print out element data
Definition: ShellNL.cpp:1037
Definition: Matrix.h:82
Base class for nine node quads.
Definition: QuadBase9N.h:49
Base class for small displacement 3D coordinate transformations.
Definition: ShellLinearCrdTransf3d.h:41
const Vector & getResistingForceIncInertia(void) const
get residual with inertia terms
Definition: ShellNL.cpp:377
void setDomain(Domain *theDomain)
set domain
Definition: ShellNL.cpp:83
void zeroLoad(void)
Zeroes the element load vector.
Definition: ShellNL.cpp:312
const Vector & getResistingForce(void) const
get residual
Definition: ShellNL.cpp:363
Base class for force deformation section models. Constitutive equations of the section.
Definition: SectionForceDeformation.h:86
int revertToLastCommit(void)
revert to last commit
Definition: ShellNL.cpp:115
int revertToStart(void)
revert to start
Definition: ShellNL.cpp:121
Communication parameters between processes.
Definition: CommParameters.h:65
ShellNL(void)
null constructor
Definition: ShellNL.cpp:62
int sendData(CommParameters &)
Send members through the channel being passed as parameter.
Definition: ShellNL.cpp:1003
int recvData(const CommParameters &)
Receives members through the channel being passed as parameter.
Definition: ShellNL.cpp:1016
================================================================================
Definition: ContinuaReprComponent.h:34
Base class for 3D coordinate transformations.
Definition: ShellCrdTransf3dBase.h:48