xc
ShellNLDKGQ.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.10 $
49 // $Date: 2014/09/30 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/ShellNLDKGQ/ShellNLDKGQ.h,v $
51 
52 //Written: Lisha Wang, Xinzheng Lu, Linlin Xie, Song Cen & Quan Gu
53 //
54 // four node flat shell element with membrane and drill DOF
55 // considering geometric nonlinear, form nonlinear shell element
56 // using updated Lagrangian formula
57 // Ref: Plate Bending Part - DKQ, thin plate element
58 // Membrane Part - GQ12, a membrane element with drilling DOF
59 //
60 
61 #ifndef ShellNLDKGQ_h
62 #define ShellNLDKGQ_h
63 
64 #include "Shell4NBase.h"
65 #include "domain/mesh/element/utils/coordTransformation/ShellNLCrdTransf3d.h"
66 
67 namespace XC {
68 
74 class ShellNLDKGQ: public Shell4NBase
75  {
76  private:
77  static ShellNLCrdTransf3d non_linear_trf;
78 
79  //last resid
80  // Vector CstrainGauss,TstrainGauss;
81  Vector CstrainGauss;
82  mutable Vector TstrainGauss; //modify for geometric nonlinearity
83 
84  //start Yuli Huang (yulihuang@gmail.com) & Xinzheng Lu (luxz@tsinghua.edu.cn)
85  void updateBasis(void) const;
86  //end Yuli Huang (yulihuang@gmail.com) & Xinzheng Lu (luxz@tsinghua.edu.cn)
87 
88  //inertia terms
89  void formInertiaTerms(int tangFlag) const;
90 
91  //form residual and tangent
92  void formResidAndTangent(int tangFlag) const;
93 
94  //compute Bdrill matrix
95  //double* computeBdrill( int node, const double shp[3][4]);
96 
97  //assemble a B matrix
98  const Matrix& assembleB( const Matrix &Bmembrane,
99  const Matrix &Bbend,
100  const Matrix &Bshear ) const;
101 
102  //compute Bmembrane matrix
103  const Matrix& computeBmembrane( int node, const double shp[3][4],
104  const double shpDrill[4][4]) const;
105 
106  //compute Bbend matrix
107  const Matrix& computeBbend( int node, const double shpBend[6][12] ) const;
108 
109  //add for geometric nonlinearity
110  const Matrix& computeBG(int node, const double shpBend[6][12]) const;
111  const Vector& computeNLdstrain(const Matrix &BG,const Vector &dispIncLocalBend) const;
112 
113  //shape function routine for four node quads
114  static void shape2d( double ss, double tt,
115  const double x[2][4],
116  double shp[3][4],
117  double &xsj ,double sx[2][2]);
118 
119  //shape function routine for membrane
120  static void shapeDrill(double ss, double tt, const double x[2][4],
121  double sx[2][2], double shpDrill[4][4]);
122  //shape function routine for bending
123  static void shapeBend(double ss, double tt, const double x[2][4], double sx[2][2], double shpBend[6][12]);
124 
125  int sendData(Communicator &);
126  int recvData(const Communicator &);
127 
128  public:
129  ShellNLDKGQ(void);
130  ShellNLDKGQ(int tag, const SectionForceDeformation *);
131  ShellNLDKGQ(int tag,
132  int node1, int node2, int node3, int node4,
133  const SectionForceDeformation &theMaterial);
134  Element *getCopy(void) const;
135 
136  //set domain because frank is a dumb ass
137  void setDomain(Domain *theDomain);
138 
139 
140  int commitState(void); //commit state
141  int revertToLastCommit(void); //revert to last commit
142  int revertToStart(void); //revert to start
143 
144  //print out element data
145  void Print(std::ostream &os, int flag) const;
146 
147  //return stiffness matrix
148  const Matrix &getTangentStiff(void) const;
149  const Matrix &getInitialStiff(void) const;
150  const Matrix &getMass(void) const;
151 
152  // public methods for element output
153  int sendSelf(Communicator &);
154  int recvSelf(const Communicator &);
155  };
156 
157 } // end of XC namespace
158 #endif
Base class for force deformation section models.
Definition: SectionForceDeformation.h:88
Float vector abstraction.
Definition: Vector.h:94
int commitState(void)
commit state
Definition: ShellNLDKGQ.cpp:99
Communication parameters between processes.
Definition: Communicator.h:66
four node flat shell element with membrane and drill DOF considering geometric nonlinear, form nonlinear shell element using updated Lagrangian formula.
Definition: ShellNLDKGQ.h:74
int revertToLastCommit(void)
revert to last commit
Definition: ShellNLDKGQ.cpp:109
Base class for the finite elements.
Definition: Element.h:112
void setDomain(Domain *theDomain)
set domain
Definition: ShellNLDKGQ.cpp:91
Base class for four node shell elements.
Definition: Shell4NBase.h:52
int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: ShellNLDKGQ.cpp:1495
Element * getCopy(void) const
Virtual constructor.
Definition: ShellNLDKGQ.cpp:87
ShellNLDKGQ(void)
Definition: ShellNLDKGQ.cpp:67
int revertToStart(void)
revert to start
Definition: ShellNLDKGQ.cpp:118
const Matrix & getTangentStiff(void) const
return stiffness matrix
Definition: ShellNLDKGQ.cpp:127
int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: ShellNLDKGQ.cpp:1480
Base class for small displacement 3D coordinate transformations.
Definition: ShellNLCrdTransf3d.h:42
const Matrix & getInitialStiff(void) const
return secant matrix
Definition: ShellNLDKGQ.cpp:136
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
const Matrix & getMass(void) const
return mass matrix
Definition: ShellNLDKGQ.cpp:491
Matrix of floats.
Definition: Matrix.h:111
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
void Print(std::ostream &os, int flag) const
print out element data
Definition: ShellNLDKGQ.cpp:1515