xc
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 Pos3dArray3d;
73 class Pos2d;
74 class Pos3d;
75 class Rect3d3dCooSys;
76 class GeomObj2d;
77 class GeomObj3d;
78 
79 namespace XC {
80 class Vector;
81 class Renderer;
82 class Info;
83 class Information;
84 class Parameter;
85 class Response;
86 class ElementalLoad;
87 class Node;
88 class NodePtrArray3d;
89 class ElemPtrArray3d;
90 class SetBase;
91 class SetEstruct;
92 class NodePtrsWithIDs;
93 class Material;
94 class DqVectors;
95 class DqMatrices;
96 class DefaultTag;
97 class GaussModel;
98 class ParticlePos3d;
99 
103 //
109 class Element: public MeshComponent
110  {
111  public:
112  static double dead_srf;
113  typedef std::vector<const Node *> NodesEdge;
114  inline static void setDeadSRF(const double &d)
116  { dead_srf= d; }
117  private:
118  int nodeIndex;
119 
120  static std::deque<Matrix> theMatrices;
121  static std::deque<Vector> theVectors1;
122  static std::deque<Vector> theVectors2;
123 
124  void compute_damping_matrix(Matrix &) const;
125  static DefaultTag defaultTag; //<! default tag for next new element.
126  protected:
127  friend class EntMdlr;
128  friend class Preprocessor;
129  virtual ElemPtrArray3d put_on_mesh(const NodePtrArray3d &,meshing_dir) const;
130  virtual ElemPtrArray3d cose(const SetEstruct &f1,const SetEstruct &f2) const;
131 
132  const Vector &getRayleighDampingForces(void) const;
133 
136  //(mutable para que getDamp pueda ser const).
137  mutable Matrix Kc;
138  //(mutable para que getDamp pueda ser const).
139  int sendData(CommParameters &cp);
140  int recvData(const CommParameters &cp);
141 
142  public:
143  Element(int tag, int classTag);
144  virtual Element *getCopy(void) const= 0;
145 
146  static DefaultTag &getDefaultTag(void);
147 
148  // methods dealing with nodes and number of external dof
150  virtual int getNumExternalNodes(void) const =0;
151  virtual int getNumEdges(void) const;
152  virtual NodePtrsWithIDs &getNodePtrs(void)= 0;
153  virtual const NodePtrsWithIDs &getNodePtrs(void) const= 0;
154  std::vector<int> getIdxNodes(void) const;
161  virtual int getNumDOF(void) const= 0;
162  virtual size_t getDimension(void) const;
163  virtual void setIdNodes(const std::vector<int> &inodes);
164  virtual void setIdNodes(const ID &inodes);
165  void setDomain(Domain *theDomain);
166 
167  // methods dealing with committed state and update
168  virtual int commitState(void);
173  virtual int revertToLastCommit(void) = 0;
174  virtual int revertToStart(void);
175  virtual int update(void);
176  virtual bool isSubdomain(void);
177 
178  // methods to return the current linearized stiffness,
179  // damping and mass matrices
180 
187  virtual const Matrix &getTangentStiff(void) const= 0;
188  virtual const Matrix &getInitialStiff(void) const= 0;
189  virtual const Matrix &getDamp(void) const;
190  virtual const Matrix &getMass(void) const;
191 
192  // methods for applying loads
193  virtual void zeroLoad(void);
194  virtual int addLoad(ElementalLoad *theLoad, double loadFactor)=0;
195  virtual int addInertiaLoadToUnbalance(const Vector &accel)=0;
196 
197  virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF) const;
198 
199  // methods for obtaining resisting force (force includes elemental loads)
200 
206  virtual const Vector &getResistingForce(void) const= 0;
207 
216  virtual const Vector &getResistingForceIncInertia(void) const;
217  const Vector &getNodeResistingComponents(const size_t &,const Vector &) const;
218  const Vector &getNodeResistingForce(const size_t &iNod) const;
219  const Vector &getNodeResistingForceIncInertia(const size_t &iNod) const;
220  const Vector &getNodeResistingForce(const Node *) const;
221  const Vector &getNodeResistingForceIncInertia(const Node *) const;
222  Vector getEquivalentStaticLoad(int mode,const double &) const;
223  Matrix getEquivalentStaticNodalLoads(int mode,const double &) const;
224 
225  // method for obtaining information specific to an element
226  virtual Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
227  virtual int getResponse(int responseID, Information &eleInformation);
228  Response *setMaterialResponse(Material *,const std::vector<std::string> &,const size_t &,Information &);
229 
230 // AddingSensitivity:BEGIN //////////////////////////////////////////
231  virtual int addInertiaLoadSensitivityToUnbalance(const Vector &accel, bool tag);
232  virtual int setParameter(const std::vector<std::string> &argv, Parameter &param);
233  int setMaterialParameter(Material *,const std::vector<std::string> &,const size_t &, Parameter &);
234  virtual int updateParameter(int parameterID, Information &info);
235  virtual int activateParameter(int parameterID);
236  virtual const Vector &getResistingForceSensitivity(int gradNumber);
237  virtual const Matrix &getInitialStiffSensitivity(int gradNumber);
238  virtual const Matrix &getDampSensitivity(int gradNumber);
239  virtual const Matrix &getMassSensitivity(int gradNumber);
240  virtual int commitSensitivity(int gradNumber, int numGrads);
241 // AddingSensitivity:END ///////////////////////////////////////////
242 
243  virtual int addResistingForceToNodalReaction(bool inclInertia);
244 
245  double MaxCooNod(int i) const;
246  double MinCooNod(int i) const;
247  const Matrix &getCooNodes(void) const;
248  virtual Matrix getLocalAxes(bool initialGeometry= true) const;
249  virtual Vector getBaseVector(size_t i,bool initialGeometry= true) const;
250  virtual Vector3d getBaseVector3d(size_t i,bool initialGeometry= true) const;
251  virtual Vector3d getIVector3d(bool initialGeometry= true) const;
252  virtual Vector3d getJVector3d(bool initialGeometry= true) const;
253  virtual Vector3d getKVector3d(bool initialGeometry= true) const;
254  virtual Rect3d3dCooSys getCooSys(bool) const;
255  Pos3d getPosNode(const size_t &i,bool initialGeometry= true) const;
256  std::list<Pos3d> getPosNodes(bool initialGeometry= true) const;
257  virtual Pos3d getCenterOfMassPosition(bool initialGeometry= true) const;
258  Vector getCenterOfMassCoordinates(bool initialGeometry= true) const;
259  Pos3dArray3d getPoints(const size_t &ni,const size_t &nj,const size_t &nk,bool initialGeometry= true);
260  bool In(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
261  bool Out(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
262  bool In(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
263  bool Out(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
264  virtual double getDist2(const Pos2d &p,bool initialGeometry= true) const;
265  virtual double getDist(const Pos2d &p,bool initialGeometry= true) const;
266  virtual double getDist2(const Pos3d &p,bool initialGeometry= true) const;
267  virtual double getDist(const Pos3d &p,bool initialGeometry= true) const;
268 
269  void resetTributaries(void) const;
270  void dumpTributaries(const std::vector<double> &) const;
271  virtual void computeTributaryLengths(bool initialGeometry= true) const;
272  virtual double getTributaryLength(const Node *) const;
273  virtual double getTributaryLengthByTag(const int &) const;
274  virtual void computeTributaryAreas(bool initialGeometry= true) const;
275  virtual double getTributaryArea(const Node *) const;
276  virtual double getTributaryAreaByTag(const int &) const;
277  virtual void computeTributaryVolumes(bool initialGeometry= true) const;
278  virtual double getTributaryVolume(const Node *) const;
279  virtual double getTributaryVolumeByTag(const int &) const;
280 
281  virtual Vector getInterpolationFactors(const ParticlePos3d &) const;
282  virtual Vector getInterpolationFactors(const Pos3d &) const;
283 
284  virtual int getVtkCellType(void) const;
285  virtual const GaussModel &getGaussModel(void) const;
286  virtual NodesEdge getNodesEdge(const size_t &) const;
287  virtual int getEdgeNodes(const Node *,const Node *) const;
288  int getEdgeNodes(const int &,const int &) const;
289  virtual ID getEdgesNode(const Node *) const;
290  std::set<int> getEdgesNodes(const NodePtrSet &) const;
291  ID getEdgesNodeByTag(const int &) const;
292  virtual ID getLocalIndexNodesEdge(const size_t &i) const;
293 
294  virtual std::set<std::string> getMaterialNames(void) const;
295  boost::python::list getMaterialNamesPy(void) const;
296 
297 
298 
299  std::set<SetBase *> get_sets(void) const;
300  void add_to_sets(std::set<SetBase *> &);
301 
302  };
303 } // end of XC namespace
304 
305 #endif
306 
virtual const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Element.cpp:265
virtual int revertToStart(void)
Reverts the element to its initial state.
Definition: Element.cpp:140
virtual void computeTributaryAreas(bool initialGeometry=true) const
Compute tributary areas for each element node.
Definition: Element.cpp:946
virtual void zeroLoad(void)
Zeroes the loads over the element.
Definition: Element.cpp:206
void resetTributaries(void) const
Resets tributary areas of connected nodes.
Definition: Element.cpp:915
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:988
std::vector< const Node * > NodesEdge
Nodes on an element edge.
Definition: Element.h:113
Float vector abstraction.
Definition: Vector.h:93
Natural coordinates of an element&#39;s particle.
Definition: ParticlePos3d.h:40
const Vector & getNodeResistingForceIncInertia(const size_t &iNod) const
Returns the generalized force (including inertia forces) of the element over the node identified by t...
Definition: Element.cpp:344
virtual int getNumEdges(void) const
Returns number of edges (it must be overloaded for elements that have nodes inside edges...
Definition: Element.cpp:111
Information about an element.
Definition: Information.h:80
virtual int getVtkCellType(void) const
Interfaz con VTK.
Definition: Element.cpp:657
virtual ID getEdgesNode(const Node *) const
Returns the edges of the element that ends in the node being passed as parameter. ...
Definition: Element.cpp:706
virtual Vector3d getBaseVector3d(size_t i, bool initialGeometry=true) const
Returns a base vector in the direction of the local i-th axis from the i-th row of the matrix returne...
Definition: Element.cpp:858
Definition: Response.h:71
Matrix Kc
pointer to hold last committed matrix if needed for rayleigh damping
Definition: Element.h:137
virtual Element * getCopy(void) const =0
Virtual constructor.
virtual size_t getDimension(void) const
Returns the element dimension (0, 1, 3 or 3).
Definition: Element.cpp:172
boost::python::list getMaterialNamesPy(void) const
Return the names of the material(s) of the element in a Python list.
Definition: Element.cpp:1067
Finite element model generation tools.
Definition: Preprocessor.h:58
virtual std::set< std::string > getMaterialNames(void) const
Return the names of the material(s) of the element.
Definition: Element.cpp:1058
virtual const Vector & getResistingForce(void) const =0
Returns the resisting force vector for the element.
const Matrix & getCooNodes(void) const
Returns the coordinates of the nodes.
Definition: Element.cpp:817
int recvData(const CommParameters &cp)
Receives object members through the channel being passed as parameter.
Definition: Element.cpp:1086
virtual Vector3d getIVector3d(bool initialGeometry=true) const
Returns a vector in the direction of the local x axis from the first row of the matrix returned by ge...
Definition: Element.cpp:866
virtual void setIdNodes(const std::vector< int > &inodes)
Set the nodes.
Definition: Element.cpp:180
void add_to_sets(std::set< SetBase *> &)
Adds the element to the sets being passed as parameters.
Definition: Element.cpp:772
virtual double getTributaryAreaByTag(const int &) const
Returns the tributary area corresponding to the node.
Definition: Element.cpp:958
virtual Vector getInterpolationFactors(const ParticlePos3d &) const
Returns interpolattion factors for a material point.
Definition: Element.cpp:637
Base class for materials.
Definition: Material.h:91
virtual int setParameter(const std::vector< std::string > &argv, Parameter &param)
Sets the value param to the parameter argv.
Definition: Element.cpp:469
Base class for nodes and elements (mesh components).
Definition: MeshComponent.h:43
virtual const Matrix & getTangentStiff(void) const =0
Return the tangent stiffness matrix.
Vector of integers.
Definition: ID.h:93
Vector getCenterOfMassCoordinates(bool initialGeometry=true) const
Returns the coordinates del centro de gravedad of the element.
Definition: Element.cpp:1032
virtual double getTributaryVolume(const Node *) const
Returns the tributary volume corresponding to the node being passed as parameter. ...
Definition: Element.cpp:975
std::set< int > getEdgesNodes(const NodePtrSet &) const
Returns the element edges that have both ends in the node set being passed as parameter.
Definition: Element.cpp:717
static double dead_srf
Stress reduction factor for foozen elements.
Definition: Element.h:112
virtual double getTributaryArea(const Node *) const
Returns the tributary area corresponding to the node.
Definition: Element.cpp:954
virtual int addResistingForceToNodalReaction(bool inclInertia)
Adds nodal reactions.
Definition: Element.cpp:540
virtual void computeTributaryLengths(bool initialGeometry=true) const
Computes the tributary lengths that corresponds to each node of the element.
Definition: Element.cpp:925
Three-dimensional array of pointers to elements.
Definition: ElemPtrArray3d.h:43
virtual double getTributaryVolumeByTag(const int &) const
Returns the tributary volume corresponding to the node which tag se pasa as parameter.
Definition: Element.cpp:980
virtual int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Element.cpp:472
Base class for the finite elements.
Definition: Element.h:109
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:997
virtual const GaussModel & getGaussModel(void) const
Returns the Gauss integration model of the element.
Definition: Element.cpp:666
RayleighDampingFactors rayFactors
Rayleigh damping factors.
Definition: Element.h:135
virtual const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: Element.cpp:231
void dumpTributaries(const std::vector< double > &) const
Adds to the tributary magnitude of each node the vector being passed as parameter.
Definition: Element.cpp:920
static void setDeadSRF(const double &d)
Assigns Stress Reduction Factor for element deactivation.
Definition: Element.h:115
virtual int commitState(void)
Commit the current element state.
Definition: Element.cpp:118
double MinCooNod(int i) const
Returns the minimum value of the i coordinate of the element nodes.
Definition: Element.cpp:813
Default tag.
Definition: DefaultTag.h:36
virtual Pos3d getCenterOfMassPosition(bool initialGeometry=true) const
Returns the coordinates del centro de gravedad of the element.
Definition: Element.cpp:1023
virtual Vector3d getJVector3d(bool initialGeometry=true) const
Returns a vector in the direction of the local y axis from the second row of the matrix returned by g...
Definition: Element.cpp:871
std::list< Pos3d > getPosNodes(bool initialGeometry=true) const
Returns the coordinates of the nodes.
Definition: Element.cpp:821
virtual int update(void)
Updates the element state.
Definition: Element.cpp:132
virtual Vector3d getKVector3d(bool initialGeometry=true) const
Returns a vector in the direction of the local z axis from the third row of the matrix returned by ge...
Definition: Element.cpp:876
virtual int activateParameter(int parameterID)
Activates the parameter identified by parameterID.
Definition: Element.cpp:475
virtual Response * setResponse(const std::vector< std::string > &argv, Information &eleInformation)
setResponse() is a method invoked to determine if the element will respond to a request for a certain...
Definition: Element.cpp:441
virtual Rect3d3dCooSys getCooSys(bool) const
Returns the element coordinate system from the matrix returned by getLocalAxes.
Definition: Element.cpp:881
Vector load
vector for applying loads
Definition: Element.h:134
virtual int revertToLastCommit(void)=0
Revert to the last commited state.
virtual double getTributaryLength(const Node *) const
Returns the tributary length corresponding to the node being passed as parameter. ...
Definition: Element.cpp:934
static DefaultTag & getDefaultTag(void)
Returns next element&#39;s tag value by default.
Definition: Element.cpp:106
double MaxCooNod(int i) const
Returns the maximum value of the i coordinate of the element nodes.
Definition: Element.cpp:809
virtual int getNumDOF(void) const =0
return the number of DOF associated with the element.
virtual int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Element.cpp:457
Base class for Gauss integration models.
Definition: GaussModel.h:40
const Vector & getNodeResistingComponents(const size_t &, const Vector &) const
Returns the generalized force of the element over the iNod-th node.
Definition: Element.cpp:324
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:370
const Vector & getNodeResistingForce(const size_t &iNod) const
Returns the generalized force of the element over the iNod-th node.
Definition: Element.cpp:336
Multiblock topology object (point, line, face, block,...).
Definition: EntMdlr.h:54
std::set< SetBase * > get_sets(void) const
Returns the sets to which the element belongs.
Definition: Element.cpp:756
Rayleigh damping factors.
Definition: RayleighDampingFactors.h:58
virtual int getEdgeNodes(const Node *, const Node *) const
Returns the edge of the element that ends in the nodes being passed as parameters.
Definition: Element.cpp:686
virtual double getTributaryLengthByTag(const int &) const
Returns the tributary length corresponding to the node which tag se pasa as parameter.
Definition: Element.cpp:939
structured set, i.
Definition: SetEstruct.h:45
ID getEdgesNodeByTag(const int &) const
Returns the edges of the element that have an en in the node with the tag being passed as parameter...
Definition: Element.cpp:737
Node pointer container for elements.
Definition: NodePtrsWithIDs.h:45
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:746
int sendData(CommParameters &cp)
Sends object members through the channel being passed as parameter.
Definition: Element.cpp:1077
virtual const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Element.cpp:252
Three-dimensional array of pointers to nodes.
Definition: NodePtrArray3d.h:50
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:836
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
virtual void computeTributaryVolumes(bool initialGeometry=true) const
Computes the tributary volumes that corresponds to each node of the element.
Definition: Element.cpp:966
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
virtual bool isSubdomain(void)
Returns true if the element is a subdomain.
Definition: Element.cpp:424
virtual Vector getBaseVector(size_t i, bool initialGeometry=true) const
Returns a base vector in the direction of the local i-th axis from the i-th row of the matrix returne...
Definition: Element.cpp:846
Definition: Parameter.h:65
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Element.cpp:188
virtual int getNumExternalNodes(void) const =0
return the number of external nodes associated with the element.
const Vector & getRayleighDampingForces(void) const
Returns element Rayleigh damping forces.
Definition: Element.cpp:388
Base class for loads over elements.
Definition: ElementalLoad.h:77
Pos3d getPosNode(const size_t &i, bool initialGeometry=true) const
Returns the position of the i-th node.
Definition: Element.cpp:901
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:101
virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF) const
Set Rayleigh damping factors.
Definition: Element.cpp:144
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:116
Mesh node.
Definition: Node.h:110
Pos3dArray3d getPoints(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:905
virtual ElemPtrArray3d put_on_mesh(const NodePtrArray3d &, meshing_dir) const
Places the element on the mesh.
Definition: Element.cpp:1043
virtual NodesEdge getNodesEdge(const size_t &) const
Returns the nodes of the element edge.
Definition: Element.cpp:675
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:381