xc
Shell4NBase.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 //Shell4NBase.h
29 
30 #ifndef Shell4NBase_h
31 #define Shell4NBase_h
32 
33 #include "domain/mesh/element/plane/QuadBase4N.h"
34 #include "domain/mesh/element/utils/physical_properties/SectionFDPhysicalProperties.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 Polygon3d;
41 
42 namespace XC {
43 
44 class ShellRawLoad;
45 class ShellUniformLoad;
46 class ShellCrdTransf3dBase;
47 class ParticlePos2d;
48 
50 //
52 class Shell4NBase: public QuadBase4N<SectionFDPhysicalProperties>
53  {
54  void free_mem(void);
55  void alloc(const ShellCrdTransf3dBase *);
56  protected:
57  mutable double xl[2][4];
58 
60  int applyLoad;
61  double appliedB[3];
62 
63  mutable Matrix Ki;
65 
66 
67  //static data
68  static Matrix stiff;
69  static Vector resid;
70  static Matrix mass;
71  static Matrix damping;
72 
73  void formInertiaTerms(int tangFlag) const;
74  virtual void formResidAndTangent(int tang_flag) const= 0;
75  static void shape2d(const double &,const double &, const double x[2][4], double shp[3][4], double &xsj, double sx[2][2]);
76  int sendCoordTransf(int posFlag,const int &,const int &,Communicator &);
77  int recvCoordTransf(int posFlag,const int &posClassTag,const int &posDbTag,const Communicator &);
78  int sendData(Communicator &);
79  int recvData(const Communicator &);
80 
81  friend class ShellCrdTransf3dBase;
82  friend class ShellLinearCrdTransf3d;
83  friend class ShellUpBasisCrdTransf3d;
84  typedef double (*pointer_to_xl)[4]; // To return xl.
85  inline const pointer_to_xl get_xl(void) const
86  { return xl; }
87  inline pointer_to_xl get_xl(void)
88  { return xl; }
89  public:
90  //null constructor
91  Shell4NBase(int classTag,const ShellCrdTransf3dBase *);
92  Shell4NBase(int tag,int classTag,const SectionForceDeformation *ptr_mat,const ShellCrdTransf3dBase *);
93  //full constructor
94  Shell4NBase(int tag,int classTag, int node1, int node2, int node3, int node4, const SectionFDPhysicalProperties &,const ShellCrdTransf3dBase *);
95  Shell4NBase(const Shell4NBase &);
97  ~Shell4NBase(void);
98 
99  //return number of dofs
100  int getNumDOF(void) const;
101 
102  int update(void);
103 
104  //return stiffness matrix
105  const Matrix &getTangentStiff(void) const;
106  const Matrix &getMass(void) const;
107 
108  const GaussModel &getGaussModel(void) const;
109 
111  Vector getInterpolationFactors(const Pos3d &) const;
114 
115  //Load definition methods.
116  const ShellRawLoad *vector3dRawLoadLocal(const std::vector<Vector> &);
117  const ShellRawLoad *vector3dRawLoadGlobal(const std::vector<Vector> &);
120  void strainLoad(const Matrix &);
121  virtual void createInertiaLoad(const Vector &);
122 
123  // methods for applying loads
124  void zeroLoad(void);
125  int addLoad(ElementalLoad *, double);
126  int addInertiaLoadToUnbalance(const Vector &accel);
127 
128  int commitState(void);
129  int revertToLastCommit(void);
130  int revertToStart(void);
131 
132  double getMeanInternalForce(const std::string &) const;
133  double getMeanInternalDeformation(const std::string &) const;
134 
135  virtual Matrix getLocalAxes(bool initialGeometry= true) const;
136  virtual ShellCrdTransf3dBase *getCoordTransf(void);
137  virtual const ShellCrdTransf3dBase *getCoordTransf(void) const;
138 
139  void computeBasis(void);
140  ParticlePos3d getLocalCoordinatesOfNode(const int &) const;
141  ParticlePos3d getNaturalCoordinates(const Pos3d &, bool initialGeometry= true) const;
142  Pos3d getCartesianCoordinates(const ParticlePos2d &,bool initialGeometry= true) const;
143  Pos3d getCartesianCoordinates(const ParticlePos3d &,bool initialGeometry= true) const;
144 
145  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
146  int getResponse(int responseID, Information &eleInformation);
147 
148  const Vector &getResistingForce(void) const;
149  const Vector &getResistingForceIncInertia(void) const;
150 
151  void Print(std::ostream &s, int flag ) const;
152  };
153 
154 } // end of XC namespace
155 #endif
Base class for force deformation section models.
Definition: SectionForceDeformation.h:88
void formInertiaTerms(int tangFlag) const
Form inertia terms.
Definition: Shell4NBase.cc:679
const ShellUniformLoad * vector3dUniformLoadGlobal(const Vector &)
Defines a load over the element from a vector in global coordinates.
Definition: Shell4NBase.cc:287
void Print(std::ostream &s, int flag) const
print out element data
Definition: Shell4NBase.cc:960
Shell4NBase & operator=(const Shell4NBase &)
Assignment operator.
Definition: Shell4NBase.cc:158
Float vector abstraction.
Definition: Vector.h:94
Natural coordinates of an element&#39;s particle.
Definition: ParticlePos3d.h:41
Plane polygon in a 3D space.
Definition: Polygon3d.h:35
Shell 3D coordinate transformations that updates vector directions (needs debugging).
Definition: ShellUpBasisCrdTransf3d.h:43
const ShellUniformLoad * vector3dUniformLoadLocal(const Vector &)
Defines a load over the element from a vector in local coordinates.
Definition: Shell4NBase.cc:248
int recvCoordTransf(int posFlag, const int &posClassTag, const int &posDbTag, const Communicator &)
Receives the coordinate transformation through the communicator argument.
Definition: Shell4NBase.cc:840
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
Base class response objects.
Definition: Response.h:81
Vector getInterpolationFactors(const ParticlePos3d &) const
Returns interpolation factors for a material point.
Definition: Shell4NBase.cc:558
int recvData(const Communicator &)
Receives members through the communicator argument.
Definition: Shell4NBase.cc:880
int revertToStart(void)
Return the initial state.
Definition: Shell4NBase.cc:549
Vector getInterpolatedDisplacements(const ParticlePos3d &) const
Returns interpolated displacements for a material point.
Definition: Shell4NBase.cc:638
Uniform load over shell elements.
Definition: ShellRawLoad.h:40
int sendCoordTransf(int posFlag, const int &, const int &, Communicator &)
Send the coordinate transformation through the communicator argument.
Definition: Shell4NBase.cc:824
double xl[2][4]
local nodal coordinates, two coordinates for each of four nodes
Definition: Shell4NBase.h:57
Shell4NBase(int classTag, const ShellCrdTransf3dBase *)
Constructor.
Definition: Shell4NBase.cc:108
ParticlePos3d getLocalCoordinatesOfNode(const int &) const
Return the local coordinates that correspond to the i-th node.
Definition: Shell4NBase.cc:630
Base class for 4 node quads.
Definition: QuadBase4N.h:45
Base class for four node shell elements.
Definition: Shell4NBase.h:52
int addLoad(ElementalLoad *, double)
Applies on the element the load being passed as parameter.
Definition: Shell4NBase.cc:391
Base class for 3D coordinate transformations.
Definition: ShellCrdTransf3dBase.h:49
int update(void)
Update state variables.
Definition: Shell4NBase.cc:346
Ingernal forces for a shell element.
Definition: FVectorShell.h:41
static void shape2d(const double &, const double &, const double x[2][4], double shp[3][4], double &xsj, double sx[2][2])
shape function routine for MITC4 elements.
Definition: Shell4NBase.cc:774
const ShellRawLoad * vector3dRawLoadGlobal(const std::vector< Vector > &)
Defines a load over the element from a vector of nodal loads in global coordinates.
Definition: Shell4NBase.cc:220
double appliedB[3]
Body forces applied with load.
Definition: Shell4NBase.h:61
Natural coordinates of an element&#39;s particle.
Definition: ParticlePos2d.h:40
const Vector & getResistingForce(void) const
get residual
Definition: Shell4NBase.cc:465
virtual void createInertiaLoad(const Vector &)
Creates the inertia load that corresponds to the acceleration argument.
Definition: Shell4NBase.cc:424
Base class for Gauss integration models.
Definition: GaussModel.h:41
int sendData(Communicator &)
Send members through the communicator argument.
Definition: Shell4NBase.cc:867
void zeroLoad(void)
Sets loads to zero.
Definition: Shell4NBase.cc:380
virtual ShellCrdTransf3dBase * getCoordTransf(void)
Returns a pointer to the coordinate transformation.
Definition: Shell4NBase.cc:898
ParticlePos3d getNaturalCoordinates(const Pos3d &, bool initialGeometry=true) const
Return the natural coordinates that correspond to the argument.
Definition: Shell4NBase.cc:634
Base class for small displacement 3D coordinate transformations.
Definition: ShellLinearCrdTransf3d.h:42
const ShellRawLoad * vector3dRawLoadLocal(const std::vector< Vector > &)
Defines a load over the element from a vector of nodal loads in local coordinates.
Definition: Shell4NBase.cc:178
Response * setResponse(const std::vector< std::string > &argv, Information &eleInformation)
Element response.
Definition: Shell4NBase.cc:906
Posición en tres dimensiones.
Definition: Pos3d.h:44
~Shell4NBase(void)
Destructor.
Definition: Shell4NBase.cc:171
const Matrix & getMass(void) const
return mass matrix
Definition: Shell4NBase.cc:365
int addInertiaLoadToUnbalance(const Vector &accel)
Add to the unbalance vector the inertial load corresponding to the acceleration argument.
Definition: Shell4NBase.cc:500
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
void strainLoad(const Matrix &)
Defines a strain load on this element.
Definition: Shell4NBase.cc:305
virtual Matrix getLocalAxes(bool initialGeometry=true) const
Returns a matrix with the axes of the element as matrix rows [[x1,y1,z1],[x2,y2,z2],...·].
Definition: Shell4NBase.cc:894
ShellCrdTransf3dBase * theCoordTransf
Coordinate transformation.
Definition: Shell4NBase.h:59
Physical properties for shells.
Definition: SectionFDPhysicalProperties.h:41
Matrix of floats.
Definition: Matrix.h:111
const Matrix & getTangentStiff(void) const
return stiffness matrix
Definition: Shell4NBase.cc:353
int commitState(void)
Consuma la coordinate transformation de acuerdo con el estado actual.
Definition: Shell4NBase.cc:533
Base class for loads over elements.
Definition: ElementalLoad.h:79
FVectorShell p0
Reactions in the basic system due to element loads.
Definition: Shell4NBase.h:64
int revertToLastCommit(void)
Returns to the last committed state.
Definition: Shell4NBase.cc:541
Uniform load over shell elements.
Definition: ShellUniformLoad.h:40
const Vector & getResistingForceIncInertia(void) const
get residual with inertia terms
Definition: Shell4NBase.cc:482
int getNumDOF(void) const
return number of dofs
Definition: Shell4NBase.cc:340
const GaussModel & getGaussModel(void) const
Return the element Gauss points.
Definition: Shell4NBase.cc:375
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Shell4NBase.cc:941
void computeBasis(void)
compute local coordinates and basis
Definition: Shell4NBase.cc:761
Pos3d getCartesianCoordinates(const ParticlePos2d &, bool initialGeometry=true) const
Return the cartesian coordinates of the argument.
Definition: Shell4NBase.cc:587