XC Open source finite element analysis program
EnhancedQuad.h
1 //----------------------------------------------------------------------------
2 // XC program; finite element analysis code
3 // for structural analysis and design.
4 //
5 // Copyright (C) Luis Claudio Pérez Tato
6 //
7 // This program derives from OpenSees <http://opensees.berkeley.edu>
8 // developed by the «Pacific earthquake engineering research center».
9 //
10 // Except for the restrictions that may arise from the copyright
11 // of the original program (see copyright_opensees.txt)
12 // XC is free software: you can redistribute it and/or modify
13 // it under the terms of the GNU General Public License as published by
14 // the Free Software Foundation, either version 3 of the License, or
15 // (at your option) any later version.
16 //
17 // This software is distributed in the hope that it will be useful, but
18 // WITHOUT ANY WARRANTY; without even the implied warranty of
19 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 // GNU General Public License for more details.
21 //
22 //
23 // You should have received a copy of the GNU General Public License
24 // along with this program.
25 // If not, see <http://www.gnu.org/licenses/>.
26 //----------------------------------------------------------------------------
27 /* ****************************************************************** **
28 ** OpenSees - Open System for Earthquake Engineering Simulation **
29 ** Pacific Earthquake Engineering Research Center **
30 ** **
31 ** **
32 ** (C) Copyright 1999, The Regents of the University of California **
33 ** All Rights Reserved. **
34 ** **
35 ** Commercial use of this program without express permission of the **
36 ** University of California, Berkeley, is strictly prohibited. See **
37 ** file 'COPYRIGHT' in main directory for information on usage and **
38 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
39 ** **
40 ** Developed by: **
41 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
42 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
43 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
44 ** **
45 ** ****************************************************************** */
46 
47 // $Revision: 1.8 $
48 // $Date: 2003/02/14 23:01:10 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/fourNodeQuad/EnhancedQuad.h,v $
50 
51 #include "domain/mesh/element/plane/QuadBase4N.h"
52 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h"
53 #include <utility/matrix/Vector.h>
54 #include <utility/matrix/Matrix.h>
55 
56 namespace XC {
58 //
61 class EnhancedQuad : public QuadBase4N<NDMaterialPhysicalProperties>
62  {
63  private:
64 
65  mutable Vector alpha;
66  mutable Matrix *Ki;
67 
68  //static data
69  static Matrix stiff;
70  static Vector resid;
71  static Matrix mass;
72  static Matrix damping;
73 
74  //quadrature data
75  static const double sg[4];
76  static const double tg[4];
77  static const double wg[4];
78 
79 
80  static double stressData[][4];
81  static double tangentData[][3][4];
82 
83 
84  //local nodal coordinates, two coordinates for each of four nodes
85  // static double xl[2][4];
86  static double xl[][4];
87 
88  //save stress and tangent data
89  static void saveData(int gp, const Vector &stress,const Matrix &tangent);
90 
91  //recover stress and tangent data
92  void getData(int gp,Vector &stress,Matrix &tangent) const;
93 
94  //compute enhanced strain B-matrices
95  const Matrix& computeBenhanced( int node, double L1, double L2, double j, const Matrix &Jinv) const;
96 
97 
98  //compute local coordinates and basis
99  void computeBasis(void) const;
100 
101  //form residual and tangent
102  void formResidAndTangent(int tang_flag) const;
103 
104  //inertia terms
105  void formInertiaTerms(int tangFlag) const;
106 
107 
108  //compute Jacobian matrix and inverse at point {L1,L2}
109  static void computeJacobian( double L1, double L2, const double x[2][4], Matrix &JJ, Matrix &JJinv );
110 
111  //compute Bbend matrix
112  const Matrix& computeB(int node, const double shp[3][4]) const;
113 
114  //Matrix transpose of a 3x2 matrix
115  static const Matrix& transpose(const Matrix &M);
116 
117  //shape function routine for four node quads
118  static void shape2d( double ss, double tt, const double x[2][4], double shp[3][4], double &xsj );
119 
120 
121  bool check_material_type(const std::string &type) const;
122  protected:
123  int sendData(CommParameters &);
124  int recvData(const CommParameters &);
125  public:
126 
127  //full constructor
128  EnhancedQuad(int tag, int nd1, int nd2, int nd3, int nd4, NDMaterial &, const std::string &);
129 
130  //null constructor
131  EnhancedQuad(void);
132  Element *getCopy(void) const;
133  //destructor
134  virtual ~EnhancedQuad(void);
135 
136  //set domain
137  void setDomain(Domain *);
138 
139  //return number of dofs
140  int getNumDOF(void) const;
141 
142  // methods dealing with state updates
143  int update(void);
144 
145  //print out element data
146  void Print(std::ostream &s, int flag);
147 
148  //return stiffness matrix
149  const Matrix &getTangentStiff(void) const;
150  const Matrix &getInitialStiff(void) const;
151  const Matrix &getMass(void) const;
152 
153  int addLoad(ElementalLoad *theLoad, double loadFactor);
154  int addInertiaLoadToUnbalance(const Vector &accel);
155 
156  //get residual and residual with inertia terms
157  const Vector &getResistingForce(void) const;
158  const Vector &getResistingForceIncInertia(void) const;
159 
160  virtual int sendSelf(CommParameters &);
161  virtual int recvSelf(const CommParameters &);
162  };
163 } // end of XC namespace
void Print(std::ostream &s, int flag)
Imprime el objeto.
Definition: EnhancedQuad.cpp:141
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: EnhancedQuad.cpp:1328
int getNumDOF(void) const
return number of dofs
Definition: EnhancedQuad.cpp:136
const Matrix & getInitialStiff(void) const
Definition: EnhancedQuad.cpp:170
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Element * getCopy(void) const
Virtual constructor.
Definition: EnhancedQuad.cpp:116
Base class for 2D and 3D materials.
Definition: NDMaterial.h:91
Definition: Vector.h:82
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: EnhancedQuad.cpp:515
Four-node quadrilateral element, which uses a bilinear isoparametric formulation with enhanced strain...
Definition: EnhancedQuad.h:61
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
int recvData(const CommParameters &)
Receives members through the channel being passed as parameter.
Definition: EnhancedQuad.cpp:1319
Base class for 4 node quads.
Definition: QuadBase4N.h:45
Definition: Matrix.h:82
int update(void)
Actualiza el estado of the element.
Definition: EnhancedQuad.cpp:1035
int sendData(CommParameters &)
Send members through the channel being passed as parameter.
Definition: EnhancedQuad.cpp:1310
void setDomain(Domain *)
Sets the domain for the element.
Definition: EnhancedQuad.cpp:129
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: EnhancedQuad.cpp:448
Communication parameters between processes.
Definition: CommParameters.h:65
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: EnhancedQuad.cpp:1342
================================================================================
Definition: ContinuaReprComponent.h:34