xc
Node.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.9 $
48 // $Date: 2005/11/23 22:48:50 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/domain/mesh/node/Node.h,v $
50 
51 
52 #ifndef Node_h
53 #define Node_h
54 
55 // Written: fmk
56 // Created: 11/96
57 // Revision: A
58 //
59 // Purpose: This file contains the class interface for Node.
60 // A Node provides the abstraction for a node in the FEM.
61 // Nodes have original position, trial displacement, velocity and
62 // acceleration and committed displacement, velocity and acceleration
63 // (the last committed() trial quantities).
64 //
65 // What: "@(#) Node.h, revA"
66 
67 #include "domain/mesh/MeshComponent.h"
68 #include "NodeDispVectors.h"
69 #include "NodeVelVectors.h"
70 #include "NodeAccelVectors.h"
71 #include "utility/matrix/Matrix.h"
72 #include <boost/python/list.hpp>
73 
74 class Pos2d;
75 class Pos3d;
76 class SlidingVectorsSystem3d;
77 class GeomObj2d;
78 class GeomObj3d;
79 
80 namespace XC {
81 class ContinuaReprComponent;
82 class Matrix;
83 class Channel;
84 class Renderer;
85 class TrfGeom;
86 class NodeLocker;
87 class DefaultTag;
88 class ElementEdges;
89 class SetBase;
90 class MeshEdge;
91 class DOF_Group;
92 class DqPtrsElem;
93 
98 //
110 class Node: public MeshComponent
111  {
112  private:
113  // private data associated with each node object
114  int numberDOF;
115  DOF_Group *theDOF_GroupPtr;
116  Vector Crd;
117 
118  NodeDispVectors disp;
119  NodeVelVectors vel;
120  NodeAccelVectors accel;
121 
122  Matrix R;
123  Matrix mass;
124  mutable Vector unbalLoad;
125  Vector unbalLoadWithInertia;
126  mutable Vector reaction;
127  double alphaM;
128  mutable double tributary;
129 
130  Matrix theEigenvectors; //Eigenvectors matrix.
131 
132  // AddingSensitivity:BEGIN /////////////////////////////////////////
133  Matrix dispSensitivity;
134  Matrix velSensitivity;
135  Matrix accSensitivity;
136  int parameterID;
137  // AddingSensitivity:END ///////////////////////////////////////////
138 
139  static std::deque<Matrix> theMatrices;
140 
141  mutable std::set<ContinuaReprComponent *> connected;
142 
143  std::set<int> freeze_constraints;
144  const ID &get_id_constraints(void) const;
145  void set_id_constraints(const ID &);
146 
147  static DefaultTag defaultTag; //<! tag for next new node.
148  protected:
149 
150  DbTagData &getDbTagData(void) const;
151  int sendData(CommParameters &);
152  int recvData(const CommParameters &);
153  public:
154  // typedefs.
155  typedef std::set<Element *> ElementPtrSet;
156  typedef std::set<const Element *> ElementConstPtrSet;
157 
158  // constructors
159  Node(int classTag);
160  Node(int tag, int classTag);
161  Node(int tag, int ndof, double Crd1);
162  Node(int tag, int ndof, double Crd1, double Crd2);
163  Node(int tag, int ndof, double Crd1, double Crd2, double Crd3);
164  Node(int tag, int ndof, const Vector &crds);
165  Node(const Node &theCopy, bool copyMass = true);
166  Node *getCopy(void) const;
167 
168  // destructor
169  virtual ~Node(void);
170 
171  static DefaultTag &getDefaultTag(void);
172 
173  // public methods dealing with the DOF at the node
174  virtual int getNumberDOF(void) const;
175  virtual void setDOF_GroupPtr(DOF_Group *theDOF_Grp);
176  virtual DOF_Group *getDOF_GroupPtr(void);
177 
178  void connect(ContinuaReprComponent *el) const;
179  void disconnect(ContinuaReprComponent *el) const;
180  ElementConstPtrSet getConnectedElements(void) const;
181  ElementPtrSet getConnectedElements(void);
182  const MeshEdge *next(const std::deque<MeshEdge> &, const std::set<const MeshEdge *> &) const;
183 
184  const bool isDead(void) const;
185  const bool isAlive(void) const;
186  const bool isFrozen(void) const;
187  const bool isFree(void) const;
188  void kill(void);
189  void alive(void);
190  void freeze_if_dead(NodeLocker *locker);
191  void melt_if_alive(NodeLocker *locker);
192 
194  void fix(const std::vector<int> &, const std::vector<double> &);
195  void fix(const ID &, const Vector &);
196 
197  // public methods for obtaining the nodal coordinates
198  virtual size_t getDim(void) const;
199  virtual const Vector &getCrds(void) const;
200  virtual Vector &getCrds(void);
201  Vector getCrds3d(void) const;
202  Pos2d getPosition2d(const Vector &) const;
203  Pos3d getPosition3d(const Vector &) const;
204  Pos2d getInitialPosition2d(void) const;
205  Pos3d getInitialPosition3d(void) const;
206  Pos2d getCurrentPosition2d(const double &factor= 1.0) const;
207  Pos3d getCurrentPosition3d(const double &factor= 1.0) const;
208  bool In(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
209  bool In(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
210  bool Out(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
211  bool Out(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
212  void setPos(const Pos3d &);
213  void Mueve(const Vector3d &desplaz);
214  void Transforma(const TrfGeom &trf);
215  double getDist2(const Pos2d &p,bool initialGeometry= true) const;
216  double getDist(const Pos2d &p,bool initialGeometry= true) const;
217  double getDist2(const Pos3d &p,bool initialGeometry= true) const;
218  double getDist(const Pos3d &p,bool initialGeometry= true) const;
219 
220  // public methods for obtaining committed and trial
221  // response quantities of the node
222  virtual const Vector &getDisp(void) const;
223  Vector getDispXYZ(void) const;
224  Vector getRotXYZ(void) const;
225  virtual const Vector &getVel(void) const;
226  Vector getVelXYZ(void) const;
227  Vector getOmegaXYZ(void) const;
228  virtual const Vector &getAccel(void) const;
229  Vector getAccelXYZ(void) const;
230  Vector getAlphaXYZ(void) const;
231  virtual const Vector &getIncrDisp(void) const;
232  virtual const Vector &getIncrDeltaDisp(void) const;
233 
234  virtual const Vector &getTrialDisp(void) const;
235  virtual const Vector &getTrialVel(void) const;
236  virtual const Vector &getTrialAccel(void) const;
237 
238  // public methods for updating the trial response quantities
239  virtual int setTrialDispComponent(double value, int dof);
240  virtual int setTrialDisp(const Vector &);
241  virtual int setTrialVel(const Vector &);
242  virtual int setTrialAccel(const Vector &);
243 
244  virtual int incrTrialDisp(const Vector &);
245  virtual int incrTrialVel(const Vector &);
246  virtual int incrTrialAccel(const Vector &);
247 
248  const NodalLoad *newLoad(const Vector &v);
249  // public methods for adding and obtaining load information
250  virtual void zeroUnbalancedLoad(void);
251  virtual int addUnbalancedLoad(const Vector &load, double fact = 1.0);
252  virtual int addInertiaLoadToUnbalance(const Vector &accel, double fact = 1.0);
253  virtual const Vector &getUnbalancedLoad(void);
254  virtual const Vector &getUnbalancedLoadIncInertia(void);
255 
256  // public methods dealing with the committed state of the node
257  virtual int commitState();
258  virtual int revertToLastCommit();
259  virtual int revertToStart();
260 
261  // public methods for dynamic analysis
262  virtual const Matrix &getMass(void);
263  virtual int setMass(const Matrix &theMass);
264  virtual int setNumColR(int numCol);
265  virtual int setR(int row, int col, double Value);
266  virtual const Vector &getRV(const Vector &V);
267 
268  virtual int setRayleighDampingFactor(double alphaM);
269  virtual const Matrix &getDamp(void);
270 
271  void addTributary(const double &) const;
272  void resetTributary(void) const;
273  const double &getTributary(void) const;
274 
275  // public methods for eigen vector
276  virtual int setNumEigenvectors(int numVectorsToStore);
277  virtual int setEigenvector(int , const Vector &);
278  inline int getNumModes(void) const
279  { return theEigenvectors.noCols(); }
280  virtual Vector getEigenvector(int ) const;
281  Vector getNormalizedEigenvector(int ) const;
282  virtual const Matrix &getEigenvectors(void);
283  Matrix getNormalizedEigenvectors(void) const;
284  Pos2d getEigenPosition2d(const double &, int) const;
285  Pos3d getEigenPosition3d(const double &, int) const;
286 
287  //Angular frequencies.
288  double getAngularFrequency(int) const;
289  Vector getAngularFrequencies(void) const;
290 
291  //Modal participation factors.
292  virtual double getModalParticipationFactor(int mode) const;
294  virtual double getModalParticipationFactor(int mode,const std::set<int> &dofs) const;
295  Vector getModalParticipationFactors(const std::set<int> &dofs) const;
296  Vector getModalParticipationFactorsForDOFs(const boost::python::list &) const;
297  //Distribution factors.
298  Vector getDistributionFactor(int mode) const;
299  Vector getDistributionFactor(int mode,const std::set<int> &dofs) const;
300  Matrix getDistributionFactors(void) const;
301 
302  //Effective modal masses
303  double getEffectiveModalMass(int mode) const;
304  Vector getEffectiveModalMasses(void) const;
305 
306  //Equivalent static load.
307  Vector getEquivalentStaticLoad(int mode,const double &) const;
308 
309  //Maximum modal values
310  Vector getMaxModalDisplacement(int mode,const double &) const;
311  Vector getMaxModalVelocity(int mode,const double &) const;
312  Vector getMaxModalAcceleration(int mode,const double &) const;
313  Vector getMaxModalDisplacementForDOFs(int mode,const double &,const std::set<int> &dofs) const;
314  Vector getMaxModalVelocityForDOFs(int mode,const double &,const std::set<int> &dofs) const;
315  Vector getMaxModalAccelerationForDOFs(int mode,const double &,const std::set<int> &dofs) const;
316  Vector getMaxModalDisplacementForDOFs(int mode,const double &,const boost::python::list &dofs) const;
317  Vector getMaxModalVelocityForDOFs(int mode,const double &,const boost::python::list &dofs) const;
318  Vector getMaxModalAccelerationForDOFs(int mode,const double &,const boost::python::list &dofs) const;
319 
320  // public methods for output
321  virtual int sendSelf(CommParameters &);
322  virtual int recvSelf(const CommParameters &);
323 
324  std::set<SetBase *> get_sets(void) const;
325  void add_to_sets(std::set<SetBase *> &);
326 
327  virtual void Print(std::ostream &s, int flag = 0);
328 
329  virtual const Vector &getReaction(void) const;
330  const Vector &getResistingForce(const ElementConstPtrSet &,const bool &) const;
331  SlidingVectorsSystem3d getResistingSlidingVectorsSystem3d(const ElementConstPtrSet &,const bool &) const;
332  virtual int addReactionForce(const Vector &, double factor);
333  virtual int resetReactionForce(bool inclInertia);
334  bool checkReactionForce(const double &) const;
335 
336  // AddingSensitivity:BEGIN /////////////////////////////////////////
337  int addInertiaLoadSensitivityToUnbalance(const Vector &accel, double fact = 1.0, bool tag=false);
338  Matrix getMassSensitivity(void);
339  virtual const Matrix &getDampSensitivity(void);
340  int getCrdsSensitivity(void);
341  int saveSensitivity(Vector *v, Vector *vdot, Vector *vdotdot, int gradNum, int numGrads);
342  double getDispSensitivity(int dof, int gradNum);
343  double getVelSensitivity(int dof, int gradNum);
344  double getAccSensitivity(int dof, int gradNum);
345  int setParameter(const std::vector<std::string> &argv, Parameter &param);
346  int updateParameter(int parameterID, Information &info);
347  int activateParameter(int parameterID);
348 
349  // AddingSensitivity:END ///////////////////////////////////////////
350 
351  };
352 
353 Pos3d pos_node(const Node &nod,bool initialGeometry= true);
354 
355 } // end of XC namespace
356 
357 #endif
358 
virtual int setTrialDisp(const Vector &)
Set the current trial displacement.
Definition: Node.cpp:878
virtual const Vector & getUnbalancedLoadIncInertia(void)
Return unbalanced load vector including inertial forces.
Definition: Node.cpp:1075
virtual size_t getDim(void) const
Return the dimension of the node vector position (1,2 or 3).
Definition: Node.cpp:487
Pos3d getInitialPosition3d(void) const
Return the initial position of the node in 3D.
Definition: Node.cpp:537
Vectors that store trial and commited values of the node acceleration.
Definition: NodeAccelVectors.h:41
SFreedom_Constraint * fix(const SFreedom_Constraint &)
Introduce en the node una constraint como la being passed as parameter.
Definition: Node.cpp:426
Vector getEffectiveModalMasses(void) const
Returns the effective modal masses.
Definition: Node.cpp:1537
Vector getMaxModalVelocity(int mode, const double &) const
Return the maximal modal velocity for the mode being passed as parameter and the acceleration corresp...
Definition: Node.cpp:1564
virtual int setMass(const Matrix &theMass)
Set the mass matrix of the node.
Definition: Node.cpp:1224
void melt_if_alive(NodeLocker *locker)
Deletes the constraint over the node DOFs previously created by the freeze method.
Definition: Node.cpp:380
Float vector abstraction.
Definition: Vector.h:93
void disconnect(ContinuaReprComponent *el) const
Removes a component (element, constraint,...) from the connected component list.
Definition: Node.cpp:297
virtual int setTrialAccel(const Vector &)
Set the current trial acceleration.
Definition: Node.cpp:896
virtual const Vector & getDisp(void) const
Returns the displacement of the node.
Definition: Node.cpp:802
Information about an element.
Definition: Information.h:80
virtual void zeroUnbalancedLoad(void)
Causes the node to zero out its unbalanced load vector.
Definition: Node.cpp:961
Vector getCrds3d(void) const
Returns the node coordinates in a 3D space.
Definition: Node.cpp:507
virtual const Vector & getTrialDisp(void) const
Returns the trial value of the displacement of the node.
Definition: Node.cpp:835
virtual Vector getEigenvector(int) const
Returns the eigenvector that corresponds to i-th mode.
Definition: Node.cpp:1367
void add_to_sets(std::set< SetBase *> &)
Adds the node to the sets being passed as parameters.
Definition: Node.cpp:1749
const NodalLoad * newLoad(const Vector &v)
Create a new load on the node and put it in the current load pattern.
Definition: Node.cpp:930
Vector getMaxModalDisplacement(int mode, const double &) const
Returns the maximal modal displacement for the node being passed as parameter and the acceleration co...
Definition: Node.cpp:1557
Vector getAccelXYZ(void) const
Returns the XYZ components of the translational acceleration of the node.
Definition: Node.cpp:719
virtual int revertToStart()
Returns to the initial state.
Definition: Node.cpp:1131
Pos2d getEigenPosition2d(const double &, int) const
Return the "modal" position of the node scaled by factor: return initPos+ factor * getNormalizedEigen...
Definition: Node.cpp:1393
bool In(const GeomObj3d &, const double &factor=1.0, const double &tol=0.0) const
Returns true if the current position of the node scaled by factor: return initPos+ factor * nodDispla...
Definition: Node.cpp:624
virtual const Matrix & getDamp(void)
Return the damping matrix of the node.
Definition: Node.cpp:1170
Vector getDispXYZ(void) const
Returns the XYZ components of node displacement.
Definition: Node.cpp:703
Node * getCopy(void) const
Virtual constructor.
Definition: Node.cpp:305
Matrix getNormalizedEigenvectors(void) const
Returns a matrix with the eigenvectors as columns.
Definition: Node.cpp:1377
const double & getTributary(void) const
Return tributary value (length, area or volume).
Definition: Node.cpp:1197
void Transforma(const TrfGeom &trf)
Applies to the node position the transformation being passed as parameter.
Definition: Node.cpp:779
bool checkReactionForce(const double &) const
Checks that reactions on the node correspond to constrained degrees of freedom.
Definition: Node.cpp:2062
Node(int classTag)
Default constructor.
Definition: Node.cpp:109
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Node.cpp:1846
Pos2d getPosition2d(const Vector &) const
Return the 2D position obtained as: initPos+ factor * v.
Definition: Node.cpp:551
void resetTributary(void) const
Zeroes tributary (length, area or volume).
Definition: Node.cpp:1193
Pos3d pos_node(const Node &nod, bool initialGeometry=true)
Return the 3D position of the node.
Definition: Node.cpp:787
virtual const Vector & getVel(void) const
Returns the velocity of the node.
Definition: Node.cpp:811
virtual int incrTrialVel(const Vector &)
Increments the current trial velocity.
Definition: Node.cpp:916
int sendData(CommParameters &)
Sends object members through the channel being passed as parameter.
Definition: Node.cpp:1638
virtual void setDOF_GroupPtr(DOF_Group *theDOF_Grp)
Sets the DOF_Group pointer.
Definition: Node.cpp:473
Base class for nodes and elements (mesh components).
Definition: MeshComponent.h:43
virtual int revertToLastCommit()
Returns to the last commited state.
Definition: Node.cpp:1116
int noCols() const
Returns the number of columns, numCols, of the Matrix.
Definition: Matrix.h:263
Vector that stores the dbTags of the class members.
Definition: DbTagData.h:43
virtual int setTrialVel(const Vector &)
Set the current trial velocity.
Definition: Node.cpp:887
virtual int setR(int row, int col, double Value)
Sets the (row,col) entry of R to be equal to Value.
Definition: Node.cpp:1257
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: Node.cpp:1659
Vector getOmegaXYZ(void) const
Returns the XYZ components of the angular velocity of the node.
Definition: Node.cpp:715
bool Out(const GeomObj3d &, const double &factor=1.0, const double &tol=0.0) const
Returns true if the current position of the node scaled by factor: return initPos+ factor * nodDispla...
Definition: Node.cpp:642
Vector of integers.
Definition: ID.h:93
DbTagData & getDbTagData(void) const
Returns a vector to store the dbTags de los miembros of the clase.
Definition: Node.cpp:1631
const bool isAlive(void) const
True if node is active.
Definition: Node.cpp:313
virtual const Vector & getUnbalancedLoad(void)
Returns unbalanced load vector.
Definition: Node.cpp:1066
void connect(ContinuaReprComponent *el) const
Inserts a component (element, constraint,...) to the connected component list.
Definition: Node.cpp:287
Vector getAngularFrequencies(void) const
Returns the angular frequencies from the modal analysis.
Definition: Node.cpp:1351
const bool isDead(void) const
True if node is inactive.
Definition: Node.cpp:309
virtual int incrTrialAccel(const Vector &)
Increments the current trial acceleration.
Definition: Node.cpp:926
virtual const Vector & getIncrDeltaDisp(void) const
Return the incremental displacement.
Definition: Node.cpp:866
virtual int commitState()
Commits the state of the node.
Definition: Node.cpp:1100
virtual int resetReactionForce(bool inclInertia)
Calculate the reactions in this node (used in Domain::calculateNodalReactions).
Definition: Node.cpp:2117
void setPos(const Pos3d &)
Sets the node position.
Definition: Node.cpp:760
virtual void Print(std::ostream &s, int flag=0)
Prints node data.
Definition: Node.cpp:1761
virtual ~Node(void)
destructor
Definition: Node.cpp:415
const bool isFree(void) const
returns true if the node has no constraints.
Definition: Node.cpp:399
virtual const Vector & getTrialAccel(void) const
Returns the trial value of the acceleration of the node.
Definition: Node.cpp:851
virtual int addInertiaLoadToUnbalance(const Vector &accel, double fact=1.0)
Adds inertial loads to the node unbalanced load vector.
Definition: Node.cpp:1006
double getDist(const Pos2d &p, bool initialGeometry=true) const
Returns the distance from the node to the point being passed as parameter.
Definition: Node.cpp:737
virtual double getModalParticipationFactor(int mode) const
Returns the modal participation factor corresponding to i-th mode.
Definition: Node.cpp:1409
Vector getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: Node.cpp:1417
Vector getMaxModalDisplacementForDOFs(int mode, const double &, const std::set< int > &dofs) const
Returns the maximal modal displacement on the DOFs and mode being passed as parameter and the acceler...
Definition: Node.cpp:1579
Vector getNormalizedEigenvector(int) const
Returns the eigenvector that corresponds to i-th mode normalized so the maximum values of the compone...
Definition: Node.cpp:1373
virtual int addUnbalancedLoad(const Vector &load, double fact=1.0)
Adds vector to the node unbalanced load.
Definition: Node.cpp:972
Vector getDistributionFactor(int mode) const
Returns the distribution factor corresponding to the i-th mode.
Definition: Node.cpp:1488
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: Node.cpp:1719
Vector getMaxModalAcceleration(int mode, const double &) const
Return the maximal modal acceleration for the mode being passed as parameter and the acceleration cor...
Definition: Node.cpp:1571
int setParameter(const std::vector< std::string > &argv, Parameter &param)
Sets the value param to the parameter argv.
Definition: Node.cpp:1801
virtual const Vector & getAccel(void) const
Returns the acceleration of the node.
Definition: Node.cpp:820
Pos2d getInitialPosition2d(void) const
Return the initial position of the node in 2D.
Definition: Node.cpp:525
Vectores to store trial and commited values («commited») on node displacements.
Definition: NodeDispVectors.h:41
Default tag.
Definition: DefaultTag.h:36
Matrix getDistributionFactors(void) const
Returns the matrix with the computed distribution factors placed by columns.
Definition: Node.cpp:1502
Pos3d getEigenPosition3d(const double &, int) const
Return the "modal" position of the node scaled by factor: return initPos+ factor * getNormalizedEigen...
Definition: Node.cpp:1401
int activateParameter(int parameterID)
Activates the parameter identified by parameterID.
Definition: Node.cpp:1874
const Vector & getResistingForce(const ElementConstPtrSet &, const bool &) const
Return the action of the elements from the set over this node.
Definition: Node.cpp:1981
virtual const Vector & getCrds(void) const
Return a const reference to the vector of nodal coordinates.
Definition: Node.cpp:495
Single freedom constraint.
Definition: SFreedom_Constraint.h:84
Vectores to store trial and commited values of node velocity.
Definition: NodeVelVectors.h:41
Base class for components used to represent the material (continuum).
Definition: ContinuaReprComponent.h:37
virtual int setNumEigenvectors(int numVectorsToStore)
Set the dimensions of the matrix that constains the eigenvectors.
Definition: Node.cpp:1298
virtual int setRayleighDampingFactor(double alphaM)
Sets the Rayleigh dumping factor.
Definition: Node.cpp:1162
A DOF_Group object is instantiated by the ConstraintHandler for every unconstrained node in the domai...
Definition: DOF_Group.h:106
virtual int addReactionForce(const Vector &, double factor)
Increments the node reaction.
Definition: Node.cpp:2040
virtual int sendSelf(CommParameters &)
Send the object through the channel being passed as parameter.
Definition: Node.cpp:1697
virtual const Matrix & getMass(void)
Return the mass matrix of the node.
Definition: Node.cpp:1158
Mesh edge.
Definition: MeshEdge.h:39
void Mueve(const Vector3d &desplaz)
Moves the node (intended only for its use from XC::Set).
Definition: Node.cpp:2130
double getAngularFrequency(int) const
Return the angular frequency corresponding to the i-th mode.
Definition: Node.cpp:1343
Geometric transformation that can be applied to the components of a set.
Definition: TrfGeom.h:49
std::set< Element * > ElementPtrSet
Container of element pointers.
Definition: Node.h:155
std::set< const Element * > ElementConstPtrSet
Container of const element pointers.
Definition: Node.h:156
const MeshEdge * next(const std::deque< MeshEdge > &, const std::set< const MeshEdge *> &) const
Returns an edge that has its origin in this node (and is not in visited).
Definition: Node.cpp:1963
Vector getMaxModalVelocityForDOFs(int mode, const double &, const std::set< int > &dofs) const
Returns the maximum modal velocity on the DOFs and mode being passed as parameter and the acceleratio...
Definition: Node.cpp:1597
virtual const Vector & getRV(const Vector &V)
This is a method provided for Element objects, the Node object returns the product of the matrix R an...
Definition: Node.cpp:1276
void freeze_if_dead(NodeLocker *locker)
Imposes zero displacement (zero value for all node DOFs).
Definition: Node.cpp:345
Vector getEquivalentStaticLoad(int mode, const double &) const
Return the equivalent static load for the mode being passed as parameter and the acceleration corresp...
Definition: Node.cpp:1548
Single freedom constraints that make part of a load pattern.
Definition: NodeLocker.h:44
virtual int incrTrialDisp(const Vector &)
Increments the current trial displacement.
Definition: Node.cpp:906
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
virtual DOF_Group * getDOF_GroupPtr(void)
Gets the DOF_Group pointer.
Definition: Node.cpp:482
virtual const Vector & getReaction(void) const
Return the node reaction.
Definition: Node.cpp:2036
virtual int setEigenvector(int, const Vector &)
Set el eigenvector eigenvector al modo mode.
Definition: Node.cpp:1318
ElementConstPtrSet getConnectedElements(void) const
Return a list of pointers to the elements that are connected with this node.
Definition: Node.cpp:1937
Vector getMaxModalAccelerationForDOFs(int mode, const double &, const std::set< int > &dofs) const
Return the maximal modal acceleration on the DOFs and mode being passed as parameter and the accelera...
Definition: Node.cpp:1614
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
Vector getModalParticipationFactorsForDOFs(const boost::python::list &) const
Returns the modal participation factors.
Definition: Node.cpp:1481
Pos3d getCurrentPosition3d(const double &factor=1.0) const
Return the current position of the node scaled by factor: return initPos+ factor * nodDisplacement...
Definition: Node.cpp:615
Vector getRotXYZ(void) const
Returns the XYZ components of node rotation.
Definition: Node.cpp:707
Definition: Parameter.h:65
const bool isFrozen(void) const
returns true if the node is frozen.
Definition: Node.cpp:395
Vector getVelXYZ(void) const
Returns the XYZ components of the translational velocity of the node.
Definition: Node.cpp:711
virtual const Vector & getTrialVel(void) const
Returns the trial value of the velocity of the node.
Definition: Node.cpp:843
virtual int getNumberDOF(void) const
Return the number of node DOFs, associated with the node.
Definition: Node.cpp:464
Vector getAlphaXYZ(void) const
Returns the XYZ components of the angular acceleration of the node.
Definition: Node.cpp:723
Mesh node.
Definition: Node.h:110
double getDist2(const Pos2d &p, bool initialGeometry=true) const
Returns the square of the distance from the node to the point being passed as parameter.
Definition: Node.cpp:730
virtual int setNumColR(int numCol)
Creates a Matrix to store the R matrix.
Definition: Node.cpp:1245
SlidingVectorsSystem3d getResistingSlidingVectorsSystem3d(const ElementConstPtrSet &, const bool &) const
Returns the sliding vector system that represents the action of the elements of the set over the node...
Definition: Node.cpp:2013
Pos3d getPosition3d(const Vector &) const
Return the 3D position obtained as: initPos+ factor * v.
Definition: Node.cpp:573
Pos2d getCurrentPosition2d(const double &factor=1.0) const
Return the current position of the node scaled by factor: return initPos+ factor * nodDisplacement...
Definition: Node.cpp:607
virtual const Matrix & getEigenvectors(void)
Returns the eigenvectors that correspond to the node.
Definition: Node.cpp:1363
virtual const Vector & getIncrDisp(void) const
Returns the displacement increment.
Definition: Node.cpp:855
Load over a node.
Definition: NodalLoad.h:82
static DefaultTag & getDefaultTag(void)
Returns a default value for node identifier.
Definition: Node.cpp:421
double getEffectiveModalMass(int mode) const
Return the effective modal mass that corresponds to i mode.
Definition: Node.cpp:1525
void addTributary(const double &) const
Adds to tributary the lentgth, area o volume being passed as parameter.
Definition: Node.cpp:1189
std::set< SetBase * > get_sets(void) const
Returns the sets that include this node.
Definition: Node.cpp:1733