xc
Node.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.9 $
49 // $Date: 2005/11/23 22:48:50 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/domain/mesh/node/Node.h,v $
51 
52 
53 #ifndef Node_h
54 #define Node_h
55 
56 // Written: fmk
57 // Created: 11/96
58 // Revision: A
59 //
60 // Purpose: This file contains the class interface for Node.
61 // A Node provides the abstraction for a node in the FEM.
62 // Nodes have original position, trial displacement, velocity and
63 // acceleration and committed displacement, velocity and acceleration
64 // (the last committed() trial quantities).
65 //
66 // What: "@(#) Node.h, revA"
67 
68 #include "domain/mesh/MeshComponent.h"
69 #include "NodeDispVectors.h"
70 #include "NodeVelVectors.h"
71 #include "NodeAccelVectors.h"
72 #include "utility/matrix/Matrix.h"
73 #include <boost/python/list.hpp>
74 
75 class Pos2d;
76 class Pos3d;
78 class GeomObj2d;
79 class GeomObj3d;
80 
81 namespace XC {
82 class ContinuaReprComponent;
83 class Matrix;
84 class Channel;
85 class Renderer;
86 class TrfGeom;
87 class NodeLocker;
88 class DefaultTag;
89 class SetBase;
90 class MeshEdge;
91 class DOF_Group;
92 class DqPtrsElem;
93 class Constraint;
94 
99 //
111 class Node: public MeshComponent
112  {
113  public:
114  // typedefs.
115  typedef std::set<Element *> ElementPtrSet;
116  typedef std::set<const Element *> ElementConstPtrSet;
117  private:
118  // private data associated with each node object
119  int numberDOF;
120  DOF_Group *theDOF_GroupPtr;
121  Vector Crd;
122 
123  NodeDispVectors disp;
124  NodeVelVectors vel;
125  NodeAccelVectors accel;
126 
127  Matrix R;
128  Matrix mass;
129  mutable Vector unbalLoad;
130  mutable Vector unbalLoadWithInertia;
131  mutable Vector reaction;
132  double alphaM;
133  mutable double tributary;
134 
135  Matrix theEigenvectors; //Eigenvectors matrix.
136 
137  // AddingSensitivity:BEGIN /////////////////////////////////////////
138  Matrix dispSensitivity;
139  Matrix velSensitivity;
140  Matrix accSensitivity;
141  int parameterID;
142  // AddingSensitivity:END ///////////////////////////////////////////
143 
144  static std::deque<Matrix> theMatrices;
145 
146  mutable std::set<ContinuaReprComponent *> connected;
147 
148  std::map<int, std::set<int> > freeze_constraints;
149  ID get_node_lockers_tags(void) const;
150  std::vector<ID> get_constraints_tags(void) const;
151  void set_id_constraints(const ID &, const std::vector<ID> &);
152  bool is_a_freeze_constraint(const ContinuaReprComponent *) const;
153 
154  Matrix get_element_stiff(const ElementConstPtrSet &,bool initial) const;
155  Matrix get_constraints_stiff(void) const;
156 
157  static DefaultTag defaultTag; //<! tag for next new node.
158  protected:
159 
160  DbTagData &getDbTagData(void) const;
161  int sendData(Communicator &);
162  int recvData(const Communicator &);
163  public:
164  // constructors
165  Node(int classTag);
166  Node(int tag, int classTag);
167  Node(int tag, int ndof, double Crd1);
168  Node(int tag, int ndof, double Crd1, double Crd2);
169  Node(int tag, int ndof, double Crd1, double Crd2, double Crd3);
170  Node(int tag, int ndof, const Vector &crds);
171  Node(const Node &theCopy, bool copyMass = true);
172  Node *getCopy(void) const;
173 
174  // destructor
175  virtual ~Node(void);
176 
177  static DefaultTag &getDefaultTag(void);
178 
179  // public methods dealing with the DOF at the node
180  virtual int getNumberDOF(void) const;
181  virtual void setDOF_GroupPtr(DOF_Group *theDOF_Grp);
182  virtual DOF_Group *getDOF_GroupPtr(void);
183  const ID &getDOFs(void) const;
184 
185  size_t getNumberOfConnectedConstraints(void) const;
186  std::set<const Constraint *>getConnectedConstraints(void) const;
187  std::set<Constraint *>getConnectedConstraints(void);
188  boost::python::list getConnectedConstraintsPy(void);
189  boost::python::list getConnectedConstraintTags(void) const;
190 
191  void connect(ContinuaReprComponent *el) const;
192  void disconnect(ContinuaReprComponent *el) const;
193  size_t getNumberOfConnectedElements(void) const;
194  ElementConstPtrSet getConnectedElements(void) const;
195  ElementPtrSet getConnectedElements(void);
196  boost::python::list getConnectedElementsPy(void);
197  boost::python::list getConnectedElementTags(void) const;
198  const MeshEdge *next(const std::deque<MeshEdge> &, const std::set<const MeshEdge *> &) const;
199 
200  const bool isDead(void) const;
201  const bool isAlive(void) const;
202  const bool isFrozen(void) const;
203  const bool isFree(void) const;
204  const bool isFixedOnlyByFreezeConstraints(void) const;
205  void kill(void);
206  void alive(void);
207  void freeze_if_dead(NodeLocker *locker);
208  void melt_if_alive(void);
209 
211  void fix(const std::vector<int> &, const std::vector<double> &);
212  void fix(const ID &, const Vector &);
213 
214  // public methods for obtaining the nodal coordinates
215  virtual size_t getDim(void) const;
216  virtual const Vector &getCrds(void) const;
217  virtual Vector &getCrds(void);
218  Vector getCrds3d(void) const;
219  Pos2d getPosition2d(const Vector &) const;
220  Pos3d getPosition3d(const Vector &) const;
221  Pos2d getInitialPosition2d(void) const;
222  Pos3d getInitialPosition3d(void) const;
223  Pos2d getCurrentPosition2d(const double &factor= 1.0) const;
224  Pos3d getCurrentPosition3d(const double &factor= 1.0) const;
225  bool In(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
226  bool In(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
227  bool Out(const GeomObj3d &,const double &factor= 1.0, const double &tol= 0.0) const;
228  bool Out(const GeomObj2d &,const double &factor= 1.0, const double &tol= 0.0) const;
229  void setPos(const Pos3d &);
230  void Move(const Vector3d &desplaz);
231  void Transform(const TrfGeom &trf);
232  double getDist2(const Pos2d &p,bool initialGeometry= true) const;
233  double getDist(const Pos2d &p,bool initialGeometry= true) const;
234  double getDist2(const Pos3d &p,bool initialGeometry= true) const;
235  double getDist(const Pos3d &p,bool initialGeometry= true) const;
236 
237  // public methods for obtaining committed and trial
238  // response quantities of the node
239  virtual const Vector &getDisp(void) const;
240  Vector getDispXYZ(void) const;
241  Vector getRotXYZ(void) const;
242  virtual const Vector &getVel(void) const;
243  Vector getVelXYZ(void) const;
244  Vector getOmegaXYZ(void) const;
245  virtual const Vector &getAccel(void) const;
246  Vector getAccelXYZ(void) const;
247  Vector getAlphaXYZ(void) const;
248  virtual const Vector &getIncrDisp(void) const;
249  virtual const Vector &getIncrDeltaDisp(void) const;
250 
251  virtual const Vector &getTrialDisp(void) const;
252  virtual const Vector &getTrialVel(void) const;
253  virtual const Vector &getTrialAccel(void) const;
254 
255  // public methods for updating the trial response quantities
256  virtual int setTrialDispComponent(double value, int dof);
257  virtual int setTrialDisp(const Vector &);
258  virtual int setTrialVel(const Vector &);
259  virtual int setTrialAccel(const Vector &);
260 
261  virtual int incrTrialDisp(const Vector &);
262  virtual int incrTrialVel(const Vector &);
263  virtual int incrTrialAccel(const Vector &);
264 
265  const NodalLoad *newLoad(const Vector &);
266  void createInertiaLoad(const Vector &);
267  // public methods for adding and obtaining load information
268  virtual void zeroUnbalancedLoad(void);
269  virtual int addUnbalancedLoad(const Vector &load, double fact = 1.0);
270  virtual int addInertiaLoadToUnbalance(const Vector &accel, double fact = 1.0);
271  virtual const Vector &getUnbalancedLoad(void) const;
272  virtual const Vector &getUnbalancedLoadIncInertia(void) const;
273 
274  // public methods dealing with the committed state of the node
275  virtual int commitState();
276  virtual int revertToLastCommit();
277  virtual int revertToStart();
278 
279  // public methods for dynamic analysis
280  virtual const Matrix &getMass(void) const;
281  double getMassComponent(const int &) const;
282  virtual int setMass(const Matrix &theMass);
283  virtual int setNumColR(int numCol);
284  virtual int setR(int row, int col, double Value);
285  virtual const Vector &getRV(const Vector &V);
286 
287  virtual int setRayleighDampingFactor(double alphaM);
288  virtual const Matrix &getDamp(void) const;
289 
290  void addTributary(const double &) const;
291  void resetTributary(void) const;
292  const double &getTributary(void) const;
293 
294  // public methods for eigen vector
295  virtual int setNumEigenvectors(int numVectorsToStore);
296  virtual int setEigenvector(int , const Vector &);
297  inline int getNumModes(void) const
298  { return theEigenvectors.noCols(); }
299  virtual Vector getEigenvector(int ) const;
300  Vector getNormalizedEigenvector(int ) const;
301  virtual const Matrix &getEigenvectors(void) const;
302  Matrix getNormalizedEigenvectors(void) const;
303  Pos2d getEigenPosition2d(const double &, int) const;
304  Pos3d getEigenPosition3d(const double &, int) const;
307  void clearEigenvectors(void);
308 
309  //Angular frequencies.
310  double getAngularFrequency(int) const;
311  Vector getAngularFrequencies(void) const;
312 
313  //Modal participation factors.
314  virtual double getModalParticipationFactor(int mode) const;
316  virtual double getModalParticipationFactor(int mode,const std::set<int> &dofs) const;
317  Vector getModalParticipationFactors(const std::set<int> &dofs) const;
318  Vector getModalParticipationFactorsForDOFs(const boost::python::list &) const;
319  //Distribution factors.
320  Vector getDistributionFactor(int mode) const;
321  Vector getDistributionFactor(int mode,const std::set<int> &dofs) const;
322  Matrix getDistributionFactors(void) const;
323 
324  //Effective modal masses
325  double getEffectiveModalMass(int mode) const;
326  Vector getEffectiveModalMasses(void) const;
327 
328  //Equivalent static load.
329  Vector getEquivalentStaticLoad(int mode,const double &) const;
330 
331  //Maximum modal values
332  Vector getMaxModalDisplacement(int mode,const double &) const;
333  Vector getMaxModalVelocity(int mode,const double &) const;
334  Vector getMaxModalAcceleration(int mode,const double &) const;
335  Vector getMaxModalDisplacementForDOFs(int mode,const double &,const std::set<int> &dofs) const;
336  Vector getMaxModalVelocityForDOFs(int mode,const double &,const std::set<int> &dofs) const;
337  Vector getMaxModalAccelerationForDOFs(int mode,const double &,const std::set<int> &dofs) const;
338  Vector getMaxModalDisplacementForDOFs(int mode,const double &,const boost::python::list &dofs) const;
339  Vector getMaxModalVelocityForDOFs(int mode,const double &,const boost::python::list &dofs) const;
340  Vector getMaxModalAccelerationForDOFs(int mode,const double &,const boost::python::list &dofs) const;
341 
342  // public methods for output
343  virtual int sendSelf(Communicator &);
344  virtual int recvSelf(const Communicator &);
345  boost::python::dict getPyDict(void) const;
346  void setPyDict(const boost::python::dict &);
347 
348  std::set<SetBase *> get_sets(void) const;
349  void add_to_sets(std::set<SetBase *> &);
350  void copySetsFrom(const Node &);
351 
352  virtual void Print(std::ostream &s, int flag = 0) const;
353 
354  Vector3d get3dForceComponents(const Vector &) const;
355  Vector3d get3dMomentComponents(const Vector &) const;
356 
357  virtual const Vector &getReaction(void) const;
358  void setReaction(const Vector &);
359  Vector3d getReactionForce3d(void) const;
360  Vector3d getReactionMoment3d(void) const;
361  const Vector &getResistingForce(const ElementConstPtrSet &,const bool &) const;
362  SlidingVectorsSystem3d getResistingSlidingVectorsSystem3d(const ElementConstPtrSet &,const bool &) const;
363  virtual int addReactionForce(const Vector &, double factor);
364  virtual int resetReactionForce(bool inclInertia);
365  bool checkReactionForce(const double &) const;
366  Matrix getTangentStiff(const ElementConstPtrSet &) const;
367  Matrix getInitialStiff(const ElementConstPtrSet &) const;
368 
369  // AddingSensitivity:BEGIN /////////////////////////////////////////
370  int addInertiaLoadSensitivityToUnbalance(const Vector &accel, double fact = 1.0, bool tag=false);
371  Matrix getMassSensitivity(void) const;
372  virtual const Matrix &getDampSensitivity(void) const;
373  int getCrdsSensitivity(void) const;
374  int saveSensitivity(Vector *v, Vector *vdot, Vector *vdotdot, int gradNum, int numGrads);
375  double getDispSensitivity(int dof, int gradNum) const;
376  double getVelSensitivity(int dof, int gradNum) const;
377  double getAccSensitivity(int dof, int gradNum) const;
378  int setParameter(const std::vector<std::string> &argv, Parameter &param);
379  int updateParameter(int parameterID, Information &info);
380  int activateParameter(int parameterID);
381 
382  // AddingSensitivity:END ///////////////////////////////////////////
383 
384  };
385 
386 Pos3d pos_node(const Node &nod,bool initialGeometry= true);
387 
388 } // end of XC namespace
389 
390 #endif
391 
virtual int setTrialDisp(const Vector &)
Set the current trial displacement.
Definition: Node.cpp:1065
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: Node.cpp:2008
Vector3d getEigenvectorDisp3dComponents(int) const
Return the eigenvector displacement components in a 3D vector.
Definition: Node.cpp:1626
virtual size_t getDim(void) const
Return the dimension of the node vector position (1,2 or 3).
Definition: Node.cpp:656
Vector3d get3dMomentComponents(const Vector &) const
Return the "moment/rotation" components in a 3d vector.
Definition: Node.cpp:2613
Pos3d getInitialPosition3d(void) const
Return the initial position of the node in 3D.
Definition: Node.cpp:706
Vectors that store trial and committed values of the node acceleration.
Definition: NodeAccelVectors.h:42
SFreedom_Constraint * fix(const SFreedom_Constraint &)
Introduce en the node una constraint como la being passed as parameter.
Definition: Node.cpp:577
void createInertiaLoad(const Vector &)
Creates the inertia load that corresponds to the acceleration argument.
Definition: Node.cpp:1143
Vector getEffectiveModalMasses(void) const
Returns the effective modal masses.
Definition: Node.cpp:1764
Vector3d get3dForceComponents(const Vector &) const
Return the "force/displacement" components in a 3d vector.
Definition: Node.cpp:2573
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:1791
virtual int setMass(const Matrix &theMass)
Set the mass matrix of the node.
Definition: Node.cpp:1430
Float vector abstraction.
Definition: Vector.h:94
boost::python::list getConnectedElementTags(void) const
Return a python list containing the tags of the elements that are connected with this node...
Definition: Node.cpp:2375
void disconnect(ContinuaReprComponent *el) const
Removes a component (element, constraint,...) from the connected component list.
Definition: Node.cpp:312
virtual int setTrialAccel(const Vector &)
Set the current trial acceleration.
Definition: Node.cpp:1083
Vector3d getReactionMoment3d(void) const
Return the reaction moment in a 3D vector.
Definition: Node.cpp:2569
virtual const Vector & getDisp(void) const
Returns the displacement of the node.
Definition: Node.cpp:989
void copySetsFrom(const Node &)
Add this element to all the sets containing the given one.
Definition: Node.cpp:2056
Information about an element.
Definition: Information.h:81
virtual void zeroUnbalancedLoad(void)
Causes the node to zero out its unbalanced load vector.
Definition: Node.cpp:1154
Communication parameters between processes.
Definition: Communicator.h:66
Vector getCrds3d(void) const
Returns the node coordinates in a 3D space.
Definition: Node.cpp:676
virtual int sendSelf(Communicator &)
Send the object through the communicator argument.
Definition: Node.cpp:1985
boost::python::list getConnectedConstraintTags(void) const
Return a Python list with the tags of the constraints connected to the node.
Definition: Node.cpp:2302
virtual const Vector & getTrialDisp(void) const
Returns the trial value of the displacement of the node.
Definition: Node.cpp:1022
virtual Vector getEigenvector(int) const
Returns the eigenvector that corresponds to i-th mode.
Definition: Node.cpp:1579
void add_to_sets(std::set< SetBase *> &)
Adds the node to the sets being passed as parameters.
Definition: Node.cpp:2040
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:1784
Vector getAccelXYZ(void) const
Returns the XYZ components of the translational acceleration of the node.
Definition: Node.cpp:906
virtual int revertToStart()
Returns to the initial state.
Definition: Node.cpp:1325
Pos2d getEigenPosition2d(const double &, int) const
Return the "modal" position of the node scaled by factor: return initPos+ factor * getNormalizedEigen...
Definition: Node.cpp:1611
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:811
Posición en dos dimensiones.
Definition: Pos2d.h:41
virtual const Matrix & getEigenvectors(void) const
Returns the eigenvectors that correspond to the node.
Definition: Node.cpp:1575
Vector getDispXYZ(void) const
Returns the XYZ components of node displacement.
Definition: Node.cpp:890
Node * getCopy(void) const
Virtual constructor.
Definition: Node.cpp:320
Matrix getNormalizedEigenvectors(void) const
Returns a matrix with the eigenvectors as columns.
Definition: Node.cpp:1589
const double & getTributary(void) const
Return tributary value (length, area or volume).
Definition: Node.cpp:1403
bool checkReactionForce(const double &) const
Checks that reactions on the node correspond to constrained degrees of freedom.
Definition: Node.cpp:2677
Node(int classTag)
Default constructor.
Definition: Node.cpp:108
void clearEigenvectors(void)
Remove the stored eigenvectors.
Definition: Node.cpp:1604
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Node.cpp:2152
Pos2d getPosition2d(const Vector &) const
Return the 2D position obtained as: initPos+ factor * v.
Definition: Node.cpp:720
void resetTributary(void) const
Zeroes tributary (length, area or volume).
Definition: Node.cpp:1399
Pos3d pos_node(const Node &nod, bool initialGeometry=true)
Return the 3D position of the node.
Definition: Node.cpp:974
virtual const Vector & getVel(void) const
Returns the velocity of the node.
Definition: Node.cpp:998
virtual int incrTrialVel(const Vector &)
Increments the current trial velocity.
Definition: Node.cpp:1103
virtual void setDOF_GroupPtr(DOF_Group *theDOF_Grp)
Sets the DOF_Group pointer.
Definition: Node.cpp:627
boost::python::list getConnectedElementsPy(void)
Return a python list of pointers to the elements that are connected with this node.
Definition: Node.cpp:2360
boost::python::dict getPyDict(void) const
Return a Python dictionary with the object members values.
Definition: Node.cpp:1914
Base class for nodes and elements (mesh components).
Definition: MeshComponent.h:44
virtual int revertToLastCommit()
Returns to the last committed state.
Definition: Node.cpp:1310
int noCols() const
Returns the number of columns, numCols, of the Matrix.
Definition: Matrix.h:273
Vector that stores the dbTags of the class members.
Definition: DbTagData.h:44
virtual int setTrialVel(const Vector &)
Set the current trial velocity.
Definition: Node.cpp:1074
virtual int setR(int row, int col, double Value)
Sets the (row,col) entry of R to be equal to Value.
Definition: Node.cpp:1464
Vector getOmegaXYZ(void) const
Returns the XYZ components of the angular velocity of the node.
Definition: Node.cpp:902
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:829
virtual const Matrix & getMass(void) const
Return the mass matrix of the node.
Definition: Node.cpp:1352
Vector of integers.
Definition: ID.h:95
DbTagData & getDbTagData(void) const
Returns a vector to store the dbTags de los miembros of the clase.
Definition: Node.cpp:1858
const bool isAlive(void) const
True if node is active.
Definition: Node.cpp:330
virtual const Matrix & getDamp(void) const
Return the damping matrix of the node.
Definition: Node.cpp:1376
void connect(ContinuaReprComponent *el) const
Inserts a component (element, constraint,...) to the connected component list.
Definition: Node.cpp:301
size_t getNumberOfConnectedConstraints(void) const
Return the nunber of connected constraints.
Definition: Node.cpp:2246
Vector getAngularFrequencies(void) const
Returns the angular frequencies from the modal analysis.
Definition: Node.cpp:1563
const bool isDead(void) const
True if node is inactive.
Definition: Node.cpp:324
double getAccSensitivity(int dof, int gradNum) const
Return acceleration sensitivity.
Definition: Node.cpp:2236
virtual int incrTrialAccel(const Vector &)
Increments the current trial acceleration.
Definition: Node.cpp:1113
virtual const Vector & getIncrDeltaDisp(void) const
Return the incremental displacement.
Definition: Node.cpp:1053
virtual int commitState()
Commits the state of the node.
Definition: Node.cpp:1294
void setReaction(const Vector &)
Set the node reaction.
Definition: Node.cpp:2548
virtual int resetReactionForce(bool inclInertia)
Calculate the reactions in this node (used in Domain::calculateNodalReactions).
Definition: Node.cpp:2739
const NodalLoad * newLoad(const Vector &)
Create a new load on the node and put it in the current load pattern.
Definition: Node.cpp:1117
void setPos(const Pos3d &)
Sets the node position.
Definition: Node.cpp:947
virtual ~Node(void)
destructor
Definition: Node.cpp:566
const bool isFree(void) const
returns true if the node has no constraints.
Definition: Node.cpp:548
void Move(const Vector3d &desplaz)
Moves the node (intended only for its use from XC::Set).
Definition: Node.cpp:2752
virtual const Vector & getTrialAccel(void) const
Returns the trial value of the acceleration of the node.
Definition: Node.cpp:1038
virtual int addInertiaLoadToUnbalance(const Vector &accel, double fact=1.0)
Adds inertial loads to the node unbalanced load vector.
Definition: Node.cpp:1199
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:924
const ID & getDOFs(void) const
Gets the equation numbers corresponding to the node DOFS.
Definition: Node.cpp:640
virtual double getModalParticipationFactor(int mode) const
Returns the modal participation factor corresponding to i-th mode.
Definition: Node.cpp:1635
Vector getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: Node.cpp:1643
Base class for the two-dimensional geometric objects.
Definition: GeomObj2d.h:37
double getMassComponent(const int &) const
Return the mass matrix component for the DOF argument.
Definition: Node.cpp:1356
void melt_if_alive(void)
Deletes the constraint over the node DOFs previously created by the freeze method.
Definition: Node.cpp:513
void Transform(const TrfGeom &trf)
Applies to the node position the transformation being passed as parameter.
Definition: Node.cpp:966
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:1806
Vector getNormalizedEigenvector(int) const
Returns the eigenvector that corresponds to i-th mode normalized so the maximum values of the compone...
Definition: Node.cpp:1585
virtual int addUnbalancedLoad(const Vector &load, double fact=1.0)
Adds vector to the node unbalanced load.
Definition: Node.cpp:1165
virtual void Print(std::ostream &s, int flag=0) const
Prints node data.
Definition: Node.cpp:2066
Vector getDistributionFactor(int mode) const
Returns the distribution factor corresponding to the i-th mode.
Definition: Node.cpp:1715
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:1798
virtual const Vector & getUnbalancedLoad(void) const
Returns unbalanced load vector.
Definition: Node.cpp:1260
int setParameter(const std::vector< std::string > &argv, Parameter &param)
Sets the value param to the parameter argv.
Definition: Node.cpp:2106
virtual const Vector & getAccel(void) const
Returns the acceleration of the node.
Definition: Node.cpp:1007
Pos2d getInitialPosition2d(void) const
Return the initial position of the node in 2D.
Definition: Node.cpp:694
Vectors to store trial and committed values («committed») on node displacements.
Definition: NodeDispVectors.h:42
Default tag.
Definition: DefaultTag.h:37
Matrix getDistributionFactors(void) const
Returns the matrix with the computed distribution factors placed by columns.
Definition: Node.cpp:1729
Pos3d getEigenPosition3d(const double &, int) const
Return the "modal" position of the node scaled by factor: return initPos+ factor * getNormalizedEigen...
Definition: Node.cpp:1619
int activateParameter(int parameterID)
Activates the parameter identified by parameterID.
Definition: Node.cpp:2180
boost::python::list getConnectedConstraintsPy(void)
Return a Python list of pointers to the constraints connected to the node.
Definition: Node.cpp:2285
const Vector & getResistingForce(const ElementConstPtrSet &, const bool &) const
Return the action of the elements from the set over this node.
Definition: Node.cpp:2406
virtual const Vector & getCrds(void) const
Return a const reference to the vector of nodal coordinates.
Definition: Node.cpp:664
Single freedom constraint.
Definition: SFreedom_Constraint.h:85
Vectors to store trial and committed values of node velocity.
Definition: NodeVelVectors.h:42
Vector3d getEigenvectorRot3dComponents(int) const
Return the eigenvector rotation components in a 3D vector.
Definition: Node.cpp:1630
double getDispSensitivity(int dof, int gradNum) const
Return displacement sensitivity.
Definition: Node.cpp:2218
Base class for components used to represent the material (continuum).
Definition: ContinuaReprComponent.h:39
virtual int setNumEigenvectors(int numVectorsToStore)
Set the dimensions of the matrix that constrains the eigenvectors.
Definition: Node.cpp:1507
virtual const Vector & getUnbalancedLoadIncInertia(void) const
Return unbalanced load vector including inertial forces.
Definition: Node.cpp:1269
virtual int setRayleighDampingFactor(double alphaM)
Sets the Rayleigh dumping factor.
Definition: Node.cpp:1368
A DOF_Group object is instantiated by the ConstraintHandler for every unconstrained node in the domai...
Definition: DOF_Group.h:107
virtual int addReactionForce(const Vector &, double factor)
Increments the node reaction.
Definition: Node.cpp:2649
const bool isFixedOnlyByFreezeConstraints(void) const
Return true if the node is fixed only by its freeze constraints (all the other constraints, if any, are dead)
Definition: Node.cpp:364
Mesh edge.
Definition: MeshEdge.h:40
double getAngularFrequency(int) const
Return the angular frequency corresponding to the i-th mode.
Definition: Node.cpp:1555
Geometric transformation that can be applied to the components of a set.
Definition: TrfGeom.h:48
std::set< Element * > ElementPtrSet
Container of element pointers.
Definition: Node.h:115
std::set< const Element * > ElementConstPtrSet
Container of const element pointers.
Definition: Node.h:116
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:2388
int addInertiaLoadSensitivityToUnbalance(const Vector &accel, double fact=1.0, bool tag=false)
??
Definition: Node.cpp:1223
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:1824
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:1484
void freeze_if_dead(NodeLocker *locker)
Imposes zero displacement (zero value for all node DOFs).
Definition: Node.cpp:394
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:1775
Posición en tres dimensiones.
Definition: Pos3d.h:44
Single freedom constraints that make part of a load pattern.
Definition: NodeLocker.h:45
Sliding vectors system en un espacio tridimensional.
Definition: SlidingVectorsSystem3d.h:39
void setPyDict(const boost::python::dict &)
Set the values of the object members from a Python dictionary.
Definition: Node.cpp:1941
Matrix getTangentStiff(const ElementConstPtrSet &) const
Return the tangent stiffness contribution of the elements argument.
Definition: Node.cpp:2531
virtual int incrTrialDisp(const Vector &)
Increments the current trial displacement.
Definition: Node.cpp:1093
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
virtual DOF_Group * getDOF_GroupPtr(void)
Gets the DOF_Group pointer.
Definition: Node.cpp:636
virtual const Vector & getReaction(void) const
Return the node reaction.
Definition: Node.cpp:2544
virtual int setEigenvector(int, const Vector &)
Set el eigenvector eigenvector al modo mode.
Definition: Node.cpp:1528
ElementConstPtrSet getConnectedElements(void) const
Return a list of pointers to the elements that are connected with this node.
Definition: Node.cpp:2332
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:1841
Matrix of floats.
Definition: Matrix.h:111
Vector getModalParticipationFactorsForDOFs(const boost::python::list &) const
Returns the modal participation factors.
Definition: Node.cpp:1708
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:793
Vector getRotXYZ(void) const
Returns the XYZ components of node rotation.
Definition: Node.cpp:894
std::set< const Constraint * > getConnectedConstraints(void) const
Return a set of const pointers the constraints connected to the node.
Definition: Node.cpp:2259
Parameter.
Definition: Parameter.h:68
const bool isFrozen(void) const
returns true if the node is frozen.
Definition: Node.cpp:544
Vector getVelXYZ(void) const
Returns the XYZ components of the translational velocity of the node.
Definition: Node.cpp:898
void alive(void)
Activates the component.
Definition: Node.cpp:558
virtual const Vector & getTrialVel(void) const
Returns the trial value of the velocity of the node.
Definition: Node.cpp:1030
virtual int getNumberDOF(void) const
Return the number of node DOFs, associated with the node.
Definition: Node.cpp:618
Vector getAlphaXYZ(void) const
Returns the XYZ components of the angular acceleration of the node.
Definition: Node.cpp:910
void kill(void)
Deactivates the component.
Definition: Node.cpp:551
size_t getNumberOfConnectedElements(void) const
Return the number of elements that are connected with this node.
Definition: Node.cpp:2319
int sendData(Communicator &)
Sends object members through the communicator argument.
Definition: Node.cpp:1865
Mesh node.
Definition: Node.h:111
Vector3d getReactionForce3d(void) const
Return the reaction force in a 3D vector.
Definition: Node.cpp:2565
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:917
virtual int setNumColR(int numCol)
Creates a Matrix to store the R matrix.
Definition: Node.cpp:1452
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:2439
Pos3d getPosition3d(const Vector &) const
Return the 3D position obtained as: initPos+ factor * v.
Definition: Node.cpp:742
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:776
virtual const Vector & getIncrDisp(void) const
Returns the displacement increment.
Definition: Node.cpp:1042
Load over a node.
Definition: NodalLoad.h:83
Matrix getInitialStiff(const ElementConstPtrSet &) const
Return the initial stiffness contribution of the elements argument.
Definition: Node.cpp:2540
static DefaultTag & getDefaultTag(void)
Returns a default value for node identifier.
Definition: Node.cpp:572
Vector en tres dimensiones.
Definition: Vector3d.h:39
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: Node.cpp:1888
double getVelSensitivity(int dof, int gradNum) const
Return velocity sensitivity.
Definition: Node.cpp:2227
double getEffectiveModalMass(int mode) const
Return the effective modal mass that corresponds to i mode.
Definition: Node.cpp:1752
void addTributary(const double &) const
Adds to tributary the lentgth, area o volume being passed as parameter.
Definition: Node.cpp:1395
Clase base para los objetos en tres dimensiones.
Definition: GeomObj3d.h:43
std::set< SetBase * > get_sets(void) const
Returns the sets that include this node.
Definition: Node.cpp:2023