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