OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
structengngmodel.C
Go to the documentation of this file.
1 /*
2  *
3  * ##### ##### ###### ###### ### ###
4  * ## ## ## ## ## ## ## ### ##
5  * ## ## ## ## #### #### ## # ##
6  * ## ## ## ## ## ## ## ##
7  * ## ## ## ## ## ## ## ##
8  * ##### ##### ## ###### ## ##
9  *
10  *
11  * OOFEM : Object Oriented Finite Element Code
12  *
13  * Copyright (C) 1993 - 2013 Borek Patzak
14  *
15  *
16  *
17  * Czech Technical University, Faculty of Civil Engineering,
18  * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * This program is distributed in the hope that it will be useful,
26  * but WITHOUT ANY WARRANTY; without even the implied warranty of
27  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28  * Lesser General Public License for more details.
29  *
30  * You should have received a copy of the GNU Lesser General Public
31  * License along with this library; if not, write to the Free Software
32  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33  */
34 
35 #include "../sm/EngineeringModels/structengngmodel.h"
36 #include "../sm/Elements/structuralelement.h"
37 #include "../sm/Elements/structuralelementevaluator.h"
38 #include "../sm/Elements/Interfaces/structuralinterfaceelement.h"
39 #include "dofmanager.h"
40 #include "dof.h"
41 #include "element.h"
42 #include "timestep.h"
43 #include "outputmanager.h"
44 #include "activebc.h"
45 #include "assemblercallback.h"
46 #include "unknownnumberingscheme.h"
47 
48 #include "../sm/Materials/structuralmaterial.h"
49 #include "../sm/CrossSections/structuralcrosssection.h"
50 
51 namespace oofem {
52 
54 {
55  //static_cast< StructuralElement & >( element ).giveInternalForcesVector(vec, tStep, 1);
56  element.giveCharacteristicVector(vec, LastEquilibratedInternalForcesVector, mode, tStep);
57 }
58 
60 {
61  StructuralElement &selem = static_cast< StructuralElement & >( element );
62 
63  vec.clear();
64  for ( auto &gp : *selem.giveDefaultIntegrationRulePtr() ) {
65  FloatArray epsilonTemperature;
66 
68  static_cast< StructuralMaterial *>( selem.giveStructuralCrossSection()->giveMaterial(gp) )->computeStressIndependentStrainVector(epsilonTemperature, gp, tStep, VM_Incremental);
69 
70  if ( epsilonTemperature.giveSize() > 0 ) {
71  FloatMatrix D, B;
72  FloatArray s;
73  double dV = selem.computeVolumeAround(gp);
74  selem.computeBmatrixAt(gp, B);
75  selem.computeConstitutiveMatrixAt(D, ElasticStiffness, gp, tStep);
76  s.beProductOf(D, epsilonTemperature);
77  vec.plusProduct(B, s, dV);
78  }
79  }
80 }
81 
83 {
84  static_cast< StructuralElement & >( element ).computeInitialStressMatrix(answer, tStep);
85 }
86 
87 
88 
90  internalVarUpdateStamp(0), internalForcesEBENorm()
91 { }
92 
93 
95 { }
96 
97 
98 void
100 //
101 // computes and prints reaction forces in all supported or restrained dofs
102 //
103 {
104  FloatArray reactions;
105  IntArray dofManMap, dofidMap, eqnMap;
106 
107  Domain *domain = this->giveDomain(di);
108 
109  // map contains corresponding dofmanager and dofs numbers corresponding to prescribed equations
110  // sorted according to dofmanger number and as a minor crit. according to dof number
111  // this is necessary for extractor, since the sorted output is expected
112  this->buildReactionTable(dofManMap, dofidMap, eqnMap, tStep, di);
113 
114  //
115  // print header
116  //
117  fprintf(out, "\n\n\tR E A C T I O N S O U T P U T:\n\t_______________________________\n\n\n");
118 
119  // compute reaction forces
120  this->computeReaction(reactions, tStep, di);
121 
122  //
123  // loop over reactions and print them
124  //
125  for ( int i = 1; i <= dofManMap.giveSize(); i++ ) {
126  if ( domain->giveOutputManager()->testDofManOutput(dofManMap.at(i), tStep) ) {
127  fprintf(out, "\tNode %8d iDof %2d reaction % .4e [bc-id: %d]\n",
128  domain->giveDofManager( dofManMap.at(i) )->giveLabel(),
129  dofidMap.at(i), reactions.at( eqnMap.at(i) ),
130  domain->giveDofManager( dofManMap.at(i) )->giveDofWithID( dofidMap.at(i) )->giveBcId() );
131  }
132  }
133 }
134 
137 }
138 
139 void
141  IntArray &eqn, TimeStep *tStep, int di)
142 {
143  // determine number of restrained dofs
144  Domain *domain = this->giveDomain(di);
146  int ndofMan = domain->giveNumberOfDofManagers();
147  int rindex, count = 0;
148 
149  // initialize corresponding dofManagers and dofs for each restrained dof
150  restrDofMans.resize(numRestrDofs);
151  restrDofs.resize(numRestrDofs);
152  eqn.resize(numRestrDofs);
153 
154  for ( int i = 1; i <= ndofMan; i++ ) {
155  DofManager *inode = domain->giveDofManager(i);
156  for ( Dof *jdof: *inode ) {
157  if ( jdof->isPrimaryDof() && ( jdof->hasBc(tStep) ) ) { // skip slave dofs
158  rindex = jdof->__givePrescribedEquationNumber();
159  if ( rindex ) {
160  count++;
161  restrDofMans.at(count) = i;
162  restrDofs.at(count) = jdof->giveDofID();
163  eqn.at(count) = rindex;
164  } else {
165  // NullDof has no equation number and no prescribed equation number
166  //_error("No prescribed equation number assigned to supported DOF");
167  }
168  }
169  }
170  }
171  // Trim to size.
172  restrDofMans.resizeWithValues(count);
173  restrDofs.resizeWithValues(count);
174  eqn.resizeWithValues(count);
175 }
176 
177 
178 void
180 {
181  FloatArray contribution;
182 
184  answer.zero();
185 
186  // Add internal forces
187  this->assembleVector( answer, tStep, LastEquilibratedInternalForceAssembler(), VM_Total,
189  // Subtract external loading
191  //this->assembleVector( answer, tStep, ExternalForceAssembler(), VM_Total,
192  // EModelDefaultPrescribedEquationNumbering(), this->giveDomain(di) );
194  this->computeExternalLoadReactionContribution(contribution, tStep, di);
195  answer.subtract(contribution);
197 }
198 
199 
200 void
202 {
204  reactions.zero();
205  this->assembleVector( reactions, tStep, ExternalForceAssembler(), VM_Total,
207 }
208 
209 
210 void
211 StructuralEngngModel :: giveInternalForces(FloatArray &answer, bool normFlag, int di, TimeStep *tStep)
212 {
213  // Simply assembles contributions from each element in domain
214  Domain *domain = this->giveDomain(di);
215  // Update solution state counter
216  tStep->incrementStateCounter();
217 
219  answer.zero();
220  this->assembleVector(answer, tStep, InternalForceAssembler(), VM_Total,
221  EModelDefaultEquationNumbering(), domain, normFlag ? & this->internalForcesEBENorm : NULL);
222 
223  // Redistributes answer so that every process have the full values on all shared equations
225 
226  // Remember last internal vars update time stamp.
228 }
229 
230 
231 void
233 {
234  this->updateInternalState(tStep);
236 }
237 
238 
239 int
241 {
242  Domain *domain = this->giveDomain(1);
243  // check for proper element type
244  for ( auto &elem : domain->giveElements() ) {
245  StructuralElement *sePtr = dynamic_cast< StructuralElement * >( elem.get() );
246  StructuralInterfaceElement *siePtr = dynamic_cast< StructuralInterfaceElement * >( elem.get() );
247  StructuralElementEvaluator *see = dynamic_cast< StructuralElementEvaluator * >( elem.get() );
248 
249  if ( sePtr == NULL && see == NULL && siePtr == NULL ) {
250  OOFEM_WARNING("Element %d has no structural support", elem->giveLabel());
251  return 0;
252  }
253  }
254 
256 
257  return 1;
258 }
259 
260 
261 void
263 {
264  if ( !this->giveDomain(1)->giveOutputManager()->testTimeStepOutput(tStep) ) {
265  return; // Do not print even Solution step header
266  }
267 
268  EngngModel :: printOutputAt(file, tStep);
269  this->printReactionForces(tStep, 1, file);
270 }
271 
272 
273 void
275 {
276  for ( auto &domain: domainList ) {
278  for ( auto &dman : domain->giveDofManagers() ) {
279  this->updateDofUnknownsDictionary(dman.get(), tStep);
280  }
281  }
282 
283  for ( auto &bc : domain->giveBcs() ) {
285 
286  if ( ( abc = dynamic_cast< ActiveBoundaryCondition * >(bc.get()) ) ) {
287  int ndman = abc->giveNumberOfInternalDofManagers();
288  for ( int j = 1; j <= ndman; j++ ) {
290  }
291  }
292  }
293 
295  for ( auto &elem : domain->giveElements() ) {
296  elem->updateInternalState(tStep);
297  }
298 
300  }
301  }
302 }
303 
304 
305 
306 #ifdef __OOFEG
307 void
309 {
310  Domain *domain = this->giveDomain(1);
311 
312  if ( type != 1 ) {
313  return;
314  }
315 
316  for ( auto &elem : domain->giveElements() ) {
317  elem->showSparseMtrxStructure(TangentStiffnessMatrix, gc, tStep);
318  }
319 }
320 #endif
321 } // end namespace oofem
This class represent a new concept on how to define elements.
int testDofManOutput(int, TimeStep *)
Tests if given dof manager is required to do its output for given time step.
The representation of EngngModel default unknown numbering.
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
Class and object Domain.
Definition: domain.h:115
virtual void vectorFromElement(FloatArray &vec, Element &element, TimeStep *tStep, ValueModeType mode) const
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
Implementation for assembling internal forces vectors in standard monolithic, nonlinear FE-problems...
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Returns number of equations for given domain in active (current time step) time step.
Definition: engngm.C:391
int giveNumberOfDofManagers() const
Returns number of dof managers in domain.
Definition: domain.h:432
std::vector< std::unique_ptr< Domain > > domainList
List of problem domains.
Definition: engngm.h:207
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
Computes constitutive matrix of receiver.
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual int checkConsistency()
Allows programmer to test some receiver&#39;s internal data, before computation begins.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual int requiresUnknownsDictionaryUpdate()
Indicates if EngngModel requires Dofs dictionaries to be updated.
Definition: engngm.h:845
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual DofManager * giveInternalDofManager(int i)
Gives an internal dof manager from receiver.
void incrementStateCounter()
Updates solution state counter.
Definition: timestep.h:190
Abstract base class for all finite elements.
Definition: element.h:145
virtual void terminate(TimeStep *tStep)
Terminates the solution of time step.
Definition: engngm.C:656
Base class for dof managers.
Definition: dofmanager.h:113
StateCounterType internalVarUpdateStamp
Contains last time stamp of internal variable update.
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual void giveInternalForces(FloatArray &answer, bool normFlag, int di, TimeStep *tStep)
Evaluates the nodal representation of internal forces by assembling contributions from individual ele...
void printReactionForces(TimeStep *tStep, int id, FILE *out)
Computes and prints reaction forces, computed from nodal internal forces.
virtual int giveNumberOfInternalDofManagers()
Gives the number of internal dof managers.
void buildReactionTable(IntArray &restrDofMans, IntArray &restrDofs, IntArray &eqn, TimeStep *tStep, int di)
Builds the reaction force table.
Abstract base class for all "structural" finite elements.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: element.h:518
StateCounterType giveSolutionStateCounter()
Returns current solution state counter.
Definition: timestep.h:188
int updateSharedDofManagers(FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
Exchanges necessary remote DofManagers data.
Definition: engngm.C:1957
virtual void terminate(TimeStep *tStep)
Terminates the solution of time step.
Implementation for assembling external forces vectors in standard monolithic FE-problems.
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void updateYourself(TimeStep *tStep)
Updates internal state after finishing time step.
Definition: engngm.C:612
void resizeWithValues(int n, int allocChunk=0)
Checks size of receiver towards requested bounds.
Definition: intarray.C:115
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
Abstract base class for all active boundary conditions.
Definition: activebc.h:63
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual void updateYourself(TimeStep *tStep)
Updates internal state after finishing time step.
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual int checkConsistency()
Allows programmer to test some receiver&#39;s internal data, before computation begins.
Definition: engngm.h:995
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void updateDofUnknownsDictionary(DofManager *, TimeStep *)
Updates necessary values in Dofs unknown dictionaries.
Definition: engngm.h:859
The representation of EngngModel default prescribed unknown numbering.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to output domain stream, for given time step.
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to output domain stream, for given time step.
Definition: engngm.C:695
virtual void showSparseMtrxStructure(int type, oofegGraphicContext &gc, TimeStep *tStep)
Shows the sparse structure of required matrix, type == 1 stiffness.
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)=0
Computes the geometrical matrix of receiver in given integration point.
Abstract base class for all "structural" constitutive models.
OutputManager * giveOutputManager()
Returns domain output manager.
Definition: domain.C:1436
FloatArray internalForcesEBENorm
Norm of nodal internal forces evaluated on element by element basis (squared)
std::vector< std::unique_ptr< Element > > & giveElements()
Definition: domain.h:279
void computeReaction(FloatArray &answer, TimeStep *tStep, int di)
Computes reaction forces.
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
Definition: element.C:580
Abstract base class representing the "problem" under consideration.
Definition: engngm.h:181
Abstract base class for all structural interface elements.
int giveSize() const
Definition: intarray.h:203
virtual void matrixFromElement(FloatMatrix &mat, Element &element, TimeStep *tStep) const
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
the oofem namespace is to define a context or scope in which all oofem names are defined.
Domain * giveDomain(int n)
Service for accessing particular problem domain.
Definition: engngm.C:1720
virtual ~StructuralEngngModel()
Destructor.
Abstract class Dof represents Degree Of Freedom in finite element mesh.
Definition: dof.h:93
DofManager * giveDofManager(int n)
Service for accessing particular domain dof manager.
Definition: domain.C:314
virtual void computeExternalLoadReactionContribution(FloatArray &reactions, TimeStep *tStep, int di)
Computes the contribution external loading to reaction forces in given domain.
void updateInternalState(TimeStep *tStep)
Updates nodal values (calls also this->updateDofUnknownsDictionary for updating dofs unknowns diction...
#define OOFEM_WARNING(...)
Definition: error.h:62
Assembles the internal forces, without updating the strain.
void assembleVector(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
Assembles characteristic vector of required type from dofManagers, element, and active boundary condi...
Definition: engngm.C:986
Class representing solution step.
Definition: timestep.h:80
virtual Material * giveMaterial(IntegrationPoint *ip)
Returns the material associated with the GP.
StructuralEngngModel(int i, EngngModel *_master=NULL)
Creates new StructuralEngngModel with number i, associated to domain d.
virtual void vectorFromElement(FloatArray &vec, Element &element, TimeStep *tStep, ValueModeType mode) const
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
Definition: floatarray.C:226

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:31 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011