OOFEM 3.0
Loading...
Searching...
No Matches
nlineardynamic.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
35#include <iostream>
36using namespace std;
37
41#include "nummet.h"
42#include "timestep.h"
43#include "metastep.h"
44#include "element.h"
45#include "dof.h"
46#include "error.h"
47#include "verbose.h"
49#include "nrsolver.h"
50#include "outputmanager.h"
51#include "datastream.h"
52#include "classfactory.h"
53#include "timer.h"
54#include "contextioerr.h"
55#include "sparsemtrx.h"
56#include "errorestimator.h"
58#include "assemblercallback.h"
59
60#ifdef __MPI_PARALLEL_MODE
61 #include "loadbalancer.h"
62 #include "problemcomm.h"
63 #include "processcomm.h"
64#endif
65
66namespace oofem {
68
69NonLinearDynamic :: NonLinearDynamic(int i, EngngModel *_master) : StructuralEngngModel(i, _master),
71{
73
74 ndomains = 1;
77}
78
79
80NonLinearDynamic :: ~NonLinearDynamic()
81{
82}
83
84
85NumericalMethod *NonLinearDynamic :: giveNumericalMethod(MetaStep *mStep)
86{
87 if ( mStep == NULL ) {
88 OOFEM_ERROR("undefined meta step");
89 }
90
91 if ( this->nMethod ) {
92 this->nMethod->reinitialize();
93 } else {
94 this->nMethod = std::make_unique<NRSolver>(this->giveDomain(1), this);
95 }
96
97 return this->nMethod.get();
98}
99
100
101void
102NonLinearDynamic :: updateAttributes(MetaStep *mStep)
103{
104 auto &ir = mStep->giveAttributesRecord();
105
106 StructuralEngngModel :: updateAttributes(mStep);
107
108 deltaT = 1.0;
110 if ( deltaT < 0. ) {
111 OOFEM_ERROR("deltaT < 0");
112 }
113
114 eta = 0;
116
117 delta = 0;
119}
120
121
122void
123NonLinearDynamic :: initializeFrom(InputRecord &ir)
124{
125 StructuralEngngModel :: initializeFrom(ir);
126
127 int val = 0;
130
131 val = 0;
134
137
138 deltaT = 1.0;
140 if ( deltaT < 0. ) {
141 OOFEM_ERROR("deltaT < 0");
142 }
143
144 val = 0;
147
148 gamma = 0.5;
149 beta = 0.25; // Default Newmark parameters.
151 OOFEM_LOG_INFO("Selecting Newmark-beta metod\n");
155 OOFEM_LOG_INFO("Selecting Two-point Backward Euler method\n");
157 OOFEM_LOG_INFO("Selecting Three-point Backward Euler metod\n");
158 } else {
159 throw ValueInputException(ir, _IFT_NonLinearDynamic_ddtScheme, "Time-stepping scheme not found!");
160 }
161
162 MANRMSteps = 0;
164
165#ifdef __MPI_PARALLEL_MODE
166 if ( isParallel() ) {
168 communicator = new NodeCommunicator(this, commBuff, this->giveRank(),
169 this->giveNumberOfProcesses());
170
172 nonlocalExt = 1;
174 this->giveNumberOfProcesses());
175 }
176 }
177#endif
178}
179
180
181double NonLinearDynamic :: giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof)
182{
183 int eq = dof->__giveEquationNumber();
184#ifdef DEBUG
185 if ( eq == 0 ) {
186 OOFEM_ERROR("invalid equation number");
187 }
188#endif
189
190 if ( tStep != this->giveCurrentStep() ) {
191 OOFEM_ERROR("unknown time step encountered");
192 }
193
194 switch ( mode ) {
195 case VM_Incremental:
196 return incrementOfDisplacement.at(eq);
197
198 case VM_Total:
199 return totalDisplacement.at(eq);
200
201 case VM_Velocity:
202 return velocityVector.at(eq);
203
204 case VM_Acceleration:
205 return accelerationVector.at(eq);
206
207 default:
208 OOFEM_ERROR("Unknown is of undefined ValueModeType for this problem");
209 }
210}
211
212
213TimeStep *NonLinearDynamic :: giveNextStep()
214{
215 int istep = giveNumberOfFirstStep();
216 int mStepNum = 1;
217 StateCounterType counter = 1;
218 double deltaTtmp = deltaT;
219 double totalTime = deltaT;
221
222 //Do not increase deltaT on microproblem
223 if ( pScale == microScale ) {
224 deltaTtmp = 0.;
225 }
226
227 if ( currentStep ) {
228 totalTime = currentStep->giveTargetTime() + deltaTtmp;
229 istep = currentStep->giveNumber() + 1;
230 counter = currentStep->giveSolutionStateCounter() + 1;
231 mStepNum = currentStep->giveMetaStepNumber();
232 td = currentStep->giveTimeDiscretization();
233
234 if ( !this->giveMetaStep(mStepNum)->isStepValid(istep) ) {
235 mStepNum++;
236 if ( mStepNum > nMetaSteps ) {
237 OOFEM_ERROR("no next step available, mStepNum=%d > nMetaSteps=%d", mStepNum, nMetaSteps);
238 }
239 }
240 }
241
242 previousStep = std :: move(currentStep);
243 currentStep = std::make_unique<TimeStep>(istep, this, mStepNum, totalTime, deltaTtmp, counter, td);
244
245 return currentStep.get();
246}
247
248
249void NonLinearDynamic :: solveYourself()
250{
251 if ( commInitFlag ) {
252 #ifdef __VERBOSE_PARALLEL
253 // Force equation numbering before setting up comm maps.
255 OOFEM_LOG_INFO("[process rank %d] neq is %d\n", this->giveRank(), neq);
256 #endif
257 if ( isParallel() ) {
258 // Set up communication patterns.
259 this->initializeCommMaps();
260 }
261 commInitFlag = 0;
262 }
263
264 StructuralEngngModel :: solveYourself();
265}
266
267
268void
269NonLinearDynamic :: solveYourselfAt(TimeStep *tStep)
270{
271 if ( commInitFlag ) {
272 #ifdef __VERBOSE_PARALLEL
273 // Force equation numbering before setting up comm maps.
275 OOFEM_LOG_INFO("[process rank %d] neq is %d\n", this->giveRank(), neq);
276 #endif
277 if ( isParallel() ) {
278 // Set up communication patterns.
279 this->initializeCommMaps();
280 }
281 commInitFlag = 0;
282 }
283
284#ifdef VERBOSE
285 OOFEM_LOG_RELEVANT( "\n\nSolving [Step number %8d, Time %15e]\n\n", tStep->giveNumber(), tStep->giveTargetTime() );
286#endif
287
288 proceedStep(1, tStep);
289}
290
291void
292NonLinearDynamic :: initializeYourself(TimeStep *tStep)
293{
294 Domain *domain = this->giveDomain(1);
296
297 if ( tStep->isTheFirstStep() && initFlag ) {
298 // Initialization
299 incrementOfDisplacement.resize(neq);
301 totalDisplacement.resize(neq);
302 totalDisplacement.zero();
303 velocityVector.resize(neq);
304 velocityVector.zero();
305 accelerationVector.resize(neq);
306 accelerationVector.zero();
307 internalForces.resize(neq);
308 internalForces.zero();
311 previousTotalDisplacement.resize(neq);
313 previousVelocityVector.resize(neq);
315 previousAccelerationVector.resize(neq);
317 previousInternalForces.resize(neq);
319
321 -deltaT, deltaT, 0);
322
323 // Considering initial conditions.
324 for ( auto &node : domain->giveDofManagers() ) {
325 for ( Dof *dof: *node ) {
326 // Ask for initial values obtained from
327 // bc (boundary conditions) and ic (initial conditions).
328 if ( !dof->isPrimaryDof() ) {
329 continue;
330 }
331
332 int jj = dof->__giveEquationNumber();
333 if ( jj ) {
334 totalDisplacement.at(jj) = dof->giveUnknown(VM_Total, stepWhenIcApply);
335 velocityVector.at(jj) = dof->giveUnknown(VM_Velocity, stepWhenIcApply);
336 accelerationVector.at(jj) = dof->giveUnknown(VM_Acceleration, stepWhenIcApply);
337 }
338 }
339 }
340 this->updateInternalRHS(internalForces, tStep, domain, &this->internalForcesEBENorm);
341 }
342}
343
344void
345NonLinearDynamic :: proceedStep(int di, TimeStep *tStep)
346{
347 // creates system of governing eq's and solves them at given time step
348 // first assemble problem at current time step
349
351
352 // Time-stepping constants
353 this->determineConstants(tStep);
354
355
356
357 if ( initFlag ) {
358 // First assemble problem at current time step.
359 // Option to take into account initial conditions.
362 // create mass matrix as compcol storage, need only mass multiplication by vector, not linear system solution.
363 massMatrix = classFactory.createSparseMtrx(SMT_CompCol);
364 }
365
367 OOFEM_ERROR("sparse matrix creation failed");
368 }
369
370 if ( nonlocalStiffnessFlag ) {
371 if ( !effectiveStiffnessMatrix->isAsymmetric() ) {
372 OOFEM_ERROR("effectiveStiffnessMatrix does not support asymmetric storage");
373 }
374 }
375
376 effectiveStiffnessMatrix->buildInternalStructure( this, di, EModelDefaultEquationNumbering() );
377 massMatrix->buildInternalStructure( this, di, EModelDefaultEquationNumbering() );
378
379 // Assemble mass matrix
380 this->assemble( *massMatrix, tStep, MassMatrixAssembler(),
382
383 // Initialize vectors
384 help.resize(neq);
385 help.zero();
386 rhs.resize(neq);
387 rhs.zero();
388 rhs2.resize(neq);
389 rhs2.zero();
390
396
397 forcesVector.resize(neq);
398 forcesVector.zero();
399
400 totIterations = 0;
401 initFlag = 0;
402 }
403
404#ifdef VERBOSE
405 OOFEM_LOG_DEBUG("Assembling load\n");
406#endif
407
408 // Assemble the external forces
409 FloatArray loadVector;
411 loadVector.zero();
412 this->assembleVector( loadVector, tStep, ExternalForceAssembler(),
413 VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(di) );
415
416 // Assembling the effective load vector
417 for ( int i = 1; i <= neq; i++ ) {
418 //help.beScaled(a2, previousVelocityVector);
419 //help.add(a3, previousAccelerationVector);
420 //help.add(a4 * eta, previousVelocityVector);
421 //help.add(a5 * eta, previousAccelerationVector);
422 //help.add(a6 * eta, previousIncrementOfDisplacement);
424 + eta * ( a4 * previousVelocityVector.at(i)
427 }
428
429 massMatrix->times(help, rhs);
430
431 if ( delta != 0 ) {
432 //help.beScaled(a4 * delta, previousVelocityVector);
433 //help.add(a5 * delta, previousAccelerationVector);
434 //help.add(a6 * delta, previousIncrementOfDisplacement);
435 for ( int i = 1; i <= neq; i++ ) {
436 help.at(i) = delta * ( a4 * previousVelocityVector.at(i)
439 }
440 this->timesMtrx(help, rhs2, TangentStiffnessMatrix, this->giveDomain(di), tStep);
441
442 help.zero();
443 rhs.add(rhs2);
444
445 //this->assembleVector(rhs, tStep, MatrixProductAssembler(TangentAssembler(), help), VM_Total,
446 // EModelDefaultEquationNumbering(), this->giveDomain(1));
447 }
448
449 rhs.add(loadVector);
450 rhs.subtract(previousInternalForces);
451
452 //
453 // Set-up numerical model.
454 //
456
457 //
458 // Call numerical model to solve problem.
459 //
460 double loadLevel = 1.0;
461
462 if ( totIterations == 0 ) {
464 }
465
466 ConvergedReason numMetStatus = nMethod->solve(*effectiveStiffnessMatrix, rhs, NULL,
468 internalForcesEBENorm, loadLevel, SparseNonLinearSystemNM :: rlm_total, currentIterations, tStep);
469 if ( numMetStatus != CR_CONVERGED ) {
470 OOFEM_ERROR("NRSolver failed to solve problem");
471 }
473 tStep->convergedReason = numMetStatus;
474
475
478 for ( int i = 1; i <= neq; i++ ) {
479 accelerationVector.at(i) = a0 * incrementOfDisplacement.at(i) - a2 *rhs.at(i) - a3 *rhs2.at(i);
480 velocityVector.at(i) = a1 * incrementOfDisplacement.at(i) - a4 *rhs.at(i) - a5 *rhs2.at(i)
482 }
484}
485
486
487void
488NonLinearDynamic :: determineConstants(TimeStep *tStep)
489{
490 TimeDiscretizationType timeDiscretization = tStep->giveTimeDiscretization();
491
492 if ( timeDiscretization == TD_Newmark ) {
493 OOFEM_LOG_DEBUG("Solving using Newmark-beta method\n");
494 } else if ( timeDiscretization == TD_TwoPointBackward ) {
495 OOFEM_LOG_DEBUG("Solving using Backward Euler method\n");
496 } else if ( timeDiscretization == TD_ThreePointBackward ) {
497 OOFEM_LOG_DEBUG("Solving using Three-point Backward Euler method\n");
498 if ( tStep->isTheFirstStep() ) {
499 timeDiscretization = TD_TwoPointBackward;
500 }
501 } else {
502 OOFEM_ERROR("Time-stepping scheme not found!\n")
503 }
504
505 deltaT = tStep->giveTimeIncrement();
506
507 double dt2 = deltaT * deltaT;
508
509 if ( timeDiscretization == TD_Newmark ) {
510 a0 = 1. / ( beta * dt2 );
511 a1 = gamma / ( beta * deltaT );
512 a2 = 1. / ( beta * deltaT );
513 a3 = 1. / ( 2. * beta ) - 1.;
514 a4 = ( gamma / beta ) - 1.;
515 a5 = deltaT / 2. * ( gamma / beta - 2. );
516 a6 = 0.;
517 } else if ( timeDiscretization == TD_TwoPointBackward ) {
518 a0 = 1. / dt2;
519 a1 = 1. / deltaT;
520 a2 = 1. / deltaT;
521 a3 = 0.;
522 a4 = 0.;
523 a5 = 0.;
524 a6 = 0.;
525 } else if ( timeDiscretization == TD_ThreePointBackward ) {
526 a0 = 2. / dt2;
527 a1 = 3. / ( 2. * deltaT );
528 a2 = 2. / deltaT;
529 a3 = 0.;
530 a4 = 0.;
531 a5 = 0.;
532 a6 = 1. / ( 2. * deltaT );
533 }
534}
535
536
537void NonLinearDynamic :: updateYourself(TimeStep *tStep)
538{
539 totIterations = 0;
540
546
547 StructuralEngngModel :: updateYourself(tStep);
548}
549
550
551void NonLinearDynamic :: updateSolution(FloatArray &solutionVector, TimeStep *tStep, Domain *d)
552{
553 // No-op; this still assumes that the solution vector is passed by reference to the NRSolver.
554 //this->field->update(VM_Total, tStep, solutionVector, EModelDefaultEquationNumbering());
555}
556
557
558void NonLinearDynamic :: updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm)
559{
560#ifdef VERBOSE
561 OOFEM_LOG_DEBUG("Updating internal RHS\n");
562#endif
563#ifdef TIME_REPORT
564 Timer timer;
565 timer.startTimer();
566#endif
567 if ( ( currentIterations != 0 ) || ( totIterations == 0 ) ) {
569
570 // Updating the residual vector @ NR-solver
571 help.beScaled(a0 + eta * a1, incrementOfDisplacement);
572
573 if (!( tStep->isTheFirstStep() && initFlag )) {
574 massMatrix->times(help, rhs2);
575
577 forcesVector.add(rhs2);
579
580 if ( delta != 0 ) {
582 this->timesMtrx(help, rhs2, TangentStiffnessMatrix, this->giveDomain(1), tStep);
583 //this->assembleVector(rhs2, tStep, MatrixProductAssembler(TangentAssembler(), help), VM_Total,
584 // EModelDefaultEquationNumbering(), this->giveDomain(1));
585 forcesVector.add(rhs2);
586 }
587 }
588 }
589#ifdef TIME_REPORT
590 timer.stopTimer();
591 OOFEM_LOG_DEBUG( "User time consumed by updating internal RHS: %.2fs\n", timer.getUtime() );
592#endif
593}
594
595
596void NonLinearDynamic :: updateMatrix(SparseMtrx &mat, TimeStep *tStep, Domain *d)
597{
598 // Prevent assembly if already assembled ( totIterations > 0 )
599 // Allow if MANRMSteps != 0
600 if ( ( totIterations == 0 ) || MANRMSteps ) {
601#ifdef VERBOSE
602 OOFEM_LOG_DEBUG("Updating effective stiffness matrix\n");
603#endif
604#ifdef TIME_REPORT
605 Timer timer;
606 timer.startTimer();
607#endif
608#if 1
609 mat.zero();
610 this->assemble(mat, tStep, EffectiveTangentAssembler(TangentStiffness, false, 1 + this->delta * a1, this->a0 + this->eta * this->a1),
612#else
613 this->assemble(mat, tStep, TangentStiffnessMatrix, EModelDefaultEquationNumbering(), d);
614 mat.times(1. + this->delta * a1);
615 mat.add(this->a0 + this->eta * this->a1, this->massMatrix);
616#endif
617#ifdef TIME_REPORT
618 timer.stopTimer();
619 OOFEM_LOG_DEBUG( "User time consumed by updating nonlinear LHS: %.2fs\n", timer.getUtime() );
620#endif
621 }
622}
623
624
625void NonLinearDynamic :: updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d)
626//
627// Updates some component, which is used by numerical method
628// to newly reached state. used mainly by numerical method
629// when new tanget stiffness is needed during finding
630// of new equlibrium stage.
631//
632{
633#ifdef TIME_REPORT
634 Timer timer;
635#endif
636
637 switch ( cmpn ) {
638 case NonLinearLhs:
639 // Prevent assembly if already assembled ( totIterations > 0 )
640 // Allow if MANRMSteps != 0
641 if ( ( totIterations == 0 ) || MANRMSteps ) {
642#ifdef VERBOSE
643 OOFEM_LOG_DEBUG("Updating effective stiffness matrix\n");
644#endif
645#ifdef TIME_REPORT
646 timer.startTimer();
647#endif
648#if 1
650 this->assemble(*effectiveStiffnessMatrix, tStep, EffectiveTangentAssembler(TangentStiffness, false, 1 + this->delta * a1, this->a0 + this->eta * this->a1),
652#else
653 this->assemble(effectiveStiffnessMatrix, tStep, TangentStiffnessMatrix,
655 effectiveStiffnessMatrix->times(1. + this->delta * a1);
656 effectiveStiffnessMatrix->add(this->a0 + this->eta * this->a1, this->massMatrix);
657#endif
658#ifdef TIME_REPORT
659 timer.stopTimer();
660 OOFEM_LOG_DEBUG( "User time consumed by updating nonlinear LHS: %.2fs\n", timer.getUtime() );
661#endif
662 }
663 break;
664 case InternalRhs:
665#ifdef VERBOSE
666 OOFEM_LOG_DEBUG("Updating internal RHS\n");
667#endif
668 {
669#ifdef TIME_REPORT
670 timer.startTimer();
671#endif
672 if ( ( currentIterations != 0 ) || ( totIterations == 0 ) ) {
674
675 // Updating the residual vector @ NR-solver
676 help.beScaled(a0 + eta * a1, incrementOfDisplacement);
677
678 massMatrix->times(help, rhs2);
679
681 forcesVector.add(rhs2);
683
684 if ( delta != 0 ) {
686 this->timesMtrx(help, rhs2, TangentStiffnessMatrix, this->giveDomain(1), tStep);
687 //this->assembleVector(rhs2, tStep, MatrixProductAssembler(TangentAssembler(), help), VM_Total,
688 // EModelDefaultEquationNumbering(), this->giveDomain(1));
689
690 forcesVector.add(rhs2);
691 }
692 }
693#ifdef TIME_REPORT
694 timer.stopTimer();
695 OOFEM_LOG_DEBUG( "User time consumed by updating internal RHS: %.2fs\n", timer.getUtime() );
696#endif
697 }
698 break;
699 default:
700 OOFEM_ERROR("Unknown Type of component.");
701 }
702}
703
704void
705NonLinearDynamic :: printOutputAt(FILE *file, TimeStep *tStep)
706{
707 if ( !this->giveDomain(1)->giveOutputManager()->testTimeStepOutput(tStep) ) {
708 return; // Do not print even Solution step header
709 }
710
711 fprintf(file, "\n\nOutput for time %.3e, solution step number %d\n", tStep->giveTargetTime(), tStep->giveNumber() );
712 fprintf(file, "Equilibrium reached in %d iterations\n\n", currentIterations);
713
714 nMethod->printState(file);
715
716 this->giveDomain(1)->giveOutputManager()->doDofManOutput(file, tStep);
717 this->giveDomain(1)->giveOutputManager()->doElementOutput(file, tStep);
718 this->printReactionForces(tStep, 1, file);
719}
720
721void
722NonLinearDynamic :: printDofOutputAt(FILE *stream, Dof *iDof, TimeStep *tStep)
723{
724 static char dofchar[] = "dva";
725 static ValueModeType dofmodes[] = {
726 VM_Total, VM_Velocity, VM_Acceleration
727 };
728
729 iDof->printMultipleOutputAt(stream, tStep, dofchar, dofmodes, 3);
730}
731
732void NonLinearDynamic :: saveContext(DataStream &stream, ContextMode mode)
733{
735
736 EngngModel :: saveContext(stream, mode);
737
738 if ( ( iores = incrementOfDisplacement.storeYourself(stream) ) != CIO_OK ) {
739 THROW_CIOERR(iores);
740 }
741
742 if ( ( iores = totalDisplacement.storeYourself(stream) ) != CIO_OK ) {
743 THROW_CIOERR(iores);
744 }
745
746 if ( ( iores = velocityVector.storeYourself(stream) ) != CIO_OK ) {
747 THROW_CIOERR(iores);
748 }
749
750 if ( ( iores = accelerationVector.storeYourself(stream) ) != CIO_OK ) {
751 THROW_CIOERR(iores);
752 }
753}
754
755
756void NonLinearDynamic :: restoreContext(DataStream &stream, ContextMode mode)
757{
759
760 EngngModel :: restoreContext(stream, mode);
761
762 if ( ( iores = incrementOfDisplacement.restoreYourself(stream) ) != CIO_OK ) {
763 THROW_CIOERR(iores);
764 }
765
766 if ( ( iores = totalDisplacement.restoreYourself(stream) ) != CIO_OK ) {
767 THROW_CIOERR(iores);
768 }
769
770 if ( ( iores = velocityVector.restoreYourself(stream) ) != CIO_OK ) {
771 THROW_CIOERR(iores);
772 }
773
774 if ( ( iores = accelerationVector.restoreYourself(stream) ) != CIO_OK ) {
775 THROW_CIOERR(iores);
776 }
777}
778
779
780void
781NonLinearDynamic :: updateDomainLinks()
782{
783 EngngModel :: updateDomainLinks();
784
785 this->giveNumericalMethod( this->giveCurrentMetaStep() )->setDomain( this->giveDomain(1) );
786#ifdef __MPI_PARALLEL_MODE
787 if ( this->giveLoadBalancer() ) {
788 this->giveLoadBalancer()->setDomain( this->giveDomain(1) );
789 }
790
791#endif
792}
793
794
795void
796NonLinearDynamic :: assemble(SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma,
797 const UnknownNumberingScheme &s, Domain *domain)
798{
799#ifdef TIME_REPORT
800 Timer timer;
801 timer.startTimer();
802#endif
803
804 EngngModel :: assemble(answer, tStep, ma, s, domain);
805
806 if ( ( nonlocalStiffnessFlag ) && dynamic_cast< const TangentAssembler* >(&ma) ) {
807 // Add nonlocal contribution.
808 for ( auto &elem : domain->giveElements() ) {
809 static_cast< NLStructuralElement * >( elem.get() )->addNonlocalStiffnessContributions(answer, s, tStep);
810 }
811
812 // Print storage statistics.
813 answer.printStatistics();
814 }
815
816#ifdef TIME_REPORT
817 timer.stopTimer();
818 OOFEM_LOG_DEBUG( "NonLinearDynamic info: user time consumed by assembly: %.2fs\n", timer.getUtime() );
819#endif
820}
821
822#ifdef __OOFEG
823void
824NonLinearDynamic :: showSparseMtrxStructure(int type, oofegGraphicContext &gc, TimeStep *tStep)
825{
826 Domain *domain = this->giveDomain(1);
827 CharType ctype;
828
829 if ( type != 1 ) {
830 return;
831 }
832
833 ctype = TangentStiffnessMatrix;
834
835 for ( auto &elem : domain->giveElements() ) {
836 elem->showSparseMtrxStructure(ctype, gc, tStep);
837 }
838
839 for ( auto &elem : domain->giveElements() ) {
840 elem->showExtendedSparseMtrxStructure(ctype, gc, tStep);
841 }
842}
843#endif
844
845
846void
847NonLinearDynamic :: timesMtrx(FloatArray &vec, FloatArray &answer, CharType type, Domain *domain, TimeStep *tStep)
848{
849 int nelem = domain->giveNumberOfElements();
850 //int neq = this->giveNumberOfDomainEquations(1, EModelDefaultEquationNumbering());
851 int jj, kk, n;
852 FloatMatrix charMtrx, R;
853 IntArray loc;
855
856 answer.resize( this->giveNumberOfDomainEquations(domain->giveNumber(), en) );
857 answer.zero();
858 for ( int i = 1; i <= nelem; i++ ) {
859 Element *element = domain->giveElement(i);
860 // Skip remote elements (these are used as mirrors of remote elements on other domains
861 // when nonlocal constitutive models are used. They introduction is necessary to
862 // allow local averaging on domains without fine grain communication between domains).
863 if ( element->giveParallelMode() == Element_remote ) {
864 continue;
865 }
866
867 element->giveLocationArray(loc, en);
868 element->giveCharacteristicMatrix(charMtrx, type, tStep);
869
870 if ( charMtrx.isNotEmpty() ) {
872 if ( element->giveRotationMatrix(R) ) {
873 charMtrx.rotatedWith(R);
874 }
875
876 //
877 // assemble it manually
878 //
879#ifdef DEBUG
880 if ( loc.giveSize() != charMtrx.giveNumberOfRows() ) {
881 OOFEM_ERROR("dimension mismatch");
882 }
883
884#endif
885
886 n = loc.giveSize();
887
888 for ( int j = 1; j <= n; j++ ) {
889 jj = loc.at(j);
890 if ( jj ) {
891 for ( int k = 1; k <= n; k++ ) {
892 kk = loc.at(k);
893 if ( kk ) {
894 answer.at(jj) += charMtrx.at(j, k) * vec.at(kk);
895 }
896 }
897 }
898 }
899 }
900 }
901
903}
904
905
906int
907NonLinearDynamic :: estimateMaxPackSize(IntArray &commMap, DataStream &buff, int packUnpackType)
908{
909 int count = 0, pcount = 0;
910 Domain *domain = this->giveDomain(1);
911
912 if ( packUnpackType == 0 ) {
913 for ( int map: commMap ) {
914 DofManager *dman = domain->giveDofManager( map );
915 for ( Dof *dof: *dman ) {
916 if ( dof->isPrimaryDof() && ( dof->__giveEquationNumber() ) ) {
917 count++;
918 } else {
919 pcount++;
920 }
921 }
922 }
923
924 return ( buff.givePackSizeOfDouble(1) * max(count, pcount) );
925 } else if ( packUnpackType == 1 ) {
926 for ( int map: commMap ) {
927 count += domain->giveElement( map )->estimatePackSize(buff);
928 }
929
930 return count;
931 }
932
933 return 0;
934}
935
936
937#ifdef __MPI_PARALLEL_MODE
939NonLinearDynamic :: giveLoadBalancer()
940{
941 if ( lb ) {
942 return lb.get();
943 }
944
945 if ( loadBalancingFlag ) {
946 lb = classFactory.createLoadBalancer( "parmetis", this->giveDomain(1) );
947 return lb.get();
948 } else {
949 return nullptr;
950 }
951}
952
953
955NonLinearDynamic :: giveLoadBalancerMonitor()
956{
957 if ( lbm ) {
958 return lbm.get();
959 }
960
961 if ( loadBalancingFlag ) {
962 lbm = classFactory.createLoadBalancerMonitor( "wallclock", this);
963 return lbm.get();
964 } else {
965 return nullptr;
966 }
967}
968#endif
969
970void
971NonLinearDynamic :: packMigratingData(TimeStep *tStep)
972{
973 Domain *domain = this->giveDomain(1);
974 int ndofman = domain->giveNumberOfDofManagers(), _eq;
975
976 for ( int idofman = 1; idofman <= ndofman; idofman++ ) {
977 DofManager *_dm = domain->giveDofManager(idofman);
978 for ( Dof *_dof: *_dm ) {
979 if ( _dof->isPrimaryDof() ) {
980 if ( ( _eq = _dof->__giveEquationNumber() ) ) {
981 // pack values in solution vectors
982 _dof->updateUnknownsDictionary( tStep, VM_Total, totalDisplacement.at(_eq) );
983 }
984 }
985 }
986 }
987}
988
989void
990NonLinearDynamic :: unpackMigratingData(TimeStep *tStep)
991{
992 Domain *domain = this->giveDomain(1);
993 int ndofman = domain->giveNumberOfDofManagers(), _eq;
994 //int myrank = this->giveRank();
995
996 // resize target arrays
998 totalDisplacement.resize(neq);
999 incrementOfDisplacement.resize(neq);
1000
1001 for ( int idofman = 1; idofman <= ndofman; idofman++ ) {
1002 DofManager *_dm = domain->giveDofManager(idofman);
1003 for ( Dof *_dof: *_dm ) {
1004 if ( _dof->isPrimaryDof() ) {
1005 if ( ( _eq = _dof->__giveEquationNumber() ) ) {
1006 // pack values in solution vectors
1007 totalDisplacement.at(_eq) = _dof->giveUnknownsDictionaryValue( tStep, VM_Total );
1008 #if 0
1009 // debug print
1010 if ( _dm->giveParallelMode() == DofManager_shared ) {
1011 fprintf(stderr, "[%d] Shared: %d(%d) -> %d\n", myrank, idofman, idof, _eq);
1012 } else {
1013 fprintf(stderr, "[%d] Local : %d(%d) -> %d\n", myrank, idofman, idof, _eq);
1014 }
1015
1016 #endif
1017 }
1018 }
1019 }
1020 }
1021
1022 this->initializeCommMaps(true);
1023 nMethod->reinitialize();
1024 // reinitialize error estimator (if any)
1025 if ( this->giveDomainErrorEstimator(1) ) {
1026 this->giveDomainErrorEstimator(1)->reinitialize();
1027 }
1028
1029 initFlag = true;
1030}
1031
1032} // end namespace oofem
#define REGISTER_EngngModel(class)
virtual int givePackSizeOfDouble(std::size_t count)=0
dofManagerParallelMode giveParallelMode() const
Definition dofmanager.h:526
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
int giveNumberOfDofManagers() const
Returns number of dof managers in domain.
Definition domain.h:461
DofManager * giveDofManager(int n)
Definition domain.C:317
Element * giveElement(int n)
Definition domain.C:165
std ::vector< std ::unique_ptr< Element > > & giveElements()
Definition domain.h:294
virtual bool giveRotationMatrix(FloatMatrix &answer)
Definition element.C:298
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Definition element.C:620
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
Definition element.C:429
int estimatePackSize(DataStream &buff)
Definition element.C:1748
elementParallelMode giveParallelMode() const
Definition element.h:1139
int giveNumberOfTimeStepWhenIcApply()
Returns the time step number, when initial conditions should apply.
Definition engngm.h:787
void initializeCommMaps(bool forceInit=false)
Definition engngm.C:2147
std::unique_ptr< LoadBalancer > lb
Load Balancer.
Definition engngm.h:304
virtual TimeStep * giveCurrentStep(bool force=false)
Definition engngm.h:717
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Definition engngm.C:452
int giveNumberOfProcesses() const
Returns the number of collaborating processes.
Definition engngm.h:1156
int giveRank() const
Returns domain rank in a group of collaborating processes (0..groupSize-1).
Definition engngm.h:1154
ProblemCommunicator * communicator
Communicator.
Definition engngm.h:315
std ::unique_ptr< TimeStep > previousStep
Previous time step.
Definition engngm.h:243
MetaStep * giveCurrentMetaStep()
Returns current meta step.
Definition engngm.C:1900
int ndomains
Number of receiver domains.
Definition engngm.h:215
ProblemCommunicator * nonlocCommunicator
NonLocal Communicator. Necessary when nonlocal constitutive models are used.
Definition engngm.h:318
Domain * giveDomain(int n)
Definition engngm.C:1936
std ::unique_ptr< TimeStep > currentStep
Current time step.
Definition engngm.h:241
int nonlocalExt
Flag indicating if nonlocal extension active, which will cause data to be sent between shared element...
Definition engngm.h:293
bool loadBalancingFlag
If set to true, load balancing is active.
Definition engngm.h:307
virtual ErrorEstimator * giveDomainErrorEstimator(int n)
Definition engngm.h:372
MetaStep * giveMetaStep(int i)
Returns the i-th meta step.
Definition engngm.h:773
EngngModelTimer timer
E-model timer.
Definition engngm.h:279
int nMetaSteps
Number of meta steps.
Definition engngm.h:235
CommunicatorBuff * commBuff
Common Communicator buffer.
Definition engngm.h:313
std::unique_ptr< LoadBalancerMonitor > lbm
Definition engngm.h:305
virtual int giveNumberOfFirstStep(bool force=false)
Definition engngm.h:763
bool isParallel() const
Returns true if receiver in parallel mode.
Definition engngm.h:1152
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
problemScale pScale
Multiscale mode.
Definition engngm.h:269
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(const FloatMatrix &r, char mode='n')
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
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
InputRecord & giveAttributesRecord()
Returns e-model attributes.
Definition metastep.h:122
TimeDiscretizationType initialTimeDiscretization
void updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm) override
FloatArray previousAccelerationVector
FloatArray previousTotalDisplacement
void determineConstants(TimeStep *tStep)
LinSystSolverType solverType
void proceedStep(int di, TimeStep *tStep)
FloatArray previousIncrementOfDisplacement
double deltaT
Intrinsic time increment.
NumericalMethod * giveNumericalMethod(MetaStep *mStep) override
Returns reference to receiver's numerical method.
std ::unique_ptr< SparseMtrx > effectiveStiffnessMatrix
LoadBalancer * giveLoadBalancer() override
std ::unique_ptr< SparseMtrx > massMatrix
std ::unique_ptr< SparseNonLinearSystemNM > nMethod
Numerical method used to solve the problem.
SparseMtrxType sparseMtrxType
void timesMtrx(FloatArray &answer, FloatArray &vec, CharType type, Domain *domain, TimeStep *tStep)
virtual void add(double x, SparseMtrx &m)
Definition sparsemtrx.h:166
virtual void zero()=0
Zeroes the receiver.
virtual void times(const FloatArray &x, FloatArray &answer) const
Definition sparsemtrx.h:137
virtual void printStatistics() const
Prints the receiver statistics (one-line) to stdout.
Definition sparsemtrx.h:273
virtual void addNonlocalStiffnessContributions(SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep)
StateCounterType internalVarUpdateStamp
FloatArray internalForcesEBENorm
Norm of nodal internal forces evaluated on element by element basis (squared).
void printReactionForces(TimeStep *tStep, int id, FILE *out)
StructuralEngngModel(int i, EngngModel *master=nullptr)
Creates new StructuralEngngModel with number i, associated to domain d.
void updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm) override
ConvergedReason convergedReason
Status of solution step (Converged,.
Definition timestep.h:120
double giveTimeIncrement()
Returns solution step associated time increment.
Definition timestep.h:168
double giveTargetTime()
Returns target time.
Definition timestep.h:164
int numberOfIterations
Number of itarations needed to achieve convergence.
Definition timestep.h:116
TimeDiscretizationType giveTimeDiscretization()
Returns time discretization.
Definition timestep.h:219
int giveNumber()
Returns receiver's number.
Definition timestep.h:144
bool isTheFirstStep()
Definition timestep.C:148
#define THROW_CIOERR(e)
#define _IFT_EngngModel_lstype
Definition engngm.h:91
#define _IFT_EngngModel_smtype
Definition engngm.h:92
#define OOFEM_ERROR(...)
Definition error.h:79
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#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.
FloatArrayF< N > max(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ DofManager_shared
Definition dofmanager.h:68
@ CBT_static
@ SMT_CompCol
Compressed column.
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.
@ NonLinearLhs
@ InternalRhs
@ microScale
Definition problemmode.h:47
#define _IFT_NonLinearDynamic_deltat
#define _IFT_NonLinearDynamic_ddtScheme
#define _IFT_NonLinearDynamic_delta
#define _IFT_NonLinearDynamic_eta
#define _IFT_NonLinearDynamic_nonlocalext
#define _IFT_NonLinearDynamic_gamma
#define _IFT_NonLinearDynamic_nonlocstiff
#define _IFT_NonLinearDynamic_beta
#define _IFT_NRSolver_manrmsteps
Definition nrsolver.h:58
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]

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