xc
ShellMITC4Base.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 //ShellMITC4Base.h
29 
30 #ifndef ShellMITC4Base_h
31 #define ShellMITC4Base_h
32 
33 #include "Shell4NBase.h"
34 
35 namespace XC {
36 
38 //
41  {
42  protected:
43  double Ktt;
44 
45  static const int ngauss= 4;
46  static const int nstress;
47  mutable std::vector<Vector> strains;
48  std::vector<Vector> persistentInitialDeformation;
49  static ShellBData BData;
50 
51 
52  void formResidAndTangent(int tang_flag) const;
53  const Matrix calculateG(void) const;
54  double *computeBdrill(int node, const double shp[3][4]) const;
55  const Matrix& assembleB(const Matrix &Bmembrane, const Matrix &Bbend, const Matrix &Bshear) const;
56  const Matrix& computeBmembrane(int node, const double shp[3][4] ) const;
57  const Matrix& computeBbend(int node, const double shp[3][4] ) const;
58  int sendData(Communicator &);
59  int recvData(const Communicator &);
60 
61  public:
62  //null constructor
63  ShellMITC4Base(int classTag,const ShellCrdTransf3dBase *);
64  ShellMITC4Base(int tag,int classTag,const SectionForceDeformation *ptr_mat,const ShellCrdTransf3dBase *);
65  //full constructor
66  ShellMITC4Base(int tag,int classTag, int node1, int node2, int node3, int node4, const SectionFDPhysicalProperties &,const ShellCrdTransf3dBase *);
67 
68  // Element birth and death stuff.
69  const std::vector<Vector> &getPersistentInitialDeformation(void) const;
71 
72  //set domain because frank is a dumb ass
73  void setDomain(Domain *theDomain);
74 
75  //return stiffness matrix
76  const Matrix &getInitialStiff(void) const;
77 
78  void alive(void);
79 
80  };
81 
82 } // end of XC namespace
83 #endif
void alive(void)
Reactivates the element.
Definition: ShellMITC4Base.cc:146
Base class for force deformation section models.
Definition: SectionForceDeformation.h:88
std::vector< Vector > persistentInitialDeformation
Persistent initial strain at element level. Used to store the deformation during the inactive phase o...
Definition: ShellMITC4Base.h:48
std::vector< Vector > strains
strains at gauss points.
Definition: ShellMITC4Base.h:47
static const int nstress
(8): three membrane, three moment, two shear
Definition: ShellMITC4Base.h:46
Communication parameters between processes.
Definition: Communicator.h:66
double * computeBdrill(int node, const double shp[3][4]) const
compute Bdrill
Definition: ShellMITC4Base.cc:772
void incrementPersistentInitialDeformationWithCurrentDeformation(void)
Increments the persistent (does not get wiped out by zeroLoad) initial deformation of the section...
Definition: ShellMITC4Base.cc:105
static ShellBData BData
B-bar data.
Definition: ShellMITC4Base.h:49
const Matrix calculateG(void) const
Computes the matrix G.
Definition: ShellMITC4Base.cc:157
Base class for four node shell elements.
Definition: Shell4NBase.h:52
Base class for 3D coordinate transformations.
Definition: ShellCrdTransf3dBase.h:49
Auxiliary data for shell elements.
Definition: ShellBData.h:43
double Ktt
drilling stiffness
Definition: ShellMITC4Base.h:43
const Matrix & computeBmembrane(int node, const double shp[3][4]) const
compute Bmembrane matrix
Definition: ShellMITC4Base.cc:820
const Matrix & getInitialStiff(void) const
return secant matrix
Definition: ShellMITC4Base.cc:202
void formResidAndTangent(int tang_flag) const
form residual and tangent
Definition: ShellMITC4Base.cc:424
ShellMITC4Base(int classTag, const ShellCrdTransf3dBase *)
Constructor.
Definition: ShellMITC4Base.cc:77
void setDomain(Domain *theDomain)
Set the element domain.
Definition: ShellMITC4Base.cc:116
int recvData(const Communicator &)
Receives members through the communicator argument.
Definition: ShellMITC4Base.cc:990
int sendData(Communicator &)
Send members through the communicator argument.
Definition: ShellMITC4Base.cc:981
const Matrix & assembleB(const Matrix &Bmembrane, const Matrix &Bbend, const Matrix &Bshear) const
assemble a B matrix
Definition: ShellMITC4Base.cc:850
static const int ngauss
Number of gauss points.
Definition: ShellMITC4Base.h:45
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Physical properties for shells.
Definition: SectionFDPhysicalProperties.h:41
Matrix of floats.
Definition: Matrix.h:111
const std::vector< Vector > & getPersistentInitialDeformation(void) const
Returns the value of the persistent (does not get wiped out by zeroLoad) initial deformation of the e...
Definition: ShellMITC4Base.cc:98
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
const Matrix & computeBbend(int node, const double shp[3][4]) const
compute Bbend matrix
Definition: ShellMITC4Base.cc:954
Base class for MIT C4 shell elements.
Definition: ShellMITC4Base.h:40