xc
ShadowSubdomain.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.5 $
48 // $Date: 2005/11/30 23:47:00 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/domain/subdomain/ShadowSubdomain.h,v $
50 
51 
52 #ifndef ShadowSubdomain_h
53 #define ShadowSubdomain_h
54 
55 // Written: fmk
56 // Revision: A
57 //
58 // Description: This file contains the class definition for ShadowSubdomain.
59 // ShadowSubdomain is a container class. The class is responsible for holding
60 // and providing access to the Elements, Nodes, LoadCases, SFreedom_Constraints
61 // and MFreedom_Constraints that have been added to the ShadowSubdomain.
62 //
63 // What: "@(#) ShadowSubdomain.h, revA"
64 
65 #include "Subdomain.h"
66 #include <utility/actor/shadow/Shadow.h>
67 #include "utility/matrix/ID.h"
68 
69 namespace XC {
71 //
73 class ShadowSubdomain: public Shadow, public Subdomain
74  {
75  private:
76  mutable ID msgData;
77  ID theElements;
78  ID theNodes;
79  ID theExternalNodes;
80  ID theLoadCases;
81 
82  TaggedObjectStorage *theShadowSPs;
83  TaggedObjectStorage *theShadowMPs;
84  TaggedObjectStorage *theShadowLPs;
85 
86  int numDOF;
87  int numElements;
88  int numNodes;
89  int numExternalNodes;
90  int numSPs;
91  int numMPs;
92  int numLoadPatterns;
93 
94  bool buildRemote;
95  bool gotRemoteData;
96 
97  FE_Element *theFEele;
98 
99  mutable Vector *theVector; // for storing residual info
100  mutable Matrix *theMatrix; // for storing tangent info
101 
102  static char *shadowSubdomainProgram;
103 
104  static int count; // MHS
105  static int numShadowSubdomains;
106  static std::deque<ShadowSubdomain *> theShadowSubdomains;
107 
108  void free_vectors(void) const;
109  void alloc_vectors(const size_t &) const;
110  void resize_vectors(const size_t &) const;
111  void free_arrays(void);
112  void alloc_arrays(const size_t &,const size_t &);
113  protected:
114  virtual int buildMap(void) const;
115  virtual int buildEleGraph(Graph &theEleGraph);
116  virtual int buildNodeGraph(Graph &theNodeGraph);
117 
118  public:
119  ShadowSubdomain(int tag, MachineBroker &theMachineBroker, FEM_ObjectBroker &theObjectBroker,DataOutputHandler::map_output_handlers *oh,CommandEntity *owr);
120  ShadowSubdomain(int tag, Channel &, FEM_ObjectBroker &,DataOutputHandler::map_output_handlers *,CommandEntity *);
121  virtual ~ShadowSubdomain(void);
122 
123  // method added for parallel domain generation
124  virtual int buildSubdomain(int numSubdomains,
125  PartitionedModelBuilder &theBuilder);
126  virtual int getRemoteData(void);
127 
128  // Methods inherited from Domain, Subdomain and Element
129  // which must be rewritten
130 
131  virtual bool addElement(Element *);
132  virtual bool addNode(Node *);
133  virtual bool addExternalNode(Node *);
136  virtual bool addLoadPattern(LoadPattern *);
137  virtual bool addNodalLoad(NodalLoad *, int loadPattern);
138  virtual bool addElementalLoad(ElementalLoad *, int loadPattern);
139  virtual bool addSFreedom_Constraint(SFreedom_Constraint *, int loadPattern);
140 
141 
142  virtual bool hasNode(int tag);
143  virtual bool hasElement(int tag);
144 
145  virtual void clearAll(void);
146  virtual bool removeElement(int tag);
147  virtual bool removeNode(int tag);
148  virtual bool removeSFreedom_Constraint(int tag);
149  virtual bool removeMFreedom_Constraint(int tag);
150  virtual bool removeLoadPattern(int tag);
151  virtual bool removeNodalLoad(int tag, int loadPattern);
152  virtual bool removeElementalLoad(int tag, int loadPattern);
153  virtual bool removeSFreedom_Constraint(int tag, int loadPattern);
154 
155  virtual ElementIter &getElements();
156  virtual NodeIter &getNodes();
157  virtual NodeIter &getInternalNodeIter(void);
158  virtual NodeIter &getExternalNodeIter(void);
159 
160  virtual Element *getElementPtr(int tag);
161  virtual Node *getNodePtr(int tag);
162 
163  virtual int getNumElements(void) const;
164  virtual int getNumNodes(void) const;
165  virtual int getNumSPs(void) const;
166  virtual int getNumMPs(void) const;
167  virtual int getNumLoadPatterns(void) const;
168 
169  virtual Graph &getElementGraph(void);
170  virtual Graph &getNodeGraph(void);
171 
172  // methods to update the domain
173  virtual void setCommitTag(int newTag);
174  virtual void setCurrentTime(double newTime);
175  virtual void setCommittedTime(double newTime);
176  virtual void applyLoad(double pseudoTime);
177  virtual void setLoadConstant(void);
178 
179  virtual int update(void);
180  virtual int update(double newTime, double dT);
181  virtual int commit(void);
182  virtual int revertToLastCommit(void);
183  virtual int revertToStart(void);
184  virtual int barrierCheckIN(void);
185  virtual int barrierCheckOUT(int);
186 
187  virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF);
188 
189  virtual int addRecorder(Recorder &theRecorder);
190  virtual int removeRecorders(void);
191 
192  virtual void wipeAnalysis(void);
193  virtual void setDomainDecompAnalysis(DomainDecompositionAnalysis &theAnalysis);
194  virtual int setAnalysisAlgorithm(EquiSolnAlgo &theAlgorithm);
195  virtual int setAnalysisIntegrator(IncrementalIntegrator &theIntegrator);
196  virtual int setAnalysisLinearSOE(LinearSOE &theSOE);
197  virtual int setAnalysisConvergenceTest(ConvergenceTest &theTest);
198  virtual void clearAnalysis(void);
199  virtual void domainChange(void);
200 
201  virtual int getNumExternalNodes(void) const;
202  virtual const ID &getExternalNodes(void) const;
203  virtual int getNumDOF(void) const;
204 
205  virtual const Matrix &getTang(void);
206  virtual const Vector &getResistingForce(void) const;
207 
208  virtual int computeTang(void);
209  virtual int computeResidual(void);
210 
211  const Vector &getLastExternalSysResponse(void);
212  virtual int computeNodalResponse(void);
213  virtual int newStep(double deltaT);
214 
215  virtual int sendSelf(CommParameters &);
216  virtual int recvSelf(const CommParameters &);
217 
218  virtual double getCost(void);
219 
220  virtual void Print(std::ostream &s, int flag =0);
221 
222  // nodal methods required in domain interface for parallel interprter
223  virtual double getNodeDisp(int nodeTag, int dof, int &errorFlag);
224  virtual int setMass(const Matrix &mass, int nodeTag);
225  };
226 } // end of XC namespace
227 
228 #endif
virtual int computeResidual(void)
The method first starts a Timer object running.
Definition: ShadowSubdomain.cpp:1121
virtual void applyLoad(double pseudoTime)
Apply the loads for the given time pseudoTime.
Definition: ShadowSubdomain.cpp:785
virtual int setRayleighDampingFactors(const RayleighDampingFactors &rF)
Set Rayleigh damping factors.
Definition: ShadowSubdomain.cpp:881
virtual double getCost(void)
Return the current value of realCost.
Definition: ShadowSubdomain.cpp:1197
Float vector abstraction.
Definition: Vector.h:93
Domain enclosed in another domain.
Definition: Subdomain.h:101
virtual bool removeElementalLoad(int tag, int loadPattern)
Removes the elemental load from container.
Definition: ShadowSubdomain.cpp:620
virtual int computeTang(void)
The method first starts a Timer object running.
Definition: ShadowSubdomain.cpp:1095
virtual bool removeLoadPattern(int tag)
Remove from domain el load pattern identified by the argument.
Definition: ShadowSubdomain.cpp:586
virtual bool addElementalLoad(ElementalLoad *, int loadPattern)
Adds a load over element to the pattern.
Definition: ShadowSubdomain.cpp:451
Linear system of equations.
Definition: LinearSOE.h:91
virtual NodeIter & getInternalNodeIter(void)
Return an iterator to the internal nodes of the subdomain, nodes that are added using the addNode() c...
Definition: ShadowSubdomain.cpp:685
virtual int buildNodeGraph(Graph &theNodeGraph)
Builds the node graph.
Definition: ShadowSubdomain.cpp:1242
virtual int commit(void)
invokes the base Domain classes commit() method.
Definition: ShadowSubdomain.cpp:921
An Recorder object is used in the program to store/restore information at each commit().
Definition: Recorder.h:86
Iterator over an element container.
Definition: ElementIter.h:73
virtual int getNumDOF(void) const
Returns the num of external dof associated with the subdomain.
Definition: ShadowSubdomain.cpp:1051
virtual int addRecorder(Recorder &theRecorder)
Adds a recorder to the model.
Definition: ShadowSubdomain.cpp:905
virtual void Print(std::ostream &s, int flag=0)
Print stuff.
Definition: ShadowSubdomain.cpp:1226
virtual const Matrix & getTang(void)
Return the Matrix obtained from invoking getTangent() on the DomainDecompositionAnalysis object...
Definition: ShadowSubdomain.cpp:1061
virtual bool removeSFreedom_Constraint(int tag)
Removes from domain the single freedom constraint identified by the argument.
Definition: ShadowSubdomain.cpp:560
FEM_ObjectBroker is is an object broker class for the finite element method.
Definition: FEM_ObjectBroker.h:145
virtual Graph & getElementGraph(void)
Builds (if necessary) the domain elements graph and returns a reference to it.
Definition: ShadowSubdomain.cpp:764
Vector of integers.
Definition: ID.h:93
virtual bool removeNode(int tag)
Remove a node from the subdomain.
Definition: ShadowSubdomain.cpp:524
const Vector & getLastExternalSysResponse(void)
Return the Vector obtained by calling getLastSysResponse() on the associated FE_Element.
Definition: ShadowSubdomain.cpp:1148
A LoadPattern object is used to to store reference loads and single point constraints and a TimeSerie...
Definition: LoadPattern.h:93
virtual int revertToLastCommit(void)
Return the domain to its last commited state.
Definition: ShadowSubdomain.cpp:932
Base class for the finite elements.
Definition: Element.h:109
virtual void clearAll(void)
Removes all components from domain (nodes, elements, loads & constraints).
Definition: ShadowSubdomain.cpp:897
virtual void setLoadConstant(void)
Set all the loads as constant.
Definition: ShadowSubdomain.cpp:831
virtual bool addElement(Element *)
Adds to the domain the element being passed as parameter.
Definition: ShadowSubdomain.cpp:267
virtual int getNumNodes(void) const
Returns the number of external and internal Nodes.
Definition: ShadowSubdomain.cpp:753
virtual int recvSelf(const CommParameters &)
Receive itself.
Definition: ShadowSubdomain.cpp:1218
Used when performing a domain decomposition analysis.
Definition: DomainDecompositionAnalysis.h:90
virtual bool addLoadPattern(LoadPattern *)
Appends the load pattern to the domain.
Definition: ShadowSubdomain.cpp:379
Local object associated with an actor.
Definition: Shadow.h:90
virtual bool removeMFreedom_Constraint(int tag)
Removes from domain the multi-freedom constraint identified by the argument.
Definition: ShadowSubdomain.cpp:573
Single freedom constraint.
Definition: SFreedom_Constraint.h:84
The Graph class provides the abstraction of a graph.
Definition: Graph.h:93
IncrementalIntegrator is an algorithmic class for setting up the finite element equations in an incre...
Definition: IncrementalIntegrator.h:96
Local representation of a remote subdomain.
Definition: ShadowSubdomain.h:73
virtual int buildEleGraph(Graph &theEleGraph)
Builds the element graph.
Definition: ShadowSubdomain.cpp:1235
virtual int update(void)
Updates the state of the domain.
Definition: ShadowSubdomain.cpp:839
Channel is an abstract base class which defines the channel interface.
Definition: Channel.h:91
EquiSolnAlgo is an abstract base class, i.e.
Definition: EquiSolnAlgo.h:88
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: ShadowSubdomain.cpp:1257
virtual bool addExternalNode(Node *)
A Method to add the node pointed to by the argument.
Definition: ShadowSubdomain.cpp:319
virtual const Vector & getResistingForce(void) const
Return the Vector obtained from invoking getCondensedRHS() on the DomainDecompositionAnalysis object...
Definition: ShadowSubdomain.cpp:1078
Rayleigh damping factors.
Definition: RayleighDampingFactors.h:58
virtual int sendSelf(CommParameters &)
Send itself.
Definition: ShadowSubdomain.cpp:1211
A MachineBroker is responsible for getting an actor process running on the parallel machine...
Definition: MachineBroker.h:76
convergence test.
Definition: ConvergenceTest.h:80
virtual ElementIter & getElements()
Returns an iterator to the element container.
Definition: ShadowSubdomain.cpp:665
virtual void setCommittedTime(double newTime)
Set the committed time to newTime.
Definition: ShadowSubdomain.cpp:821
virtual int computeNodalResponse(void)
Set the nodal responses for the nodes in the subdomain.
Definition: ShadowSubdomain.cpp:1158
virtual Graph & getNodeGraph(void)
Builds (if necessary) the domain node graph and returns a reference to it.
Definition: ShadowSubdomain.cpp:774
virtual void setDomainDecompAnalysis(DomainDecompositionAnalysis &theAnalysis)
Sets the corresponding DomainDecompositionAnalysis object to be {theAnalysis}.
Definition: ShadowSubdomain.cpp:958
Finite element as seen by analysis.
Definition: FE_Element.h:107
virtual int revertToStart(void)
Return the domain to its initial state and triggers the "restart" method for all the recorders...
Definition: ShadowSubdomain.cpp:943
The PartitionedModelBuilder class is an abstract class.
Definition: PartitionedModelBuilder.h:86
virtual void setCommitTag(int newTag)
Set the committed tag to newTag.
Definition: ShadowSubdomain.cpp:800
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
virtual bool addNodalLoad(NodalLoad *, int loadPattern)
Appends a nodal load to the pattern being passed as parameter.
Definition: ShadowSubdomain.cpp:426
virtual int getNumExternalNodes(void) const
Returns the number of external nodes that have been successfully added to the subdomain as external n...
Definition: ShadowSubdomain.cpp:1032
Iterator over the nodes.
Definition: NodeIter.h:73
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
Multi-freedom constraint.
Definition: MFreedom_Constraint.h:110
virtual bool addNode(Node *)
Method to add a node to the subdomain.
Definition: ShadowSubdomain.cpp:299
virtual const ID & getExternalNodes(void) const
Returns an ID containing the tags of all nodes added to the subdomain as external nodes and have yet ...
Definition: ShadowSubdomain.cpp:1042
virtual int setMass(const Matrix &mass, int nodeTag)
Set the mass matrix for the node identified by the argument.
Definition: ShadowSubdomain.cpp:1273
virtual bool removeNodalLoad(int tag, int loadPattern)
Removes the nodal load from container.
Definition: ShadowSubdomain.cpp:600
Base class for loads over elements.
Definition: ElementalLoad.h:77
virtual bool removeElement(int tag)
Remove the element identified by the argument.
Definition: ShadowSubdomain.cpp:493
Mesh node.
Definition: Node.h:110
virtual bool addMFreedom_Constraint(MFreedom_Constraint *)
Adds to the domain a multi-freedom constraint.
Definition: ShadowSubdomain.cpp:361
virtual int getNumElements(void) const
Return the number of elements.
Definition: ShadowSubdomain.cpp:747
virtual int removeRecorders(void)
Remove the recorders.
Definition: ShadowSubdomain.cpp:914
virtual bool addSFreedom_Constraint(SFreedom_Constraint *)
Adds a single freedom constraint to the domain.
Definition: ShadowSubdomain.cpp:342
virtual void setCurrentTime(double newTime)
Set the current time to newTime.
Definition: ShadowSubdomain.cpp:808
Load over a node.
Definition: NodalLoad.h:82
virtual void domainChange(void)
Sets a flag indicating that the integer returned in the next call to hasDomainChanged() must be incre...
Definition: ShadowSubdomain.cpp:1018
virtual NodeIter & getNodes()
Return an iter to all nodes that have been added to the subdomain.
Definition: ShadowSubdomain.cpp:675