OOFEM 3.0
Loading...
Searching...
No Matches
diidynamic.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 - 2025 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
38#include "timestep.h"
39#include "element.h"
40#include "dofmanager.h"
41#include "dof.h"
42#include "contextioerr.h"
43#include "datastream.h"
44#include "verbose.h"
46#include "classfactory.h"
47
48namespace oofem {
50
60
61DIIDynamic :: ~DIIDynamic()
62{
63}
64
65NumericalMethod *DIIDynamic :: giveNumericalMethod(MetaStep *mStep)
66// Only one has reason for DIIDynamic
67// - SolutionOfLinearEquations
68{
69 if ( nMethod ) {
70 return nMethod.get();
71 }
72
73 nMethod = classFactory.createSparseLinSolver(solverType, this->giveDomain(1), this);
74 if ( !nMethod ) {
75 OOFEM_ERROR("linear solver creation failed for lstype %d", solverType);
76 }
77
78 return nMethod.get();
79}
80
81
82void
83DIIDynamic :: initializeFrom(InputRecord &ir)
84{
85 StructuralEngngModel :: initializeFrom(ir);
86
87 int val = 0;
90
91 val = 0;
94
95 val = 0;
97
98 eta = 0;
100
101 delta = 0;
103
105
106 gamma = 0.5;
107 beta = 0.25; // Default Newmark parameters.
108 theta = 1.37; // Default Wilson parameter.
110 OOFEM_LOG_INFO("Selecting Newmark-beta metod\n");
114 OOFEM_LOG_INFO("Selecting Two-point Backward Euler method\n");
116 OOFEM_LOG_INFO("Selecting Three-point Backward Euler metod\n");
117 } else if ( initialTimeDiscretization == TD_Wilson ) {
118 OOFEM_LOG_INFO("Selecting Wilson-theta metod\n");
120 if ( theta < 1.37 ) {
121 OOFEM_WARNING("Found theta < 1.37. Performing correction, theta = 1.37");
122 theta = 1.37;
123 }
124 } else {
125 throw ValueInputException(ir, _IFT_DIIDynamic_ddtScheme, "Time-stepping scheme not found!");
126 }
127
129}
130
131
132double DIIDynamic :: giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof)
133// Returns unknown quantity displacement, velocity or acceleration of equation eq.
134// This function translates this request to numerical method language.
135{
136 int eq = dof->__giveEquationNumber();
137#ifdef DEBUG
138 if ( eq == 0 ) {
139 OOFEM_ERROR("invalid equation number");
140 }
141#endif
142
143 if ( tStep != this->giveCurrentStep() ) {
144 OOFEM_ERROR("unknown time step encountered");
145 }
146
147 switch ( mode ) {
148 case VM_Total:
149 return displacementVector.at(eq);
150
151 case VM_Velocity:
152 return velocityVector.at(eq);
153
154 case VM_Acceleration:
155 return accelerationVector.at(eq);
156
157 default:
158 OOFEM_ERROR("Unknown is of undefined ValueModeType for this problem");
159 }
160
161 // return 0.0;
162}
163
164
165TimeStep *DIIDynamic :: giveNextStep()
166{
167 int istep = giveNumberOfFirstStep();
168 double totalTime = deltaT;
169 StateCounterType counter = 1;
171
172 if ( currentStep ) {
173 totalTime = currentStep->giveTargetTime() + deltaT;
174 istep = currentStep->giveNumber() + 1;
175 counter = currentStep->giveSolutionStateCounter() + 1;
176 td = currentStep->giveTimeDiscretization();
177 if ( currentStep->isTheFirstStep() &&
180 }
181 }
182
183 previousStep = std :: move(currentStep);
184
185 currentStep = std::make_unique<TimeStep>(istep, this, 1, totalTime, deltaT, counter, td);
186
187 return currentStep.get();
188}
189
190void DIIDynamic :: initializeYourself(TimeStep *tStep)
191{
192 Domain *domain = this->giveDomain(1);
194
195 if ( tStep->isTheFirstStep() ) {
197 -deltaT, deltaT, 0);
198 //
199 // Determine initial displacements, velocities and accelerations.
200 //
201 loadVector.resize(neq);
202 loadVector.zero();
203 displacementVector.resize(neq);
204 displacementVector.zero();
205 velocityVector.resize(neq);
206 velocityVector.zero();
207 accelerationVector.resize(neq);
208 accelerationVector.zero();
211
212 for ( auto &node : domain->giveDofManagers()) {
213 for ( Dof *iDof: *node ) {
214 //
215 // Ask for initial values obtained from boundary conditions and initial conditions.
216 //
217 if ( !iDof->isPrimaryDof() ) {
218 continue;
219 }
220
221 int jj = iDof->__giveEquationNumber();
222 if ( jj ) {
223 displacementVector.at(jj) = iDof->giveUnknown(VM_Total, stepWhenIcApply);
224 velocityVector.at(jj) = iDof->giveUnknown(VM_Velocity, stepWhenIcApply);
225 accelerationVector.at(jj) = iDof->giveUnknown(VM_Acceleration, stepWhenIcApply);
226 }
227 }
228 }
229
230 delete stepWhenIcApply;
231 } // End of initialization.
232}
233
234void DIIDynamic :: solveYourself()
235{
236 StructuralEngngModel :: solveYourself();
237}
238
239void DIIDynamic :: solveYourselfAt(TimeStep *tStep)
240{
241 // Determine the constants.
242 this->determineConstants(tStep);
243
244 Domain *domain = this->giveDomain(1);
246
247
248 if ( initFlag ) {
249#ifdef VERBOSE
250 OOFEM_LOG_DEBUG("Assembling stiffness matrix\n");
251#endif
252 stiffnessMatrix = classFactory.createSparseMtrx(sparseMtrxType);
253 if ( !stiffnessMatrix ) {
254 OOFEM_ERROR("sparse matrix creation failed");
255 }
256
257 stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
258
259 this->assemble(*stiffnessMatrix, tStep, EffectiveTangentAssembler(TangentStiffness, false, 1 + this->delta * a1, this->a0 + this->eta * this->a1),
261
262 help.resize(neq);
263 help.zero();
264
269
270 initFlag = 0;
271 }
272
273 if ( ( tStep->givePreviousStep() != NULL ) && ( tStep->giveTimeDiscretization() != tStep->givePreviousStep()->giveTimeDiscretization() ) ) {
274#ifdef VERBOSE
275 OOFEM_LOG_DEBUG("Assembling stiffness matrix\n");
276#endif
277 stiffnessMatrix->zero();
278 this->assemble(*stiffnessMatrix, tStep, EffectiveTangentAssembler(TangentStiffness, false, 1 + this->delta * a1, this->a0 + this->eta * this->a1),
280 }
281
282#ifdef VERBOSE
283 OOFEM_LOG_DEBUG("Assembling load\n");
284#endif
285
286 this->assembleLoadVector(loadVector, domain, VM_Total, tStep);
287
288 //
289 // contribution of internal, viscous and inertia forces due to prescribed (dirichlet) BCs
290 //
291 {
292 FloatArray bcLoadVector;
293 this->assembleDirichletBcRhsVector(bcLoadVector, domain, tStep);
294 loadVector.subtract(bcLoadVector);
295 }
296
297
298 //
299 // Assemble effective load vector (rhs).
300 //
301
302 for ( int i = 1; i <= neq; i++ ) {
303 help.at(i) = a0 * previousDisplacementVector.at(i)
306 + eta * ( a1 * previousDisplacementVector.at(i)
307 + a4 * previousVelocityVector.at(i)
310 }
311
312 this->timesMtrx(help, rhs, MassMatrix, domain, tStep);
313 //this->assembleVector(help, tStep, MatrixProductAssembler(MassMatrix(), rhs), VM_Total,
314 // EModelDefaultEquationNumbering(), this->giveDomain(1));
315
316 help.zero();
317
318 if ( delta != 0 ) {
319 for ( int i = 1; i <= neq; i++ ) {
320 help.at(i) = delta * ( a1 * previousDisplacementVector.at(i)
321 + a4 * previousVelocityVector.at(i)
324 }
325 this->timesMtrx(help, rhs2, TangentStiffnessMatrix, domain, tStep);
326 //this->assembleVector(help, tStep, MatrixProductAssembler(TangentAssembler(), rhs2), VM_Total,
327 // EModelDefaultEquationNumbering(), this->giveDomain(1));
328
329 help.zero();
330 for ( int i = 1; i <= neq; i++ ) {
331 rhs.at(i) += rhs2.at(i);
332 }
333 }
334
335 if ( tStep->giveTimeDiscretization() == TD_Wilson ) {
336 for ( int i = 1; i <= neq; i++ ) {
337 rhs.at(i) += previousLoadVector.at(i)
338 + theta * ( loadVector.at(i) - previousLoadVector.at(i) );
339 }
340 } else {
341 for ( int i = 1; i <= neq; i++ ) {
342 rhs.at(i) += loadVector.at(i);
343 }
344 }
345
346 //
347 // Call numerical model to solve arised problem.
348 //
349#ifdef VERBOSE
350 OOFEM_LOG_RELEVANT( "\n\nSolving [step number %8d, time %15e]\n", tStep->giveNumber(), tStep->giveTargetTime() );
351#endif
352
354
355 if ( tStep->giveTimeDiscretization() == TD_Wilson ) {
356 OOFEM_LOG_INFO("TD_Wilson: Updating acceleration, velocity and displacement.\n");
357 for ( int i = 1; i <= neq; i++ ) {
361
364
367 + a10 * ( accelerationVector.at(i) + 2 * previousAccelerationVector.at(i) );
368 }
369 } else if ( ( tStep->giveTimeDiscretization() == TD_ThreePointBackward ) ||
371 for ( int i = 1; i <= neq; i++ ) {
373 - a2 *previousVelocityVector.at(i);
374
377 }
378 } else {
379 for ( int i = 1; i <= neq; i++ ) {
383
386 + a7 *accelerationVector.at(i);
387 }
388 }
389}
390
391
392void DIIDynamic :: updateYourself(TimeStep *tStep)
393{
399
400 StructuralEngngModel :: updateYourself(tStep);
401}
402
403
404void
405DIIDynamic :: printDofOutputAt(FILE *stream, Dof *iDof, TimeStep *tStep)
406{
407 static char dofchar[] = "dva";
408 static ValueModeType dofmodes[] = {
409 VM_Total, VM_Velocity, VM_Acceleration
410 };
411
412 iDof->printMultipleOutputAt(stream, tStep, dofchar, dofmodes, 3);
413}
414
415
416void
417DIIDynamic :: timesMtrx(FloatArray &vec, FloatArray &answer, CharType type, Domain *domain, TimeStep *tStep)
418{
419 int nelem = domain->giveNumberOfElements();
421 int i, j, k, jj, kk, n;
422 FloatMatrix charMtrx, R;
423 IntArray loc;
424 Element *element;
426
427 answer.resize(neq);
428 answer.zero();
429
430 for ( i = 1; i <= nelem; i++ ) {
431 element = domain->giveElement(i);
432 // Skip remote elements.
433 if ( element->giveParallelMode() == Element_remote ) {
434 continue;
435 }
436
437 element->giveLocationArray(loc, en);
438 element->giveCharacteristicMatrix(charMtrx, type, tStep);
439 if ( charMtrx.isNotEmpty() ) {
441 if ( element->giveRotationMatrix(R) ) {
442 charMtrx.rotatedWith(R);
443 }
444
445
446#ifdef DEBUG
447 if ( loc.giveSize() != charMtrx.giveNumberOfRows() ) {
448 OOFEM_ERROR("dimension mismatch");
449 }
450
451#endif
452
453 n = loc.giveSize();
454 for ( j = 1; j <= n; j++ ) {
455 jj = loc.at(j);
456 if ( jj ) {
457 for ( k = 1; k <= n; k++ ) {
458 kk = loc.at(k);
459 if ( kk ) {
460 answer.at(jj) += charMtrx.at(j, k) * vec.at(kk);
461 }
462 }
463 }
464 }
465 }
466 }
467}
468
469
470void
471DIIDynamic :: assembleLoadVector(FloatArray &_loadVector, Domain *domain, ValueModeType mode, TimeStep *tStep)
472{
474 _loadVector.zero();
475
476 this->assembleVector(_loadVector, tStep, ExternalForceAssembler(), mode,
479}
480
481
482void
483DIIDynamic :: assembleDirichletBcRhsVector (FloatArray& answer, Domain* d, TimeStep *tStep)
484{
485 IntArray loc, dofids;
486 int nelem = d->giveNumberOfElements();
487 FloatArray rp, charVec;
488 FloatMatrix k,m, R;
489 FloatMatrix capacity;
490
492 answer.zero();
493
494 for ( int ielem = 1; ielem <= nelem; ielem++ ) {
495 Element *element = d->giveElement(ielem);
496
497 element->giveElementDofIDMask(dofids);
499 element->computeVectorOfPrescribed(dofids, VM_Total, tStep, rp);
500 int rotFlag = element->giveRotationMatrix(R);
501
502 if ( rp.containsOnlyZeroes() ) {
503 continue;
504 } else {
505 element->giveCharacteristicMatrix(k, TangentStiffnessMatrix, tStep);
506 charVec.beProductOf(k, rp);
507 if (rotFlag) charVec.rotatedWith(R, 't');
508 answer.assemble(charVec, loc);
509 }
510
511 // add acceleration forces due to BCs
512 element->computeVectorOfPrescribed(dofids, VM_Acceleration, tStep, rp);
513 if ( rp.containsOnlyZeroes() ) {
514 continue;
515 } else {
516 element->giveCharacteristicMatrix(m, MassMatrix, tStep);
517 charVec.beProductOf(m, rp);
518 if (rotFlag) charVec.rotatedWith(R, 't');
519 answer.assemble(charVec, loc);
520 }
521
522 // add dumping forces due to BCs
523 element->computeVectorOfPrescribed(dofids, VM_Velocity, tStep, rp);
524 if ( rp.containsOnlyZeroes() ) {
525 continue;
526 } else {
527 if (!k.isNotEmpty()) element->giveCharacteristicMatrix(k, TangentStiffnessMatrix, tStep);
528 if (!m.isNotEmpty()) element->giveCharacteristicMatrix(m, MassMatrix, tStep);
529 k.times(this->delta);
530 m.times(this->eta);
531 m.add(k);
532 charVec.beProductOf(m, rp);
533 if (rotFlag) charVec.rotatedWith(R, 't');
534 answer.assemble(charVec, loc);
535 }
536
537
538
539
540 } // end element loop
541
542}
543
544
545
546void
547DIIDynamic :: determineConstants(TimeStep *tStep)
548{
549 if ( this->giveCurrentStep()->isTheFirstStep() && initialTimeDiscretization == TD_ThreePointBackward ) {
550 this->giveCurrentStep()->setTimeDiscretization( TD_TwoPointBackward );
551 }
552
553 deltaT = tStep->giveTimeIncrement();
554
555 if ( tStep->giveTimeDiscretization() == TD_Newmark ) {
556 double dt2 = deltaT * deltaT;
557 a0 = 1 / ( beta * dt2 );
558 a1 = gamma / ( beta * deltaT );
559 a2 = 1 / ( beta * deltaT );
560 a3 = 1 / ( 2 * beta ) - 1;
561 a4 = ( gamma / beta ) - 1;
562 a5 = deltaT / 2 * ( gamma / beta - 2 );
563 a6 = deltaT * ( 1 - gamma );
564 a7 = deltaT * gamma;
565 a8 = dt2 * ( ( 1 / 2 ) - beta );
566 a9 = dt2 * beta;
567 a10 = 0;
568 a11 = 0;
569 } else if ( tStep->giveTimeDiscretization() == TD_TwoPointBackward ) {
570 double dt2 = deltaT * deltaT;
571 a0 = 1 / dt2;
572 a1 = 1 / deltaT;
573 a2 = 1 / deltaT;
574 a3 = 0;
575 a4 = 0;
576 a5 = 0;
577 a6 = 0;
578 a7 = 0;
579 a8 = 0;
580 a9 = 0;
581 a10 = 0;
582 a11 = 0;
583 } else if ( tStep->giveTimeDiscretization() == TD_ThreePointBackward ) {
584 double dt2 = deltaT * deltaT;
585 a0 = 2 / dt2;
586 a1 = 3 / ( 2 * deltaT );
587 a2 = 2 / deltaT;
588 a3 = 0;
589 a4 = 0;
590 a5 = 0;
591 a6 = 0;
592 a7 = 0;
593 a8 = 0;
594 a9 = 0;
595 a10 = 0;
596 a11 = 1 / ( 2 * deltaT );
597 } else if ( tStep->giveTimeDiscretization() == TD_Wilson ) {
598 a0 = 6 / ( ( theta * deltaT ) * ( theta * deltaT ) );
599 a1 = 3 / ( theta * deltaT );
600 a2 = 2 * a1;
601 a3 = 2;
602 a4 = 2;
603 a5 = theta * deltaT / 2;
604 a6 = a0 / theta;
605 a7 = a2 / theta;
606 a8 = 1 - 3 / theta;
607 a9 = deltaT / 2;
608 a10 = deltaT * deltaT / 6;
609 a11 = 0;
610 } else {
611 OOFEM_ERROR("Time-stepping scheme not found!");
612 }
613}
614
615
616void DIIDynamic :: saveContext(DataStream &stream, ContextMode mode)
617{
619
620 EngngModel :: saveContext(stream, mode);
621
622 if ( ( iores = displacementVector.storeYourself(stream) ) != CIO_OK ) {
623 THROW_CIOERR(iores);
624 }
625
626 if ( ( iores = velocityVector.storeYourself(stream) ) != CIO_OK ) {
627 THROW_CIOERR(iores);
628 }
629
630 if ( ( iores = accelerationVector.storeYourself(stream) ) != CIO_OK ) {
631 THROW_CIOERR(iores);
632 }
633
634 if ( ( iores = loadVector.storeYourself(stream) ) != CIO_OK ) {
635 THROW_CIOERR(iores);
636 }
637
638 if ( ( iores = previousIncrementOfDisplacement.storeYourself(stream) ) != CIO_OK ) {
639 THROW_CIOERR(iores);
640 }
641}
642
643void DIIDynamic :: restoreContext(DataStream &stream, ContextMode mode)
644{
646
647 EngngModel :: restoreContext(stream, mode);
648
649 if ( ( iores = displacementVector.restoreYourself(stream) ) != CIO_OK ) {
650 THROW_CIOERR(iores);
651 }
652
653 if ( ( iores = velocityVector.restoreYourself(stream) ) != CIO_OK ) {
654 THROW_CIOERR(iores);
655 }
656
657 if ( ( iores = accelerationVector.restoreYourself(stream) ) != CIO_OK ) {
658 THROW_CIOERR(iores);
659 }
660
661 if ( ( iores = loadVector.restoreYourself(stream) ) != CIO_OK ) {
662 THROW_CIOERR(iores);
663 }
664
665 if ( ( iores = previousIncrementOfDisplacement.restoreYourself(stream) ) != CIO_OK ) {
666 THROW_CIOERR(iores);
667 }
668}
669} // end namespace oofem
#define REGISTER_EngngModel(class)
void timesMtrx(FloatArray &answer, FloatArray &vec, CharType type, Domain *domain, TimeStep *tStep)
Definition diidynamic.C:417
FloatArray previousDisplacementVector
Definition diidynamic.h:80
void determineConstants(TimeStep *tStep)
Definition diidynamic.C:547
FloatArray rhs
Definition diidynamic.h:78
void assembleDirichletBcRhsVector(FloatArray &answer, Domain *d, TimeStep *tStep)
Definition diidynamic.C:483
std ::unique_ptr< SparseLinearSystemNM > nMethod
Numerical method used to solve the problem.
Definition diidynamic.h:95
FloatArray previousVelocityVector
Definition diidynamic.h:80
LinSystSolverType solverType
Definition diidynamic.h:88
void assembleLoadVector(FloatArray &_loadVector, Domain *domain, ValueModeType mode, TimeStep *tStep)
Definition diidynamic.C:471
TimeDiscretizationType initialTimeDiscretization
Definition diidynamic.h:92
FloatArray velocityVector
Definition diidynamic.h:79
FloatArray accelerationVector
Definition diidynamic.h:79
FloatArray previousAccelerationVector
Definition diidynamic.h:80
FloatArray previousIncrementOfDisplacement
Definition diidynamic.h:81
FloatArray previousLoadVector
Definition diidynamic.h:78
FloatArray displacementVector
Definition diidynamic.h:79
FloatArray help
Definition diidynamic.h:82
FloatArray loadVector
Definition diidynamic.h:78
FloatArray rhs2
Definition diidynamic.h:78
std ::unique_ptr< SparseMtrx > stiffnessMatrix
Definition diidynamic.h:77
SparseMtrxType sparseMtrxType
Definition diidynamic.h:89
virtual void printMultipleOutputAt(FILE *File, TimeStep *tStep, char *ch, ValueModeType *mode, int nite)
Definition dof.C:92
virtual int __giveEquationNumber() const =0
int giveNumber()
Returns domain number.
Definition domain.h:281
int giveNumberOfElements() const
Returns number of elements in domain.
Definition domain.h:463
std ::vector< std ::unique_ptr< DofManager > > & giveDofManagers()
Definition domain.h:427
Element * giveElement(int n)
Definition domain.C:165
virtual void giveElementDofIDMask(IntArray &answer) const
Definition element.h:510
virtual bool giveRotationMatrix(FloatMatrix &answer)
Definition element.C:298
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Definition element.C:620
void computeVectorOfPrescribed(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:212
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
Definition element.C:429
elementParallelMode giveParallelMode() const
Definition element.h:1139
int giveNumberOfTimeStepWhenIcApply()
Returns the time step number, when initial conditions should apply.
Definition engngm.h:787
virtual TimeStep * giveCurrentStep(bool force=false)
Definition engngm.h:717
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Definition engngm.C:452
std ::unique_ptr< TimeStep > previousStep
Previous time step.
Definition engngm.h:243
int ndomains
Number of receiver domains.
Definition engngm.h:215
Domain * giveDomain(int n)
Definition engngm.C:1936
std ::unique_ptr< TimeStep > currentStep
Current time step.
Definition engngm.h:241
virtual int giveNumberOfFirstStep(bool force=false)
Definition engngm.h:763
int updateSharedDofManagers(FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
Definition engngm.C:2162
std ::unique_ptr< TimeStep > stepWhenIcApply
Solution step when IC (initial conditions) apply.
Definition engngm.h:239
void assembleVector(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
Definition engngm.C:1106
bool containsOnlyZeroes() const
Definition floatarray.C:671
void assemble(const FloatArray &fe, const IntArray &loc)
Definition floatarray.C:616
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void rotatedWith(FloatMatrix &r, char mode)
Definition floatarray.C:814
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void times(double f)
void rotatedWith(const FloatMatrix &r, char mode='n')
void add(const FloatMatrix &a)
bool isNotEmpty() const
Tests for empty matrix.
int giveNumberOfRows() const
Returns number of rows of receiver.
double at(std::size_t i, std::size_t j) const
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
StructuralEngngModel(int i, EngngModel *master=nullptr)
Creates new StructuralEngngModel with number i, associated to domain d.
double giveTimeIncrement()
Returns solution step associated time increment.
Definition timestep.h:168
double giveTargetTime()
Returns target time.
Definition timestep.h:164
TimeDiscretizationType giveTimeDiscretization()
Returns time discretization.
Definition timestep.h:219
int giveNumber()
Returns receiver's number.
Definition timestep.h:144
TimeStep * givePreviousStep()
Returns pointer to previous solution step.
Definition timestep.C:132
bool isTheFirstStep()
Definition timestep.C:148
#define THROW_CIOERR(e)
#define _IFT_DIIDynamic_eta
Definition diidynamic.h:50
#define _IFT_DIIDynamic_theta
Definition diidynamic.h:52
#define _IFT_DIIDynamic_gamma
Definition diidynamic.h:48
#define _IFT_DIIDynamic_deltat
Definition diidynamic.h:46
#define _IFT_DIIDynamic_delta
Definition diidynamic.h:51
#define _IFT_DIIDynamic_ddtScheme
Definition diidynamic.h:47
#define _IFT_DIIDynamic_beta
Definition diidynamic.h:49
#define _IFT_EngngModel_lstype
Definition engngm.h:91
#define _IFT_EngngModel_smtype
Definition engngm.h:92
#define OOFEM_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
#define OOFEM_LOG_INFO(...)
Definition logger.h:143
#define OOFEM_LOG_RELEVANT(...)
Definition logger.h:142
#define OOFEM_LOG_DEBUG(...)
Definition logger.h:144
long ContextMode
Definition contextmode.h:43
@ Element_remote
Element in active domain is only mirror of some remote element.
Definition element.h:89
FloatArrayF< N > assemble(const FloatArrayF< M > &x, int const (&c)[M])
Assemble components into zero matrix.
long StateCounterType
StateCounterType type used to indicate solution state.
ClassFactory & classFactory
TimeDiscretizationType
Time discretization used by transient solvers.
@ TD_TwoPointBackward
Two-point Backward Euler method.
@ TD_Newmark
Newmark-beta method.
@ TD_ThreePointBackward
Three-point Backward Euler method.
@ TD_Wilson
Wilson-theta method.

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011