XC Open source finite element analysis program
Element.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 /* ****************************************************************** **
28 ** OpenSees - Open System for Earthquake Engineering Simulation **
29 ** Pacific Earthquake Engineering Research Center **
30 ** **
31 ** **
32 ** (C) Copyright 1999, The Regents of the University of California **
33 ** All Rights Reserved. **
34 ** **
35 ** Commercial use of this program without express permission of the **
36 ** University of California, Berkeley, is strictly prohibited. See **
37 ** file 'COPYRIGHT' in main directory for information on usage and **
38 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
39 ** **
40 ** Developed by: **
41 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
42 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
43 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
44 ** **
45 ** ****************************************************************** */
46 
47 // $Revision: 1.12 $
48 // $Date: 2005/02/17 22:29:54 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/Element.h,v $
50 
51 
52 // Written: fmk
53 // Created: 11/96
54 // Revision: A
55 //
56 // Description: This file contains the class definition for Element.
57 // Element is an abstract base class and thus no objects of it's type
58 // can be instantiated. It has pure virtual functions which must be
59 // implemented in it's derived classes.
60 //
61 // What: "@(#) Element.h, revA"
62 
63 #ifndef Element_h
64 #define Element_h
65 
66 #include "domain/mesh/MeshComponent.h"
67 #include "preprocessor/MeshingParams.h"
68 #include "domain/mesh/element/utils/RayleighDampingFactors.h"
69 #include "utility/matrix/Matrix.h"
70 #include "domain/mesh/node/NodeTopology.h"
71 
72 class TritrizPos3d;
73 class Pos2d;
74 class Pos3d;
75 
76 namespace XC {
77 class Vector;
78 class Renderer;
79 class Info;
80 class Information;
81 class Parameter;
82 class Response;
83 class ElementalLoad;
84 class Node;
85 class TritrizPtrNod;
86 class TritrizPtrElem;
87 class SetBase;
88 class SetEstruct;
89 class NodePtrsWithIDs;
90 class Material;
91 class DqVectors;
92 class DqMatrices;
93 class DefaultTag;
94 class GaussModel;
95 class MEDGaussModel;
96 class ParticlePos3d;
97 
99 //
101 //
104 class Element: public MeshComponent
105  {
106  public:
107  static double dead_srf;
108  typedef std::vector<const Node *> NodesEdge;
109  inline static void setDeadSRF(const double &d)
111  { dead_srf= d; }
112  private:
113  int nodeIndex;
114 
115  static std::deque<Matrix> theMatrices;
116  static std::deque<Vector> theVectors1;
117  static std::deque<Vector> theVectors2;
118 
119  void compute_damping_matrix(Matrix &) const;
120  static DefaultTag defaultTag; //<! default tag for next new element.
121  protected:
122  friend class EntMdlr;
123  friend class Preprocessor;
124  virtual TritrizPtrElem put_on_mesh(const TritrizPtrNod &,meshing_dir) const;
125  virtual TritrizPtrElem cose(const SetEstruct &f1,const SetEstruct &f2) const;
126 
127  const Vector &getRayleighDampingForces(void) const;
128 
131  //(mutable para que getDamp pueda ser const).
132  mutable Matrix Kc;
133  //(mutable para que getDamp pueda ser const).
134  int sendData(CommParameters &cp);
135  int recvData(const CommParameters &cp);
136 
137  public:
138  Element(int tag, int classTag);
139  virtual Element *getCopy(void) const= 0;
140 
141  static DefaultTag &getDefaultTag(void);
142 
143  // methods dealing with nodes and number of external dof
144  virtual int getNumExternalNodes(void) const =0;
145  virtual int getNumEdges(void) const;
146  virtual NodePtrsWithIDs &getNodePtrs(void)= 0;
147  virtual const NodePtrsWithIDs &getNodePtrs(void) const= 0;
148  std::vector<int> getIdxNodes(void) const;
149  virtual int getNumDOF(void) const= 0;
150  virtual size_t getDimension(void) const;
151  virtual void setIdNodos(const std::vector<int> &inodos);
152  virtual void setIdNodos(const ID &inodos);
153  void setDomain(Domain *theDomain);
154 
155  // methods dealing with committed state and update
156  virtual int commitState(void);
157  virtual int revertToLastCommit(void) = 0;
158  virtual int revertToStart(void);
159  virtual int update(void);
160  virtual bool isSubdomain(void);
161 
162  // methods to return the current linearized stiffness,
163  // damping and mass matrices
164  virtual const Matrix &getTangentStiff(void) const= 0;
165  virtual const Matrix &getInitialStiff(void) const= 0;
166  virtual const Matrix &getDamp(void) const;
167  virtual const Matrix &getMass(void) const;
168 
169  // methods for applying loads
170  virtual void zeroLoad(void);
171  virtual int addLoad(ElementalLoad *theLoad, double loadFactor)=0;
172  virtual int addInertiaLoadToUnbalance(const Vector &accel)=0;
173 
174  virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF) const;
175 
176  // methods for obtaining resisting force (force includes elemental loads)
177  virtual const Vector &getResistingForce(void) const= 0;
178  virtual const Vector &getResistingForceIncInertia(void) const;
179  const Vector &getNodeResistingComponents(const size_t &,const Vector &) const;
180  const Vector &getNodeResistingForce(const size_t &iNod) const;
181  const Vector &getNodeResistingForceIncInertia(const size_t &iNod) const;
182  const Vector &getNodeResistingForce(const Node *) const;
183  const Vector &getNodeResistingForceIncInertia(const Node *) const;
184  Vector getEquivalentStaticLoad(int mode,const double &) const;
185  Matrix getEquivalentStaticNodalLoads(int mode,const double &) const;
186 
187  // method for obtaining information specific to an element
188  virtual Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
189  virtual int getResponse(int responseID, Information &eleInformation);
190  Response *setMaterialResponse(Material *,const std::vector<std::string> &,const size_t &,Information &);
191 
192 // AddingSensitivity:BEGIN //////////////////////////////////////////
193  virtual int addInertiaLoadSensitivityToUnbalance(const Vector &accel, bool tag);
194  virtual int setParameter(const std::vector<std::string> &argv, Parameter &param);
195  int setMaterialParameter(Material *,const std::vector<std::string> &,const size_t &, Parameter &);
196  virtual int updateParameter(int parameterID, Information &info);
197  virtual int activateParameter(int parameterID);
198  virtual const Vector &getResistingForceSensitivity(int gradNumber);
199  virtual const Matrix &getInitialStiffSensitivity(int gradNumber);
200  virtual const Matrix &getDampSensitivity(int gradNumber);
201  virtual const Matrix &getMassSensitivity(int gradNumber);
202  virtual int commitSensitivity(int gradNumber, int numGrads);
203 // AddingSensitivity:END ///////////////////////////////////////////
204 
205  virtual int addResistingForceToNodalReaction(bool inclInertia);
206 
207  double MaxCooNod(int i) const;
208  double MinCooNod(int i) const;
209  const Matrix &getCooNodos(void) const;
210  virtual Matrix getLocalAxes(bool initialGeometry= true) const;
211  Pos3d getPosNodo(const size_t &i,bool initialGeometry= true) const;
212  std::list<Pos3d> getPosNodos(bool initialGeometry= true) const;
213  virtual Pos3d getPosCdg(bool initialGeometry= true) const;
214  Vector getCooCdg(bool initialGeometry= true) const;
215  TritrizPos3d getPuntos(const size_t &ni,const size_t &nj,const size_t &nk,bool initialGeometry= true);
216  virtual double getDist2(const Pos2d &p,bool initialGeometry= true) const;
217  virtual double getDist(const Pos2d &p,bool initialGeometry= true) const;
218  virtual double getDist2(const Pos3d &p,bool initialGeometry= true) const;
219  virtual double getDist(const Pos3d &p,bool initialGeometry= true) const;
220 
221  void resetTributarias(void) const;
222  void vuelcaTributarias(const std::vector<double> &) const;
223  virtual void calculaLongsTributarias(bool initialGeometry= true) const;
224  virtual double getLongTributaria(const Node *) const;
225  virtual double getLongTributariaByTag(const int &) const;
226  virtual void calculaAreasTributarias(bool initialGeometry= true) const;
227  virtual double getAreaTributaria(const Node *) const;
228  virtual double getAreaTributariaByTag(const int &) const;
229  virtual void calculaVolsTributarios(bool initialGeometry= true) const;
230  virtual double getVolTributario(const Node *) const;
231  virtual double getVolTributarioByTag(const int &) const;
232 
233  virtual Vector getInterpolationFactors(const ParticlePos3d &) const;
234  virtual Vector getInterpolationFactors(const Pos3d &) const;
235 
236  virtual int getVtkCellType(void) const;
237  virtual int getMEDCellType(void) const;
238  virtual const GaussModel &getGaussModel(void) const;
239  MEDGaussModel getMEDGaussModel(void) const;
240  virtual NodesEdge getNodesEdge(const size_t &) const;
241  virtual int getEdgeNodes(const Node *,const Node *) const;
242  int getEdgeNodes(const int &,const int &) const;
243  virtual ID getEdgesNode(const Node *) const;
244  std::set<int> getEdgesNodes(const NodePtrSet &) const;
245  ID getEdgesNodeByTag(const int &) const;
246  virtual ID getLocalIndexNodesEdge(const size_t &i) const;
247 
248  std::set<SetBase *> get_sets(void) const;
249  void add_to_sets(std::set<SetBase *> &);
250 
251  };
252 } // end of XC namespace
253 
254 #endif
255 
void resetTributarias(void) const
Resets tributary areas of connected nodes.
Definition: Element.cpp:795
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: Element.cpp:770
virtual Element * getCopy(void) const =0
Virtual constructor.
Vector load
vector for applying loads
Definition: Element.h:129
virtual Pos3d getPosCdg(bool initialGeometry=true) const
Returns the coordinates del centro de gravedad of the element.
Definition: Element.cpp:906
virtual int commitState(void)
Consuma el estado of the element.
Definition: Element.cpp:111
Vector getEquivalentStaticLoad(int mode, const double &) const
Returns the equivalent static load for the mode being passed as parameter and the acceleration that c...
Definition: Element.cpp:335
"Tritriz" of pointers to elements.
Definition: TritrizPtrElem.h:43
Multiblock topology object (point, line, face, block,...).
Definition: EntMdlr.h:53
virtual double getLongTributaria(const Node *) const
Returns the tributary length corresponding to the node being passed as parameter. ...
Definition: Element.cpp:814
std::vector< const Node * > NodesEdge
Definition: Element.h:108
Base class for materials.
Definition: Material.h:85
Element(int tag, int classTag)
Constructor that takes the element&#39;s unique tag and the number of external nodes for the element...
Definition: Element.cpp:97
virtual double getVolTributario(const Node *) const
Returns the tributary volume corresponding to the node being passed as parameter. ...
Definition: Element.cpp:858
Natural coordinates of an element&#39;s particle.
Definition: ParticlePos3d.h:40
std::set< int > getEdgesNodes(const NodePtrSet &) const
Returns the bordes of the element que tienen ambos extremos en el node set being passed as parameter...
Definition: Element.cpp:661
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Base class for Gauss integration models.
Definition: GaussModel.h:40
"Tritriz" of pointers to elements.
Definition: TritrizPtrNod.h:51
virtual int addResistingForceToNodalReaction(bool inclInertia)
Adds nodal reactions.
Definition: Element.cpp:476
Definition: Vector.h:82
const Vector & getNodeResistingForce(const size_t &iNod) const
Returns the generalized force of the element over the iNod-th node.
Definition: Element.cpp:301
static double dead_srf
Stress reduction factor for foozen elements.
Definition: Element.h:107
virtual const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Element.cpp:217
virtual double getAreaTributariaByTag(const int &) const
Returns the área tributaria corresponding to the node cuyo tag se pasa as parameter.
Definition: Element.cpp:841
static void setDeadSRF(const double &d)
Assigns Stress Reduction Factor for element deactivation.
Definition: Element.h:110
Mesh node.
Definition: Node.h:99
std::list< Pos3d > getPosNodos(bool initialGeometry=true) const
Returns the coordinates of the nodes.
Definition: Element.cpp:765
Information about an element.
Definition: Information.h:80
virtual double getDist2(const Pos2d &p, bool initialGeometry=true) const
Returns the squared distance from the element to the point being passed as parameter.
Definition: Element.cpp:871
Information about Gauss model.
Definition: MEDGaussModel.h:38
virtual Vector getInterpolationFactors(const ParticlePos3d &) const
Returns interpolattion factors for a material point.
Definition: Element.cpp:572
Base class for nodes and elements (mesh components).
Definition: MeshComponent.h:41
virtual int getVtkCellType(void) const
Interfaz con VTK.
Definition: Element.cpp:590
virtual void calculaVolsTributarios(bool initialGeometry=true) const
Calcula los volúmenes tributarios that corresponds to cada nodo of the element.
Definition: Element.cpp:849
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
MEDGaussModel getMEDGaussModel(void) const
Returns the Gauss integration model of the element for MED library.
Definition: Element.cpp:611
virtual NodesEdge getNodesEdge(const size_t &) const
Returns the nodos del borde (o arista) of the element which index is being passed as parameter...
Definition: Element.cpp:619
const Vector & getNodeResistingComponents(const size_t &, const Vector &) const
Returns the generalized force of the element over the iNod-th node.
Definition: Element.cpp:289
Matrix Kc
pointer to hold last committed matrix if needed for rayleigh damping
Definition: Element.h:132
const Vector & getNodeResistingForceIncInertia(const size_t &iNod) const
Returns the fuerza generalizada (incluyendo fuerzas de inercia) of the element sobre el nodo which in...
Definition: Element.cpp:309
Default tag.
Definition: DefaultTag.h:36
virtual double getVolTributarioByTag(const int &) const
Returns the tributary volume corresponding to the node cuyo tag se pasa as parameter.
Definition: Element.cpp:863
Finite element model generation tools.
Definition: Preprocessor.h:58
virtual void calculaAreasTributarias(bool initialGeometry=true) const
Calcula las áreas tributarias that corresponds to cada nodo of the element.
Definition: Element.cpp:827
virtual double getDist(const Pos2d &p, bool initialGeometry=true) const
Returns the the distance from the element to the point being passed as parameter. ...
Definition: Element.cpp:880
double MinCooNod(int i) const
Returns the minimum value de la coordenada i of the nodes of the element.
Definition: Element.cpp:757
RayleighDampingFactors rayFactors
Rayleigh damping factors.
Definition: Element.h:130
virtual void zeroLoad(void)
Anula el load vector aplicadas of the element.
Definition: Element.cpp:185
Node pointer container for elements.
Definition: NodePtrsWithIDs.h:45
Definition: ID.h:77
structured set, i. e. a set that can return a pointer a to a node or an element from its indices i...
Definition: SetEstruct.h:45
virtual double getLongTributariaByTag(const int &) const
Returns the tributary length corresponding to the node cuyo tag se pasa as parameter.
Definition: Element.cpp:819
Matrix getEquivalentStaticNodalLoads(int mode, const double &) const
Returns the equivalent static load on each node for the mode being passed as parameter and the corres...
Definition: Element.cpp:346
virtual int getMEDCellType(void) const
Interfaz con el formato MED de Salome.
Definition: Element.cpp:597
double MaxCooNod(int i) const
Returns the valor máximo de la coordenada i of the nodes of the element.
Definition: Element.cpp:753
Definition: Matrix.h:82
virtual ID getEdgesNode(const Node *) const
Returns the bordes of the element que tienen por extremo el nodo being passed as parameter.
Definition: Element.cpp:650
Rayleigh damping factors.
Definition: RayleighDampingFactors.h:45
virtual TritrizPtrElem put_on_mesh(const TritrizPtrNod &, meshing_dir) const
Places the element on the mesh.
Definition: Element.cpp:926
virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF) const
Asigna valores a los coeficientes de amortiguamiento de Rayleigh.
Definition: Element.cpp:126
virtual const GaussModel & getGaussModel(void) const
Returns the Gauss integration model of the element.
Definition: Element.cpp:604
TritrizPos3d getPuntos(const size_t &ni, const size_t &nj, const size_t &nk, bool initialGeometry=true)
Returns a grid of points distributed along the line.
Definition: Element.cpp:786
void add_to_sets(std::set< SetBase * > &)
Adds the element to the sets being passed as parameters.
Definition: Element.cpp:716
virtual size_t getDimension(void) const
Returns the element dimension (0, 1, 3 or 3).
Definition: Element.cpp:154
virtual ID getLocalIndexNodesEdge(const size_t &i) const
Returns the local indexes of the element nodes that lie over the i-th edge.
Definition: Element.cpp:690
Definition: Parameter.h:65
void vuelcaTributarias(const std::vector< double > &) const
Adds to the tributary magnitude of each node the vector being passed as parameter.
Definition: Element.cpp:800
int sendData(CommParameters &cp)
Sends object members through the channel being passed as parameter.
Definition: Element.cpp:939
const Vector & getRayleighDampingForces(void) const
Returns element Rayleigh damping forces.
Definition: Element.cpp:353
virtual const Matrix & getDamp(void) const
Returns the matriz de amortiguamiento.
Definition: Element.cpp:203
virtual double getAreaTributaria(const Node *) const
Returns the área tributaria corresponding to the node being passed as parameter.
Definition: Element.cpp:836
virtual const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: Element.cpp:230
virtual int getNumEdges(void) const
Returns number of edges (it must be overloaded for elements that have nodes inside edges...
Definition: Element.cpp:107
virtual int update(void)
Actualiza el estado of the element.
Definition: Element.cpp:119
Vector getCooCdg(bool initialGeometry=true) const
Returns the coordinates del centro de gravedad of the element.
Definition: Element.cpp:915
ID getEdgesNodeByTag(const int &) const
Returns the bordes of the element que tienen por extremo el nodo cuyo tag is being passed as paramete...
Definition: Element.cpp:681
Communication parameters between processes.
Definition: CommParameters.h:65
================================================================================
Definition: ContinuaReprComponent.h:34
int recvData(const CommParameters &cp)
Receives object members through the channel being passed as parameter.
Definition: Element.cpp:948
virtual void calculaLongsTributarias(bool initialGeometry=true) const
Computes the tributary lengths that corresponds to each node of the element.
Definition: Element.cpp:805
std::set< SetBase * > get_sets(void) const
Returns the sets to which the element belongs.
Definition: Element.cpp:700
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Element.cpp:170
Pos3d getPosNodo(const size_t &i, bool initialGeometry=true) const
Returns the position of the i-th node.
Definition: Element.cpp:782
virtual void setIdNodos(const std::vector< int > &inodos)
Asigna los nodos.
Definition: Element.cpp:162
const Matrix & getCooNodos(void) const
Returns the coordinates of the nodes.
Definition: Element.cpp:761
Definition: Response.h:71
static DefaultTag & getDefaultTag(void)
Returns next element&#39;s tag value by default.
Definition: Element.cpp:102
virtual int getEdgeNodes(const Node *, const Node *) const
Returns the borde (o arista) of the element que tiene por extremos los nodos being passed as paramete...
Definition: Element.cpp:630