XC Open source finite element analysis program
Domain.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.16 $
48 // $Date: 2005/11/30 23:32:32 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/domain/domain/Domain.h,v $
50 
51 // Written: fmk
52 // Created: Fri Sep 20 15:27:47: 1996
53 // Revision: A
54 //
55 // Description: This file contains the class definition for Domain.
56 // Domain is a container class. The class is responsible for holding
57 // and providing access to the Elements, Nodes, SFreedom_Constraints
58 // MFreedom_Constraints, and LoadPatterns.
59 //
60 // What: "@(#) Domain.h, revA"
61 
62 #ifndef Domain_h
63 #define Domain_h
64 
65 #include "utility/recorder/ObjWithRecorders.h"
66 #include "PseudoTimeTracker.h"
67 #include "../mesh/Mesh.h"
68 #include "../constraints/ConstrContainer.h"
69 #include "utility/matrix/Vector.h"
70 
71 namespace XC {
72 class Element;
73 class Node;
74 class Preprocessor;
75 
76 class ElementIter;
77 class NodeIter;
78 
79 class SingleDomEleIter;
80 class SingleDomNodIter;
81 
82 class LoadCombination;
83 
84 class MeshRegion;
85 class DqMeshRegion;
86 class Recorder;
87 class Graph;
88 class NodeGraph;
89 class ElementGraph;
90 class FEM_ObjectBroker;
91 class RayleighDampingFactors;
92 
94 
96 //
99  {
100  private:
101  PseudoTimeTracker timeTracker;
102  std::string CallbackCommit;
103 
104  int dbTag;
105  int currentGeoTag;
106  bool hasDomainChangedFlag;
107  int commitTag;
108  Mesh mesh;
109  ConstrContainer constraints;
110  Vector theEigenvalues;
111  Vector modalParticipationFactors;
112  DqMeshRegion *theRegions;
113  std::string nmbCombActual;
114 
115  int lastChannel;
116  int lastGeoSendTag;
117 
118  protected:
119  virtual int buildEleGraph(Graph &theEleGraph);
120  virtual int buildNodeGraph(Graph &theNodeGraph);
121  inline virtual Domain *get_domain_ptr(void)
122  { return this; }
123 
124  void libera(void);
125  DbTagData &getDbTagData(void) const;
126  int sendData(CommParameters &cp);
127  int recvData(const CommParameters &cp);
128  public:
129  Domain(EntCmd *owr,DataOutputHandler::map_output_handlers *oh);
130  Domain(EntCmd *owr,int numNods, int numElements, int numSPs, int numMPs,int numLPatterns,int numNLockers,DataOutputHandler::map_output_handlers *oh);
131 
132  virtual ~Domain(void);
133 
134  // methods to populate a domain
135  virtual bool addElement(Element *);
136  virtual bool addNode(Node *);
140  virtual bool addLoadPattern(LoadPattern *);
141  virtual bool addNodeLocker(NodeLocker *);
142  virtual bool addLoadCombination(LoadCombination *);
143 
144  void setNodeReactionException(const int &);
145  void checkNodalReactions(const double &);
146 
147  // methods to add components to a LoadPattern object
148  virtual bool addSFreedom_Constraint(SFreedom_Constraint *, int loadPatternTag);
149  virtual bool addNodalLoad(NodalLoad *, int loadPatternTag);
150  virtual bool addElementalLoad(ElementalLoad *, int loadPatternTag);
151 
152  // methods to remove the components
153  virtual void clearAll(void);
154  virtual bool removeElement(int tag);
155  virtual bool removeNode(int tag);
156  virtual bool removeSFreedom_Constraint(int theNode, int theDOF, int loadPatternTag);
157  virtual bool removeSFreedom_Constraint(int tag);
158  virtual bool removeMFreedom_Constraint(int tag);
159  virtual bool removeMRMFreedom_Constraint(int tag);
160  virtual bool removeLoadPattern(int loadTag);
161  virtual bool removeNodeLocker(int nlTag);
162  bool removeLoadPattern(LoadPattern *lp);
163  bool removeNodeLocker(NodeLocker *lp);
165  void removeLPs(void);
166  void removeNLs(void);
167 
168  virtual bool removeNodalLoad(int tag, int loadPattern);
169  virtual bool removeElementalLoad(int tag, int loadPattern);
170  virtual bool removeSFreedom_Constraint(int tag, int loadPattern);
171 
172  virtual void clearDOF_GroupPtr(void);
173 
174  // methods to access the components of a domain
175  virtual ElementIter &getElements(void);
176  virtual NodeIter &getNodes(void);
177  virtual Mesh &getMesh(void);
178  virtual const Mesh &getMesh(void) const;
179  virtual ConstrContainer &getConstraints(void);
180  virtual const ConstrContainer &getConstraints(void) const;
181 
182  const std::string &getNombreCombActual(void) const;
183 
184  bool existElement(int tag);
185  virtual Element *getElement(int tag);
186  virtual const Element *getElement(int tag) const;
187  bool existNode(int tag);
188  virtual Node *getNode(int tag);
189  virtual const Node *getNode(int tag) const;
190 
191 
192  // methods to query the state of the domain
193  inline const PseudoTimeTracker &getTimeTracker(void) const
194  { return timeTracker; }
195  inline int getCurrentGeoTag(void) const
196  { return currentGeoTag; }
197  virtual int getCommitTag(void) const;
198  virtual int getNumElements(void) const;
199  virtual int getNumNodes(void) const;
200  virtual const Vector &getPhysicalBounds(void);
201 
202 
203  // methods to get element and node graphs
204  virtual Graph &getElementGraph(void);
205  virtual Graph &getNodeGraph(void);
206 
207  // methods to update the domain
208  virtual void setCommitTag(int newTag);
209  virtual void setCurrentTime(double newTime);
210  virtual void setCommittedTime(double newTime);
211  virtual void setTime(double newTime);
212  virtual void applyLoad(double pseudoTime);
213  virtual void setLoadConstant(void);
214  virtual int initialize(void);
215  virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF);
216 
217  virtual int commit(void);
218  virtual int revertToLastCommit(void);
219  virtual int revertToStart(void);
220  virtual int update(void);
221  virtual int update(double newTime, double dT);
222  virtual int newStep(double dT);
223 
224  void resetLoadCase(void);
225 
226  // methods for eigenvalue analysis
227  int getNumModes(void) const;
228  virtual int setEigenvalues(const Vector &);
229  virtual const double &getEigenvalue(int) const;
230  double getAngularFrequency(int) const;
231  double getPeriodo(int) const;
232  double getFrecuencia(int) const;
233  virtual const Vector &getEigenvalues(void) const;
234  Vector getAngularFrequencies(void) const;
235  Vector getPeriodos(void) const;
236  Vector getFrecuencias(void) const;
237  virtual int setModalParticipationFactors(const Vector &);
238  virtual const double &getModalParticipationFactor(int mode) const;
239  virtual const Vector &getModalParticipationFactors(void) const;
240  const double getEffectiveModalMass(int mode) const;
241  Vector getEffectiveModalMasses(void) const;
242  double getTotalMass(void) const;
243 
244  // methods for other objects to determine if model has changed
245  virtual void domainChange(void);
246  virtual int hasDomainChanged(void);
247  virtual void setDomainChangeStamp(int newStamp);
248 
249  virtual int addRegion(MeshRegion &theRegion);
250  virtual MeshRegion *getRegion(int region);
251 
252  virtual void Print(std::ostream &s, int flag =0);
253  friend std::ostream &operator<<(std::ostream &s, Domain &M);
254 
255  virtual int sendSelf(CommParameters &);
256  virtual int recvSelf(const CommParameters &);
257  friend int sendDomain(Domain &, int posDbTag,DbTagData &,CommParameters &cp);
258  friend int receiveDomain(Domain &, int posDbTag,DbTagData &,const CommParameters &cp);
259 
260  const Preprocessor *GetPreprocessor(void) const;
262 
263  // nodal methods required in domain interface for parallel interprter
264  virtual double getNodeDisp(int nodeTag, int dof, int &errorFlag);
265  virtual int setMass(const Matrix &mass, int nodeTag);
266 
267  virtual int calculateNodalReactions(bool inclInertia,const double &);
268  virtual int addRecorder(Recorder &theRecorder);
269 
270  static void setDeadSRF(const double &);
271  };
272 
273 int sendDomain(Domain &,int posDbTag,DbTagData &,CommParameters &cp);
274 int receiveDomain(Domain &, int posDbTag,DbTagData &,CommParameters &cp);
275 
276 
277 } // end of XC namespace
278 
279 #endif
280 
281 
A LoadPattern object is used to to store reference loads and single point constraints and a TimeSerie...
Definition: LoadPattern.h:87
virtual int setModalParticipationFactors(const Vector &)
Sets the values of the modal participation factors.
Definition: Domain.cpp:758
virtual bool addSFreedom_Constraint(SFreedom_Constraint *)
Adds to the domain una constraint monopunto.
Definition: Domain.cpp:192
double getAngularFrequency(int) const
Return the angular frequency of the i-th mode.
Definition: Domain.cpp:708
virtual void domainChange(void)
Establece que the model ha cambiado.
Definition: Domain.cpp:799
Single freedom constraints that make part of a load pattern.
Definition: NodeLocker.h:44
virtual const double & getEigenvalue(int) const
Returns the eigenvalue of the i-th mode.
Definition: Domain.cpp:704
virtual bool addMRMFreedom_Constraint(MRMFreedom_Constraint *)
Adds to the domain una constraint multi retained node.
Definition: Domain.cpp:216
void setNodeReactionException(const int &)
Asigna la excepción para comprobación de reacciones (ver Domain::checkNodalReactions).
Definition: Domain.cpp:993
int sendData(CommParameters &cp)
Send data through the channel being passed as parameter.
Definition: Domain.cpp:885
virtual const Vector & getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: Domain.cpp:769
virtual Graph & getNodeGraph(void)
Builds (if necessary) the domain nodes graph and returns a reference to it.
Definition: Domain.cpp:568
virtual int update(void)
Updates the state of the mesh.
Definition: Domain.cpp:673
virtual Mesh & getMesh(void)
Returns a reference to the domain mesh.
Definition: Domain.cpp:499
static void setDeadSRF(const double &)
Assigns Stress Reduction Factor for element deactivation.
Definition: Domain.cpp:179
virtual ConstrContainer & getConstraints(void)
Returns domain constraints.
Definition: Domain.cpp:507
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: Domain.cpp:943
int getNumModes(void) const
Returns the number of computed eigenvalues.
Definition: Domain.cpp:754
virtual int revertToLastCommit(void)
Returns the domain to its last commited state.
Definition: Domain.cpp:632
virtual bool removeNodeLocker(int nlTag)
Elimina del domain el.
Definition: Domain.cpp:383
Finite element mesh.
Definition: Mesh.h:64
void resetLoadCase(void)
Prepares the domain to solve for a new load pattern.
Definition: Domain.cpp:170
Definition: MeshRegion.h:74
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Vector getPeriodos(void) const
Returns a vector with the computed periods (for each mode).
Definition: Domain.cpp:734
double getFrecuencia(int) const
Return the frequency of the i-th mode.
Definition: Domain.cpp:716
double getPeriodo(int) const
Returns the period of the i-th mode.
Definition: Domain.cpp:712
Definition: Vector.h:82
virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF)
Asigna valores a los coeficientes de amortiguamiento de Rayleigh.
Definition: Domain.cpp:606
virtual void clearDOF_GroupPtr(void)
Clears the pointers to DOF groups.
Definition: Domain.cpp:483
virtual bool addNodalLoad(NodalLoad *, int loadPatternTag)
Appends a nodal load to the pattern being passed as parameter.
Definition: Domain.cpp:244
virtual MeshRegion * getRegion(int region)
Returns a pointer to the la región cuyo tag being passed as parameter.
Definition: Domain.cpp:860
Vector que almacena los dbTags de los miembros de la clase.
Definition: DbTagData.h:43
void removeNLs(void)
Elimina del domain todos los bloqueos de nodos.
Definition: Domain.cpp:449
Iterator over an element container.
Definition: ElementIter.h:73
virtual bool removeNode(int tag)
Elimina the node identified by the tag being passed as parameter.
Definition: Domain.cpp:276
int recvData(const CommParameters &cp)
Receive data through the channel being passed as parameter.
Definition: Domain.cpp:908
virtual bool addLoadPattern(LoadPattern *)
Adds al modelo la hipótesis simple being passed as parameter.
Definition: Domain.cpp:316
Vector getFrecuencias(void) const
Returns a vector with the computed frequencies (for each mode).
Definition: Domain.cpp:744
const Preprocessor * GetPreprocessor(void) const
Returns (if possible) a pointer to the preprocessor.
Definition: Domain.cpp:1008
virtual bool removeElementalLoad(int tag, int loadPattern)
Removes from domain the elemental load being passed as parameter.
Definition: Domain.cpp:468
virtual const Vector & getPhysicalBounds(void)
Returns the BND of the model.
Definition: Domain.cpp:558
virtual bool addMFreedom_Constraint(MFreedom_Constraint *)
Adds to the domain una constraint multipunto.
Definition: Domain.cpp:204
Mesh node.
Definition: Node.h:99
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: Domain.cpp:962
virtual int addRecorder(Recorder &theRecorder)
Adds a recorder to the model.
Definition: Domain.cpp:838
virtual void setLoadConstant(void)
Set all the loads as constant.
Definition: Domain.cpp:598
virtual bool addNodeLocker(NodeLocker *)
Adds al modelo.
Definition: Domain.cpp:335
Objeto capaz de gestionar Recorders.
Definition: ObjWithRecorders.h:44
Load pattern combination (1.5*PP+1.0*CP+1.6*SC ...).
Definition: LoadCombination.h:43
virtual int commit(void)
Commits domain state and triggers "record" method for all defined recorders.
Definition: Domain.cpp:611
virtual int setMass(const Matrix &mass, int nodeTag)
Asigna la matriz de masas al nudo cuyo tag being passed as parameter.
Definition: Domain.cpp:989
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
virtual ~Domain(void)
Destructor.
Definition: Domain.cpp:166
virtual bool removeElement(int tag)
Clears the element identified by the tag being passed as parameter.
Definition: Domain.cpp:272
Finite element model generation tools.
Definition: Preprocessor.h:58
Iterador sobre el conteedor de nodos.
Definition: NodeIter.h:73
virtual int getNumElements(void) const
Returns the number of elements.
Definition: Domain.cpp:550
virtual int calculateNodalReactions(bool inclInertia, const double &)
Calculate nodal reaction forces and moments.
Definition: Domain.cpp:1000
virtual int addRegion(MeshRegion &theRegion)
Adds a region.
Definition: Domain.cpp:851
virtual void clearAll(void)
Removes all components from domain (nodes, elements, loads & constraints). GENERAL NOTE ON REMOVAL OF...
Definition: Domain.cpp:133
virtual Node * getNode(int tag)
Returns a pointer to the nodo cuyo tag being passed as parameter.
Definition: Domain.cpp:539
virtual const double & getModalParticipationFactor(int mode) const
Returns the modal participation factor of the i-th mode.
Definition: Domain.cpp:765
Vector getAngularFrequencies(void) const
Returns a vector with the computed angular frequencies (for each mode).
Definition: Domain.cpp:724
DbTagData & getDbTagData(void) const
Returns a vector para almacenar los dbTags de los miembros de la clase.
Definition: Domain.cpp:878
An Recorder object is used in the program to store/restore information at each commit().
Definition: Recorder.h:79
virtual int buildNodeGraph(Graph &theNodeGraph)
Builds the node graph.
Definition: Domain.cpp:873
Vector de celdas.
Definition: DqMeshRegion.h:40
Base class for distributed processing.
Definition: DistributedBase.h:40
virtual bool addLoadCombination(LoadCombination *)
Adds to the domain the load combination being passed as parameter.
Definition: Domain.cpp:347
virtual int initialize(void)
Inicializa.
Definition: Domain.cpp:602
virtual bool removeNodalLoad(int tag, int loadPattern)
Removes from domain the nodal load being passed as parameter.
Definition: Domain.cpp:461
Multiple retained nodes constraint.
Definition: MRMFreedom_Constraint.h:58
bool existNode(int tag)
Returns true if the mesh has a node with the tag being passed as parameter.
Definition: Domain.cpp:534
const double getEffectiveModalMass(int mode) const
Return the effective modal mass of the i-th mode.
Definition: Domain.cpp:773
bool existElement(int tag)
Returns true if the element identified by the tag being passed as parameter already exists en el doma...
Definition: Domain.cpp:519
Definition: Matrix.h:82
virtual bool addElement(Element *)
Adds to the domain the element being passed as parameter.
Definition: Domain.cpp:184
virtual bool removeMFreedom_Constraint(int tag)
Elimina del domain la constraint multipunto cuyo tag being passed as parameter.
Definition: Domain.cpp:297
virtual ElementIter & getElements(void)
Returns an iterator to the element container.
Definition: Domain.cpp:487
Rayleigh damping factors.
Definition: RayleighDampingFactors.h:45
virtual double getNodeDisp(int nodeTag, int dof, int &errorFlag)
Returns the value of dof component of displacement for the node with the tag being passed as paramete...
Definition: Domain.cpp:985
virtual int buildEleGraph(Graph &theEleGraph)
Builds the element graph.
Definition: Domain.cpp:869
The Graph class provides the abstraction of a graph, a collection of vertices and edges...
Definition: Graph.h:84
void removeLoadCombination(LoadCombination *comb)
Removes from the domain the load combination being passed as parameter.
Definition: Domain.cpp:422
virtual bool addNode(Node *)
Adds to the domain el nodo being passed as parameter.
Definition: Domain.cpp:188
virtual int revertToStart(void)
Returns the domain to its initial state and triggers the "restart" method for all the recorders...
Definition: Domain.cpp:650
virtual int setEigenvalues(const Vector &)
Sets eigenvalues.
Definition: Domain.cpp:693
Multi-freedom constraint. Objectt of this class store the information for a multifreedom constraint...
Definition: MFreedom_Constraint.h:84
void removeLPs(void)
Elimina del domain todos los load patterns.
Definition: Domain.cpp:439
virtual NodeIter & getNodes(void)
Returns an iterator a los nodos del domain.
Definition: Domain.cpp:491
virtual int hasDomainChanged(void)
Returns true if the modelo ha cambiado.
Definition: Domain.cpp:803
Single freedom constraint.
Definition: SFreedom_Constraint.h:79
Vector getEffectiveModalMasses(void) const
Returns the effective modal masses for each mode.
Definition: Domain.cpp:777
Communication parameters between processes.
Definition: CommParameters.h:65
Load over a node.
Definition: NodalLoad.h:76
virtual bool removeMRMFreedom_Constraint(int tag)
Elimina del domain la constraint multi retained node cuyo tag being passed as parameter.
Definition: Domain.cpp:307
virtual int getNumNodes(void) const
Returns the número de nodos.
Definition: Domain.cpp:554
virtual const Vector & getEigenvalues(void) const
Returns the eigenvalues vector.
Definition: Domain.cpp:720
virtual bool addElementalLoad(ElementalLoad *, int loadPatternTag)
Adds to the caso being passed as parameter una load over elements.
Definition: Domain.cpp:256
================================================================================
Definition: ContinuaReprComponent.h:34
virtual Graph & getElementGraph(void)
Builds (if necessary) the domain elements graph and returns a reference to it.
Definition: Domain.cpp:563
Constraint (essential and natural boundary conditions) container.
Definition: ConstrContainer.h:61
virtual bool removeLoadPattern(int loadTag)
Elimina del domain el load pattern cuyo tag being passed as parameter.
Definition: Domain.cpp:367
Registro del tiempo.
Definition: PseudoTimeTracker.h:40
virtual void Print(std::ostream &s, int flag=0)
Imprime el domain.
Definition: Domain.cpp:821
virtual Element * getElement(int tag)
Returns a pointer to the element identified by the tag being passed as parameter. ...
Definition: Domain.cpp:524
double getTotalMass(void) const
Return the total effective modal mass.
Definition: Domain.cpp:784
Domain(EntCmd *owr, DataOutputHandler::map_output_handlers *oh)
Constructor.
Definition: Domain.cpp:115
const std::string & getNombreCombActual(void) const
Returns the name of the current load combination.
Definition: Domain.cpp:417