xc
Element.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 /* ****************************************************************** **
29 ** OpenSees - Open System for Earthquake Engineering Simulation **
30 ** Pacific Earthquake Engineering Research Center **
31 ** **
32 ** **
33 ** (C) Copyright 1999, The Regents of the University of California **
34 ** All Rights Reserved. **
35 ** **
36 ** Commercial use of this program without express permission of the **
37 ** University of California, Berkeley, is strictly prohibited. See **
38 ** file 'COPYRIGHT' in main directory for information on usage and **
39 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
40 ** **
41 ** Developed by: **
42 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
43 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
44 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
45 ** **
46 ** ****************************************************************** */
47 
48 // $Revision: 1.12 $
49 // $Date: 2005/02/17 22:29:54 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/Element.h,v $
51 
52 
53 // Written: fmk
54 // Created: 11/96
55 // Revision: A
56 //
57 // Description: This file contains the class definition for Element.
58 // Element is an abstract base class and thus no objects of it's type
59 // can be instantiated. It has pure virtual functions which must be
60 // implemented in it's derived classes.
61 //
62 // What: "@(#) Element.h, revA"
63 
64 #ifndef Element_h
65 #define Element_h
66 
67 #include "domain/mesh/MeshComponent.h"
68 #include "preprocessor/MeshingParams.h"
69 #include "domain/mesh/element/utils/RayleighDampingFactors.h"
70 #include "utility/matrix/Matrix.h"
71 #include "domain/mesh/node/NodeTopology.h"
72 #include "utility/matrices/3d_arrays/BoolArray3d.h"
73 
74 class Pos3dArray3d;
75 class Pos2d;
76 class Pos3d;
77 class Rect3d3dCooSys;
78 class Rect2d2dCooSys;
79 class GeomObj2d;
80 class GeomObj3d;
81 
82 namespace XC {
83 class Vector;
84 class Renderer;
85 class Info;
86 class Information;
87 class Parameter;
88 class Response;
89 class ElementalLoad;
90 class Node;
91 class NodePtrArray3d;
92 class ElemPtrArray3d;
93 class SetBase;
94 class SetEstruct;
95 class NodePtrsWithIDs;
96 class Material;
97 class DqVectors;
98 class DqMatrices;
99 class DefaultTag;
100 class GaussModel;
101 class ParticlePos3d;
102 
106 //
112 class Element: public MeshComponent
113  {
114  public:
115  static double dead_srf;
116  typedef std::vector<const Node *> NodesEdge;
117  inline static void setDeadSRF(const double &d)
119  { dead_srf= d; }
120  private:
121  int nodeIndex;
122 
123  static std::deque<Matrix> theMatrices;
124  static std::deque<Vector> theVectors1;
125  static std::deque<Vector> theVectors2;
126 
127  void compute_damping_matrix(Matrix &) const;
128  static DefaultTag defaultTag; //<! default tag for next new element.
129  protected:
130  friend class EntMdlr;
131  friend class Preprocessor;
132  virtual ElemPtrArray3d put_on_mesh(const NodePtrArray3d &,meshing_dir) const;
133  virtual ElemPtrArray3d sew(const SetEstruct &f1,const SetEstruct &f2) const;
134 
135  const Vector &getRayleighDampingForces(void) const;
136 
139  //(mutable to allow getDamp being const).
140  mutable Matrix Kc;
141  //(mutable to allow getDamp being const).
142 
143  const Material *get_material_ptr(const std::string &) const;
144 
145  int sendData(Communicator &comm);
146  int recvData(const Communicator &comm);
147 
148  public:
149  Element(int tag, int classTag);
150  virtual Element *getCopy(void) const= 0;
151 
152  static DefaultTag &getDefaultTag(void);
153 
154  // methods dealing with nodes and number of external dof
156  virtual int getNumExternalNodes(void) const =0;
157  virtual int getNumEdges(void) const;
158  virtual BoolArray3d getNodePattern(void) const;
159  virtual NodePtrsWithIDs &getNodePtrs(void)= 0;
160  virtual const NodePtrsWithIDs &getNodePtrs(void) const= 0;
161  std::vector<int> getIdxNodes(void) const;
168  virtual int getNumDOF(void) const= 0;
169  virtual size_t getDimension(void) const;
170  virtual double getLength(bool initialGeometry= true) const;
171  virtual double getArea(bool initialGeometry= true) const;
172  virtual double getVolume(bool initialGeometry= true) const;
173  virtual void setIdNodes(const std::vector<int> &inodes);
174  virtual void setIdNodes(const ID &inodes);
175  virtual void setIdNode(const int &i, const int &inode);
176  int find(const Node *) const;
177  void replaceNode(Node *, Node *);
178  void setDomain(Domain *theDomain);
179 
180 
181 
182  // methods dealing with committed state and update
183  virtual int commitState(void);
188  virtual int revertToLastCommit(void) = 0;
189  virtual int revertToStart(void);
190  virtual int update(void);
191  virtual bool isSubdomain(void);
192 
193  // methods to return the current linearized stiffness,
194  // damping and mass matrices
195 
202  virtual const Matrix &getTangentStiff(void) const= 0;
203  virtual const Matrix &getInitialStiff(void) const= 0;
204  virtual const Matrix &getDamp(void) const;
205  virtual const Matrix &getMass(void) const;
206  virtual Matrix getMass(const Node *) const;
207  Matrix getTotalMass(void) const;
208  double getTotalMassComponent(const int &) const;
209  Matrix getNodeMatrixComponents(const Node *,const Matrix &) const;
210 
211  // methods for applying loads
212  virtual void zeroLoad(void);
213  virtual void createInertiaLoad(const Vector &);
214  virtual int addLoad(ElementalLoad *, double loadFactor)=0;
215  virtual int addInertiaLoadToUnbalance(const Vector &accel)=0;
216 
217  virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF) const;
218 
219  // methods for obtaining resisting force (force includes elemental loads)
220 
226  virtual const Vector &getResistingForce(void) const= 0;
227 
236  virtual const Vector &getResistingForceIncInertia(void) const;
237  const Vector &getNodeResistingComponents(const size_t &,const Vector &) const;
238  const Vector &getNodeResistingForce(const size_t &iNod) const;
239  const Vector &getNodeResistingForceIncInertia(const size_t &iNod) const;
240  const Vector &getNodeResistingForce(const Node *) const;
241  const Vector &getNodeResistingForceIncInertia(const Node *) const;
242  const Vector &getNodeResistingForceByTag(const int &) const;
243  const Vector &getNodeResistingForceIncInertiaByTag(const int &) const;
244  Vector getEquivalentStaticLoad(int mode,const double &) const;
245  Matrix getEquivalentStaticNodalLoads(int mode,const double &) const;
246 
247  // method for obtaining information specific to an element
248  virtual Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
249  virtual int getResponse(int responseID, Information &eleInformation);
250  Response *setMaterialResponse(Material *,const std::vector<std::string> &,const size_t &,Information &);
251 
252 // AddingSensitivity:BEGIN //////////////////////////////////////////
253  virtual int addInertiaLoadSensitivityToUnbalance(const Vector &accel, bool tag);
254  virtual int setParameter(const std::vector<std::string> &argv, Parameter &param);
255  int setMaterialParameter(Material *,const std::vector<std::string> &,const size_t &, Parameter &);
256  virtual int updateParameter(int parameterID, Information &info);
257  virtual int activateParameter(int parameterID);
258  virtual const Vector &getResistingForceSensitivity(int gradNumber);
259  virtual const Matrix &getInitialStiffSensitivity(int gradNumber);
260  virtual const Matrix &getDampSensitivity(int gradNumber);
261  virtual const Matrix &getMassSensitivity(int gradNumber);
262  virtual int commitSensitivity(int gradNumber, int numGrads);
263 // AddingSensitivity:END ///////////////////////////////////////////
264 
265  virtual int addResistingForceToNodalReaction(bool inclInertia);
266 
267  double MaxCooNod(int i) const;
268  double MinCooNod(int i) const;
269  const Matrix &getCooNodes(void) const;
270  virtual Matrix getLocalAxes(bool initialGeometry= true) const;
271  virtual Vector getBaseVector(size_t i,bool initialGeometry= true) const;
272  virtual Vector3d getBaseVector3d(size_t i,bool initialGeometry= true) const;
273  virtual Vector3d getIVector3d(bool initialGeometry= true) const;
274  virtual Vector3d getJVector3d(bool initialGeometry= true) const;
275  virtual Vector3d getKVector3d(bool initialGeometry= true) const;
276  virtual Rect3d3dCooSys getCooSys(bool) const;
277  virtual Rect2d2dCooSys getCooSys2d(bool) const;
278  Pos3d getPosNode(const size_t &i,bool initialGeometry= true) const;
279  std::deque<Pos3d> getPosNodes(bool initialGeometry= true) const;
280  virtual Pos3d getCenterOfMassPosition(bool initialGeometry= true) const;
281  Vector getCenterOfMassCoordinates(bool initialGeometry= true) const;
282  Pos3dArray3d getPoints(const size_t &ni,const size_t &nj,const size_t &nk,bool initialGeometry= true);
283  bool In(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
284  bool Out(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
285  bool In(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
286  bool Out(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
287  bool Crosses(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
288  bool Crosses(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
289  virtual double getDist2(const Pos2d &p,bool initialGeometry= true) const;
290  virtual double getDist(const Pos2d &p,bool initialGeometry= true) const;
291  virtual double getDist2(const Pos3d &p,bool initialGeometry= true) const;
292  virtual double getDist(const Pos3d &p,bool initialGeometry= true) const;
293  virtual Pos2d getProjection(const Pos2d &p,bool initialGeometry= true) const;
294  virtual Pos3d getProjection(const Pos3d &p,bool initialGeometry= true) const;
295 
296  void resetTributaries(void) const;
297  void dumpTributaries(const std::vector<double> &) const;
298  virtual void computeTributaryLengths(bool initialGeometry= true) const;
299  virtual double getTributaryLength(const Node *) const;
300  virtual double getTributaryLengthByTag(const int &) const;
301  virtual void computeTributaryAreas(bool initialGeometry= true) const;
302  virtual double getTributaryArea(const Node *) const;
303  virtual double getTributaryAreaByTag(const int &) const;
304  virtual void computeTributaryVolumes(bool initialGeometry= true) const;
305  virtual double getTributaryVolume(const Node *) const;
306  virtual double getTributaryVolumeByTag(const int &) const;
307 
308  virtual ParticlePos3d getNaturalCoordinates(const Pos3d &, bool initialGeometry= true) const;
309  virtual Vector getInterpolationFactors(const ParticlePos3d &) const;
310  virtual Vector getInterpolationFactors(const Pos3d &) const;
311 
312  virtual int getVtkCellType(void) const;
313  virtual const GaussModel &getGaussModel(void) const;
314  virtual NodesEdge getNodesEdge(const size_t &) const;
315  virtual int getEdgeNodes(const Node *,const Node *) const;
316  int getEdgeNodes(const int &,const int &) const;
317  virtual ID getEdgesNode(const Node *) const;
318  std::set<int> getEdgesNodes(const NodePtrSet &) const;
319  ID getEdgesNodeByTag(const int &) const;
320  virtual ID getLocalIndexNodesEdge(const size_t &i) const;
321 
322  virtual std::set<std::string> getMaterialNames(void) const;
323  boost::python::list getMaterialNamesPy(void) const;
324 
325  virtual boost::python::list getValuesAtNodes(const std::string &, bool silent= false) const;
326 
327  std::set<SetBase *> get_sets(void) const;
328  void add_to_sets(std::set<SetBase *> &);
329  void copySetsFrom(const Element &);
330 
331  boost::python::dict getPyDict(void) const;
332  void setPyDict(const boost::python::dict &);
333  };
334 } // end of XC namespace
335 
336 #endif
337 
virtual const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Element.cpp:434
const Material * get_material_ptr(const std::string &) const
Return a pointer to the material that corresponds to the name.
Definition: Element.cpp:1024
virtual int revertToStart(void)
Reverts the element to its initial state.
Definition: Element.cpp:148
virtual void computeTributaryAreas(bool initialGeometry=true) const
Compute tributary areas for each element node.
Definition: Element.cpp:1268
virtual void zeroLoad(void)
Zeroes the loads over the element.
Definition: Element.cpp:318
void resetTributaries(void) const
Resets tributary areas of connected nodes.
Definition: Element.cpp:1235
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:1314
int find(const Node *) const
Returns the index of the node whose pointer is being passed as parameter.
Definition: Element.cpp:256
std::vector< const Node * > NodesEdge
Nodes on an element edge.
Definition: Element.h:116
void setPyDict(const boost::python::dict &)
Set the values of the object members from a Python dictionary.
Definition: Element.cpp:1516
Float vector abstraction.
Definition: Vector.h:94
Natural coordinates of an element&#39;s particle.
Definition: ParticlePos3d.h:41
int sendData(Communicator &comm)
Sends object members through the communicator argument.
Definition: Element.cpp:1489
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:539
Base class for grids of bool in 3D (used to express if something exists or not in a (layer...
Definition: BoolArray3d.h:34
virtual int getNumEdges(void) const
Returns number of edges (it must be overloaded for elements that have nodes inside edges...
Definition: Element.cpp:119
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
virtual int getVtkCellType(void) const
Interfaz con VTK.
Definition: Element.cpp:883
virtual ID getEdgesNode(const Node *) const
Returns the edges of the element that ends in the node being passed as parameter. ...
Definition: Element.cpp:932
virtual ParticlePos3d getNaturalCoordinates(const Pos3d &, bool initialGeometry=true) const
Return the natural coordinates that correspond to the given position.
Definition: Element.cpp:852
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:1151
Base class response objects.
Definition: Response.h:81
Matrix Kc
pointer to hold last committed matrix if needed for rayleigh damping
Definition: Element.h:140
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:183
std::deque< Pos3d > getPosNodes(bool initialGeometry=true) const
Returns the coordinates of the nodes.
Definition: Element.cpp:1089
bool Crosses(const GeomObj3d &, const double &factor=1.0, const double &tol=0.0) const
Return true if the element cross (i.e.
Definition: Element.cpp:1110
boost::python::list getMaterialNamesPy(void) const
Return the names of the material(s) of the element in a Python list.
Definition: Element.cpp:1445
Finite element model generation tools.
Definition: Preprocessor.h:59
virtual std::set< std::string > getMaterialNames(void) const
Return the names of the material(s) of the element.
Definition: Element.cpp:1435
virtual double getLength(bool initialGeometry=true) const
Return the element length.
Definition: Element.cpp:193
virtual const Vector & getResistingForce(void) const =0
Returns the resisting force vector for the element.
Posición en dos dimensiones.
Definition: Pos2d.h:41
const Matrix & getCooNodes(void) const
Returns the coordinates of the nodes.
Definition: Element.cpp:1084
int recvData(const Communicator &comm)
Receives object members through the communicator argument.
Definition: Element.cpp:1498
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:1160
virtual void setIdNodes(const std::vector< int > &inodes)
Set the nodes.
Definition: Element.cpp:222
void add_to_sets(std::set< SetBase *> &)
Adds the element to the sets being passed as parameters.
Definition: Element.cpp:999
Three-dimensional coordinate system defined in a three-dimensional space.
Definition: Rect3d3dCooSys.h:33
virtual double getTributaryAreaByTag(const int &) const
Returns the tributary area corresponding to the node.
Definition: Element.cpp:1280
virtual Vector getInterpolationFactors(const ParticlePos3d &) const
Returns interpolation factors for a material point.
Definition: Element.cpp:863
const Vector & getNodeResistingForceIncInertiaByTag(const int &) const
Returns the generalized force (including inertia forces) of the element over the node identified by t...
Definition: Element.cpp:574
Base class for materials.
Definition: Material.h:93
virtual int setParameter(const std::vector< std::string > &argv, Parameter &param)
Sets the value param to the parameter argv.
Definition: Element.cpp:682
Base class for nodes and elements (mesh components).
Definition: MeshComponent.h:44
virtual const Matrix & getTangentStiff(void) const =0
Return the tangent stiffness matrix.
int setMaterialParameter(Material *, const std::vector< std::string > &, const size_t &, Parameter &)
Set the value of a parameter of the material.
Definition: Element.cpp:1059
Matrix getTotalMass(void) const
Returns the sum of the mass matrices corresponding to the nodes.
Definition: Element.cpp:404
virtual void createInertiaLoad(const Vector &)
Creates the inertia load that corresponds to the acceleration argument.
Definition: Element.cpp:307
virtual void setIdNode(const int &i, const int &inode)
Set the i-th node.
Definition: Element.cpp:244
Vector of integers.
Definition: ID.h:95
Vector getCenterOfMassCoordinates(bool initialGeometry=true) const
Returns the coordinates of the center of mass of the element.
Definition: Element.cpp:1394
virtual double getTributaryVolume(const Node *) const
Returns the tributary volume corresponding to the node being passed as parameter. ...
Definition: Element.cpp:1298
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:943
static double dead_srf
Stress reduction factor for foozen elements.
Definition: Element.h:115
virtual double getTributaryArea(const Node *) const
Returns the tributary area corresponding to the node.
Definition: Element.cpp:1276
virtual int addResistingForceToNodalReaction(bool inclInertia)
Adds nodal reactions.
Definition: Element.cpp:753
virtual void computeTributaryLengths(bool initialGeometry=true) const
Computes the tributary lengths that corresponds to each node of the element.
Definition: Element.cpp:1246
Three-dimensional array of pointers to elements.
Definition: ElemPtrArray3d.h:47
virtual double getTributaryVolumeByTag(const int &) const
Returns the tributary volume corresponding to the node which tag se pasa as parameter.
Definition: Element.cpp:1303
virtual int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Element.cpp:685
Base class for the finite elements.
Definition: Element.h:112
virtual Pos2d getProjection(const Pos2d &p, bool initialGeometry=true) const
Return the projection of the argument on the element.
Definition: Element.cpp:1360
void replaceNode(Node *, Node *)
Replace the old node by the new one.
Definition: Element.cpp:262
const Vector & getNodeResistingForceByTag(const int &) const
Returns the generalized force of the element over the node identified by the given integer...
Definition: Element.cpp:565
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:1327
Base class for the two-dimensional geometric objects.
Definition: GeomObj2d.h:37
virtual const GaussModel & getGaussModel(void) const
Returns the Gauss integration model of the element.
Definition: Element.cpp:892
RayleighDampingFactors rayFactors
Rayleigh damping factors.
Definition: Element.h:138
virtual const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: Element.cpp:343
void dumpTributaries(const std::vector< double > &) const
Adds to the tributary magnitude of each node the vector being passed as parameter.
Definition: Element.cpp:1240
static void setDeadSRF(const double &d)
Assigns Stress Reduction Factor for element deactivation.
Definition: Element.h:118
virtual int commitState(void)
Commit the current element state.
Definition: Element.cpp:126
double MinCooNod(int i) const
Returns the minimum value of the i coordinate of the element nodes.
Definition: Element.cpp:1080
Default tag.
Definition: DefaultTag.h:37
virtual Pos3d getCenterOfMassPosition(bool initialGeometry=true) const
Returns the coordinates of the center of gravity of the element.
Definition: Element.cpp:1383
Matrix getNodeMatrixComponents(const Node *, const Matrix &) const
Returns the components of the matrix relative to the node.
Definition: Element.cpp:495
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:1166
boost::python::dict getPyDict(void) const
Return a Python dictionary with the object members values.
Definition: Element.cpp:1507
virtual int update(void)
Updates the element state.
Definition: Element.cpp:140
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:1172
virtual int activateParameter(int parameterID)
Activates the parameter identified by parameterID.
Definition: Element.cpp:688
double getTotalMassComponent(const int &) const
Return the mass matrix component for the DOF argument.
Definition: Element.cpp:420
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:654
virtual Rect3d3dCooSys getCooSys(bool) const
Returns the element coordinate system from the matrix returned by getLocalAxes.
Definition: Element.cpp:1178
Vector load
vector for applied nodal loads.
Definition: Element.h:137
virtual int revertToLastCommit(void)=0
Revert to the last committed state.
virtual double getTributaryLength(const Node *) const
Returns the tributary length corresponding to the node being passed as parameter. ...
Definition: Element.cpp:1255
static DefaultTag & getDefaultTag(void)
Returns next element&#39;s tag value by default.
Definition: Element.cpp:114
double MaxCooNod(int i) const
Returns the maximum value of the i coordinate of the element nodes.
Definition: Element.cpp:1076
bool Out(const GeomObj3d &, const double &factor=1.0, const double &tol=0.0) const
Return true if the element is outside the given object.
Definition: Element.cpp:1097
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:670
Base class for Gauss integration models.
Definition: GaussModel.h:41
const Vector & getNodeResistingComponents(const size_t &, const Vector &) const
Returns the generalized force of the element over the iNod-th node.
Definition: Element.cpp:519
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:583
const Vector & getNodeResistingForce(const size_t &iNod) const
Returns the generalized force of the element over the iNod-th node.
Definition: Element.cpp:531
Multiblock topology object (point, line, face, block,...).
Definition: EntMdlr.h:55
std::set< SetBase * > get_sets(void) const
Returns the sets to which the element belongs.
Definition: Element.cpp:982
Rayleigh damping factors.
Definition: RayleighDampingFactors.h:59
Twp-dimensional rectangular coordinate system defined in a two-dimensional space. ...
Definition: Rect2d2dCooSys.h:33
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:912
virtual double getTributaryLengthByTag(const int &) const
Returns the tributary length corresponding to the node which tag se pasa as parameter.
Definition: Element.cpp:1260
structured set, i.
Definition: SetEstruct.h:47
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:963
Position array in a three-dimensional space.
Definition: Pos3dArray3d.h:37
virtual BoolArray3d getNodePattern(void) const
Return a grid of booleans, one for each of the element nodes.
Definition: Element.cpp:1407
Node pointer container for elements.
Definition: NodePtrsWithIDs.h:46
Posición en tres dimensiones.
Definition: Pos3d.h:44
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:972
virtual const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Element.cpp:364
Three-dimensional array of pointers to nodes.
Definition: NodePtrArray3d.h:51
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: Element.cpp:1121
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
virtual boost::python::list getValuesAtNodes(const std::string &, bool silent=false) const
Return a python list with the values of the argument property at element nodes.
Definition: Element.cpp:1463
virtual void computeTributaryVolumes(bool initialGeometry=true) const
Computes the tributary volumes that corresponds to each node of the element.
Definition: Element.cpp:1289
Matrix of floats.
Definition: Matrix.h:111
virtual Rect2d2dCooSys getCooSys2d(bool) const
Returns the element coordinate system from the matrix returned by getLocalAxes.
Definition: Element.cpp:1199
virtual bool isSubdomain(void)
Returns true if the element is a subdomain.
Definition: Element.cpp:637
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:1134
Parameter.
Definition: Parameter.h:68
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Element.cpp:291
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:601
Base class for loads over elements.
Definition: ElementalLoad.h:79
Pos3d getPosNode(const size_t &i, bool initialGeometry=true) const
Returns the position of the i-th node.
Definition: Element.cpp:1220
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:104
virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF) const
Set Rayleigh damping factors.
Definition: Element.cpp:155
void copySetsFrom(const Element &)
Add this element to all the sets containing the given one.
Definition: Element.cpp:1015
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
Mesh node.
Definition: Node.h:111
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 element.
Definition: Element.cpp:1225
virtual ElemPtrArray3d put_on_mesh(const NodePtrArray3d &, meshing_dir) const
Places the element on the mesh.
Definition: Element.cpp:1418
bool In(const GeomObj3d &, const double &factor=1.0, const double &tol=0.0) const
Return true if the element is inside the given object.
Definition: Element.cpp:1093
Vector en tres dimensiones.
Definition: Vector3d.h:39
virtual NodesEdge getNodesEdge(const size_t &) const
Returns the nodes of the element edge.
Definition: Element.cpp:901
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:594
std::vector< int > getIdxNodes(void) const
Return the indexes of the nodes (used when creating VTK meshes).
Definition: Element.cpp:1072
virtual double getVolume(bool initialGeometry=true) const
Return the element volume.
Definition: Element.cpp:213
virtual double getArea(bool initialGeometry=true) const
Return the element area.
Definition: Element.cpp:203
Clase base para los objetos en tres dimensiones.
Definition: GeomObj3d.h:43