XC Open source finite element analysis program
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 SVD3d;
77 
78 namespace XC {
79 class ContinuaReprComponent;
80 class Matrix;
81 class Channel;
82 class Renderer;
83 class TrfGeom;
84 class NodeLocker;
85 class DefaultTag;
86 class ElementEdges;
87 class SetBase;
88 class MeshEdge;
89 class DOF_Group;
90 class DqPtrsElem;
91 
96 //
99 class Node: public MeshComponent
100  {
101  private:
102  // private data associated with each node object
103  int numberDOF;
104  DOF_Group *theDOF_GroupPtr;
105  Vector Crd;
106 
107  NodeDispVectors disp;
108  NodeVelVectors vel;
109  NodeAccelVectors accel;
110 
111  Matrix R;
112  Matrix mass;
113  mutable Vector unbalLoad;
114  Vector unbalLoadWithInertia;
115  mutable Vector reaction;
116  double alphaM;
117  mutable double tributaria;
118 
119  Matrix theEigenvectors; //Eigenvectors matrix.
120 
121  // AddingSensitivity:BEGIN /////////////////////////////////////////
122  Matrix dispSensitivity;
123  Matrix velSensitivity;
124  Matrix accSensitivity;
125  int parameterID;
126  // AddingSensitivity:END ///////////////////////////////////////////
127 
128  static std::deque<Matrix> theMatrices;
129 
130  mutable std::set<ContinuaReprComponent *> connected;
131 
132  std::set<int> coacciones_freeze;
133  const ID &get_id_coacciones(void) const;
134  void set_id_coacciones(const ID &);
135 
136  static DefaultTag defaultTag; //<! tag for next new node.
137  protected:
138 
139  DbTagData &getDbTagData(void) const;
140  int sendData(CommParameters &);
141  int recvData(const CommParameters &);
142  public:
143  // typedefs.
144  typedef std::set<Element *> ElementPtrSet;
145  typedef std::set<const Element *> ElementConstPtrSet;
146 
147  // constructors
148  Node(int classTag);
149  Node(int tag, int classTag);
150  Node(int tag, int ndof, double Crd1);
151  Node(int tag, int ndof, double Crd1, double Crd2);
152  Node(int tag, int ndof, double Crd1, double Crd2, double Crd3);
153  Node(int tag, int ndof, const Vector &crds);
154  Node(const Node &theCopy, bool copyMass = true);
155  Node *getCopy(void) const;
156 
157  // destructor
158  virtual ~Node(void);
159 
160  static DefaultTag &getDefaultTag(void);
161 
162  // public methods dealing with the DOF at the node
163  virtual int getNumberDOF(void) const;
164  virtual void setDOF_GroupPtr(DOF_Group *theDOF_Grp);
165  virtual DOF_Group *getDOF_GroupPtr(void);
166 
167  void connect(ContinuaReprComponent *el) const;
168  void disconnect(ContinuaReprComponent *el) const;
169  ElementConstPtrSet getConnectedElements(void) const;
170  ElementPtrSet getConnectedElements(void);
171  const MeshEdge *next(const std::deque<MeshEdge> &, const std::set<const MeshEdge *> &) const;
172 
173  const bool isDead(void) const;
174  const bool isAlive(void) const;
175  const bool isFrozen(void) const;
176  const bool isFree(void) const;
177  void kill(void);
178  void alive(void);
179  void freeze_if_dead(NodeLocker *locker);
180  void melt_if_alive(NodeLocker *locker);
181 
183  void fix(const std::vector<int> &, const std::vector<double> &);
184  void fix(const ID &, const Vector &);
185 
186  // public methods for obtaining the nodal coordinates
187  virtual size_t getDim(void) const;
188  virtual const Vector &getCrds(void) const;
189  virtual Vector &getCrds(void);
190  Vector getCrds3d(void) const;
191  Pos2d getPosInicial2d(void) const;
192  Pos3d getPosInicial3d(void) const;
193  Pos2d getPosFinal2d(void) const;
194  Pos3d getPosFinal3d(void) const;
195  void setPos(const Pos3d &);
196  void Mueve(const Vector3d &desplaz);
197  void Transforma(const TrfGeom &trf);
198  double getDist2(const Pos2d &p,bool initialGeometry= true) const;
199  double getDist(const Pos2d &p,bool initialGeometry= true) const;
200  double getDist2(const Pos3d &p,bool initialGeometry= true) const;
201  double getDist(const Pos3d &p,bool initialGeometry= true) const;
202 
203  // public methods for obtaining committed and trial
204  // response quantities of the node
205  virtual const Vector &getDisp(void) const;
206  Vector getDispXYZ(void) const;
207  Vector getRotXYZ(void) const;
208  virtual const Vector &getVel(void) const;
209  Vector getVelXYZ(void) const;
210  Vector getOmegaXYZ(void) const;
211  virtual const Vector &getAccel(void) const;
212  Vector getAccelXYZ(void) const;
213  Vector getAlphaXYZ(void) const;
214  virtual const Vector &getIncrDisp(void) const;
215  virtual const Vector &getIncrDeltaDisp(void) const;
216 
217  virtual const Vector &getTrialDisp(void) const;
218  virtual const Vector &getTrialVel(void) const;
219  virtual const Vector &getTrialAccel(void) const;
220 
221  // public methods for updating the trial response quantities
222  virtual int setTrialDispComponent(double value, int dof);
223  virtual int setTrialDisp(const Vector &);
224  virtual int setTrialVel(const Vector &);
225  virtual int setTrialAccel(const Vector &);
226 
227  virtual int incrTrialDisp(const Vector &);
228  virtual int incrTrialVel(const Vector &);
229  virtual int incrTrialAccel(const Vector &);
230 
231  // public methods for adding and obtaining load information
232  virtual void zeroUnbalancedLoad(void);
233  virtual int addUnbalancedLoad(const Vector &load, double fact = 1.0);
234  virtual int addInertiaLoadToUnbalance(const Vector &accel, double fact = 1.0);
235  virtual const Vector &getUnbalancedLoad(void);
236  virtual const Vector &getUnbalancedLoadIncInertia(void);
237 
238  // public methods dealing with the committed state of the node
239  virtual int commitState();
240  virtual int revertToLastCommit();
241  virtual int revertToStart();
242 
243  // public methods for dynamic analysis
244  virtual const Matrix &getMass(void);
245  virtual int setMass(const Matrix &theMass);
246  virtual int setNumColR(int numCol);
247  virtual int setR(int row, int col, double Value);
248  virtual const Vector &getRV(const Vector &V);
249 
250  virtual int setRayleighDampingFactor(double alphaM);
251  virtual const Matrix &getDamp(void);
252 
253  void addTributaria(const double &) const;
254  void resetTributaria(void) const;
255  const double &getTributaria(void) const;
256 
257  // public methods for eigen vector
258  virtual int setNumEigenvectors(int numVectorsToStore);
259  virtual int setEigenvector(int mode, const Vector &eigenVector);
260  inline int getNumModes(void) const
261  { return theEigenvectors.noCols(); }
262  virtual Vector getEigenvector(int mode) const;
263  Vector getNormalizedEigenvector(int mode) const;
264  virtual const Matrix &getEigenvectors(void);
265  Matrix getNormalizedEigenvectors(void) const;
266 
267  //Angular frequencies.
268  double getAngularFrequency(int) const;
269  Vector getAngularFrequencies(void) const;
270 
271  //Modal participation factors.
272  virtual double getModalParticipationFactor(int mode) const;
274  virtual double getModalParticipationFactor(int mode,const std::set<int> &gdls) const;
275  Vector getModalParticipationFactors(const std::set<int> &gdls) const;
276  Vector getModalParticipationFactorsForGdls(const boost::python::list &) const;
277  //Distribution factors.
278  Vector getDistributionFactor(int mode) const;
279  Vector getDistributionFactor(int mode,const std::set<int> &gdls) const;
280  Matrix getDistributionFactors(void) const;
281 
282  //Effective modal masses
283  double getEffectiveModalMass(int mode) const;
284  Vector getEffectiveModalMasses(void) const;
285 
286  //Equivalent static load.
287  Vector getEquivalentStaticLoad(int mode,const double &) const;
288 
289  //Maximum modal values
290  Vector getMaxModalDisplacement(int mode,const double &) const;
291  Vector getMaxModalVelocity(int mode,const double &) const;
292  Vector getMaxModalAcceleration(int mode,const double &) const;
293  Vector getMaxModalDisplacementForGdls(int mode,const double &,const std::set<int> &gdls) const;
294  Vector getMaxModalVelocityForGdls(int mode,const double &,const std::set<int> &gdls) const;
295  Vector getMaxModalAccelerationForGdls(int mode,const double &,const std::set<int> &gdls) const;
296  Vector getMaxModalDisplacementForGdls(int mode,const double &,const boost::python::list &gdls) const;
297  Vector getMaxModalVelocityForGdls(int mode,const double &,const boost::python::list &gdls) const;
298  Vector getMaxModalAccelerationForGdls(int mode,const double &,const boost::python::list &gdls) const;
299 
300  // public methods for output
301  virtual int sendSelf(CommParameters &);
302  virtual int recvSelf(const CommParameters &);
303 
304  std::set<SetBase *> get_sets(void) const;
305  void add_to_sets(std::set<SetBase *> &);
306 
307  virtual void Print(std::ostream &s, int flag = 0);
308 
309  virtual const Vector &getReaction(void) const;
310  const Vector &getResistingForce(const ElementConstPtrSet &,const bool &) const;
311  SVD3d getResistingSVD3d(const ElementConstPtrSet &,const bool &) const;
312  virtual int addReactionForce(const Vector &, double factor);
313  virtual int resetReactionForce(bool inclInertia);
314  void checkReactionForce(const double &);
315 
316  // AddingSensitivity:BEGIN /////////////////////////////////////////
317  int addInertiaLoadSensitivityToUnbalance(const Vector &accel, double fact = 1.0, bool tag=false);
318  Matrix getMassSensitivity(void);
319  virtual const Matrix &getDampSensitivity(void);
320  int getCrdsSensitivity(void);
321  int saveSensitivity(Vector *v, Vector *vdot, Vector *vdotdot, int gradNum, int numGrads);
322  double getDispSensitivity(int dof, int gradNum);
323  double getVelSensitivity(int dof, int gradNum);
324  double getAccSensitivity(int dof, int gradNum);
325  int setParameter(const std::vector<std::string> &argv, Parameter &param);
326  int updateParameter(int parameterID, Information &info);
327  int activateParameter(int parameterID);
328 
329  // AddingSensitivity:END ///////////////////////////////////////////
330 
331  };
332 
333 Pos3d pos_node(const Node &nod,bool initialGeometry= true);
334 
335 } // end of XC namespace
336 
337 #endif
338 
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:616
std::set< SetBase * > get_sets(void) const
Returns the sets a los que pertenece este nodo.
Definition: Node.cpp:1401
virtual void Print(std::ostream &s, int flag=0)
Imprime el objeto.
Definition: Node.cpp:1426
Single freedom constraints that make part of a load pattern.
Definition: NodeLocker.h:44
int sendData(CommParameters &)
Sends object members through the channel being passed as parameter.
Definition: Node.cpp:1331
static DefaultTag & getDefaultTag(void)
Returns a default value for node identifier.
Definition: Node.cpp:369
Mesh edge.
Definition: MeshEdge.h:39
virtual const Vector & getUnbalancedLoad(void)
Returns unbalanced load vector.
Definition: Node.cpp:824
Pos2d getPosFinal2d(void) const
Return the final position of the node.
Definition: Node.cpp:480
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:1626
void add_to_sets(std::set< SetBase * > &)
Adds the node to the sets being passed as parameters.
Definition: Node.cpp:1417
Vector getCrds3d(void) const
Returns the coordenadas del nodo en un espacio 3d.
Definition: Node.cpp:436
Vectores to store trial and commited values of node velocity.
Definition: NodeVelVectors.h:41
A DOF_Group object is instantiated by the ConstraintHandler for every unconstrained node in the domai...
Definition: DOF_Group.h:94
void resetTributaria(void) const
Anula la magnitud tributaria (longitud, área o volumen).
Definition: Node.cpp:924
Pos3d pos_node(const Node &nod, bool initialGeometry=true)
Return the 3D position of the node.
Definition: Node.cpp:666
Vector getVelXYZ(void) const
Returns the XYZ components of the translational velocity of the node.
Definition: Node.cpp:590
void Mueve(const Vector3d &desplaz)
Moves the node (intended only for its use from XC::Set).
Definition: Node.cpp:1751
virtual int getNumberDOF(void) const
Return the number of node DOFs.
Definition: Node.cpp:410
virtual void setDOF_GroupPtr(DOF_Group *theDOF_Grp)
Sets the DOF_Group pointer.
Definition: Node.cpp:414
Definition: Vector.h:82
Vector que almacena los dbTags de los miembros de la clase.
Definition: DbTagData.h:43
Vectors that store trial and commited values of the node acceleration.
Definition: NodeAccelVectors.h:41
Vector getMaxModalAccelerationForGdls(int mode, const double &, const std::set< int > &gdls) const
Return the maximal modal acceleration on the DOFs and mode being passed as parameter and the accelera...
Definition: Node.cpp:1307
virtual const Vector & getIncrDisp(void) const
Returns the displacement increment.
Definition: Node.cpp:707
virtual DOF_Group * getDOF_GroupPtr(void)
Gets the DOF_Group pointer.
Definition: Node.cpp:419
Pos3d getPosInicial3d(void) const
Return the intial position of the node in 3D.
Definition: Node.cpp:466
Vector getOmegaXYZ(void) const
Returns the XYZ components of the angular velocity of the node.
Definition: Node.cpp:594
void connect(ContinuaReprComponent *el) const
Inserts a component (element, constraint,...) to the connected component list.
Definition: Node.cpp:245
virtual size_t getDim(void) const
Return the dimension of the node vector position (1,2 or 3).
Definition: Node.cpp:424
double getEffectiveModalMass(int mode) const
Return the masa modal efectiva correspondiente al modo i.
Definition: Node.cpp:1218
Base class for components used to represent the material (continuum).
Definition: ContinuaReprComponent.h:37
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: Node.cpp:1351
void freeze_if_dead(NodeLocker *locker)
Imposes zero displacement (zero value for all node DOFs).
Definition: Node.cpp:299
virtual int setRayleighDampingFactor(double alphaM)
Returns the coeficiente de amortiguamiento de Rayleigh.
Definition: Node.cpp:893
virtual int setMass(const Matrix &theMass)
Asigna la matriz de masas to the node.
Definition: Node.cpp:950
Mesh node.
Definition: Node.h:99
std::set< const Element * > ElementConstPtrSet
Container of const element pointers.
Definition: Node.h:145
Information about an element.
Definition: Information.h:80
virtual int setEigenvector(int mode, const Vector &eigenVector)
Asigna el autovector eigenvector al modo mode.
Definition: Node.cpp:1032
Vector getAccelXYZ(void) const
Returns the XYZ components of the translational acceleration of the node.
Definition: Node.cpp:598
Pos2d getPosInicial2d(void) const
Return the intial position of the node in 2D.
Definition: Node.cpp:454
const bool isFrozen(void) const
returns true if the node is frozen.
Definition: Node.cpp:349
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:1257
SFreedom_Constraint * fix(const SFreedom_Constraint &)
Introduce en el nodo una constraint como la being passed as parameter.
Definition: Node.cpp:374
Base class for nodes and elements (mesh components).
Definition: MeshComponent.h:41
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:1250
void disconnect(ContinuaReprComponent *el) const
Removes a component (element, constraint,...) from the connected component list.
Definition: Node.cpp:254
virtual const Vector & getDisp(void) const
Returns the displacement of the node.
Definition: Node.cpp:676
virtual int sendSelf(CommParameters &)
Envia el objeto through the channel being passed as parameter.
Definition: Node.cpp:1375
double getAngularFrequency(int) const
Return the angular frequency corresponding to the i-th mode.
Definition: Node.cpp:1052
virtual const Matrix & getDamp(void)
Return the matriz de amortiguamiento del nodo.
Definition: Node.cpp:901
Default tag.
Definition: DefaultTag.h:36
Pos3d getPosFinal3d(void) const
Return the final position of the node.
Definition: Node.cpp:503
Vector getMaxModalDisplacementForGdls(int mode, const double &, const std::set< int > &gdls) const
Returns the maximal modal displacement on the DOFs and mode being passed as parameter and the acceler...
Definition: Node.cpp:1272
virtual const Matrix & getMass(void)
Return the matriz de masas del nodo.
Definition: Node.cpp:889
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:1241
Vector getAlphaXYZ(void) const
Returns the XYZ components of the angular acceleration of the node.
Definition: Node.cpp:602
Definition: ID.h:77
virtual int commitState()
Commits the state of the node.
Definition: Node.cpp:848
virtual const Vector & getUnbalancedLoadIncInertia(void)
Return unbalanced load vector including inertial forces.
Definition: Node.cpp:829
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:609
virtual int revertToLastCommit()
Returns to the last commited state.
Definition: Node.cpp:858
const double & getTributaria(void) const
Return the magnitud tributaria (longitud, área o volumen).
Definition: Node.cpp:928
virtual int revertToStart()
Returns to the initial state.
Definition: Node.cpp:869
Vector getDispXYZ(void) const
Returns the XYZ components of node displacement.
Definition: Node.cpp:582
Definition: TrfGeom.h:49
Vectores to store trial and commited values («commited») on node displacements.
Definition: NodeDispVectors.h:41
Vector getEffectiveModalMasses(void) const
Returns the effective modal mases.
Definition: Node.cpp:1230
Definition: Matrix.h:82
void setPos(const Pos3d &)
Sets the node position.
Definition: Node.cpp:639
virtual int setNumEigenvectors(int numVectorsToStore)
Dimensiona la matriz que contiene a los eigenvectors.
Definition: Node.cpp:1014
virtual const Vector & getAccel(void) const
Returns the acceleration of the node.
Definition: Node.cpp:684
Vector getModalParticipationFactorsForGdls(const boost::python::list &) const
Returns the modal participation factors. If gdls argument is not empty it "projects" the i-th mode ov...
Definition: Node.cpp:1174
virtual const Vector & getReaction(void) const
Return the node reaction.
Definition: Node.cpp:1697
virtual const Vector & getVel(void) const
Returns the velocity of the node.
Definition: Node.cpp:680
void addTributaria(const double &) const
Adds to the magnitud tributaria la longitud, área o volumen being passed as parameter.
Definition: Node.cpp:920
ElementConstPtrSet getConnectedElements(void) const
Return a list of pointers to the elements that are connected with this node.
Definition: Node.cpp:1600
virtual int addUnbalancedLoad(const Vector &load, double fact=1.0)
Adds vector to unbalanced load.
Definition: Node.cpp:740
Vector getNormalizedEigenvector(int mode) const
Returns the autovector correspondiente al modo i normalizado de modo que la componente máxima valga 1...
Definition: Node.cpp:1081
void checkReactionForce(const double &)
Checks reaction on the node.
Definition: Node.cpp:1721
virtual int addInertiaLoadToUnbalance(const Vector &accel, double fact=1.0)
Adds inertial loads to unbalanced load vector.
Definition: Node.cpp:768
const bool isAlive(void) const
True if node is active.
Definition: Node.cpp:270
virtual ~Node(void)
destructor
Definition: Node.cpp:363
virtual const Matrix & getEigenvectors(void)
Returns the eigenvectors correspondientes to the node.
Definition: Node.cpp:1072
Vector getRotXYZ(void) const
Returns the XYZ components of node rotation.
Definition: Node.cpp:586
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:1264
Definition: Parameter.h:65
Vector getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: Node.cpp:1110
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: Node.cpp:1388
virtual const Vector & getTrialDisp(void) const
Returns the trial value of the displacement of the node.
Definition: Node.cpp:695
virtual const Vector & getTrialAccel(void) const
Returns the trial value of the acceleration of the node.
Definition: Node.cpp:703
virtual const Vector & getCrds(void) const
return the vector of nodal coordinates
Definition: Node.cpp:428
void Transforma(const TrfGeom &trf)
Applies to the node position the transformation being passed as parameter.
Definition: Node.cpp:658
Vector getMaxModalVelocityForGdls(int mode, const double &, const std::set< int > &gdls) const
Returns the maximum modal velocity on the DOFs and mode being passed as parameter and the acceleratio...
Definition: Node.cpp:1290
Single freedom constraint.
Definition: SFreedom_Constraint.h:79
Vector getAngularFrequencies(void) const
Returns the angular frequencies from the modal analysis.
Definition: Node.cpp:1060
Matrix getDistributionFactors(void) const
Returns the matrix with the computed distribution factors placed by columns.
Definition: Node.cpp:1195
SVD3d getResistingSVD3d(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:1675
Communication parameters between processes.
Definition: CommParameters.h:65
virtual const Vector & getTrialVel(void) const
Returns the trial value of the velocity of the node.
Definition: Node.cpp:699
std::set< Element * > ElementPtrSet
Container of element pointers.
Definition: Node.h:144
Matrix getNormalizedEigenvectors(void) const
Returns a matriz con los eigenvectors normalizados colocados por columnas (norma_infinito).
Definition: Node.cpp:1086
const bool isDead(void) const
True if node is inactive.
Definition: Node.cpp:266
virtual double getModalParticipationFactor(int mode) const
Returns the modal participation factor corresponding to i-th mode.
Definition: Node.cpp:1102
void melt_if_alive(NodeLocker *locker)
Deletes the constraint over the node DOFs previously created by the freeze method.
Definition: Node.cpp:334
================================================================================
Definition: ContinuaReprComponent.h:34
Vector getDistributionFactor(int mode) const
Returns the distribution factor of the i-th mode.
Definition: Node.cpp:1181
Node(int classTag)
Default constructor.
Definition: Node.cpp:102
virtual Vector getEigenvector(int mode) const
Returns the autovector correspondiente al modo i.
Definition: Node.cpp:1076
const bool isFree(void) const
returns true if the node has no constraints.
Definition: Node.cpp:353
Node * getCopy(void) const
Virtual constructor.
Definition: Node.cpp:262
virtual int addReactionForce(const Vector &, double factor)
Increments the node reaction.
Definition: Node.cpp:1701
DbTagData & getDbTagData(void) const
Returns a vector para almacenar los dbTags de los miembros de la clase.
Definition: Node.cpp:1324
virtual int resetReactionForce(bool inclInertia)
Calculate the reactions in this node (used in en Domain::calculateNodalReactions).
Definition: Node.cpp:1738
const Vector & getResistingForce(const ElementConstPtrSet &, const bool &) const
Return the action of the elements from the set over this node.
Definition: Node.cpp:1644