xc
EnhancedQuad.h
1 // -*-c++-*-
2 //----------------------------------------------------------------------------
3 // XC program; finite element analysis code
4 // for structural analysis and design.
5 //
6 // Copyright (C) Luis C. Pérez Tato
7 //
8 // This program derives from OpenSees <http://opensees.berkeley.edu>
9 // developed by the «Pacific earthquake engineering research center».
10 //
11 // Except for the restrictions that may arise from the copyright
12 // of the original program (see copyright_opensees.txt)
13 // XC is free software: you can redistribute it and/or modify
14 // it under the terms of the GNU General Public License as published by
15 // the Free Software Foundation, either version 3 of the License, or
16 // (at your option) any later version.
17 //
18 // This software is distributed in the hope that it will be useful, but
19 // WITHOUT ANY WARRANTY; without even the implied warranty of
20 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 // GNU General Public License for more details.
22 //
23 //
24 // You should have received a copy of the GNU General Public License
25 // along with this program.
26 // If not, see <http://www.gnu.org/licenses/>.
27 //----------------------------------------------------------------------------
28 /* ****************************************************************** **
29 ** OpenSees - Open System for Earthquake Engineering Simulation **
30 ** Pacific Earthquake Engineering Research Center **
31 ** **
32 ** **
33 ** (C) Copyright 1999, The Regents of the University of California **
34 ** All Rights Reserved. **
35 ** **
36 ** Commercial use of this program without express permission of the **
37 ** University of California, Berkeley, is strictly prohibited. See **
38 ** file 'COPYRIGHT' in main directory for information on usage and **
39 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
40 ** **
41 ** Developed by: **
42 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
43 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
44 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
45 ** **
46 ** ****************************************************************** */
47 
48 // $Revision: 1.8 $
49 // $Date: 2003/02/14 23:01:10 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/fourNodeQuad/EnhancedQuad.h,v $
51 
52 #include "domain/mesh/element/plane/SolidMech4N.h"
53 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h"
54 #include <utility/matrix/Vector.h>
55 #include <utility/matrix/Matrix.h>
56 
57 namespace XC {
59 //
63  {
64  private:
65 
66  mutable Vector alpha;
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  protected:
122  int sendData(Communicator &);
123  int recvData(const Communicator &);
124  public:
125  EnhancedQuad(int tag= 0,const NDMaterial *ptr_mat= nullptr);
126  //full constructor
127  EnhancedQuad(int tag, int nd1, int nd2, int nd3, int nd4, NDMaterial &, const std::string &);
128 
129  Element *getCopy(void) const;
130  //destructor
131  virtual ~EnhancedQuad(void);
132 
133  void setDomain(Domain *); //set domain
134 
135  //return number of dofs
136  int getNumDOF(void) const;
137 
138  // methods dealing with state updates
139  int update(void);
140 
141  //print out element data
142  void Print(std::ostream &s, int flag) const;
143 
144  //return stiffness matrix
145  const Matrix &getTangentStiff(void) const;
146  const Matrix &getInitialStiff(void) const;
147  const Matrix &getMass(void) const;
148 
149  void alive(void);
150  int addLoad(ElementalLoad *theLoad, double loadFactor);
151  int addInertiaLoadToUnbalance(const Vector &accel);
152 
153  //get residual and residual with inertia terms
154  const Vector &getResistingForce(void) const;
155  const Vector &getResistingForceIncInertia(void) const;
156 
157  int setParameter(const std::vector<std::string> &argv, Parameter &param);
158  int updateParameter(int parameterID, Information &info);
159 
160  virtual int sendSelf(Communicator &);
161  virtual int recvSelf(const Communicator &);
162  };
163 } // end of XC namespace
const Matrix & getMass(void) const
return mass matrix
Definition: EnhancedQuad.cpp:442
void setDomain(Domain *)
Set domain.
Definition: EnhancedQuad.cpp:125
int update(void)
Updates the element state.
Definition: EnhancedQuad.cpp:1018
Float vector abstraction.
Definition: Vector.h:94
int getNumDOF(void) const
return number of dofs
Definition: EnhancedQuad.cpp:132
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
void Print(std::ostream &s, int flag) const
Print stuff.
Definition: EnhancedQuad.cpp:137
void alive(void)
Reactivates the element.
Definition: EnhancedQuad.cpp:453
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: EnhancedQuad.cpp:154
virtual int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: EnhancedQuad.cpp:1336
int sendData(Communicator &)
Send members through the communicator argument.
Definition: EnhancedQuad.cpp:1280
Element * getCopy(void) const
Virtual constructor.
Definition: EnhancedQuad.cpp:116
const Vector & getResistingForce(void) const
Get resisting force (residual vector).
Definition: EnhancedQuad.cpp:507
Four-node quadrilateral element, which uses a bilinear isoparametric formulation with enhanced strain...
Definition: EnhancedQuad.h:62
EnhancedQuad(int tag=0, const NDMaterial *ptr_mat=nullptr)
Constructor.
Definition: EnhancedQuad.cpp:90
Base class for the finite elements.
Definition: Element.h:112
int updateParameter(int parameterID, Information &info)
Update parameter value.
Definition: EnhancedQuad.cpp:1316
int recvData(const Communicator &)
Receives members through the communicator argument.
Definition: EnhancedQuad.cpp:1288
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: EnhancedQuad.cpp:1351
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
Four node quadrilateral element for two-dimensional problems.
Definition: SolidMech4N.h:76
const Matrix & getInitialStiff(void) const
Return secant matrix.
Definition: EnhancedQuad.cpp:166
Parameter.
Definition: Parameter.h:68
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
const Vector & getResistingForceIncInertia(void) const
Get resisting force with inertia terms.
Definition: EnhancedQuad.cpp:521
int setParameter(const std::vector< std::string > &argv, Parameter &param)
Set parameter value.
Definition: EnhancedQuad.cpp:1296