OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
linearstability.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 // please activate or de-activate next line
36 //#define LIN_STAB_COMPATIBILITY_MODE
37 
38 #include "../sm/EngineeringModels/linearstability.h"
39 #include "timestep.h"
40 #include "element.h"
41 #include "contextioerr.h"
42 #include "floatmatrix.h"
43 #include "verbose.h"
44 #include "floatarray.h"
45 #include "classfactory.h"
46 #include "datastream.h"
47 #include "exportmodulemanager.h"
48 #include "dofmanager.h"
49 #include "dof.h"
50 #include "unknownnumberingscheme.h"
51 #include "outputmanager.h"
52 
53 #ifdef __OOFEG
54  #include "oofeggraphiccontext.h"
55 #endif
56 
57 namespace oofem {
58 REGISTER_EngngModel(LinearStability);
59 
61 {
62  if ( !nMethod ) {
63  nMethod.reset( classFactory.createGeneralizedEigenValueSolver(solverType, this->giveDomain(1), this) );
64  if ( !nMethod ) {
65  OOFEM_ERROR("solver creation failed");
66  }
67  }
68 
69  return nMethod.get();
70 }
71 
73 {
74  if ( !nMethodLS ) {
75  nMethodLS.reset( classFactory.createSparseLinSolver(ST_Direct, this->giveDomain(1), this) );
76  if ( !nMethodLS ) {
77  OOFEM_ERROR("solver creation failed");
78  }
79  }
80 
81  return nMethodLS.get();
82 }
83 
86 {
87  IRResultType result; // Required by IR_GIVE_FIELD macro
88 
89  //StructuralEngngModel::instanciateFrom(ir);
91  // numberOfSteps set artifficially to numberOfRequiredEigenValues
92  // in order to allow
93  // use restoreContext function for different eigenValues
95 
97  if ( rtolv < 1.e-12 ) {
98  rtolv = 1.e-12;
99  }
100 
101  if ( rtolv > 0.01 ) {
102  rtolv = 0.01;
103  }
104 
105  int val = 0;
108 
109 
110  nMetaSteps = 0;
111 
113 
114  if(suppressOutput) {
115  printf("Suppressing output.\n");
116  }
117  else {
118 
119  if ( ( outputStream = fopen(this->dataOutputFileName.c_str(), "w") ) == NULL ) {
120  OOFEM_ERROR("Can't open output file %s", this->dataOutputFileName.c_str());
121  }
122 
123  fprintf(outputStream, "%s", PRG_HEADER);
124  fprintf(outputStream, "\nStarting analysis on: %s\n", ctime(& this->startTime) );
125  fprintf(outputStream, "%s\n", simulationDescription.c_str());
126  }
127 
128  return IRRT_OK;
129 }
130 
131 
133 // returns unknown quantity like displacement, eigen value.
134 {
135  int eq = dof->__giveEquationNumber();
136 #ifdef DEBUG
137  if ( eq == 0 ) {
138  OOFEM_ERROR("invalid equation number");
139  }
140 #endif
141 
142  int activeVector = ( int ) tStep->giveTargetTime();
143  switch ( mode ) {
144  case VM_Total: // EigenVector
145  if ( activeVector ) {
146  return eigVec.at(eq, activeVector);
147  }
148 
149  return displacementVector.at(eq);
150 
151  default:
152  OOFEM_ERROR("Unknown is of undefined type for this problem");
153  }
154 
155  return 0.;
156 }
157 
159 {
160  int istep = giveNumberOfFirstStep();
161  StateCounterType counter = 1;
162 
163  if ( currentStep ) {
164  istep = currentStep->giveNumber() + 1;
165  counter = currentStep->giveSolutionStateCounter() + 1;
166  }
167 
168  previousStep = std :: move(currentStep);
169  currentStep.reset( new TimeStep(istep, this, 1, 0., 0., counter) );
170 
171  return currentStep.get();
172 }
173 
175 {
177  // update state according to new meta step
178  this->giveNextStep();
179  this->updateAttributes( this->giveCurrentMetaStep() );
180  this->solveYourselfAt( this->giveCurrentStep() );
181  this->terminate( this->giveCurrentStep() );
182 }
183 
184 
186 {
187  //
188  // creates system of governing eq's and solves them at given time step
189  //
190  this->giveNumericalMethod( this->giveMetaStep( tStep->giveMetaStepNumber() ) );
192 
193  // first assemble problem at current time step
194 
195  if ( tStep->giveNumber() == 1 ) {
196  //
197  // first step - solve linear static problem
198  //
200  stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
201 
202  //
203  // allocate space for displacement Vector
204  //
206  //
207  // allocate space for load vector
208  //
210  }
211 
212 #ifndef LIN_STAB_COMPATIBILITY_MODE
213  #ifdef VERBOSE
214  OOFEM_LOG_INFO("Assembling stiffness matrix\n");
215  #endif
216  stiffnessMatrix->zero();
217  this->assemble( *stiffnessMatrix, tStep, TangentAssembler(TangentStiffness),
219 #endif
220 
221 
222 #ifdef VERBOSE
223  OOFEM_LOG_INFO("Assembling load\n");
224 #endif
225 
227  loadVector.zero();
228 
229  // Internal forces first, negated;
230  this->assembleVector( loadVector, tStep, InternalForceAssembler(), VM_Total,
233 
234  this->assembleVector( loadVector, tStep, ExternalForceAssembler(), VM_Total,
237 
238  //
239  // call numerical model to solve problem
240  //
241 #ifdef VERBOSE
242  OOFEM_LOG_INFO("Solving linear static problem\n");
243 #endif
244 
246  // terminate linear static computation (necessary, in order to compute stresses in elements).
247  this->terminateLinStatic( this->giveCurrentStep() );
248  /*
249  * Normal forces already known, proceed with linear stability
250  */
251 
252  stiffnessMatrix->zero();
253  if ( !initialStressMatrix ) {
254  initialStressMatrix.reset( stiffnessMatrix->GiveCopy() );
255  } else {
256  initialStressMatrix->zero();
257  }
258 
259 #ifdef VERBOSE
260  OOFEM_LOG_INFO("Assembling stiffness matrix\n");
261 #endif
262  this->assemble( *stiffnessMatrix, tStep, TangentAssembler(TangentStiffness),
264 #ifdef VERBOSE
265  OOFEM_LOG_INFO("Assembling initial stress matrix\n");
266 #endif
269  initialStressMatrix->times(-1.0);
270 
271  // stiffnessMatrix->printYourself();
272  // initialStressMatrix->printYourself();
273 
274  //
275  // create resulting objects eigVec and eigVal
276  //
279  eigVec.zero();
280  eigVal.zero();
281 
282  //
283  // call numerical model to solve arised problem
284  //
285 #ifdef VERBOSE
286  OOFEM_LOG_INFO("Solving ...\n");
287 #endif
288 
290  // compute eigen frequencies
291  //for (i = 1; i<= numberOfRequiredEigenValues; i++)
292  // eigVal.at(i) = sqrt(eigVal.at(i));
293 }
294 
296 { }
297 
298 void
300 {
301  Domain *domain = this->giveDomain(1);
302  tStep->setTime(0.);
303 
305  for ( auto &dman : domain->giveDofManagers() ) {
306  this->updateDofUnknownsDictionary(dman.get(), tStep);
307  }
308  }
309 
310  for ( auto &dman : domain->giveDofManagers() ) {
311  dman->updateYourself(tStep);
312  }
313 
314 # ifdef VERBOSE
315  VERBOSE_PRINT0("Updated nodes & sides ", domain->giveNumberOfDofManagers())
316 # endif
317 
318 
319  for ( auto &elem : domain->giveElements() ) {
320  elem->updateInternalState(tStep);
321  elem->updateYourself(tStep);
322  }
323 
324 # ifdef VERBOSE
325  VERBOSE_PRINT0("Updated Elements ", domain->giveNumberOfElements())
326 # endif
327 
328 #if 0
329  // save context if required
330  // default - save only if ALWAYS is set ( see cltypes.h )
331 
332  if ((domain->giveContextOutputMode() == COM_Always ) ||
333  (domain->giveContextOutputMode() == COM_Required )) {
334  this->saveContext(NULL);
335  } else if (domain->giveContextOutputMode() == COM_UserDefined ) {
336  if (tStep->giveNumber()%domain->giveContextOutputStep() == 0)
337  this->saveContext(NULL);
338  }
339 #endif
340 }
341 
342 
344 {
345  if ( !suppressOutput ) {
346  this->printOutputAt(this->giveOutputStream(), tStep);
347  fflush( this->giveOutputStream() );
348  }
349 
350  Domain *domain = this->giveDomain(1);
351  // i = 0 represents the linear solution, which is followed by the eigen vectors starting at i = 1
352  for ( int i = 0; i <= numberOfRequiredEigenValues; i++ ) {
353  tStep->setTime( ( double ) i );
354 
355  if ( this->requiresUnknownsDictionaryUpdate() ) {
356  for ( auto &dman : domain->giveDofManagers() ) {
357  this->updateDofUnknownsDictionary(dman.get(), tStep);
358  }
359  }
360 
361  for ( auto &dman : domain->giveDofManagers() ) {
362  dman->updateYourself(tStep);
363  }
364 
365  tStep->setNumber(i);
367  }
368 }
369 
370 
372 {
373  Domain *domain = this->giveDomain(1);
374  if ( !domain->giveOutputManager()->testTimeStepOutput(tStep) ) {
375  return;
376  }
377 
378  fprintf(file, "\nLinear Stability:");
379  fprintf(file, "\nEigen Values are:\n-----------------\n");
380 
381  for ( int i = 1; i <= numberOfRequiredEigenValues; i++ ) {
382  fprintf(file, "%15.8e ", eigVal.at(i) );
383  if ( ( i % 5 ) == 0 ) {
384  fprintf(file, "\n");
385  }
386  }
387 
388  fprintf(file, "\n\n");
389 
390  for ( int i = 0; i <= numberOfRequiredEigenValues; i++ ) {
391  if ( i == 0 ) {
392  fprintf(file, "\nLinear solution\n\n");
393  } else {
394  fprintf(file, "\nEigen vector no. %d, correposnding eigen value is %15.8e\n\n", i, eigVal.at(i));
395  }
396  tStep->setTime( ( double ) i );
397 
398  if ( this->requiresUnknownsDictionaryUpdate() ) {
399  for ( auto &dman : domain->giveDofManagers() ) {
400  this->updateDofUnknownsDictionary(dman.get(), tStep);
401  }
402  }
403 
404  for ( auto &dman : domain->giveDofManagers() ) {
405  dman->updateYourself(tStep);
406  dman->printOutputAt(file, tStep);
407  }
408 
409  tStep->setNumber(i);
410 
411  if ( i == 0 ) {
412  for ( auto &elem : domain->giveElements() ) {
413  elem->printOutputAt(file, tStep);
414  }
415  this->printReactionForces(tStep, 1., file);
416  }
417  }
418 }
419 
420 
422 {
423  this->giveCurrentStep()->setTime( ( double ) activeVector );
424 }
425 
426 
428 {
429  contextIOResultType iores;
430 
431  if ( ( iores = StructuralEngngModel :: saveContext(stream, mode) ) != CIO_OK ) {
432  THROW_CIOERR(iores);
433  }
434 
435  if ( ( iores = displacementVector.storeYourself(stream) ) != CIO_OK ) {
436  THROW_CIOERR(iores);
437  }
438 
439  if ( ( iores = eigVal.storeYourself(stream) ) != CIO_OK ) {
440  THROW_CIOERR(iores);
441  }
442 
443  if ( ( iores = eigVec.storeYourself(stream) ) != CIO_OK ) {
444  THROW_CIOERR(iores);
445  }
446 
447  return CIO_OK;
448 }
449 
450 
452 {
453  contextIOResultType iores;
454 
455  if ( ( iores = StructuralEngngModel :: restoreContext(stream, mode) ) != CIO_OK ) {
456  THROW_CIOERR(iores);
457  }
458 
459  if ( ( iores = displacementVector.restoreYourself(stream) ) != CIO_OK ) {
460  THROW_CIOERR(iores);
461  }
462 
463  if ( ( iores = eigVal.restoreYourself(stream) ) != CIO_OK ) {
464  THROW_CIOERR(iores);
465  }
466 
467  if ( ( iores = eigVec.restoreYourself(stream) ) != CIO_OK ) {
468  THROW_CIOERR(iores);
469  }
470 
471  return CIO_OK;
472 }
473 
474 } // end namespace oofem
EngngModelTimer timer
E-model timer.
Definition: engngm.h:267
GenEigvalSolverType
Types of general eigenvalue solvers.
virtual void updateAttributes(MetaStep *mStep)
Update receiver attributes according to step metaStep attributes.
Definition: engngm.C:589
The representation of EngngModel default unknown numbering.
FILE * outputStream
Output stream.
Definition: engngm.h:242
std::unique_ptr< TimeStep > currentStep
Current time step.
Definition: engngm.h:231
std::string simulationDescription
Definition: engngm.h:316
Class and object Domain.
Definition: domain.h:115
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
std::unique_ptr< SparseGeneralEigenValueSystemNM > nMethod
std::string dataOutputFileName
Path to output stream.
Definition: engngm.h:238
int giveNumberOfDofManagers() const
Returns number of dof managers in domain.
Definition: domain.h:432
Class representing meta step.
Definition: metastep.h:62
virtual void updateYourself(TimeStep *tStep)
Updates internal state after finishing time step.
virtual void doStepOutput(TimeStep *tStep)
Prints the ouput of the solution step (using virtual this->printOutputAtservice) to the stream detemi...
SparseMtrx * createSparseMtrx(SparseMtrxType type)
Creates new instance of sparse matrix corresponding to given keyword.
Definition: classfactory.C:105
void startTimer(EngngModelTimerType t)
Definition: timer.h:128
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
contextIOResultType storeYourself(DataStream &stream) const
Definition: floatarray.C:872
std::unique_ptr< TimeStep > previousStep
Previous time step.
Definition: engngm.h:233
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
MetaStep * giveMetaStep(int i)
Returns the i-th meta step.
Definition: engngm.C:1765
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
long StateCounterType
StateCounterType type used to indicate solution state.
This base class is an abstraction for all numerical methods solving sparse linear system of equations...
virtual void setActiveVector(int i)
Only relevant for eigen value analysis. Otherwise does noting.
This base class is an abstraction for numerical algorithm.
Definition: nummet.h:80
Implementation for assembling tangent matrices in standard monolithic FE-problems.
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
time_t startTime
Solution start time.
Definition: engngm.h:259
double giveTargetTime()
Returns target time.
Definition: timestep.h:146
void setTime(double newt)
Sets target and intrinsic time to be equal.
Definition: timestep.h:154
void doOutput(TimeStep *tStep, bool substepFlag=false)
Writes the output.
std::unique_ptr< SparseLinearSystemNM > nMethodLS
Numerical method used to solve the static problem.
std::vector< std::unique_ptr< DofManager > > & giveDofManagers()
Definition: domain.h:400
virtual int __giveEquationNumber() const =0
Returns equation number of receiver, usually assigned by emodel.
int giveNumberOfElements() const
Returns number of elements in domain.
Definition: domain.h:434
REGISTER_EngngModel(ProblemSequence)
virtual NumericalMethod * giveNumericalMethod(MetaStep *mStep)
Returns reference to receiver&#39;s numerical method.
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
void printReactionForces(TimeStep *tStep, int id, FILE *out)
Computes and prints reaction forces, computed from nodal internal forces.
contextIOResultType storeYourself(DataStream &stream) const
Definition: floatmatrix.C:1852
If required (for backtracking computation).
#define VERBOSE_PRINT0(str, number)
Definition: verbose.h:56
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain)
Assembles characteristic matrix of required type into given sparse matrix.
Definition: engngm.C:803
MetaStep * giveCurrentMetaStep()
Returns current meta step.
Definition: engngm.C:1684
#define OOFEM_LOG_INFO(...)
Definition: logger.h:127
int giveNumber()
Returns receiver&#39;s number.
Definition: timestep.h:129
FILE * giveOutputStream()
Returns file descriptor of output file.
Definition: engngm.C:1791
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode)
Restores the state of model from output stream.
#define OOFEM_ERROR(...)
Definition: error.h:61
int numberOfSteps
Total number of time steps.
Definition: engngm.h:209
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.
GenEigvalSolverType solverType
Numerical method used to solve the problem.
#define _IFT_LinearStability_rtolv
std::unique_ptr< SparseMtrx > stiffnessMatrix
SparseLinearSystemNM * giveNumericalMethodForLinStaticProblem(TimeStep *tStep)
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatarray.C:895
#define _IFT_EngngModel_suppressOutput
Definition: engngm.h:84
int testTimeStepOutput(TimeStep *)
Tests if given time step output is required.
int nMetaSteps
Number of meta steps.
Definition: engngm.h:225
virtual int giveNumberOfFirstStep(bool force=false)
Returns number of first time step used by receiver.
Definition: engngm.h:730
Enable for post-processing.
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
ExportModuleManager * exportModuleManager
Export module manager.
Definition: engngm.h:250
#define _IFT_LinearStability_stype
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description in input reader.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode)
Stores the state of model to output stream.
Definition: engngm.C:1527
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatmatrix.C:1877
Initializes the variable VERBOSE, in order to get a few intermediate messages on screen: beginning an...
virtual double giveUnknownComponent(ValueModeType type, TimeStep *tStep, Domain *d, Dof *dof)
Returns requested unknown.
SparseLinearSystemNM * createSparseLinSolver(LinSystSolverType st, Domain *d, EngngModel *m)
Creates new instance of SparseLinearSystemNM corresponding to given type.
Definition: classfactory.C:120
Input attribute of domain (each n-th step).
virtual TimeStep * giveNextStep()
Returns next time step (next to current step) of receiver.
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual void updateDofUnknownsDictionary(DofManager *, TimeStep *)
Updates necessary values in Dofs unknown dictionaries.
Definition: engngm.h:859
int giveMetaStepNumber()
Returns receiver&#39;s meta step number.
Definition: timestep.h:135
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
Class representing the general Input Record.
Definition: inputrecord.h:101
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
ClassFactory & classFactory
Definition: classfactory.C:59
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
OutputManager * giveOutputManager()
Returns domain output manager.
Definition: domain.C:1436
bool suppressOutput
Flag for suppressing output to file.
Definition: engngm.h:314
Callback class for assembling initial stress matrices.
std::vector< std::unique_ptr< Element > > & giveElements()
Definition: domain.h:279
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode)
Stores the state of model to output stream.
SparseGeneralEigenValueSystemNM * createGeneralizedEigenValueSolver(GenEigvalSolverType name, Domain *d, EngngModel *m)
Definition: classfactory.C:411
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to output domain stream, for given time step.
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
virtual TimeStep * giveCurrentStep(bool force=false)
Returns current time step.
Definition: engngm.h:683
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode)
Restores the state of model from output stream.
Definition: engngm.C:1592
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
void terminateLinStatic(TimeStep *tStep)
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
Abstract class Dof represents Degree Of Freedom in finite element mesh.
Definition: dof.h:93
Symmetric skyline.
void negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
virtual void solveYourselfAt(TimeStep *tStep)
Solves problem for given time step.
virtual void solveYourself()
Starts solution process.
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
void setNumber(int i)
Set receiver&#39;s number.
Definition: timestep.h:131
#define _IFT_LinearStability_nroot
std::unique_ptr< SparseMtrx > initialStressMatrix
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

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:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011