xc
Domain.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.16 $
49 // $Date: 2005/11/30 23:32:32 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/domain/domain/Domain.h,v $
51 
52 // Written: fmk
53 // Created: Fri Sep 20 15:27:47: 1996
54 // Revision: A
55 //
56 // Description: This file contains the class definition for Domain.
57 // Domain is a container class. The class is responsible for holding
58 // and providing access to the Elements, Nodes, SFreedom_Constraints
59 // MFreedom_Constraints, and LoadPatterns.
60 //
61 // What: "@(#) Domain.h, revA"
62 
63 #ifndef Domain_h
64 #define Domain_h
65 
66 #include "utility/recorder/ObjWithRecorders.h"
67 #include "PseudoTimeTracker.h"
68 #include "../mesh/Mesh.h"
69 #include "../constraints/ConstrContainer.h"
70 #include "utility/matrix/Vector.h"
71 #include "domain/mesh/region/DqMeshRegion.h"
72 
73 namespace XC {
74 class Element;
75 class Node;
76 class Preprocessor;
77 
78 class ElementIter;
79 class NodeIter;
80 
81 class SingleDomParamIter;
82 
83 class LoadCombination;
84 
85 class MeshRegion;
86 class DqMeshRegion;
87 class Recorder;
88 class Graph;
89 class NodeGraph;
90 class ElementGraph;
91 class FEM_ObjectBroker;
92 class RayleighDampingFactors;
93 
95 //
96 //
118  {
119  private:
120  PseudoTimeTracker timeTracker;
121  std::string callbackCommit;
122 
123  int dbTag;
124  int currentGeoTag;
125  bool hasDomainChangedFlag;
126  int commitTag;
127  Mesh mesh;
128  ConstrContainer constraints;
129  Vector theEigenvalues;
130  Vector modalParticipationFactors;
131  DqMeshRegion theRegions;
132  std::deque<std::string> activeCombinations;
133 
134  TaggedObjectStorage *theParameters;
135  SingleDomParamIter *theParamIter;
136  std::deque<int> paramIndex;
137 
138  int lastChannel;
139  int lastGeoSendTag;
140 
141  void alloc_containers(void);
142  void alloc_iters(void);
143  bool check_containers(void) const;
144  protected:
145  virtual int buildEleGraph(Graph &theEleGraph);
146  virtual int buildNodeGraph(Graph &theNodeGraph);
147  inline virtual Domain *get_domain_ptr(void)
148  { return this; }
149 
150  void free_mem(void);
151  DbTagData &getDbTagData(void) const;
152  int sendData(Communicator &comm);
153  int recvData(const Communicator &comm);
154  public:
155  Domain(CommandEntity *owr,DataOutputHandler::map_output_handlers *oh);
156  Domain(CommandEntity *owr,int numNods, int numElements, int numSPs, int numMPs,int numLPatterns,int numNLockers,DataOutputHandler::map_output_handlers *oh);
157 
158  virtual ~Domain(void);
159 
160  // methods to populate a domain
161  virtual bool addElement(Element *);
162  virtual bool addNode(Node *);
166  virtual bool addLoadPattern(LoadPattern *);
167  virtual bool isLoadPatternActive(const LoadPattern *) const;
168  virtual bool addNodeLocker(NodeLocker *);
169  virtual bool addLoadCombination(LoadCombination *);
170  virtual bool addParameter(Parameter *);
171 
172  void setNodeReactionException(const int &);
173  bool checkNodalReactions(const double &);
174 
175  // methods to add components to a LoadPattern object
176  virtual bool addSFreedom_Constraint(SFreedom_Constraint *, int loadPatternTag);
177  virtual bool addNodalLoad(NodalLoad *, int loadPatternTag);
178  virtual bool addElementalLoad(ElementalLoad *, int loadPatternTag);
179 
180  // methods to remove the components
181  virtual void clearAll(void);
182  virtual bool removeElement(int tag);
183  bool remove(Element *);
184  virtual bool removeNode(int tag);
185  bool remove(Node *);
186  virtual bool removeSFreedom_Constraint(int theNode, int theDOF, int loadPatternTag);
187  virtual bool removeSFreedom_Constraint(int tag);
188  virtual bool removeMFreedom_Constraint(int tag);
189  virtual bool removeMRMFreedom_Constraint(int tag);
190  bool remove(Constraint *);
191  virtual bool removeLoadPattern(int loadTag);
192  virtual bool removeNodeLocker(int nlTag);
193  bool removeLoadPattern(LoadPattern *lp);
194  bool removeAllLoadPatterns(void);
195  bool removeNodeLocker(NodeLocker *lp);
197  void removeAllLoadCombinations(void);
199  void removeLPs(void);
200  void removeNLs(void);
201 
202  virtual bool removeNodalLoad(int tag, int loadPattern);
203  virtual bool removeElementalLoad(int tag, int loadPattern);
204  virtual bool removeSFreedom_Constraint(int tag, int loadPattern);
205 
206  virtual void clearDOF_GroupPtr(void);
207 
208  // methods to access the components of a domain
209  virtual ElementIter &getElements(void);
210  virtual NodeIter &getNodes(void);
211  virtual Mesh &getMesh(void);
212  virtual const Mesh &getMesh(void) const;
213  virtual ConstrContainer &getConstraints(void);
214  virtual const ConstrContainer &getConstraints(void) const;
215 
216  std::string getCurrentCombinationName(void) const;
217  std::string getCurrentLoadCaseDescription(void) const;
218 
219  bool existElement(int tag);
220  virtual Element *getElement(int tag);
221  virtual const Element *getElement(int tag) const;
222  bool existNode(int tag);
223  virtual Node *getNode(int tag);
224  virtual const Node *getNode(int tag) const;
225  bool existSFreedom_Constraint(int tag);
226  virtual Constraint *getSFreedom_Constraint(int tag);
227  virtual const Constraint *getSFreedom_Constraint(int tag) const;
228  bool existMFreedom_Constraint(int tag);
229  virtual Constraint *getMFreedom_Constraint(int tag);
230  virtual const Constraint *getMFreedom_Constraint(int tag) const;
231  bool existMRMFreedom_Constraint(int tag);
232  virtual Constraint *getMRMFreedom_Constraint(int tag);
233  virtual const Constraint *getMRMFreedom_Constraint(int tag) const;
234 
235 
236  virtual Parameter *getParameter(int tag);
237  virtual const Parameter *getParameter(int tag) const;
238 
239  // methods to query the state of the domain
241  inline const PseudoTimeTracker &getTimeTracker(void) const
242  { return timeTracker; }
244  inline double getCommittedTime(void) const
245  { return timeTracker.getCommittedTime(); }
247  inline double getCurrentTime(void) const
248  { return timeTracker.getCurrentTime(); }
249 
250 
251  inline int getCurrentGeoTag(void) const
252  { return currentGeoTag; }
253  virtual int getCommitTag(void) const;
254  virtual int getNumElements(void) const;
255  virtual int getNumNodes(void) const;
256  virtual const Vector &getPhysicalBounds(void);
257 
258 
259  // methods to get element and node graphs
260  virtual Graph &getElementGraph(void);
261  virtual Graph &getNodeGraph(void);
262 
263  // methods to update the domain
264  virtual void setCommitTag(int newTag);
265  virtual void setCurrentTime(double newTime);
266  virtual void setCommittedTime(double newTime);
267  virtual void setTime(double newTime);
268  virtual void applyLoad(double pseudoTime);
269  virtual void setLoadConstant(void);
270  virtual int initialize(void);
271  virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF);
272 
273  virtual int commit(void);
274  virtual int revertToLastCommit(void);
275  virtual int revertToStart(void);
276  virtual int update(void);
277  virtual int update(double newTime, double dT);
278  virtual int newStep(double dT);
279 
280  void resetLoadCase(void);
281 
282  // methods for eigenvalue analysis
283  int getNumModes(void) const;
284  virtual int setEigenvalues(const Vector &);
285  virtual const double &getEigenvalue(int) const;
286  double getAngularFrequency(int) const;
287  double getPeriod(int) const;
288  double getFrequency(int) const;
289  virtual const Vector &getEigenvalues(void) const;
290  boost::python::list getEigenvaluesPy(void) const;
291  void clearEigenvectors(void);
292  void clearEigenvalues(void);
293  void clearEigendata(void);
294  Vector getAngularFrequencies(void) const;
295  Vector getPeriods(void) const;
296  Vector getFrequencies(void) const;
297  virtual int setModalParticipationFactors(const Vector &);
298  virtual const double &getModalParticipationFactor(int mode) const;
299  virtual const Vector &getModalParticipationFactors(void) const;
300  const double getEffectiveModalMass(int mode) const;
301  Vector getEffectiveModalMasses(void) const;
302  double getTotalEffectiveModalMass(void) const;
303 
304  // mass distribution
305  Matrix getTotalMass(void) const;
306  double getTotalMassComponent(const int &) const;
307 
308  // methods for other objects to determine if model has changed
309  virtual void domainChange(void);
310  virtual int hasDomainChanged(void);
311  virtual void setDomainChangeStamp(int newStamp);
312 
313  virtual int addRegion(MeshRegion &theRegion);
314  virtual MeshRegion *getRegion(int region);
315 
316  virtual void Print(std::ostream &s, int flag =0) const;
317  friend std::ostream &operator<<(std::ostream &, const Domain &);
318 
319  virtual int sendSelf(Communicator &);
320  virtual int recvSelf(const Communicator &);
321  friend int sendDomain(Domain &, int posDbTag,DbTagData &,Communicator &comm);
322  friend int receiveDomain(Domain &, int posDbTag,DbTagData &,const Communicator &comm);
323 
324  boost::python::dict getPyDict(void) const;
325  void setPyDict(const boost::python::dict &);
326 
327  const Preprocessor *getPreprocessor(void) const;
329 
330  // nodal methods required in domain interface for parallel interprter
331  virtual double getNodeDisp(int nodeTag, int dof, int &errorFlag);
332  virtual int setMass(const Matrix &mass, int nodeTag);
333 
334  virtual int calculateNodalReactions(bool inclInertia,const double &);
335  virtual int addRecorder(Recorder &theRecorder);
336 
337  static void setDeadSRF(const double &);
338  };
339 
340 std::ostream &operator<<(std::ostream &, const Domain &);
341 int sendDomain(Domain &,int posDbTag,DbTagData &,Communicator &comm);
342 int receiveDomain(Domain &, int posDbTag,DbTagData &,const Communicator &comm);
343 
344 
345 } // end of XC namespace
346 
347 #endif
348 
349 
virtual double getNodeDisp(int nodeTag, int dof, int &errorFlag)
Return the value of the dof component of displacement for the node with the tag being passed as param...
Definition: Domain.cpp:1565
virtual int calculateNodalReactions(bool inclInertia, const double &)
Calculate nodal reaction forces and moments.
Definition: Domain.cpp:1581
virtual void clearDOF_GroupPtr(void)
Clears the pointers to DOF groups.
Definition: Domain.cpp:790
virtual void domainChange(void)
Sets a flag indicating that the integer returned in the next call to hasDomainChanged() must be incre...
Definition: Domain.cpp:1288
virtual bool addNode(Node *)
Adds to the domain the node being passed as parameter.
Definition: Domain.cpp:251
virtual const Vector & getPhysicalBounds(void)
Definition: Domain.cpp:959
DbTagData & getDbTagData(void) const
Returns a vector to store the dbTags de los miembros of the clase.
Definition: Domain.cpp:1389
virtual int addRegion(MeshRegion &theRegion)
Adds a region.
Definition: Domain.cpp:1357
virtual bool addSFreedom_Constraint(SFreedom_Constraint *)
Adds a single freedom constraint to the domain.
Definition: Domain.cpp:261
Float vector abstraction.
Definition: Vector.h:94
virtual Parameter * getParameter(int tag)
Return a pointer to the parameter identified by the argument.
Definition: Domain.cpp:920
virtual int setMass(const Matrix &mass, int nodeTag)
Set the mass matrix for the node identified by the argument.
Definition: Domain.cpp:1569
virtual int setModalParticipationFactors(const Vector &)
Sets the values of the modal participation factors.
Definition: Domain.cpp:1231
virtual bool addNodeLocker(NodeLocker *)
Appends the node locker object to the domain.
Definition: Domain.cpp:539
virtual Mesh & getMesh(void)
Returns a reference to the domain mesh.
Definition: Domain.cpp:806
virtual void clearAll(void)
Removes all components from domain (nodes, elements, loads & constraints).
Definition: Domain.cpp:201
virtual int revertToStart(void)
Return the domain to its initial state and triggers the "restart" method for all the recorders...
Definition: Domain.cpp:1097
void removeNLs(void)
Remove all node lockers from domain.
Definition: Domain.cpp:756
virtual bool removeNodalLoad(int tag, int loadPattern)
Removes from domain the nodal load being passed as parameter.
Definition: Domain.cpp:768
Communication parameters between processes.
Definition: Communicator.h:66
virtual bool removeMRMFreedom_Constraint(int tag)
Removes from domain the multi-freedom multi-retained node constraint identified by the argument...
Definition: Domain.cpp:464
int getNumModes(void) const
Return the number of computed eigenvalues.
Definition: Domain.cpp:1227
virtual const double & getModalParticipationFactor(int mode) const
Return the modal participation factor of the i-th mode.
Definition: Domain.cpp:1238
virtual Graph & getNodeGraph(void)
Builds (if necessary) the domain node graph and returns a reference to it.
Definition: Domain.cpp:983
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: Domain.cpp:1538
virtual bool addMFreedom_Constraint(MFreedom_Constraint *)
Adds to the domain a multi-freedom constraint.
Definition: Domain.cpp:284
Finite element model generation tools.
Definition: Preprocessor.h:59
virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF)
Set Rayleigh damping factors.
Definition: Domain.cpp:1041
An Recorder object is used in the program to store/restore information at each commit().
Definition: Recorder.h:87
virtual int revertToLastCommit(void)
Return the domain to its last committed state.
Definition: Domain.cpp:1079
Matrix getTotalMass(void) const
Return the total mass matrix.
Definition: Domain.cpp:1268
Iterator over an element container.
Definition: ElementIter.h:74
virtual ElementIter & getElements(void)
Returns an iterator to the element container.
Definition: Domain.cpp:794
void clearEigenvectors(void)
Remove the stored eigenvectors.
Definition: Domain.cpp:1182
double getPeriod(int) const
Return the period of the i-th mode.
Definition: Domain.cpp:1163
virtual bool removeLoadPattern(int loadTag)
Remove from domain el load pattern identified by the argument.
Definition: Domain.cpp:628
boost::python::dict getPyDict(void) const
Return a Python dictionary with the object members values.
Definition: Domain.cpp:1456
static void setDeadSRF(const double &)
Assigns Stress Reduction Factor for element deactivation.
Definition: Domain.cpp:242
virtual int buildNodeGraph(Graph &theNodeGraph)
Builds the node graph.
Definition: Domain.cpp:1384
Vector that stores the dbTags of the class members.
Definition: DbTagData.h:44
Base class for distributed processing.
Definition: DistributedBase.h:41
virtual Node * getNode(int tag)
Return a pointer to the node identified by the argument.
Definition: Domain.cpp:853
virtual void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: Domain.cpp:1319
int recvData(const Communicator &comm)
Receive data through the communicator argument.
Definition: Domain.cpp:1420
virtual int addRecorder(Recorder &theRecorder)
Adds a recorder to the model.
Definition: Domain.cpp:1342
virtual bool removeElementalLoad(int tag, int loadPattern)
Removes from domain the elemental load being passed as parameter.
Definition: Domain.cpp:775
virtual int hasDomainChanged(void)
Returns true if the model has changed.
Definition: Domain.cpp:1298
void removeAllLoadCombinations(void)
Remove all the load combinations currently in activeCombinations.
Definition: Domain.cpp:734
A load pattern is the spatial distribution as well as its variation in time of a specific set of forc...
Definition: LoadPattern.h:97
virtual int commit(void)
Commits domain state and triggers "record" method for all defined recorders.
Definition: Domain.cpp:1054
virtual ~Domain(void)
Destructor.
Definition: Domain.cpp:229
bool existElement(int tag)
Returns true if the element identified by the tag being passed as parameter already exists in the dom...
Definition: Domain.cpp:826
void setNodeReactionException(const int &)
Sets the exception for node reactions checking (see Domain::checkNodalReactions). ...
Definition: Domain.cpp:1573
const Preprocessor * getPreprocessor(void) const
Returns (if possible) a pointer to the preprocessor.
Definition: Domain.cpp:1589
const double getEffectiveModalMass(int mode) const
Return the effective modal mass of the i-th mode.
Definition: Domain.cpp:1246
virtual bool addElement(Element *)
Adds to the domain the element being passed as parameter.
Definition: Domain.cpp:247
virtual int initialize(void)
Initialize mesh.
Definition: Domain.cpp:1037
std::string getCurrentCombinationName(void) const
Return the name of the current load combination.
Definition: Domain.cpp:699
virtual int setEigenvalues(const Vector &)
Sets eigenvalues.
Definition: Domain.cpp:1144
void clearEigendata(void)
Remove the stored eigenvalues and eigenvectors.
Definition: Domain.cpp:1190
Base class for the finite elements.
Definition: Element.h:112
void clearEigenvalues(void)
Remove the stored eigenvalues.
Definition: Domain.cpp:1186
std::string getCurrentLoadCaseDescription(void) const
Return the name of the current load case.
Definition: Domain.cpp:716
virtual bool addNodalLoad(NodalLoad *, int loadPatternTag)
Appends a nodal load to the pattern being passed as parameter.
Definition: Domain.cpp:339
bool existMRMFreedom_Constraint(int tag)
Return true if the MRMFreedom constraint exists.
Definition: Domain.cpp:901
virtual const double & getEigenvalue(int) const
Return the eigenvalue of the i-th mode.
Definition: Domain.cpp:1155
double getTotalMassComponent(const int &) const
Return the total mass matrix component for the DOF argument.
Definition: Domain.cpp:1272
void resetLoadCase(void)
Prepares the domain to solve for a new load pattern.
Definition: Domain.cpp:233
virtual bool isLoadPatternActive(const LoadPattern *) const
Return true if the load pattern is already added to the domain.
Definition: Domain.cpp:534
virtual void setLoadConstant(void)
Set all the loads as constant.
Definition: Domain.cpp:1033
bool removeAllLoadPatterns(void)
Remove from all load patterns from domain.
Definition: Domain.cpp:669
bool existMFreedom_Constraint(int tag)
Return true if the MFreedom constraint exists.
Definition: Domain.cpp:883
virtual bool addLoadPattern(LoadPattern *)
Appends the load pattern to the domain.
Definition: Domain.cpp:508
void removeAllLoadsAndCombinations(void)
Remove all the load patterns and load combinations currently in this domain.
Definition: Domain.cpp:739
virtual void setCommitTag(int newTag)
Set the committed tag to newTag.
Definition: Domain.cpp:987
Vector getPeriods(void) const
Returns a vector with the computed periods (for each mode).
Definition: Domain.cpp:1207
virtual int getNumNodes(void) const
Return the number of nodes.
Definition: Domain.cpp:949
Vector getEffectiveModalMasses(void) const
Return the effective modal masses for each mode.
Definition: Domain.cpp:1250
Registro del tiempo.
Definition: PseudoTimeTracker.h:41
virtual MeshRegion * getRegion(int region)
Returns a pointer to the region identified by the argument.
Definition: Domain.cpp:1364
virtual bool removeElement(int tag)
Remove the element identified by the argument.
Definition: Domain.cpp:385
Single freedom constraint.
Definition: SFreedom_Constraint.h:85
Domain(CommandEntity *owr, DataOutputHandler::map_output_handlers *oh)
Constructor.
Definition: Domain.cpp:153
virtual Constraint * getSFreedom_Constraint(int tag)
Return a pointer to the SFreedom constraint identified by the argument.
Definition: Domain.cpp:871
Objet that can execute python scripts.
Definition: CommandEntity.h:40
Definition: SingleDomParamIter.h:38
virtual Element * getElement(int tag)
Return a pointer to the element identified by the argument.
Definition: Domain.cpp:833
The Graph class provides the abstraction of a graph.
Definition: Graph.h:94
virtual Constraint * getMRMFreedom_Constraint(int tag)
Return a pointer to the MRMFreedom constraint identified by the argument.
Definition: Domain.cpp:907
virtual bool addLoadCombination(LoadCombination *)
Adds to the domain the load combination being passed as parameter.
Definition: Domain.cpp:551
virtual int update(void)
Updates the state of the domain.
Definition: Domain.cpp:1124
virtual bool removeNodeLocker(int nlTag)
Remove from domain el.
Definition: Domain.cpp:644
Object that can manage Recorders.
Definition: ObjWithRecorders.h:42
virtual ConstrContainer & getConstraints(void)
Returns domain constraints.
Definition: Domain.cpp:814
const PseudoTimeTracker & getTimeTracker(void) const
Return a constant reference to the internal time tracker.
Definition: Domain.h:241
virtual const Vector & getModalParticipationFactors(void) const
Return the modal participation factors.
Definition: Domain.cpp:1242
double getCurrentTime(void) const
Return the current value fo the pseudo-time.
Definition: Domain.h:247
virtual bool removeNode(int tag)
Remove the node identified by the argument.
Definition: Domain.cpp:395
Rayleigh damping factors.
Definition: RayleighDampingFactors.h:59
bool checkNodalReactions(const double &)
Check that al free nodes have zero reaction.
Definition: Domain.cpp:1577
double getTotalEffectiveModalMass(void) const
Return the total effective modal mass.
Definition: Domain.cpp:1257
bool existSFreedom_Constraint(int tag)
Return true if the SFreedom constraint exists.
Definition: Domain.cpp:865
virtual bool removeMFreedom_Constraint(int tag)
Removes from domain the multi-freedom constraint identified by the argument.
Definition: Domain.cpp:451
virtual bool addElementalLoad(ElementalLoad *, int loadPatternTag)
Adds a load over element to the pattern.
Definition: Domain.cpp:364
virtual Constraint * getMFreedom_Constraint(int tag)
Return a pointer to the MFreedom constraint identified by the argument.
Definition: Domain.cpp:889
void removeLPs(void)
Remove all the load patterns from this domain.
Definition: Domain.cpp:746
Multiple retained nodes constraint.
Definition: MRMFreedom_Constraint.h:59
bool existNode(int tag)
Return true if the mesh has a node with this tag.
Definition: Domain.cpp:847
virtual void setCurrentTime(double newTime)
Set the current time to newTime.
Definition: Domain.cpp:991
Single freedom constraints that make part of a load pattern.
Definition: NodeLocker.h:45
virtual void setDomainChangeStamp(int newStamp)
Set the domain stamp to be newStamp.
Definition: Domain.cpp:1277
Vector getAngularFrequencies(void) const
Returns a vector with the computed angular frequencies (for each mode).
Definition: Domain.cpp:1197
void setPyDict(const boost::python::dict &)
Set the values of the object members from a Python dictionary.
Definition: Domain.cpp:1486
double getCommittedTime(void) const
Return the committed value of the pseudo-time.
Definition: Domain.h:244
Vector getFrequencies(void) const
Returns a vector with the computed frequencies (for each mode).
Definition: Domain.cpp:1217
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
double getAngularFrequency(int) const
Return the angular frequency of the i-th mode.
Definition: Domain.cpp:1159
Iterator over the nodes.
Definition: NodeIter.h:74
boost::python::list getEigenvaluesPy(void) const
Returns a Python list with the computed eigenvalues for each mode.
Definition: Domain.cpp:1175
Matrix of floats.
Definition: Matrix.h:111
virtual int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: Domain.cpp:1517
Multi-freedom constraint.
Definition: MFreedom_Constraint.h:113
virtual int getNumElements(void) const
Return the number of elements.
Definition: Domain.cpp:945
Dinamically allocated mesh region deque container.
Definition: DqMeshRegion.h:44
Parameter.
Definition: Parameter.h:68
double getFrequency(int) const
Return the frequency of the i-th mode.
Definition: Domain.cpp:1167
Base class for model constraints.
Definition: Constraint.h:48
virtual void setCommittedTime(double newTime)
Set the committed time to newTime.
Definition: Domain.cpp:995
Base class for loads over elements.
Definition: ElementalLoad.h:79
Base class for load pattern combinations (1.5*selfWeight+1.0*permanentLoad+1.6*trafficLoad ...
Definition: LoadCombination.h:45
virtual void applyLoad(double pseudoTime)
Apply the loads for the given time pseudoTime.
Definition: Domain.cpp:1015
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
Mesh node.
Definition: Node.h:111
Constraint (essential and natural boundary conditions) container.
Definition: ConstrContainer.h:64
Nodes and elements of a mesh region.
Definition: MeshRegion.h:78
virtual bool removeSFreedom_Constraint(int theNode, int theDOF, int loadPatternTag)
Remove the single freedom constraint from the load pattern identified by the argument.
Definition: Domain.cpp:419
Load over a node.
Definition: NodalLoad.h:83
virtual const Vector & getEigenvalues(void) const
Return the eigenvalues vector.
Definition: Domain.cpp:1171
Finite element mesh.
Definition: Mesh.h:65
int sendData(Communicator &comm)
Send data through the communicator argument.
Definition: Domain.cpp:1396
virtual Graph & getElementGraph(void)
Builds (if necessary) the domain elements graph and returns a reference to it.
Definition: Domain.cpp:971
void removeLoadCombination(LoadCombination *comb)
Removes from the domain the load combination being passed as parameter.
Definition: Domain.cpp:721
virtual bool addMRMFreedom_Constraint(MRMFreedom_Constraint *)
Adds to the domain a multi-freedom multi-retained node constraint.
Definition: Domain.cpp:296
virtual NodeIter & getNodes(void)
Returns an iterator to the node container.
Definition: Domain.cpp:798
virtual int buildEleGraph(Graph &theEleGraph)
Builds the element graph.
Definition: Domain.cpp:1375