OOFEM 3.0
Loading...
Searching...
No Matches
adaptnlinearstatic.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
36#include "mathfem.h"
37#include "verbose.h"
38#include "timer.h"
39#include "metastep.h"
40#include "timestep.h"
41#include "nummet.h"
42#include "element.h"
43#include "node.h"
44#include "domain.h"
45#include "datareader.h"
46#include "oofemtxtdatareader.h"
47#include "remeshingcrit.h"
48#include "mesherinterface.h"
49#include "dof.h"
51#include "errorestimator.h"
52#include "classfactory.h"
53#include "datastream.h"
54#include "contextioerr.h"
55#include "oofem_terminate.h"
57
58#ifdef __MPI_PARALLEL_MODE
59 #include "parallelcontext.h"
60 #include "loadbalancer.h"
61#endif
62
63#ifdef __OOFEG
64 #include "oofeggraphiccontext.h"
65#endif
66
67#include <cstdlib>
68
69namespace oofem {
71
72AdaptiveNonLinearStatic :: AdaptiveNonLinearStatic(int i, EngngModel *_master) : NonLinearStatic(i, _master),
74{
75 meshPackage = MPT_T3D;
77
78#ifdef __MPI_PARALLEL_MODE
79 this->preMappingLoadBalancingFlag = false;
80#endif
81}
82
83
84AdaptiveNonLinearStatic :: ~AdaptiveNonLinearStatic()
85{ }
86
87
88void
89AdaptiveNonLinearStatic :: initializeFrom(InputRecord &ir)
90{
91 NonLinearStatic :: initializeFrom(ir);
92
93 int _val;
94
95 int meshPackageId = 0;
97 meshPackage = ( MeshPackageType ) meshPackageId;
98
101 _val = 0;
104
105
106 // check if error estimator initioalized
107 if (this->defaultErrEstimator == NULL) {
108 OOFEM_ERROR ("AdaptiveNonLinearStatic :: initializeFrom: Error estimator not defined [eetype missing]");
109 }
110}
111
112void
113AdaptiveNonLinearStatic :: solveYourselfAt(TimeStep *tStep)
114{
115 proceedStep(1, tStep);
116 this->updateYourself(tStep);
117
118#ifdef __OOFEG
119 ESIEventLoop( YES, const_cast< char * >("AdaptiveNonLinearStatic: Solution finished; Press Ctrl-p to continue") );
120#endif
121
122 //this->terminate( this->giveCurrentStep() );
123
124#ifdef __MPI_PARALLEL_MODE
126 this->balanceLoad( this->giveCurrentStep() );
127 }
128
129#endif
130
131 // evaluate error of the reached solution
132 this->defaultErrEstimator->estimateError( equilibratedEM, this->giveCurrentStep() );
133 //this->defaultErrEstimator->estimateError( temporaryEM, this->giveCurrentStep() );
134 this->defaultErrEstimator->giveRemeshingCrit()->estimateMeshDensities( this->giveCurrentStep() );
135 RemeshingStrategy strategy = this->defaultErrEstimator->giveRemeshingCrit()->giveRemeshingStrategy( this->giveCurrentStep() );
136
137 // if ((strategy == RemeshingFromCurrentState_RS) && (this->giveDomain(1)->giveSerialNumber() == 0))
138 // strategy = RemeshingFromPreviousState_RS;
139
140 if ( strategy == NoRemeshing_RS ) {
141 //
142 } else if ( ( strategy == RemeshingFromCurrentState_RS ) || ( strategy == RemeshingFromPreviousState_RS ) ) {
143
144 this->terminate( this->giveCurrentStep() ); // make output
145
146 // do remeshing
147 auto mesher = classFactory.createMesherInterface( meshPackage, this->giveDomain(1) );
148
149 Domain *newDomain;
150 MesherInterface :: returnCode result = mesher->createMesh(this->giveCurrentStep(), 1,
151 this->giveDomain(1)->giveSerialNumber() + 1, & newDomain);
152
153 if ( result == MesherInterface :: MI_OK ) {
154 this->initFlag = 1;
155 this->adaptiveRemap(newDomain);
156 } else if ( result == MesherInterface :: MI_NEEDS_EXTERNAL_ACTION ) {
157 if ( strategy == RemeshingFromCurrentState_RS ) {
158 // ensure the updating the step
160 //this->terminate (this->giveCurrentStep());
161 } else {
162 // save previous step (because update not called)
163 }
164
165 this->terminateAnalysis();
166 throw OOFEM_Terminate();
167 } else {
168 OOFEM_ERROR("createMesh failed");
169 }
170 }
171}
172
173void
174AdaptiveNonLinearStatic :: updateYourself(TimeStep *tStep)
175{
176 if ( timeStepLoadLevels.isEmpty() ) {
177 timeStepLoadLevels.resize( this->giveNumberOfSteps() );
178 }
179
180 // in case of adaptive restart from given timestep
181 // a time step with same number and incremented version will be generated.
182 // Then load level reached on new discretization will overwrite the old one obtained for
183 // old discretization for this step. But this is consistent, since when initialLoadVector
184 // is requested to be recovered the reference load vectors are assembled
185 // on actual discretization.
186 timeStepLoadLevels.at( tStep->giveNumber() ) = loadLevel;
187
188 NonLinearStatic :: updateYourself(tStep);
189}
190
191
192double AdaptiveNonLinearStatic :: giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof)
193{
194 int eq = dof->__giveEquationNumber();
195#ifdef DEBUG
196 if ( eq == 0 ) {
197 OOFEM_ERROR("invalid equation number");
198 }
199#endif
200
201 if ( tStep != this->giveCurrentStep() ) {
202 OOFEM_ERROR("unknown time step encountered");
203 }
204
205 if ( d->giveNumber() == 2 ) {
206 switch ( mode ) {
207 case VM_Incremental:
208 if ( d2_incrementOfDisplacement.isNotEmpty() ) {
209 return d2_incrementOfDisplacement.at(eq);
210 } else {
211 return 0.;
212 }
213
214 case VM_Total:
215 if ( d2_totalDisplacement.isNotEmpty() ) {
216 return d2_totalDisplacement.at(eq);
217 } else {
218 return 0.;
219 }
220
221 default:
222 OOFEM_ERROR("Unknown is of undefined ValueModeType for this problem");
223 }
224 } else {
225 return NonLinearStatic :: giveUnknownComponent(mode, tStep, d, dof);
226 }
227
228 // return 0.0;
229}
230
231
232int
233AdaptiveNonLinearStatic :: initializeAdaptiveFrom(EngngModel *sourceProblem)
234{
235 int result = 1;
236
237 // measure time consumed by mapping
238 Timer timer;
239 double mc1, mc2, mc3;
240 timer.startTimer();
241
242 if ( dynamic_cast< AdaptiveNonLinearStatic * >(sourceProblem) ) {
243 OOFEM_ERROR("source problem must also be AdaptiveNonlinearStatic.");
244 }
245
246 this->currentStep = std::make_unique<TimeStep>( * ( sourceProblem->giveCurrentStep() ) );
247 if ( sourceProblem->givePreviousStep() ) {
248 this->previousStep = std::make_unique<TimeStep>( * ( sourceProblem->givePreviousStep() ) );
249 }
250
251 // map primary unknowns
253
256 totalDisplacement.zero();
258
259 result &= mapper.mapAndUpdate( totalDisplacement, VM_Total,
260 sourceProblem->giveDomain(1), this->giveDomain(1), sourceProblem->giveCurrentStep() );
261
262 result &= mapper.mapAndUpdate( incrementOfDisplacement, VM_Incremental,
263 sourceProblem->giveDomain(1), this->giveDomain(1), sourceProblem->giveCurrentStep() );
264
265 timer.stopTimer();
266 mc1 = timer.getUtime();
267 timer.startTimer();
268
269 // map internal ip state
270 for ( auto &e : this->giveDomain(1)->giveElements() ) {
271 result &= e->adaptiveMap( sourceProblem->giveDomain(1), sourceProblem->giveCurrentStep() );
272 }
273
274 timer.stopTimer();
275 mc2 = timer.getUtime();
276 timer.startTimer();
277
278 // computes the stresses and calls updateYourself to mapped state
279 for ( auto &e : this->giveDomain(1)->giveElements() ) {
280 result &= e->adaptiveUpdate(currentStep.get());
281 }
282
283 // finish mapping process
284 for ( auto &e : this->giveDomain(1)->giveElements() ) {
285 result &= e->adaptiveFinish(currentStep.get());
286 }
287
288
289 // increment time step if mapped state will be considered as new solution stepL
290 /*
291 * this->giveNextStep();
292 * if (equilibrateMappedConfigurationFlag) {
293 * // we need to assemble the load vector in same time as the restarted step,
294 * // so new time step is generated with same intrincic time as has the
295 * // previous step if equilibrateMappedConfigurationFlag is set.
296 * // this allows to equlibrate the previously reached state
297 * TimeStep* cts = this->giveCurrentStep();
298 * cts->setTime(cts->giveTime()-cts->giveTimeIncrement());
299 * cts = this->givePreviousStep();
300 * cts->setTime(cts->giveTime()-cts->giveTimeIncrement());
301 * }
302 *
303 *
304 * if (this->giveCurrentStep()->giveNumber() ==
305 * this->giveMetaStep(this->giveCurrentStep()->giveMetaStepNumber())->giveFirstStepNumber()) {
306 * this->updateAttributes (this->giveCurrentStep());
307 * }
308 */
309 this->updateAttributes( this->giveCurrentMetaStep() );
310
311 // increment solution state counter - not needed, IPs are updated by adaptiveUpdate previously called
312 // and there is no change in primary vars.
313 // this->giveCurrentStep()->incrementStateCounter();
314
315 // assemble new initial load for new discretization
317 static_cast< AdaptiveNonLinearStatic * >(sourceProblem), 1, this->giveCurrentStep() );
318 // assemble new total load for new discretization
319 // this->assembleCurrentTotalLoadVector (totalLoadVector, totalLoadVectorOfPrescribed, this->giveCurrentStep());
320 // set bcloadVector to zero (no increment within same step)
321
322 timer.stopTimer();
323 mc3 = timer.getUtime();
324
325 // compute processor time used by the program
326 OOFEM_LOG_INFO("user time consumed by primary mapping: %.2fs\n", mc1);
327 OOFEM_LOG_INFO("user time consumed by ip mapping: %.2fs\n", mc2);
328 OOFEM_LOG_INFO("user time consumed by ip update: %.2fs\n", mc3);
329 OOFEM_LOG_INFO("user time consumed by mapping: %.2fs\n", mc1 + mc2 + mc3);
330
331 //
332
333 //
334 // bring mapped configuration into equilibrium
335 //
337 // use secant stiffness to resrore equilibrium
338 NonLinearStatic_stiffnessMode oldStiffMode = this->stiffMode;
340
341
342
343 if ( initFlag ) {
344 stiffnessMatrix = classFactory.createSparseMtrx(sparseMtrxType);
345 if ( !stiffnessMatrix ) {
346 OOFEM_ERROR("sparse matrix creation failed");
347 }
348
349 if ( nonlocalStiffnessFlag ) {
350 if ( !stiffnessMatrix->isAsymmetric() ) {
351 OOFEM_ERROR("stiffnessMatrix does not support asymmetric storage");
352 }
353 }
354
355 stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
356 stiffnessMatrix->zero(); // zero stiffness matrix
357 this->assemble( *stiffnessMatrix, this->giveCurrentStep(), TangentAssembler(SecantStiffness),
359 initFlag = 0;
360 }
361
362 // updateYourself() not necessary - the adaptiveUpdate previously called does the job
363 //this->updateYourself(this->giveCurrentStep());
364#ifdef VERBOSE
365 OOFEM_LOG_INFO( "Equilibrating mapped configuration [step number %5d.%d]\n",
366 this->giveCurrentStep()->giveNumber(), this->giveCurrentStep()->giveVersion() );
367#endif
368
369 //double deltaL = nMethod->giveUnknownComponent (StepLength, 0);
370 double deltaL = nMethod->giveCurrentStepLength();
372 refLoadInputMode, this->giveDomain(1), this->giveCurrentStep() );
373 //
374 // call numerical model to solve arised problem
375 //
376#ifdef VERBOSE
377 OOFEM_LOG_RELEVANT( "Solving [step number %5d.%d]\n",
378 this->giveCurrentStep()->giveNumber(), this->giveCurrentStep()->giveVersion() );
379#endif
380
381 //nMethod -> solveYourselfAt(this->giveCurrentStep()) ;
382 nMethod->setStepLength(deltaL / 5.0);
383 if ( initialLoadVector.isNotEmpty() ) {
387 } else {
391 }
392
393
394 loadVector.zero();
395
396 this->updateYourself( this->giveCurrentStep() );
397 this->terminate( this->giveCurrentStep() );
398 this->updateLoadVectors( this->giveCurrentStep() );
399
400 // restore old step length
401 nMethod->setStepLength(deltaL);
402
403 stiffMode = oldStiffMode;
404 } else {
405 // comment this, if output for mapped configuration (not equilibrated) not wanted
406 this->printOutputAt( this->giveOutputStream(), this->giveCurrentStep() );
407 }
408
409 return result;
410}
411
412
413int
414AdaptiveNonLinearStatic :: initializeAdaptive(int tStepNumber)
415{
416 try {
417 FileDataStream stream(this->giveContextFileName(tStepNumber, 0), false);
418 this->restoreContext(stream, CM_State);
419 } catch(ContextIOERR & c) {
420 c.print();
421 exit(1);
422 }
423
424 this->initStepIncrements();
425
426 int sernum = this->giveDomain(1)->giveSerialNumber();
427 OOFEM_LOG_INFO("restoring domain %d.%d\n", 1, sernum + 1);
428 Domain *dNew = new Domain(2, sernum + 1, this);
429 OOFEMTXTDataReader domainDr(this->giveDomainFileName(1, sernum + 1));
430 if ( !dNew->instanciateYourself(domainDr) ) {
431 OOFEM_ERROR("domain Instanciation failed");
432 }
433
434 // remap solution to new domain
435 return this->adaptiveRemap(dNew);
436}
437
438
439int
440AdaptiveNonLinearStatic :: adaptiveRemap(Domain *dNew)
441{
442 int result = 1;
443
444 this->initStepIncrements();
445
446 this->ndomains = 2;
447 this->domainNeqs.resize(2);
448 this->domainPrescribedNeqs.resize(2);
449 this->domainNeqs.at(2) = 0;
450 this->domainPrescribedNeqs.at(2) = 0;
451 this->domainList.emplace(domainList.begin() + 1, dNew);
452 this->parallelContextList.emplace(parallelContextList.begin() + 1, this);
453
454 // init equation numbering
455 //this->forceEquationNumbering(2);
457
458 // measure time consumed by mapping
459 Timer timer;
460 double mc1, mc2, mc3;
461
462 timer.startTimer();
463
464 // map primary unknowns
466
471
472 result &= mapper.mapAndUpdate( d2_totalDisplacement, VM_Total,
473 this->giveDomain(1), this->giveDomain(2), this->giveCurrentStep() );
474
475 result &= mapper.mapAndUpdate( d2_incrementOfDisplacement, VM_Incremental,
476 this->giveDomain(1), this->giveDomain(2), this->giveCurrentStep() );
477
478 timer.stopTimer();
479 mc1 = timer.getUtime();
480 timer.startTimer();
481
482 // map internal ip state
483 for ( auto &e : this->giveDomain(2)->giveElements() ) {
484 /* HUHU CHEATING */
485
486 if ( e->giveParallelMode() == Element_remote ) {
487 continue;
488 }
489
490 result &= e->adaptiveMap( this->giveDomain(1), this->giveCurrentStep() );
491 }
492
493 /* replace domains */
494 OOFEM_LOG_DEBUG("deleting old domain\n");
495 //delete domainList->at(1);
496 //domainList->put(1, dNew);
497 //dNew->setNumber(1);
498 //domainList->put(2, NULL);
499
500 //domainList[0] = std :: move(domainList[1]);
501 domainList.erase(domainList.begin());
502 domainList[0]->setNumber(1);
503
505
506 // keep equation numbering of new domain
507 this->numberOfEquations = this->domainNeqs.at(1) = this->domainNeqs.at(2);
510
511 // update solution
514
515
516 this->ndomains = 1;
517 // init equation numbering
518 // this->forceEquationNumbering();
519 this->updateDomainLinks();
520
521#ifdef __MPI_PARALLEL_MODE
522 if ( isParallel() ) {
523 // set up communication patterns
524 this->initializeCommMaps(true);
526 }
527
528#endif
529
530 timer.stopTimer();
531 mc2 = timer.getUtime();
532 timer.startTimer();
533
534 // computes the stresses and calls updateYourself to mapped state
535 for ( auto &e : this->giveDomain(1)->giveElements() ) {
536 /* HUHU CHEATING */
537 if ( e->giveParallelMode() == Element_remote ) {
538 continue;
539 }
540
541 result &= e->adaptiveUpdate( this->giveCurrentStep() );
542 }
543
544 // finish mapping process
545 for ( auto &e : this->giveDomain(1)->giveElements() ) {
546 /* HUHU CHEATING */
547 if ( e->giveParallelMode() == Element_remote ) {
548 continue;
549 }
550
551 result &= e->adaptiveFinish( this->giveCurrentStep() );
552 }
553
554 nMethod->reinitialize();
555
556
557 // increment time step if mapped state will be considered as new solution stepL
558 // this->giveNextStep();
560 // we need to assemble the load vector in same time as the restarted step,
561 // so new time step is generated with same intrincic time as has the
562 // previous step if equilibrateMappedConfigurationFlag is set.
563 // this allows to equlibrate the previously reached state
564 TimeStep *cts = this->giveCurrentStep();
565 // increment version of solution step
566 cts->incrementVersion();
567
568 //cts->setTime(cts->giveTime()-cts->giveTimeIncrement());
569 //cts = this->givePreviousStep();
570 //cts->setTime(cts->giveTime()-cts->giveTimeIncrement());
571 }
572
573 if ( this->giveCurrentStep()->giveNumber() == this->giveCurrentMetaStep()->giveFirstStepNumber() ) {
574 this->updateAttributes( this->giveCurrentMetaStep() );
575 }
576
577 // increment solution state counter - not needed, IPs are updated by adaptiveUpdate previously called
578 // and there is no change in primary vars.
579 // this->giveCurrentStep()->incrementStateCounter();
580
581 // assemble new initial load for new discretization
583 this, 1, this->giveCurrentStep() );
586 this->giveCurrentStep() );
587
588 // assemble new total load for new discretization
589 // this->assembleCurrentTotalLoadVector (totalLoadVector, totalLoadVectorOfPrescribed, this->giveCurrentStep());
590 // set bcloadVector to zero (no increment within same step)
591
592 timer.stopTimer();
593 mc3 = timer.getUtime();
594
595 // compute processor time used by the program
596 OOFEM_LOG_INFO("user time consumed by primary mapping: %.2fs\n", mc1);
597 OOFEM_LOG_INFO("user time consumed by ip mapping: %.2fs\n", mc2);
598 OOFEM_LOG_INFO("user time consumed by ip update: %.2fs\n", mc3);
599 OOFEM_LOG_INFO("user time consumed by mapping: %.2fs\n", mc1 + mc2 + mc3);
600
601 //
602
603 /********
604 * #if 0
605 * {
606 *
607 * // evaluate the force error of mapped configuration
608 * this->updateComponent (this->giveCurrentStep(), InternalRhs);
609 * FloatArray rhs;
610 *
611 * loadVector.resize (this->numberOfEquations);
612 * loadVector.zero();
613 * this->assemble (loadVector, this->giveCurrentStep(), ExternalForcesVector_Total, this->giveDomain(1)) ;
614 * this->assemble (loadVector, this->giveCurrentStep(), ExternalForcesVector_Total, this->giveDomain(1));
615 *
616 * rhs = loadVector;
617 * rhs.times(loadLevel);
618 * if (initialLoadVector.isNotEmpty()) rhs.add(initialLoadVector);
619 * rhs.subtract(internalForces);
620 *
621 * //
622 * // compute forceError
623 * //
624 * // err is relative error of unbalanced forces
625 * double RR, RR0, forceErr = dotProduct (rhs.givePointer(),rhs.givePointer(),rhs.giveSize());
626 * if (initialLoadVector.isNotEmpty())
627 * RR0 = dotProduct (initialLoadVector.givePointer(), initialLoadVector.givePointer(), initialLoadVector.giveSize());
628 * else
629 * RR0 = 0.0;
630 * RR = dotProduct(loadVector.givePointer(),loadVector.givePointer(),loadVector.giveSize());
631 * // we compute a relative error norm
632 * if ((RR0 + RR * loadLevel * loadLevel) < calm_SMALL_NUM) forceErr = 0.;
633 * else forceErr = sqrt (forceErr / (RR0+RR * loadLevel * loadLevel));
634 *
635 * printf ("Relative Force Error of Mapped Configuration is %-15e\n", forceErr);
636 *
637 * }
638 **#endif
639 *************/
640
641#ifdef __OOFEG
642 ESIEventLoop( YES, const_cast< char * >("AdaptiveRemap: Press Ctrl-p to continue") );
643#endif
644
645 //
646 // bring mapped configuration into equilibrium
647 //
649 // use secant stiffness to resrore equilibrium
650 NonLinearStatic_stiffnessMode oldStiffMode = this->stiffMode;
652
653
654
655 if ( initFlag ) {
656 if ( !stiffnessMatrix ) {
657 stiffnessMatrix = classFactory.createSparseMtrx(sparseMtrxType);
658 if ( !stiffnessMatrix ) {
659 OOFEM_ERROR("sparse matrix creation failed");
660 }
661 }
662
663 if ( nonlocalStiffnessFlag ) {
664 if ( !stiffnessMatrix->isAsymmetric() ) {
665 OOFEM_ERROR("stiffnessMatrix does not support asymmetric storage");
666 }
667 }
668
669 stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
670 stiffnessMatrix->zero(); // zero stiffness matrix
671 this->assemble( *stiffnessMatrix, this->giveCurrentStep(), TangentAssembler(SecantStiffness),
673
674 initFlag = 0;
675 }
676
677 // updateYourself() not necessary - the adaptiveUpdate previously called does the job
678 //this->updateYourself(this->giveCurrentStep());
679#ifdef VERBOSE
680 OOFEM_LOG_INFO( "Equilibrating mapped configuration [step number %5d.%d]\n",
681 this->giveCurrentStep()->giveNumber(), this->giveCurrentStep()->giveVersion() );
682#endif
683
684 //double deltaL = nMethod->giveUnknownComponent (StepLength, 0);
685 double deltaL = nMethod->giveCurrentStepLength();
686 //
687 // call numerical model to solve arised problem
688 //
689#ifdef VERBOSE
690 OOFEM_LOG_RELEVANT( "Solving [step number %5d.%d]\n",
691 this->giveCurrentStep()->giveNumber(), this->giveCurrentStep()->giveVersion() );
692#endif
693
694 //nMethod -> solveYourselfAt(this->giveCurrentStep()) ;
695 nMethod->setStepLength(deltaL / 5.0);
696 if ( initialLoadVector.isNotEmpty() ) {
700 } else {
704 }
705
706
707 loadVector.zero();
708
709 this->updateYourself( this->giveCurrentStep() );
710 this->terminate( this->giveCurrentStep() );
711 // this->updateLoadVectors (this->giveCurrentStep()); // already in terminate
712
713 // restore old step length
714 nMethod->setStepLength(deltaL);
715
716 stiffMode = oldStiffMode;
717 } else {
718 // comment this, if output for mapped configuration (not equilibrated) not wanted
719 this->printOutputAt( this->giveOutputStream(), this->giveCurrentStep() );
720 }
721
722 return result;
723}
724
725
726void
727AdaptiveNonLinearStatic :: saveContext(DataStream &stream, ContextMode mode)
728{
729 NonLinearStatic :: saveContext(stream, mode);
730
732 if ( ( iores = timeStepLoadLevels.storeYourself(stream) ) != CIO_OK ) {
733 THROW_CIOERR(iores);
734 }
735}
736
737void
738AdaptiveNonLinearStatic :: restoreContext(DataStream &stream, ContextMode mode)
739{
740 NonLinearStatic :: restoreContext(stream, mode);
741
743 if ( ( iores = timeStepLoadLevels.restoreYourself(stream) ) != CIO_OK ) {
744 THROW_CIOERR(iores);
745 }
746}
747
748
749void
750AdaptiveNonLinearStatic :: updateDomainLinks()
751{
752 NonLinearStatic :: updateDomainLinks();
753 this->defaultErrEstimator->setDomain( this->giveDomain(1) );
754}
755
756
757void
758AdaptiveNonLinearStatic :: assembleInitialLoadVector(FloatArray &loadVector, FloatArray &loadVectorOfPrescribed,
759 AdaptiveNonLinearStatic *sourceProblem, int domainIndx,
760 TimeStep *tStep)
761{
762 int mStepNum = tStep->giveMetaStepNumber();
763 int hasfixed, mode;
764 FloatArray _incrementalLoadVector, _incrementalLoadVectorOfPrescribed;
765 SparseNonLinearSystemNM :: referenceLoadInputModeType rlm;
766 //Domain* sourceDomain = sourceProblem->giveDomain(domainIndx);
767
769 loadVectorOfPrescribed.resize( this->giveNumberOfDomainEquations( domainIndx, EModelDefaultPrescribedEquationNumbering() ) );
770 loadVector.zero();
771 loadVectorOfPrescribed.zero();
772 _incrementalLoadVector.resize( this->giveNumberOfDomainEquations( domainIndx, EModelDefaultEquationNumbering() ) );
773 _incrementalLoadVectorOfPrescribed.resize( this->giveNumberOfDomainEquations( domainIndx, EModelDefaultPrescribedEquationNumbering() ) );
774 _incrementalLoadVector.zero();
775 _incrementalLoadVectorOfPrescribed.zero();
776
777 for ( int imstep = 1; imstep < mStepNum; imstep++ ) {
778 auto iMStep = this->giveMetaStep(imstep);
779 auto &ir = iMStep->giveAttributesRecord();
780 //hasfixed = ir.hasField("fixload");
781 hasfixed = 1;
782 if ( hasfixed ) {
783 // test for control mode
784 // here the algorithm works only for direct load control.
785 // Direct displacement control requires to know the quasi-rections, and the controlled nodes
786 // should have corresponding node on new mesh -> not supported
787 // Indirect control -> the load level from prevous steps is required, currently nt supported.
788
789 // additional problem: direct load control supports the reduction of step legth if convergence fails
790 // if this happens, this implementation does not work correctly.
791 // But there is NO WAY HOW TO TEST IF THIS HAPPEN
792
793 mode = 0;
795
796 // check if displacement control takes place
797 if ( ir.hasField(_IFT_AdaptiveNonLinearStatic_ddm) ) {
798 OOFEM_ERROR("fixload recovery not supported for direct displacement control");
799 }
800
801 int firststep = iMStep->giveFirstStepNumber();
802 int laststep = iMStep->giveLastStepNumber();
803
804 int _val = 0;
806 rlm = ( SparseNonLinearSystemNM :: referenceLoadInputModeType ) _val;
807
808 if ( mode == ( int ) nls_directControl ) { // and only load control
809 for ( int istep = firststep; istep <= laststep; istep++ ) {
810 // bad practise here
812 TimeStep *old = new TimeStep(istep, this, imstep, istep - 1.0, deltaT, 0);
813 this->assembleIncrementalReferenceLoadVectors(_incrementalLoadVector, _incrementalLoadVectorOfPrescribed,
814 rlm, this->giveDomain(domainIndx), old);
815
816 _incrementalLoadVector.times( sourceProblem->giveTimeStepLoadLevel(istep) );
817 loadVector.add(_incrementalLoadVector);
818 loadVectorOfPrescribed.add(_incrementalLoadVectorOfPrescribed);
819 }
820 } else if ( mode == ( int ) nls_indirectControl ) {
821 // bad practise here
822 if ( !ir.hasField(_IFT_NonLinearStatic_donotfixload) ) {
823 TimeStep *old = new TimeStep(firststep, this, imstep, firststep - 1.0, deltaT, 0);
824 this->assembleIncrementalReferenceLoadVectors(_incrementalLoadVector, _incrementalLoadVectorOfPrescribed,
825 rlm, this->giveDomain(domainIndx), old);
826
827 _incrementalLoadVector.times( sourceProblem->giveTimeStepLoadLevel(laststep) );
828 loadVector.add(_incrementalLoadVector);
829 loadVectorOfPrescribed.add(_incrementalLoadVectorOfPrescribed);
830 }
831 } else {
832 OOFEM_ERROR("fixload recovery not supported");
833 }
834 }
835 } // end loop over meta-steps
836
837 /* if direct control; add to initial load also previous steps in same metestep */
838 auto iMStep = this->giveMetaStep(mStepNum);
839 auto &ir = iMStep->giveAttributesRecord();
840 mode = 0;
842 int firststep = iMStep->giveFirstStepNumber();
843 int laststep = tStep->giveNumber();
844 int _val = 0;
846 rlm = ( SparseNonLinearSystemNM :: referenceLoadInputModeType ) _val;
847
848 if ( mode == ( int ) nls_directControl ) { // and only load control
849 for ( int istep = firststep; istep <= laststep; istep++ ) {
850 // bad practise here
851 TimeStep *old = new TimeStep(istep, this, mStepNum, istep - 1.0, deltaT, 0);
852 this->assembleIncrementalReferenceLoadVectors(_incrementalLoadVector, _incrementalLoadVectorOfPrescribed,
853 rlm, this->giveDomain(domainIndx), old);
854
855 _incrementalLoadVector.times( sourceProblem->giveTimeStepLoadLevel(istep) );
856 loadVector.add(_incrementalLoadVector);
857 loadVectorOfPrescribed.add(_incrementalLoadVectorOfPrescribed);
858 }
859 }
860}
861
862/*
863 * void
864 * AdaptiveNonLinearStatic::assembleCurrentTotalLoadVector (FloatArray& loadVector,
865 * FloatArray& loadVectorOfPrescribed,
866 * AdaptiveNonLinearStatic* sourceProblem, int domainIndx,
867 * TimeStep* tStep)
868 * {
869 * int mStepNum = tStep->giveMetaStepNumber() ;
870 * int mode;
871 * InputRecord* ir;
872 * MetaStep* mStep = sourceProblem->giveMetaStep(mStepNum);
873 * FloatArray _incrementalLoadVector, _incrementalLoadVectorOfPrescribed;
874 * SparseNonLinearSystemNM::referenceLoadInputModeType rlm;
875 * //Domain* sourceDomain = sourceProblem->giveDomain(domainIndx);
876 *
877 * loadVector.resize(this->giveNumberOfDomainEquations());
878 * loadVectorOfPrescribed.resize(this->giveNumberOfPrescribedDomainEquations());
879 * loadVector.zero(); loadVectorOfPrescribed.zero();
880 * _incrementalLoadVector.resize(this->giveNumberOfDomainEquations());
881 * _incrementalLoadVectorOfPrescribed.resize(this->giveNumberOfPrescribedDomainEquations());
882 * _incrementalLoadVector.zero(); _incrementalLoadVectorOfPrescribed.zero();
883 *
884 * // ASK CURRENT MSTEP FOR ITS RECORD
885 * ir = mStep->giveAttributesRecord();
886 *
887 * mode = 0;
888 * IR_GIVE_OPTIONAL_FIELD (ir, mode, _IFT_AdaptiveNonLinearStatic_controlmode, "controlmode");
889 *
890 * // check if displacement control takes place
891 * if (ir->hasField(_IFT_AdaptiveNonLinearStatic_ddm, "ddm"))
892 * OOFEM_ERROR("fixload recovery not supported for direct displacement control");
893 * int _val = 0;
894 * IR_GIVE_OPTIONAL_FIELD (ir, _val, _IFT_AdaptiveNonLinearStatic_refloadmode, "refloadmode");
895 *
896 * int firststep = mStep->giveFirstStepNumber();
897 * int laststep = tStep->giveNumber()-1;
898 *
899 * rlm = (SparseNonLinearSystemNM::referenceLoadInputModeType) _val;
900 *
901 * if (mode == (int)nls_directControl) { // and only load control
902 * for (int istep = firststep; istep<=laststep; istep++) {
903 * // bad practise here
904 * TimeStep* old = new TimeStep (istep, this, mStepNum, istep-1.0, deltaT, 0);
905 * this->assembleIncrementalReferenceLoadVectors (_incrementalLoadVector, _incrementalLoadVectorOfPrescribed,
906 * rlm, this->giveDomain(domainIndx), old);
907 *
908 * _incrementalLoadVector.times(sourceProblem->giveTimeStepLoadLevel(istep));
909 * loadVector.add(_incrementalLoadVector);
910 * loadVectorOfPrescribed.add(_incrementalLoadVectorOfPrescribed);
911 * }
912 * } else if (mode == (int)nls_indirectControl) {
913 * // bad practise here
914 * TimeStep* old = new TimeStep (firststep, this, mStepNum, firststep-1.0, deltaT, 0);
915 * this->assembleIncrementalReferenceLoadVectors (_incrementalLoadVector, _incrementalLoadVectorOfPrescribed,
916 * rlm, this->giveDomain(domainIndx), old);
917 *
918 * _incrementalLoadVector.times(sourceProblem->giveTimeStepLoadLevel(laststep));
919 * loadVector.add(_incrementalLoadVector);
920 * loadVectorOfPrescribed.add(_incrementalLoadVectorOfPrescribed);
921 * } else {
922 * OOFEM_ERROR("totalload recovery not supported");
923 * }
924 *
925 *
926 * }
927 */
928
929
930double
931AdaptiveNonLinearStatic :: giveTimeStepLoadLevel(int istep)
932{
933 if ( ( istep >= this->giveNumberOfFirstStep() ) && ( istep <= this->giveNumberOfSteps() ) ) {
934 return timeStepLoadLevels.at(istep);
935 } else {
936 OOFEM_ERROR("solution step out of range");
937 }
938
939 // return 0.0; // to make compiler happy
940}
941
942
943#ifdef __MPI_PARALLEL_MODE
945AdaptiveNonLinearStatic :: giveLoadBalancer()
946{
947 if ( lb ) {
948 return lb.get();
949 }
950
952 lb = classFactory.createLoadBalancer( "parmetis", this->giveDomain(1) );
953 return lb.get();
954 } else {
955 return nullptr;
956 }
957}
959AdaptiveNonLinearStatic :: giveLoadBalancerMonitor()
960{
961 if ( lbm ) {
962 return lbm.get();
963 }
964
966 lbm = classFactory.createLoadBalancerMonitor( "wallclock", this);
967 return lbm.get();
968 } else {
969 return nullptr;
970 }
971}
972#endif
973} // end namespace oofem
#define _IFT_AdaptiveNonLinearStatic_ddm
#define _IFT_AdaptiveNonLinearStatic_preMappingLoadBalancingFlag
#define _IFT_AdaptiveNonLinearStatic_controlmode
#define _IFT_AdaptiveNonLinearStatic_meshpackage
#define _IFT_AdaptiveNonLinearStatic_equilmc
#define _IFT_AdaptiveNonLinearStatic_refloadmode
#define REGISTER_EngngModel(class)
void assembleInitialLoadVector(FloatArray &loadVector, FloatArray &loadVectorOfPrescribed, AdaptiveNonLinearStatic *sourceProblem, int domainIndx, TimeStep *tStep)
virtual double giveTimeStepLoadLevel(int istep)
virtual int adaptiveRemap(Domain *dNew)
void restoreContext(DataStream &stream, ContextMode mode) override
AdaptiveNonLinearStatic(int i, EngngModel *master=nullptr)
int equilibrateMappedConfigurationFlag
Flag indication whether to restore equilibrium after adaptive remapping.
void updateYourself(TimeStep *tStep) override
virtual int __giveEquationNumber() const =0
int instanciateYourself(DataReader &dr)
Definition domain.C:480
int giveNumber()
Returns domain number.
Definition domain.h:281
int mapAndUpdate(FloatArray &answer, ValueModeType mode, Domain *oldd, Domain *newd, TimeStep *tStep) override
virtual void initStepIncrements()
Definition engngm.C:1702
void initializeCommMaps(bool forceInit=false)
Definition engngm.C:2147
std::unique_ptr< LoadBalancer > lb
Load Balancer.
Definition engngm.h:304
std ::vector< std ::unique_ptr< Domain > > domainList
List of problem domains.
Definition engngm.h:217
void setContextOutputMode(ContextOutputMode contextMode)
Definition engngm.h:411
int numberOfEquations
Total number of equation in current time step.
Definition engngm.h:221
virtual void balanceLoad(TimeStep *tStep)
Definition engngm.C:2238
virtual TimeStep * giveCurrentStep(bool force=false)
Definition engngm.h:717
void terminateAnalysis()
Definition engngm.C:2022
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Definition engngm.C:452
std ::unique_ptr< TimeStep > previousStep
Previous time step.
Definition engngm.h:243
std::string giveContextFileName(int tStepNumber, int stepVersion) const
Definition engngm.C:1907
int equationNumberingCompleted
Equation numbering completed flag.
Definition engngm.h:233
MetaStep * giveCurrentMetaStep()
Returns current meta step.
Definition engngm.C:1900
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
bool loadBalancingFlag
If set to true, load balancing is active.
Definition engngm.h:307
FILE * giveOutputStream()
Returns file descriptor of output file.
Definition engngm.C:1994
int giveNumberOfSteps(bool force=false)
Definition engngm.h:777
MetaStep * giveMetaStep(int i)
Returns the i-th meta step.
Definition engngm.h:773
EngngModelTimer timer
E-model timer.
Definition engngm.h:279
@ RemoteElementExchangeTag
Definition engngm.h:321
std ::vector< ParallelContext > parallelContextList
List where parallel contexts are stored.
Definition engngm.h:323
int exchangeRemoteElementData(int ExchangeTag)
Definition engngm.C:2200
virtual TimeStep * givePreviousStep(bool force=false)
Definition engngm.h:738
std::string giveDomainFileName(int domainNum, int domainSerNum) const
Definition engngm.C:1917
virtual int forceEquationNumbering(int i)
Definition engngm.C:469
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
IntArray domainPrescribedNeqs
Number of prescribed equations per domain.
Definition engngm.h:227
std::unique_ptr< ErrorEstimator > defaultErrEstimator
Error estimator. Useful for adaptivity, or simply printing errors output.
Definition engngm.h:285
int numberOfPrescribedEquations
Total number or prescribed equations in current time step.
Definition engngm.h:223
IntArray domainNeqs
Number of equations per domain.
Definition engngm.h:225
void resize(Index s)
Definition floatarray.C:94
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void add(const FloatArray &src)
Definition floatarray.C:218
void times(double s)
Definition floatarray.C:834
SparseMtrxType sparseMtrxType
FloatArray loadVector
std ::unique_ptr< SparseMtrx > stiffnessMatrix
SparseNonLinearSystemNM::referenceLoadInputModeType refLoadInputMode
FloatArray initialLoadVectorOfPrescribed
A load vector which does not scale for prescribed DOFs.
void printOutputAt(FILE *file, TimeStep *tStep) override
FloatArray incrementalLoadVector
FloatArray incrementalLoadVectorOfPrescribed
Incremental Load Vector for prescribed DOFs.
double deltaT
Intrinsic time increment.
void terminate(TimeStep *tStep) override
NonLinearStatic(int i, EngngModel *master=nullptr)
void proceedStep(int di, TimeStep *tStep)
void assembleIncrementalReferenceLoadVectors(FloatArray &_incrementalLoadVector, FloatArray &_incrementalLoadVectorOfPrescribed, SparseNonLinearSystemNM ::referenceLoadInputModeType _refMode, Domain *sourceDomain, TimeStep *tStep)
virtual void updateLoadVectors(TimeStep *tStep)
void updateAttributes(MetaStep *mStep) override
NonLinearStatic_stiffnessMode stiffMode
ConvergedReason numMetStatus
FloatArray incrementOfDisplacement
SparseNonLinearSystemNM * nMethod
Numerical method used to solve the problem.
FloatArray initialLoadVector
A load vector already applied, which does not scales.
OOFEM terminate exception class.
FloatArray internalForcesEBENorm
Norm of nodal internal forces evaluated on element by element basis (squared).
void incrementVersion()
Increments receiver's version.
Definition timestep.h:215
int giveMetaStepNumber()
Returns receiver's meta step number.
Definition timestep.h:150
int giveNumber()
Returns receiver's number.
Definition timestep.h:144
#define THROW_CIOERR(e)
#define CM_State
Definition contextmode.h:46
#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
@ equilibratedEM
@ 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.
@ COM_Always
Enable for post-processing.
NonLinearStatic_stiffnessMode
Type determining the stiffness mode.
@ nls_secantStiffness
The secant stiffness is used and updated whenever requested.
RemeshingStrategy
Type representing the remeshing strategy.
@ RemeshingFromPreviousState_RS
@ RemeshingFromCurrentState_RS
@ NoRemeshing_RS
ClassFactory & classFactory
@ nls_indirectControl
A generalized norm of displacement and loading vectors is controlled. In current implementation,...
@ nls_directControl
Describes the direct control where load or displacement (or both) are controlled.
#define _IFT_NonLinearStatic_donotfixload

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