XC Open source finite element analysis program
ShellMITC4Base.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 //ShellMITC4Base.h
28 
29 #ifndef ShellMITC4Base_h
30 #define ShellMITC4Base_h
31 
32 #include "domain/mesh/element/plane/QuadBase4N.h"
33 #include "domain/mesh/element/utils/physical_properties/SectionFDPhysicalProperties.h"
34 #include "ShellCrdTransf3dBase.h"
35 #include "ShellBData.h"
36 #include <utility/matrix/Vector.h>
37 #include <utility/matrix/Matrix.h>
38 #include "domain/mesh/element/utils/fvectors/FVectorShell.h"
39 
40 class Poligono3d;
41 
42 namespace XC {
43 
44 class ShellUniformLoad;
45 
47 //
49 class ShellMITC4Base : public QuadBase4N<SectionFDPhysicalProperties>
50  {
51  protected:
52  double Ktt;
53  mutable double xl[2][4];
54 
56 
58 
59  mutable Matrix Ki;
60 
61  std::vector<Vector> inicDisp;
62 
63  //static data
64  static Matrix stiff;
65  static Vector resid;
66  static Matrix mass;
67  static Matrix damping;
68 
69  static ShellBData BData;
70 
71  void libera(void);
72  void alloc(const ShellCrdTransf3dBase *);
73 
74  void setupInicDisp(void);
75  void capturaInicDisp(void);
76  void zeroInicDisp(void);
77 
78  void formInertiaTerms(int tangFlag) const;
79  void formResidAndTangent(int tang_flag) const;
80  const Matrix calculateG(void) const;
81  double *computeBdrill(int node, const double shp[3][4]) const;
82  const Matrix& assembleB(const Matrix &Bmembrane, const Matrix &Bbend, const Matrix &Bshear) const;
83  const Matrix& computeBmembrane(int node, const double shp[3][4] ) const;
84  const Matrix& computeBbend(int node, const double shp[3][4] ) const;
85  static void shape2d(const double &,const double &, const double x[2][4], double shp[3][4], double &xsj);
86  int sendCoordTransf(int posFlag,const int &,const int &,CommParameters &);
87  int recvCoordTransf(int posFlag,const int &posClassTag,const int &posDbTag,const CommParameters &);
88  int sendData(CommParameters &);
89  int recvData(const CommParameters &);
90 
91  public:
92  //null constructor
93  ShellMITC4Base(int classTag,const ShellCrdTransf3dBase *);
94  ShellMITC4Base(int tag,int classTag,const SectionForceDeformation *ptr_mat,const ShellCrdTransf3dBase *);
95  //full constructor
96  ShellMITC4Base(int tag,int classTag, int node1, int node2, int node3, int node4, const SectionFDPhysicalProperties &,const ShellCrdTransf3dBase *);
97  ShellMITC4Base(const ShellMITC4Base &otro);
99  ~ShellMITC4Base(void);
100 
101  //set domain because frank is a dumb ass
102  void setDomain(Domain *theDomain);
103 
104  //return number of dofs
105  int getNumDOF(void) const;
106 
107  int update(void);
108 
109  //return stiffness matrix
110  const Matrix &getTangentStiff(void) const;
111  const Matrix &getInitialStiff(void) const;
112  const Matrix &getMass(void) const;
113 
114  const GaussModel &getGaussModel(void) const;
115 
117  Vector getInterpolationFactors(const Pos3d &) const;
119  Vector getInterpolatedDisplacements(const Pos3d &) const;
120 
121  //Load definition methods.
124  void strainLoad(const Matrix &);
125 
126  // methods for applying loads
127  void zeroLoad(void);
128  int addLoad(ElementalLoad *theLoad, double loadFactor);
129  int addInertiaLoadToUnbalance(const Vector &accel);
130 
131  int commitState(void);
132  int revertToLastCommit(void);
133  int revertToStart(void);
134  void alive(void);
135 
136  const Vector &getResistingForce(void) const;
137  const Vector &getResistingForceIncInertia(void) const;
138  double getMeanInternalForce(const std::string &) const;
139  double getMeanInternalDeformation(const std::string &) const;
140 
141  virtual Matrix getLocalAxes(bool initialGeometry= true) const;
142  virtual ShellCrdTransf3dBase *getCoordTransf(void);
143  virtual const ShellCrdTransf3dBase *getCoordTransf(void) const;
144 
145  void computeBasis(void);
146  ParticlePos3d getLocalCoordinatesOfNode(const int &) const;
147  ParticlePos3d getNaturalCoordinates(const Pos3d &) const;
148 
149  void Print(std::ostream &s, int flag );
150  };
151 
152 } // end of XC namespace
153 #endif
Base class for MIT C4 shell elements.
Definition: ShellMITC4Base.h:49
void Print(std::ostream &s, int flag)
print out element data
Definition: ShellMITC4Base.cc:298
void alive(void)
Reactivates the element.
Definition: ShellMITC4Base.cc:282
const Vector & getResistingForceIncInertia(void) const
get residual with inertia terms
Definition: ShellMITC4Base.cc:720
static void shape2d(const double &, const double &, const double x[2][4], double shp[3][4], double &xsj)
shape function routine for MITC4 elements.
Definition: ShellMITC4Base.cc:1452
Natural coordinates of an element&#39;s particle.
Definition: ParticlePos3d.h:40
void setDomain(Domain *theDomain)
set domain
Definition: ShellMITC4Base.cc:229
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
double Ktt
drilling stiffness
Definition: ShellMITC4Base.h:52
int commitState(void)
Consuma la coordinate transformation de acuerdo con el estado actual.
Definition: ShellMITC4Base.cc:741
Base class for Gauss integration models.
Definition: GaussModel.h:40
Definition: Vector.h:82
std::vector< Vector > inicDisp
Initial displacements.
Definition: ShellMITC4Base.h:61
const ShellUniformLoad * vector3dUniformLoadLocal(const Vector &)
Defines a load over the element from a vector in local coordinates.
Definition: ShellMITC4Base.cc:143
static ShellBData BData
B-bar data.
Definition: ShellMITC4Base.h:69
const Vector & getResistingForce(void) const
get residual
Definition: ShellMITC4Base.cc:703
int sendData(CommParameters &)
Send members through the channel being passed as parameter.
Definition: ShellMITC4Base.cc:1544
void formInertiaTerms(int tangFlag) const
form inertia terms
Definition: ShellMITC4Base.cc:843
int addLoad(ElementalLoad *theLoad, double loadFactor)
Applies on the element the load being passed as parameter.
Definition: ShellMITC4Base.cc:653
const Matrix & computeBbend(int node, const double shp[3][4]) const
compute Bbend matrix
Definition: ShellMITC4Base.cc:1420
FVectorShell p0
Reactions in the basic system due to element loads.
Definition: ShellMITC4Base.h:57
int revertToStart(void)
Returns the initial state.
Definition: ShellMITC4Base.cc:757
void computeBasis(void)
compute local coordinates and basis
Definition: ShellMITC4Base.cc:1231
Auxiliary data for shell elements.
Definition: ShellBData.h:42
const ShellUniformLoad * vector3dUniformLoadGlobal(const Vector &)
Defines a load over the element from a vector in global coordinates.
Definition: ShellMITC4Base.cc:179
const GaussModel & getGaussModel(void) const
Returns the element Gauss points.
Definition: ShellMITC4Base.cc:642
double * computeBdrill(int node, const double shp[3][4]) const
compute Bdrill
Definition: ShellMITC4Base.cc:1238
int recvCoordTransf(int posFlag, const int &posClassTag, const int &posDbTag, const CommParameters &)
Recibe la coordinate transformation through the channel being passed as parameter.
Definition: ShellMITC4Base.cc:1519
void formResidAndTangent(int tang_flag) const
form residual and tangent
Definition: ShellMITC4Base.cc:927
const Matrix & getInitialStiff(void) const
return secant matrix
Definition: ShellMITC4Base.cc:408
~ShellMITC4Base(void)
Destructor.
Definition: ShellMITC4Base.cc:139
Base class for loads over elements.
Definition: ElementalLoad.h:73
virtual ShellCrdTransf3dBase * getCoordTransf(void)
Returns a pointer to the coordinate transformation.
Definition: ShellMITC4Base.cc:1575
const Matrix & getTangentStiff(void) const
return stiffness matrix
Definition: ShellMITC4Base.cc:351
const Matrix & getMass(void) const
return mass matrix
Definition: ShellMITC4Base.cc:632
Ingernal forces for a shell element.
Definition: FVectorShell.h:40
ShellCrdTransf3dBase * theCoordTransf
Coordinate transformation.
Definition: ShellMITC4Base.h:55
Vector getInterpolatedDisplacements(const ParticlePos3d &) const
Returns interpolated displacements for a material point.
Definition: ShellMITC4Base.cc:802
int sendCoordTransf(int posFlag, const int &, const int &, CommParameters &)
Envía la coordinate transformation through the channel being passed as parameter. ...
Definition: ShellMITC4Base.cc:1503
Uniform load over shell elements.
Definition: ShellUniformLoad.h:39
Base class for 4 node quads.
Definition: QuadBase4N.h:45
Definition: Matrix.h:82
int revertToLastCommit(void)
Returns to the last commited state.
Definition: ShellMITC4Base.cc:749
ShellMITC4Base(int classTag, const ShellCrdTransf3dBase *)
Constructor.
Definition: ShellMITC4Base.cc:101
int getNumDOF(void) const
return number of dofs
Definition: ShellMITC4Base.cc:278
int recvData(const CommParameters &)
Receives members through the channel being passed as parameter.
Definition: ShellMITC4Base.cc:1557
Physical properties for shells.
Definition: SectionFDPhysicalProperties.h:40
void zeroLoad(void)
Zeroes the element load vector.
Definition: ShellMITC4Base.cc:646
Base class for force deformation section models. Constitutive equations of the section.
Definition: SectionForceDeformation.h:86
Vector getInterpolationFactors(const ParticlePos3d &) const
Returns interpolattion factors for a material point.
Definition: ShellMITC4Base.cc:766
const Matrix & computeBmembrane(int node, const double shp[3][4]) const
compute Bmembrane matrix
Definition: ShellMITC4Base.cc:1286
const Matrix & assembleB(const Matrix &Bmembrane, const Matrix &Bbend, const Matrix &Bshear) const
assemble a B matrix
Definition: ShellMITC4Base.cc:1316
Communication parameters between processes.
Definition: CommParameters.h:65
virtual Matrix getLocalAxes(bool initialGeometry=true) const
Returs a matrix with the axes of the element as matrix rows [[x1,y1,z1],[x2,y2,z2],...·].
Definition: ShellMITC4Base.cc:1571
ShellMITC4Base & operator=(const ShellMITC4Base &otro)
Assignment operator.
Definition: ShellMITC4Base.cc:129
int update(void)
Actualiza los valores de las variables de estado.
Definition: ShellMITC4Base.cc:344
const Matrix calculateG(void) const
Calcula la matriz G.
Definition: ShellMITC4Base.cc:363
================================================================================
Definition: ContinuaReprComponent.h:34
Base class for 3D coordinate transformations.
Definition: ShellCrdTransf3dBase.h:48
double xl[2][4]
local nodal coordinates, two coordinates for each of four nodes
Definition: ShellMITC4Base.h:53