XC Open source finite element analysis program
NineNodeMixedQuad.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.5 $
48 // $Date: 2003/02/14 23:01:12 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/fourNodeQuad/NineNodeMixedQuad.h,v $
50 
51 // Ed "C++" Love
52 //
53 // Constant Presssure/Volume Four Node Quadrilateral
54 // Plane Strain (NOT PLANE STRESS)
55 //
56 
57 #include <domain/mesh/element/ElemWithMaterial.h>
58 #include <utility/matrix/Vector.h>
59 #include <utility/matrix/Matrix.h>
60 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h"
61 
62 namespace XC {
63 
64 class NDMaterial;
65 
67 //
69 class NineNodeMixedQuad : public ElemWithMaterial<9,NDMaterialPhysicalProperties>
70  {
71  private:
72  static double xl[][9];
73  mutable Matrix *Ki;
74 
75  //static data
76  static Matrix stiff;
77  static Vector resid;
78  static Matrix mass;
79  static Matrix damping;
80 
81 
82  //quadrature data
83  static double root06;
84  static double sg[3];
85  static double wg[3];
86 
87  //node information
88 
89 
90  //form residual and tangent
91  void formResidAndTangent(int tang_flag) const;
92  //inertia terms
93  void formInertiaTerms(int tangFlag) const;
94 
95  const Matrix& computeBbar( int node,
96  const double natCoor[2],
97  const double shp[3][9],
98  double shpBar[3][9][3] ) const;
99 
100  static void shape2dNine( double coor[2], const double x[2][9], double shp[3][9], double &xsj );
101  void computeBasis(void) const;
102  static double shape1d( int code, int node, double xi );
103 
104  protected:
105  int sendData(CommParameters &);
106  int recvData(const CommParameters &);
107  public :
108 
109  //null constructor
110  NineNodeMixedQuad(void);
111 
112  //full constructor
113  NineNodeMixedQuad(int tag,
114  int node1, int node2, int node3, int node4, int node5, int node6, int node7, int node8, int node9,
115  NDMaterial &theMaterial);
116  Element *getCopy(void) const;
117  //destructor
118  virtual ~NineNodeMixedQuad(void);
119 
120  //set domain
121  void setDomain(Domain *);
122 
123  //return number of dofs
124  int getNumDOF(void) const;
125 
126  //print out element data
127  void Print(std::ostream &s, int flag);
128 
129  //return stiffness matrix
130  const Matrix &getTangentStiff(void) const;
131  const Matrix &getInitialStiff(void) const;
132  const Matrix &getMass(void) const;
133 
134  int addLoad(ElementalLoad *theLoad, double loadFactor);
135  int addInertiaLoadToUnbalance(const Vector &accel);
136 
137  //get residual
138  const Vector &getResistingForce(void) const;
139 
140  //get residual with inertia terms
141  const Vector &getResistingForceIncInertia(void) const;
142 
143  virtual int sendSelf(CommParameters &);
144  virtual int recvSelf(const CommParameters &);
145  };
146 } // end of XC namespace
const Matrix & getMass(void) const
return mass matrix
Definition: NineNodeMixedQuad.cpp:367
Element with material.
Definition: ElemWithMaterial.h:40
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Base class for 2D and 3D materials.
Definition: NDMaterial.h:91
Definition: Vector.h:82
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: NineNodeMixedQuad.cpp:1113
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
void setDomain(Domain *)
Sets the domain for the element.
Definition: NineNodeMixedQuad.cpp:111
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: NineNodeMixedQuad.cpp:1102
Nine node quad.
Definition: NineNodeMixedQuad.h:69
Element * getCopy(void) const
Virtual constructor.
Definition: NineNodeMixedQuad.cpp:100
Definition: Matrix.h:82
virtual ~NineNodeMixedQuad(void)
destructor
Definition: NineNodeMixedQuad.cpp:104
void Print(std::ostream &s, int flag)
Imprime el objeto.
Definition: NineNodeMixedQuad.cpp:123
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: NineNodeMixedQuad.cpp:437
Communication parameters between processes.
Definition: CommParameters.h:65
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: NineNodeMixedQuad.cpp:1124
================================================================================
Definition: ContinuaReprComponent.h:34