OOFEM 3.0
Loading...
Searching...
No Matches
supg.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 "supg.h"
36#include "sparsemtrx.h"
37#include "nrsolver.h"
38#include "timestep.h"
39#include "element.h"
40#include "dofmanager.h"
41#include "dof.h"
42#include "initialcondition.h"
43#include "maskedprimaryfield.h"
44#include "verbose.h"
45#include "supgelement.h"
46#include "classfactory.h"
47#include "mathfem.h"
49#include "leplic.h"
50#include "levelsetpcs.h"
51#include "datastream.h"
52#include "function.h"
53#include "contextioerr.h"
54#include "timer.h"
56
57namespace oofem {
58/* define if implicit interface update required */
59//#define SUPG_IMPLICIT_INTERFACE
60
62
63
64SUPGInternalForceAssembler :: SUPGInternalForceAssembler(double l, double d, double u) :
66{}
67
68void SUPGInternalForceAssembler :: vectorFromElement(FloatArray &vec, Element &element, TimeStep *tStep, ValueModeType mode) const
69{
70 SUPGElement *eptr = static_cast< SUPGElement * >( &element );
71 FloatMatrix m1;
72 FloatArray vacc, vtot, ptot, h;
73 IntArray vloc, ploc;
74
75 int size = eptr->computeNumberOfDofs();
76 eptr->giveLocalVelocityDofMap(vloc);
77 eptr->giveLocalPressureDofMap(ploc);
78 vec.resize(size);
79 vec.zero();
80
81 eptr->computeVectorOfVelocities(VM_Acceleration, tStep, vacc);
82 eptr->computeVectorOfVelocities(VM_Total, tStep, vtot);
83 eptr->computeVectorOfPressures(VM_Total, tStep, ptot);
84
85 // MB contributions:
86 // add (M+M_delta)*a
87 eptr->computeAccelerationTerm_MB(m1, tStep);
88 h.beProductOf(m1, vacc);
89 vec.assemble(h, vloc);
90 // add advection terms N+N_delta
91 eptr->computeAdvectionTerm_MB(h, tStep);
92 vec.assemble(h, vloc);
93 // add diffusion terms (K+K_delta)h
94 eptr->computeDiffusionTerm_MB(h, tStep);
95 vec.assemble(h, vloc);
96 // add lsic stabilization term
97 eptr->computeLSICStabilizationTerm_MB(m1, tStep);
98 m1.times( lscale / ( dscale * uscale * uscale ) );
99 h.beProductOf(m1, vtot);
100 vec.assemble(h, vloc);
101 eptr->computeBCLhsTerm_MB(m1, tStep);
102 if ( m1.isNotEmpty() ) {
103 h.beProductOf(m1, vtot);
104 vec.assemble(h, vloc);
105 }
106
107 // add pressure term
108 eptr->computePressureTerm_MB(m1, tStep);
109 h.beProductOf(m1, ptot); // term due to prescribed pressure
110 vec.assemble(h, vloc);
111 eptr->computeBCLhsPressureTerm_MB(m1, tStep);
112 if ( m1.isNotEmpty() ) {
113 h.beProductOf(m1, ptot);
114 vec.assemble(h, vloc);
115 }
116
117 // MC contributions:
118 // G^T term - linear advection term
119 eptr->computeLinearAdvectionTerm_MC(m1, tStep);
120 //m1.times(uscale/lscale);
121 m1.times( 1. / ( dscale * uscale ) );
122 eptr->computeVectorOfVelocities(VM_Total, tStep, vtot);
123 h.beProductOf(m1, vtot); // term due to prescribed velocity
124 vec.assemble(h, ploc);
125 // Diffusion term
126 eptr->computeDiffusionTerm_MC(h, tStep);
127 vec.assemble(h, ploc);
128 // AccelerationTerm
129 eptr->computeAccelerationTerm_MC(m1, tStep);
130 h.beProductOf(m1, vacc);
131 vec.assemble(h, ploc);
132 eptr->computeBCLhsPressureTerm_MC(m1, tStep);
133 if ( m1.isNotEmpty() ) {
134 h.beProductOf(m1, vacc);
135 vec.assemble(h, ploc);
136 }
137
138 // advection N term (nonlinear)
139 eptr->computeAdvectionTerm_MC(h, tStep);
140 vec.assemble(h, ploc);
141 // pressure term
142 eptr->computePressureTerm_MC(m1, tStep);
143 h.beProductOf(m1, ptot);// term due to prescribed pressure
144 vec.assemble(h, ploc);
145};
146
147SUPGTangentAssembler :: SUPGTangentAssembler(MatResponseMode m, double l, double d, double u, double a) :
148 MatrixAssembler(), rmode(m), lscale(l), dscale(d), uscale(u), alpha(a)
149{}
150
151
152void SUPGTangentAssembler :: matrixFromElement(FloatMatrix &answer, Element &el, TimeStep *tStep) const
153{
154 SUPGElement *element = static_cast< SUPGElement * >( &el );
155 IntArray vloc, ploc;
156 FloatMatrix h;
157 int size = element->computeNumberOfDofs();
158 element->giveLocalVelocityDofMap(vloc);
159 element->giveLocalPressureDofMap(ploc);
160 answer.resize(size, size);
161 answer.zero();
162
163 element->computeAccelerationTerm_MB(h, tStep);
164 answer.assemble(h, vloc);
165 element->computeAdvectionDerivativeTerm_MB(h, tStep);
166 h.times( alpha * tStep->giveTimeIncrement() );
167 answer.assemble(h, vloc);
168 element->computeDiffusionDerivativeTerm_MB(h, rmode, tStep);
169
170 h.times( alpha * tStep->giveTimeIncrement() );
171 answer.assemble(h, vloc);
172 element->computePressureTerm_MB(h, tStep);
173 answer.assemble(h, vloc, ploc);
174 element->computeLSICStabilizationTerm_MB(h, tStep);
175 h.times( alpha * tStep->giveTimeIncrement() * lscale / ( dscale * uscale * uscale ) );
176 answer.assemble(h, vloc);
177 element->computeBCLhsTerm_MB(h, tStep);
178 if ( h.isNotEmpty() ) {
179 h.times( alpha * tStep->giveTimeIncrement() );
180 answer.assemble(h, vloc);
181 }
182
183 element->computeBCLhsPressureTerm_MB(h, tStep);
184 if ( h.isNotEmpty() ) {
185 answer.assemble(h, vloc, ploc);
186 }
187
188 // conservation eq part
189 element->computeLinearAdvectionTerm_MC(h, tStep);
190 h.times( alpha * tStep->giveTimeIncrement() * 1.0 / ( dscale * uscale ) );
191 answer.assemble(h, ploc, vloc);
192 element->computeAdvectionDerivativeTerm_MC(h, tStep);
193 h.times( alpha * tStep->giveTimeIncrement() );
194 answer.assemble(h, ploc, vloc);
195 element->computeAccelerationTerm_MC(h, tStep);
196 answer.assemble(h, ploc, vloc);
197 element->computeBCLhsPressureTerm_MC(h, tStep);
198 if ( h.isNotEmpty() ) {
199 answer.assemble(h, ploc, vloc);
200 }
201
202 element->computeDiffusionDerivativeTerm_MC(h, tStep);
203 h.times( alpha * tStep->giveTimeIncrement() );
204 answer.assemble(h, ploc, vloc);
205 element->computePressureTerm_MC(h, tStep);
206 answer.assemble(h, ploc);
207}
208
209
210SUPG :: SUPG(int i, EngngModel * _master) : FluidModel(i, _master), accelerationVector()
211{
212 ndomains = 1;
213}
214
215
216NumericalMethod *SUPG :: giveNumericalMethod(MetaStep *mStep)
217{
218 if ( !this->nMethod ) {
219 this->nMethod = classFactory.createSparseLinSolver(solverType, this->giveDomain(1), this);
220 if ( !this->nMethod ) {
221 OOFEM_ERROR("linear solver creation failed for lstype %d", solverType);
222 }
223 }
224 return this->nMethod.get();
225}
226
227void
228SUPG :: initializeFrom(InputRecord &ir)
229{
230 FluidModel :: initializeFrom(ir);
231
233 atolv = 1.e-15;
235
236
238
239 maxiter = 200;
241
242 int val = 0;
245
246 val = 0;
249
251 deltaTF = 0;
253
255
256 alpha = 0.5;
258
259 val = 0;
261 equationScalingFlag = val > 0;
262 if ( equationScalingFlag ) {
266 double vref = 1.0; // reference viscosity
267 Re = dscale * uscale * lscale / vref;
268 } else {
269 lscale = uscale = dscale = 1.0;
270 Re = 1.0;
271 }
272
274 VelocityPressureField = std::make_unique<DofDistributedPrimaryField>(this, 1, FT_VelocityPressure, 1);
275 } else {
276 VelocityPressureField = std::make_unique<PrimaryField>(this, 1, FT_VelocityPressure, 1);
277 }
278
279 val = 0;
281 if ( val == 1 ) {
282 this->materialInterface = std::make_unique<LEPlic>( 1, this->giveDomain(1) );
283 this->materialInterface->initializeFrom(ir);
284 // export velocity field
285 FieldManager *fm = this->giveContext()->giveFieldManager();
286 IntArray mask = {V_u, V_v, V_w};
287
288 std :: shared_ptr< Field > _velocityField = std::make_shared<MaskedPrimaryField>( FT_Velocity, this->VelocityPressureField.get(), mask );
289 fm->registerField(_velocityField, FT_Velocity);
290
291 //fsflag = 0;
292 //IR_GIVE_OPTIONAL_FIELD (ir, fsflag, _IFT_SUPG_fsflag, "fsflag");
293 } else if ( val == 2 ) {
294 // positive coefficient scheme level set alg
295 this->materialInterface = std::make_unique<LevelSetPCS>( 1, this->giveDomain(1) );
296 this->materialInterface->initializeFrom(ir);
297 }
298}
299
300
301double
302SUPG :: giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof)
303// returns unknown quantity like displaacement, velocity of equation eq
304// This function translates this request to numerical method language
305{
306 if ( this->requiresUnknownsDictionaryUpdate() ) {
307 int hash = this->giveUnknownDictHashIndx(mode, tStep);
308 if ( dof->giveUnknowns()->includes(hash) ) {
309 return dof->giveUnknowns()->at(hash);
310 } else {
311 OOFEM_ERROR("Dof unknowns dictionary does not contain unknown of value mode (%s)", __ValueModeTypeToString(mode));
312 }
313 } else {
314 int eq = dof->__giveEquationNumber();
315 if ( eq == 0 ) {
316 OOFEM_ERROR("invalid equation number");
317 }
318
319 if ( mode == VM_Acceleration ) {
320 return accelerationVector.at(eq);
321 /*
322 * if (tStep->isTheCurrentTimeStep()) {
323 * return accelerationVector.at(eq);
324 * } else if (tStep == givePreviousStep()) {
325 * return previousAccelerationVector.at(eq);
326 * } else OOFEM_ERROR("VM_Acceleration history past previous step not supported");
327 */
328 } else if ( mode == VM_Velocity ) {
329 return VelocityPressureField->giveUnknownValue(dof, VM_Total, tStep);
330 } else {
331 return VelocityPressureField->giveUnknownValue(dof, mode, tStep);
332 }
333 }
334
335 return 0;
336}
337
338
339void
340SUPG :: updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d)
341{
342 // update element stabilization
343 for ( auto &elem : d->giveElements() ) {
344 static_cast< FMElement * >( elem.get() )->updateStabilizationCoeffs(tStep);
345 }
346
347 if ( cmpn == InternalRhs ) {
348 this->internalForces.zero();
352 return;
353 } else if ( cmpn == NonLinearLhs ) {
354 this->lhs->zero();
355 //if ( 1 ) { //if ((nite > 5)) // && (rnorm < 1.e4))
356 this->assemble( *lhs, tStep, SUPGTangentAssembler(TangentStiffness, lscale, dscale, uscale, alpha),
358 // } else {
359 // this->assemble( lhs, tStep, SUPGTangentAssembler(SecantStiffness),
360 // EModelDefaultEquationNumbering(), d );
361 // }
362 return;
363 } else {
364 OOFEM_ERROR("Unknown component");
365 }
366}
367
368void SUPG :: updateSolution(FloatArray &solutionVector, TimeStep *tStep, Domain *d)
369{
370 this->VelocityPressureField->update(VM_Total, tStep, solutionVector, EModelDefaultEquationNumbering());
371 // update element stabilization
372 for ( auto &elem : d->giveElements() ) {
373 static_cast< FMElement * >( elem.get() )->updateStabilizationCoeffs(tStep);
374 }
375}
376
377
378void SUPG :: updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *enorm)
379{
380 answer.zero();
381 this->assembleVector(answer, tStep, SUPGInternalForceAssembler(lscale, dscale, uscale), VM_Total,
384}
385
386
387void SUPG :: updateMatrix(SparseMtrx &mat, TimeStep *tStep, Domain *d)
388{
389 mat.zero();
390 //if ( 1 ) { //if ((nite > 5)) // && (rnorm < 1.e4))
391 this->assemble( mat, tStep, SUPGTangentAssembler(TangentStiffness, lscale, dscale, uscale, alpha),
393 // } else {
394 // this->assemble( lhs, tStep, SUPGTangentAssembler(SecantStiffness),
395 // EModelDefaultEquationNumbering(), d );
396 // }
397}
398
399
400double
401SUPG :: giveReynoldsNumber()
402{
403 if ( equationScalingFlag ) {
404 return this->Re;
405 } else {
406 return 1.0;
407 }
408}
409
410
411TimeStep *
412SUPG :: giveSolutionStepWhenIcApply(bool force)
413{
414 if ( master && (!force)) {
415 return master->giveSolutionStepWhenIcApply();
416 } else {
417 if ( !stepWhenIcApply ) {
418 double dt = deltaT / this->giveVariableScale(VST_Time);
419
420 stepWhenIcApply = std::make_unique<TimeStep>(giveNumberOfTimeStepWhenIcApply(), this, 0, 0.0, dt, 0);
421 }
422
423 return stepWhenIcApply.get();
424 }
425}
426
427TimeStep *
428SUPG :: giveNextStep()
429{
430 double dt = deltaT;
431
432 Domain *domain = this->giveDomain(1);
433
434 if ( !currentStep ) {
435 // first step -> generate initial step
436 currentStep = std::make_unique<TimeStep>( *giveSolutionStepWhenIcApply() );
437 }
438
439 previousStep = std :: move(currentStep);
440
441 if ( deltaTF ) {
442 dt *= domain->giveFunction(deltaTF)->evaluateAtTime(previousStep->giveNumber() + 1);
443 }
444
445 // check for critical time step
446 for ( auto &elem : domain->giveElements() ) {
447 dt = min( dt, static_cast< SUPGElement & >( *elem ).computeCriticalTimeStep(previousStep.get()) );
448 }
449
450 if ( materialInterface ) {
451 dt = min( dt, materialInterface->computeCriticalTimeStep(previousStep.get()) );
452 }
453
454 // dt *= 0.6;
455 dt /= this->giveVariableScale(VST_Time);
456
457 currentStep = std::make_unique<TimeStep>(*previousStep, dt);
458
459 OOFEM_LOG_INFO( "SolutionStep %d : t = %e, dt = %e\n", currentStep->giveNumber(),
460 currentStep->giveTargetTime() * this->giveVariableScale(VST_Time), dt * this->giveVariableScale(VST_Time) );
461 return currentStep.get();
462}
463
464
465void
466SUPG :: solveYourselfAt(TimeStep *tStep)
467{
469 FloatArray *solutionVector = NULL, *prevSolutionVector = NULL;
470 FloatArray externalForces(neq);
472 this->internalForces.resize(neq);
473
474 if ( tStep->isTheFirstStep() ) {
475 TimeStep *_stepWhenIcApply = tStep->givePreviousStep();
476 if ( materialInterface ) {
477 materialInterface->initialize();
478 }
479
480 this->applyIC(_stepWhenIcApply);
481 //if (this->fsflag) this->updateDofManActivityMap(tStep);
482 }
483
485 VelocityPressureField->advanceSolution(tStep);
486 solutionVector = VelocityPressureField->giveSolutionVector(tStep);
487 prevSolutionVector = VelocityPressureField->giveSolutionVector( tStep->givePreviousStep() );
488 }
489
490 if ( initFlag ) {
492 //previousAccelerationVector.resize(neq);
493 accelerationVector.resize(neq);
494 solutionVector->resize(neq);
495 prevSolutionVector->resize(neq);
496 }
497
498 incrementalSolutionVector.resize(neq);
499
500 lhs = classFactory.createSparseMtrx(sparseMtrxType);
501 if ( !lhs ) {
502 OOFEM_ERROR("sparse matrix creation failed");
503 }
504
505 lhs->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
506
507 if ( materialInterface ) {
509 }
510
511 initFlag = 0;
512 } else if ( requiresUnknownsDictionaryUpdate() ) {
513 // rebuild lhs structure and resize solution vector
514 incrementalSolutionVector.resize(neq);
515 lhs->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
516 }
517
518
520 *solutionVector = *prevSolutionVector;
521 }
522
523 //previousAccelerationVector=accelerationVector;
524
525 // evaluate element supg and sppg stabilization coeffs
527
528 //
529 // predictor
530 //
531
534 } else {
535 this->updateSolutionVectors_predictor(* solutionVector, accelerationVector, tStep);
536 }
537
538 if ( tStep->giveNumber() != 1 ) {
539 if ( materialInterface ) {
540 //if (this->fsflag) updateDofManVals(tStep);
541#ifdef SUPG_IMPLICIT_INTERFACE
542 #ifdef TIME_REPORT
543 Timer timer;
544 timer.startTimer();
545 #endif
546 materialInterface->updatePosition( this->giveCurrentStep() );
548 #ifdef TIME_REPORT
549 timer.stopTimer();
550 OOFEM_LOG_INFO("SUPG info: user time consumed by updating interfaces: %.2fs\n", timer.getUtime();
551 );
552 #endif
553#else
554 //updateElementsForNewInterfacePosition (tStep);
555#endif
556 //if (this->fsflag) this->updateDofManActivityMap(tStep);
557 }
558 }
559
560 // evaluate element supg and sppg stabilization coeffs
561 // this -> evaluateElementStabilizationCoeffs (tStep);
562
563 //
564 // assemble rhs (residual)
565 //
566 externalForces.zero();
567 this->assembleVector( externalForces, tStep, ExternalForceAssembler(), VM_Total,
570
571#if 0
575 double loadLevel;
576 int currentIterations=0;
577 this->updateInternalRHS( this->internalForces, tStep, this->giveDomain(1), &this->eNorm );
578 ConvergedReason status = this->nMethod->solve(*this->lhs,
579 externalForces,
580 NULL,
581 solutionVector,
583 this->internalForces,
584 this->eNorm,
585 loadLevel, // Only relevant for incrementalBCLoadVector?
586 SparseNonLinearSystemNM :: rlm_total,
587 currentIterations,
588 tStep);
589
590 if ( status != CR_CONVERGED) {
591 OOFEM_ERROR("No success in solving problem at time step", tStep->giveNumber());
592 }
593
594#else
595 int nite = 0;
596 double _absErrResid, err, avn = 0.0, aivn = 0.0, rnorm;
597 int jj, nnodes = this->giveDomain(1)->giveNumberOfDofManagers();
598 FloatArray rhs, iForces(neq);
599
600
601 // algoritmic rhs part (assembled by e-model (in giveCharComponent service) from various element contribs)
602 iForces.zero();
603 this->assembleVector( iForces, tStep, SUPGInternalForceAssembler(lscale, dscale, uscale), VM_Total,
606
607 rhs.beDifferenceOf(externalForces, iForces);
608
609 //
610 // corrector
611 //
612 //OOFEM_LOG_INFO("Iteration IncrSolErr RelResidErr AbsResidErr\n_______________________________________________________\n");
613 OOFEM_LOG_INFO("Iteration IncrSolErr RelResidErr (Rel_MB ,Rel_MC ) AbsResidErr\n__________________________________________________________________________________\n");
614 do {
615 nite++;
616 //
617 // Assemble lhs
618 //
619 //if (nite == 1)
620 {
621 // momentum balance part
622 lhs->zero();
623 if ( 1 ) { //if ((nite > 5)) // && (rnorm < 1.e4))
624 this->assemble( *lhs, tStep, SUPGTangentAssembler(TangentStiffness, lscale, dscale, uscale, alpha),
626 } else {
627 this->assemble( *lhs, tStep, SUPGTangentAssembler(SecantStiffness, lscale, dscale, uscale, alpha),
629 }
630 }
631 //if (this->fsflag) this->imposeAmbientPressureInOuterNodes(lhs,&rhs,tStep);
632
633#if 1
635#else
636 {
637 auto lhs_copy = lhs->clone();
638
640
641 // check solver
642 FloatArray __rhs(neq);
643 double __absErr = 0., __relErr = 0.;
644 lhs_copy->times(incrementalSolutionVector, __rhs);
645 for ( int i = 1; i <= neq; i++ ) {
646 if ( fabs( __rhs.at(i) - rhs.at(i) ) > __absErr ) {
647 __absErr = fabs( __rhs.at(i) - rhs.at(i) );
648 }
649
650 if ( fabs( ( __rhs.at(i) - rhs.at(i) ) / rhs.at(i) ) > __relErr ) {
651 __relErr = fabs( ( __rhs.at(i) - rhs.at(i) ) / rhs.at(i) );
652 }
653 }
654
655 OOFEM_LOG_INFO("SUPG: solver check: absoluteError %e, relativeError %e\n", __absErr, __relErr);
656 }
657#endif
658
659
660
663 } else {
664 //update
666 avn = accelerationVector.computeSquaredNorm();
667 aivn = incrementalSolutionVector.computeSquaredNorm();
668 } // end update
669
670#if 0
671 #ifdef SUPG_IMPLICIT_INTERFACE
672 if ( materialInterface ) {
673 #ifdef TIME_REPORT
674 Timer timer;
675 timer.startTimer();
676 #endif
677 //if (this->fsflag) updateDofManVals(tStep);
678 materialInterface->updatePosition( this->giveCurrentStep() );
680 #ifdef TIME_REPORT
681 timer.stopTimer();
682 OOFEM_LOG_INFO( "SUPG info: user time consumed by updating interfaces: %.2fs\n", timer.getUtime() );
683 #endif
684 //if (this->fsflag) this->updateDofManActivityMap(tStep);
685 }
686
687 #endif
688#else
689 if ( materialInterface ) {
690 //if (this->fsflag) updateDofManVals(tStep);
691 //if (this->fsflag) this->updateDofManActivityMap(tStep);
692 }
693
694#endif
695
696 // check convergence and repeat iteration if desired
697 if ( avn < 1.e-12 ) {
698 err = aivn;
699 } else {
700 err = sqrt(aivn / avn);
701 }
702
703 //
704 // assemble rhs (residual)
705 //
706 iForces.zero();
707 this->assembleVector( iForces, tStep, SUPGInternalForceAssembler(lscale, dscale, uscale), VM_Total,
710 rhs.beDifferenceOf(externalForces, iForces);
711
712 // check convergence and repeat iteration if desired
713 double rnorm_mb = 0.0, rnorm_mc = 0.0;
714 for ( int i = 1; i <= nnodes; i++ ) {
715 DofManager *inode = this->giveDomain(1)->giveDofManager(i);
716 for ( Dof *dof: *inode ) {
717 if ( ( jj = dof->__giveEquationNumber() ) ) {
718 DofIDItem type = dof->giveDofID();
719 double val = rhs.at(jj);
720 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) {
721 rnorm_mb += val * val;
722 } else {
723 rnorm_mc += val * val;
724 }
725 }
726 }
727 }
728
729 rnorm_mb = sqrt(rnorm_mb);
730 rnorm_mc = sqrt(rnorm_mc);
731
732 rnorm = rhs.computeNorm();
733
734 _absErrResid = 0.0;
735 for ( int i = 1; i <= neq; i++ ) {
736 _absErrResid = max( _absErrResid, fabs( rhs.at(i) ) );
737 }
738
739 //if ( requiresUnknownsDictionaryUpdate() ) {
740 // OOFEM_LOG_INFO("%-10d n/a %-15e %-15e\n", nite, rnorm, _absErrResid);
741 //} else {
742 OOFEM_LOG_INFO("%-10d %-15e %-15e (%-10e,%-10e) %-15e\n", nite, err, rnorm, rnorm_mb, rnorm_mc, _absErrResid);
743 //}
744
745 if ( 0 ) {
746 // evaluate element supg and sppg stabilization coeffs
748 //
749 // assemble rhs (residual)
750 //
751 iForces.zero();
752 this->assembleVector( iForces, tStep, SUPGInternalForceAssembler(lscale, dscale, uscale), VM_Total,
755 rhs.beDifferenceOf(externalForces, iForces);
756 }
757 } while ( ( rnorm > rtolv ) && ( _absErrResid > atolv ) && ( nite <= maxiter ) );
758
759 if ( nite <= maxiter ) {
760 OOFEM_LOG_INFO("SUPG info: number of iterations: %d\n", nite);
761 status = CR_CONVERGED;
762 } else {
763 OOFEM_WARNING("Convergence not reached, number of iterations: %d\n", nite);
764 if ( stopmaxiter ) {
765 exit(1);
766 }
767 status = CR_DIVERGED_ITS;
768 }
769#endif
770
771#ifndef SUPG_IMPLICIT_INTERFACE
772 if ( materialInterface ) {
773 #ifdef TIME_REPORT
774 Timer _timer;
775 _timer.startTimer();
776 #endif
777 materialInterface->updatePosition( this->giveCurrentStep() );
779 #ifdef TIME_REPORT
780 _timer.stopTimer();
781 OOFEM_LOG_INFO( "SUPG info: user time consumed by updating interfaces: %.2fs\n", _timer.getUtime() );
782 #endif
783 //if (this->fsflag) this->updateDofManActivityMap(tStep);
784 }
785
786#endif
787
788 // update solution state counter
789 tStep->incrementStateCounter();
790 tStep->convergedReason = status;
791}
792
793
794void
795SUPG :: updateYourself(TimeStep *tStep)
796{
797 this->updateInternalState(tStep);
798 EngngModel :: updateYourself(tStep);
799 if ( materialInterface ) {
800 materialInterface->updateYourself(tStep);
801 }
802
803 //previousSolutionVector = solutionVector;
804}
805
806
807void
808SUPG :: updateInternalState(TimeStep *tStep)
809{
810 for ( auto &domain: domainList ) {
812 for ( auto &dman : domain->giveDofManagers() ) {
813 this->updateDofUnknownsDictionary(dman.get(), tStep);
814 }
815 }
816
817 for ( auto &elem : domain->giveElements() ) {
818 elem->updateInternalState(tStep);
819 }
820 }
821}
822
823
824void
825SUPG :: saveContext(DataStream &stream, ContextMode mode)
826{
827 EngngModel :: saveContext(stream, mode);
828
829 VelocityPressureField->saveContext(stream);
830
832 if ( ( iores = accelerationVector.storeYourself(stream) ) != CIO_OK ) {
833 THROW_CIOERR(iores);
834 }
835
836 if ( materialInterface ) {
837 materialInterface->saveContext(stream, mode);
838 }
839}
840
841
842void
843SUPG :: restoreContext(DataStream &stream, ContextMode mode)
844{
845 EngngModel :: restoreContext(stream, mode);
846
847 VelocityPressureField->restoreContext(stream);
848
850 if ( ( iores = accelerationVector.restoreYourself(stream) ) != CIO_OK ) {
851 THROW_CIOERR(iores);
852 }
853
854 if ( materialInterface ) {
855 materialInterface->restoreContext(stream, mode);
856 }
857}
858
859
860int
861SUPG :: checkConsistency()
862{
863 // check internal consistency
864 // if success returns nonzero
865 Domain *domain = this->giveDomain(1);
866
867 // check for proper element type
868 for ( auto &elem : domain->giveElements() ) {
869 if ( dynamic_cast< SUPGElement * >( elem.get() ) == NULL ) {
870 OOFEM_WARNING("Element %d has no SUPG base", elem->giveLabel());
871 return 0;
872 }
873 }
874
875 int ret = EngngModel :: checkConsistency();
876 if ( ret == 0 ) {
877 return 0;
878 }
879
880
881 // scale boundary and initial conditions
882 if ( equationScalingFlag ) {
883 for ( auto &bc : domain->giveBcs() ) {
884 if ( bc->giveBCValType() == VelocityBVT ) {
885 bc->scale(1. / uscale);
886 } else if ( bc->giveBCValType() == PressureBVT ) {
887 bc->scale( 1. / this->giveVariableScale(VST_Pressure) );
888 } else if ( bc->giveBCValType() == ForceLoadBVT ) {
889 bc->scale( 1. / this->giveVariableScale(VST_Force) );
890 } else {
891 OOFEM_ERROR("unknown bc/ic type");
892 }
893 }
894
895 for ( auto &ic : domain->giveIcs() ) {
896 if ( ic->giveICValType() == VelocityBVT ) {
897 ic->scale(VM_Total, 1. / uscale);
898 } else if ( ic->giveICValType() == PressureBVT ) {
899 ic->scale( VM_Total, 1. / this->giveVariableScale(VST_Pressure) );
900 } else {
901 OOFEM_ERROR("unknown bc/ic type");
902 }
903 }
904 }
905
906 return 1;
907}
908
909
910void
911SUPG :: updateDomainLinks()
912{
913 EngngModel :: updateDomainLinks();
914 this->giveNumericalMethod( this->giveCurrentMetaStep() )->reinitialize();
915 //this->giveNumericalMethod( this->giveCurrentMetaStep() )->setDomain( this->giveDomain(1) );
916}
917
918void
919SUPG :: printDofOutputAt(FILE *stream, Dof *iDof, TimeStep *tStep)
920{
921 double pscale = ( dscale * uscale * uscale );
922
923 DofIDItem type = iDof->giveDofID();
924 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) {
925 iDof->printSingleOutputAt(stream, tStep, 'd', VM_Total, uscale);
926 } else if ( type == P_f ) {
927 iDof->printSingleOutputAt(stream, tStep, 'd', VM_Total, pscale);
928 } else {
929 OOFEM_ERROR("unsupported dof type");
930 }
931}
932
933void
934SUPG :: applyIC(TimeStep *_stepWhenIcApply)
935{
936 Domain *domain = this->giveDomain(1);
938 FloatArray *vp_vector = NULL;
939
940#ifdef VERBOSE
941 OOFEM_LOG_INFO("Applying initial conditions\n");
942#endif
943 int nman = domain->giveNumberOfDofManagers();
944
945 accelerationVector.resize(neq);
946 accelerationVector.zero();
947
949 VelocityPressureField->advanceSolution(_stepWhenIcApply);
950 vp_vector = VelocityPressureField->giveSolutionVector(_stepWhenIcApply);
951 vp_vector->resize(neq);
952 vp_vector->zero();
953
954 for ( int j = 1; j <= nman; j++ ) {
955 DofManager *node = domain->giveDofManager(j);
956
957 for ( Dof *dof: *node ) {
958 // ask for initial values obtained from
959 // bc (boundary conditions) and ic (initial conditions)
960 if ( !dof->isPrimaryDof() ) {
961 continue;
962 }
963
964 int jj = dof->__giveEquationNumber();
965 DofIDItem type = dof->giveDofID();
966 if ( jj ) {
967 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) {
968 vp_vector->at(jj) = dof->giveUnknown(VM_Total, _stepWhenIcApply);
969 } else {
970 vp_vector->at(jj) = dof->giveUnknown(VM_Total, _stepWhenIcApply);
971 }
972 }
973 }
974 }
975 }
976
977 //ICs and BCs are scaled already in CheckConsistency
978 //vp_vector->times(1./this->giveVariableScale(VST_Velocity));
979
980 // update element state according to given ic
981
982 //this->initElementsForNewStep (stepWhenIcApply);
983 if ( materialInterface ) {
984 this->updateElementsForNewInterfacePosition(_stepWhenIcApply);
985 }
986
987 for ( auto &elem : domain->giveElements() ) {
988 SUPGElement *element = static_cast< SUPGElement * >( elem.get() );
989 element->updateInternalState(_stepWhenIcApply);
990 element->updateYourself(_stepWhenIcApply);
991 }
992}
993
994double
995SUPG :: giveVariableScale(VarScaleType varID)
996{
997 if ( varID == VST_Length ) {
998 return this->lscale;
999 } else if ( varID == VST_Velocity ) {
1000 return this->uscale;
1001 } else if ( varID == VST_Density ) {
1002 return this->dscale;
1003 } else if ( varID == VST_Time ) {
1004 return ( lscale / uscale );
1005 } else if ( varID == VST_Pressure ) {
1006 return ( dscale * uscale * uscale );
1007 } else if ( varID == VST_Force ) {
1008 return ( uscale * uscale / lscale );
1009 } else if ( varID == VST_Viscosity ) {
1010 return 1.0;
1011 } else {
1012 OOFEM_ERROR("unknown variable type");
1013 }
1014
1015 //return 0.0;
1016}
1017
1018void
1019SUPG :: evaluateElementStabilizationCoeffs(TimeStep *tStep)
1020{
1021 Domain *domain = this->giveDomain(1);
1022
1023 for ( auto &elem : domain->giveElements() ) {
1024 FMElement *ePtr = static_cast< FMElement * >( elem.get() );
1025 ePtr->updateStabilizationCoeffs(tStep);
1026 }
1027}
1028
1029void
1030SUPG :: updateElementsForNewInterfacePosition(TimeStep *tStep)
1031{
1032 Domain *domain = this->giveDomain(1);
1033
1034 OOFEM_LOG_DEBUG("SUPG :: updateElements - updating elements for interface position");
1035
1036
1037 for ( auto &elem : domain->giveElements() ) {
1038 SUPGElement *ePtr = static_cast< SUPGElement * >( elem.get() );
1040 }
1041}
1042
1043
1044void
1045SUPG :: updateDofUnknownsDictionary(DofManager *inode, TimeStep *tStep)
1046{
1047 // update DOF unknowns dictionary, where
1048 // unknowns are hold instead of keeping them in global unknowns
1049 // vectors in engng instancies
1050 // this is necessary, because during solution equation numbers for
1051 // particular DOFs may changed, and it is necesary to keep them
1052 // in DOF level.
1053
1054 // empty -> all done in updateDofUnknownsDictionary_corrector (TimeStep* tStep);
1055}
1056
1057
1058void
1059SUPG :: updateDofUnknownsDictionary_predictor(TimeStep *tStep)
1060{
1061 double dt = tStep->giveTimeIncrement();
1062 Domain *domain = this->giveDomain(1);
1063
1064 int nnodes = domain->giveNumberOfDofManagers();
1065
1067 for ( int j = 1; j <= nnodes; j++ ) {
1068 DofManager *inode = domain->giveDofManager(j);
1069 for ( Dof *dof: *inode ) {
1070 double val, prev_val, accel;
1071 DofIDItem type = dof->giveDofID();
1072 if ( !dof->hasBc(tStep) ) {
1073 prev_val = dof->giveUnknown( VM_Total, tStep->givePreviousStep() );
1074 accel = dof->giveUnknown( VM_Acceleration, tStep->givePreviousStep() );
1075 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) {
1076 val = prev_val + dt * accel;
1077 } else {
1078 val = prev_val;
1079 }
1080
1081 dof->updateUnknownsDictionary(tStep, VM_Total, val); // velocity
1082 dof->updateUnknownsDictionary(tStep, VM_Acceleration, accel); // acceleration
1083 } else {
1084 val = dof->giveBcValue(VM_Total, tStep);
1085 dof->updateUnknownsDictionary(tStep, VM_Total, val); // velocity
1086 dof->updateUnknownsDictionary(tStep, VM_Acceleration, 0.0); // acceleration
1087 }
1088 }
1089 }
1090 }
1091}
1092
1093
1094void
1095SUPG :: updateDofUnknownsDictionary_corrector(TimeStep *tStep)
1096{
1097 double dt = tStep->giveTimeIncrement();
1098 Domain *domain = this->giveDomain(1);
1099 int nnodes = domain->giveNumberOfDofManagers();
1100
1102 for ( int j = 1; j <= nnodes; j++ ) {
1103 DofManager *inode = domain->giveDofManager(j);
1104 for ( Dof *iDof: *inode ) {
1105 DofIDItem type = iDof->giveDofID();
1106 if ( !iDof->hasBc(tStep) ) {
1107 double val, prev_val;
1108 prev_val = iDof->giveUnknown(VM_Total, tStep);
1109 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) {
1110 val = prev_val + dt *alpha *incrementalSolutionVector.at( iDof->__giveEquationNumber() );
1111 } else {
1112 val = prev_val + incrementalSolutionVector.at( iDof->__giveEquationNumber() );
1113 }
1114
1115 iDof->updateUnknownsDictionary(tStep, VM_Total, val); // velocity
1116
1117 prev_val = iDof->giveUnknown(VM_Acceleration, tStep);
1118 val = prev_val + incrementalSolutionVector.at( iDof->__giveEquationNumber() );
1119 iDof->updateUnknownsDictionary(tStep, VM_Acceleration, val); // acceleration
1120 }
1121 }
1122 }
1123 }
1124}
1125
1126int
1127SUPG :: giveUnknownDictHashIndx(ValueModeType mode, TimeStep *tStep)
1128{
1129 if ( ( tStep == this->giveCurrentStep() ) || ( tStep == this->givePreviousStep() ) ) {
1130 return ( tStep->giveNumber() % 2 ) * 100 + mode;
1131 } else {
1132 OOFEM_ERROR("unsupported solution step");
1133 }
1134}
1135
1136void
1137SUPG :: updateSolutionVectors_predictor(FloatArray &solutionVector, FloatArray &_accelerationVector, TimeStep *tStep)
1138{
1139 double dt = tStep->giveTimeIncrement();
1140 Domain *domain = this->giveDomain(1);
1141
1142 for ( auto &node : domain->giveDofManagers() ) {
1143 for ( Dof *iDof: *node ) {
1144 if ( !iDof->isPrimaryDof() ) {
1145 continue;
1146 }
1147
1148 int jj = iDof->__giveEquationNumber();
1149 DofIDItem type = iDof->giveDofID();
1150
1151 if ( jj ) {
1152 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) { // v = v + (dt)*a
1153 solutionVector.at(jj) += dt * _accelerationVector.at(jj);
1154 }
1155 }
1156 }
1157 }
1158
1159 for ( auto &elem : domain->giveElements() ) {
1160 int ndofman = elem->giveNumberOfInternalDofManagers();
1161 for ( int ji = 1; ji <= ndofman; ji++ ) {
1162 DofManager *dofman = elem->giveInternalDofManager(ji);
1163 for ( Dof *iDof: *dofman ) {
1164 if ( !iDof->isPrimaryDof() ) {
1165 continue;
1166 }
1167
1168 int jj = iDof->__giveEquationNumber();
1169 DofIDItem type = iDof->giveDofID();
1170
1171 if ( jj ) {
1172 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) { // v = v + (dt*alpha)*da
1173 solutionVector.at(jj) += dt * _accelerationVector.at(jj);
1174 }
1175 }
1176 }
1177 } // end loop over elem internal dofmans
1178 } // end loop over elems
1179}
1180
1181
1182void
1183SUPG :: updateSolutionVectors(FloatArray &solutionVector, FloatArray &_accelerationVector, FloatArray &incrementalSolutionVector, TimeStep *tStep)
1184{
1185 double dt = tStep->giveTimeIncrement();
1186 Domain *domain = this->giveDomain(1);
1187
1188 _accelerationVector.add(incrementalSolutionVector);
1189
1190 for ( auto &node : domain->giveDofManagers() ) {
1191
1192 for ( Dof *iDof: *node ) {
1193 if ( !iDof->isPrimaryDof() ) {
1194 continue;
1195 }
1196
1197 int jj = iDof->__giveEquationNumber();
1198 DofIDItem type = iDof->giveDofID();
1199
1200 if ( jj ) {
1201 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) { // v = v + (dt*alpha)*da
1202 solutionVector.at(jj) += dt * alpha * incrementalSolutionVector.at(jj);
1203 } else {
1204 solutionVector.at(jj) += incrementalSolutionVector.at(jj); // p = p + dp
1205 }
1206 }
1207 }
1208 } // end loop over dnam
1209
1210 for ( auto &elem : domain->giveElements() ) {
1211 int ndofman = elem->giveNumberOfInternalDofManagers();
1212 for ( int ji = 1; ji <= ndofman; ji++ ) {
1213 DofManager *dofman = elem->giveInternalDofManager(ji);
1214 for ( Dof *iDof: *dofman ) {
1215 if ( !iDof->isPrimaryDof() ) {
1216 continue;
1217 }
1218
1219 int jj = iDof->__giveEquationNumber();
1220 DofIDItem type = iDof->giveDofID();
1221
1222 if ( jj ) {
1223 if ( ( type == V_u ) || ( type == V_v ) || ( type == V_w ) ) { // v = v + (dt*alpha)*da
1224 solutionVector.at(jj) += dt * alpha * incrementalSolutionVector.at(jj);
1225 } else {
1226 solutionVector.at(jj) += incrementalSolutionVector.at(jj); // p = p + dp
1227 }
1228 }
1229 }
1230 } // end loop over elem internal dofmans
1231 } // end loop over elems
1232}
1233
1234
1235#define __VOF_TRESHOLD 1.e-5
1236
1237#define __NODE_OUT 0
1238#define __NODE_IN 1
1239#define __NODE_OUT_ACTIVE 2
1240#define __NODE_INTERPOL_CANDIDATE 3
1241
1242/*
1243 * void
1244 * SUPG::updateDofManActivityMap (TimeStep* tStep)
1245 * {
1246 * // loop over inactive nodes and test if some of its shared element has become active (temp_vof>0)
1247 * // then initialize its velocity and pressure
1248 *
1249 * Domain* domain = this->giveDomain(1);
1250 * Element* ielem;
1251 * int nnodes = domain->giveNumberOfDofManagers();
1252 * int nelems = domain->giveNumberOfElements();
1253 * int i, inode, ie;
1254 * //int inodes, _code;
1255 * LEPlicElementInterface *interface;
1256 * IntArray dofMask(2);
1257 * FloatArray vv;
1258 *
1259 * dofMask.at(1) = V_u;
1260 * dofMask.at(2) = V_v;
1261 *
1262 * printf ("SUPG::updateDofManActivityMap called\n");
1263 *
1264 * IntArray __old_DofManActivityMask (__DofManActivityMask);
1265 * __DofManActivityMask.resize(nnodes);
1266 * for (i=1; i<=nnodes;i++) __DofManActivityMask.at(i) = __NODE_OUT;
1267 *
1268 *
1269 * for (ie=1; ie<= nelems; ie++) {
1270 * ielem = domain->giveElement(ie);
1271 * if ((interface = (LEPlicElementInterface*) (ielem->giveInterface(LEPlicElementInterfaceType)))) {
1272 * if (interface->giveVolumeFraction() > 0.9999) {
1273 * // mark all its nodes as active
1274 * nnodes=ielem->giveNumberOfNodes();
1275 * for (i=1; i<=nnodes;i++) this->__DofManActivityMask.at(ielem->giveNode(i)->giveNumber())=__NODE_IN;
1276 * } else if (interface->giveVolumeFraction() > 0.0) {
1277 * nnodes=ielem->giveNumberOfNodes();
1278 * for (i=1; i<=nnodes;i++) {
1279 * inode = ielem->giveNode(i)->giveNumber();
1280 * //if (__DofManActivityMask.at(inode) == __NODE_OUT) {
1281 * FloatArray normal; double p,x,y;
1282 * interface->giveTempInterfaceNormal(normal);
1283 * p = interface->giveTempLineConstant();
1284 *
1285 * x=domain->giveNode(inode)->giveCoordinate(1);
1286 * y=domain->giveNode(inode)->giveCoordinate(2);
1287 *
1288 * if ((normal.at(1)*x+normal.at(2)*y+p) >= 0.) {
1289 * __DofManActivityMask.at(inode) = __NODE_IN;
1290 * } else {
1291 * if (__DofManActivityMask.at(inode) == __NODE_OUT) {
1292 * if (interface->giveVolumeFraction() < 0.5) __DofManActivityMask.at(inode) = __NODE_INTERPOL_CANDIDATE;
1293 * else __DofManActivityMask.at(inode) = __NODE_OUT_ACTIVE;
1294 * }
1295 * }
1296 * //}
1297 * }
1298 * }
1299 * }
1300 * }
1301 *
1302 * }
1303 *
1304 *
1305 * void
1306 * SUPG:: updateDofManVals (TimeStep* tStep)
1307 * {
1308 * Domain* domain = this->giveDomain(1);
1309 * int nnodes = domain->giveNumberOfDofManagers();
1310 * int inode;
1311 * const IntArray* ct;
1312 * IntArray dofMask(2);
1313 * FloatArray vv;
1314 * Dof* iDof;
1315 *
1316 * dofMask.at(1) = V_u;
1317 * dofMask.at(2) = V_v;
1318 *
1319 * printf ("SUPG::updateDofManVals called\n");
1320 * for (inode=1; inode <= nnodes; inode++) {
1321 * if (__DofManActivityMask.at(inode) == __NODE_INTERPOL_CANDIDATE) {
1322 * // loop over shared elements
1323 * ct = domain->giveConnectivityTable()->giveDofManConnectivityArray(inode);
1324 * // now we have to find active dofmans in the neighborhood to interpolate velocity
1325 * double vx=0.,vy=0.;
1326 * int _j, _k, _d, _cnt=0;
1327 * Element* _je;
1328 * for (_j=1; _j<=ct->giveSize(); _j++) {
1329 * _je=domain->giveElement(ct->at(_j));
1330 * for (_k=1; _k<=_je->giveNumberOfNodes(); _k++) {
1331 * int _code = __DofManActivityMask.at(_je->giveNode(_k)->giveNumber());
1332 * if ((_code == __NODE_IN) || (_code == __NODE_OUT_ACTIVE)) { // active node (not candidate)
1333 * _je->giveNode(_k)->giveUnknownVector(vv, dofMask, VM_Total, tStep);
1334 * vx += vv.at(1); vy += vv.at(2);
1335 * _cnt++;
1336 * }
1337 * }
1338 * }
1339 * if (_cnt) printf ("%d ", inode);
1340 * vx /= _cnt; vy /= _cnt; // average
1341 * for (Dof *iDof: domain->giveDofManager(inode)) {
1342 * DofIDItem type = iDof->giveDofID();
1343 * if (!iDof->hasBc(tStep)) {
1344 * if (type == V_u)
1345 * iDof->updateUnknownsDictionary (tStep, VM_Total, vx); // velocity
1346 * else if (type == V_v)
1347 * iDof->updateUnknownsDictionary (tStep, VM_Total, vy); // velocity
1348 * else if (type == V_w) {
1349 * OOFEM_ERROR("3d not yet supported");
1350 * } else if (type == P_f) {
1351 * //iDof->updateUnknownsDictionary (tStep, VM_Total, 0.0);
1352 * } else {_error2 ("unknown DOF type encountered (%s)", __DofIDItemToString (type));}
1353 * }
1354 * }
1355 * }
1356 * }
1357 * }
1358 *
1359 * void
1360 * SUPG::imposeAmbientPressureInOuterNodes(SparseMtrx* lhs, FloatArray* rhs, TimeStep* tStep)
1361 * {
1362 * Domain* domain = this->giveDomain(1);
1363 * int nnodes = domain->giveNumberOfDofManagers();
1364 * int inode, _code;
1365 * auto pdof;
1366 *
1367 * IntArray inMask(nnodes); inMask.zero();
1368 * for (inode=1;inode<=nnodes;inode++) {
1369 * _code = __DofManActivityMask.at(inode);
1370 * if ((_code == __NODE_OUT_ACTIVE) || (_code == __NODE_INTERPOL_CANDIDATE) || (_code == __NODE_OUT)) {
1371 * if ((pdof=domain->giveDofManager(inode)->findDofWithDofId(P_f)) != domain->giveDofManager(inode)->end() ) {
1372 * int _eq = (*pdof)->giveEquationNumber();
1373 * if (_eq) {
1374 * if (fabs(lhs->at(_eq,_eq)) > 0.0) {
1375 * lhs->at(_eq,_eq) *= 1.e6;
1376 * rhs->at(_eq) = 0.0;
1377 * } else {
1378 * lhs->at(_eq,_eq) = 1.0;
1379 * rhs->at(_eq) = 0.0;
1380 * }
1381 * }
1382 * }
1383 * }
1384 * }
1385 * }
1386 *
1387 * void
1388 * SUPG::__debug (TimeStep* tStep)
1389 * {
1390 * int i, in, id, nincr = 1000;
1391 * int neq = this -> giveNumberOfDomainEquations (1, EModelDefaultEquationNumbering());
1392 * IntArray loc;
1393 * SUPGElement* element = (SUPGElement*) giveDomain(1)->giveElement(1);
1394 * FloatArray vincr(6), F(6), fi(6), fprev(6), fapprox(6), *solutionVector;
1395 * FloatMatrix d(6,6);
1396 *
1397 * vincr.at(1) = 1.e-2;
1398 * vincr.at(2) = 1.e-7;
1399 * vincr.at(3) = 2.e-2;
1400 * element -> giveLocationArray (loc);
1401 * VelocityPressureField->advanceSolution(tStep);
1402 * solutionVector = VelocityPressureField->giveSolutionVector(tStep);
1403 * solutionVector->resize(neq);
1404 * // restrict vincr
1405 * for (i=1; i<=6; i++)
1406 * if (loc.at(i) <= 0) vincr.at(i) = 0.0;
1407 *
1408 *
1409 * for (in=1; in<=nincr; in++) {
1410 * for (i=1; i<=6; i++)
1411 * if ((id = loc.at(i))) solutionVector->at(id) += vincr.at(i);
1412 * // compute F from derivative
1413 * EngngModel::giveElementCharacteristicMatrix(d, 1, TangentDiffusionDerivativeTerm_MB, tStep, giveDomain(1));
1414 * fi.beProductOf (d, vincr);
1415 * fapprox = fprev;
1416 * fapprox.add(fi);
1417 * // compute true F
1418 * this->giveElementCharacteristicVector(F, 1, DiffusionTerm_MB, VM_Total, tStep, giveDomain(1));
1419 * // print
1420 * printf ("%d %e %e %e %e %e %e %e %e %e %e %e %e\n", in, F.at(1), F.at(2), F.at(3), F.at(4), F.at(5), F.at(6),
1421 * fapprox.at(1), fapprox.at(2), fapprox.at(3), fapprox.at(4), fapprox.at(5), fapprox.at(6));
1422 * // if (in==1) fapprox = F;
1423 * fprev = fapprox;
1424 * }
1425 *
1426 * }
1427 *
1428 */
1429} // end namespace oofem
#define REGISTER_EngngModel(class)
DofIDItem giveDofID() const
Definition dof.h:276
virtual int __giveEquationNumber() const =0
virtual void giveUnknowns(FloatArray &masterUnknowns, ValueModeType mode, TimeStep *tStep)
Definition dof.C:159
virtual void printSingleOutputAt(FILE *file, TimeStep *tStep, char ch, ValueModeType mode, double scale=1.0)
Definition dof.C:76
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
Function * giveFunction(int n)
Definition domain.C:271
std ::vector< std ::unique_ptr< GeneralBoundaryCondition > > & giveBcs()
Definition domain.h:349
DofManager * giveDofManager(int n)
Definition domain.C:317
std ::vector< std ::unique_ptr< Element > > & giveElements()
Definition domain.h:294
std ::vector< std ::unique_ptr< InitialCondition > > & giveIcs()
Definition domain.h:356
virtual void updateYourself(TimeStep *tStep)
Definition element.C:824
virtual int computeNumberOfDofs()
Definition element.h:436
int giveNumberOfTimeStepWhenIcApply()
Returns the time step number, when initial conditions should apply.
Definition engngm.h:787
std ::vector< std ::unique_ptr< Domain > > domainList
List of problem domains.
Definition engngm.h:217
virtual TimeStep * giveCurrentStep(bool force=false)
Definition engngm.h:717
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Definition engngm.C:452
EngngModelContext * giveContext()
Context requesting service.
Definition engngm.h:1174
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
Domain * giveDomain(int n)
Definition engngm.C:1936
std ::unique_ptr< TimeStep > currentStep
Current time step.
Definition engngm.h:241
EngngModelTimer timer
E-model timer.
Definition engngm.h:279
@ InternalForcesExchangeTag
Definition engngm.h:321
virtual TimeStep * givePreviousStep(bool force=false)
Definition engngm.h:738
EngngModel * master
Master e-model; if defined receiver is in maintained (slave) mode.
Definition engngm.h:274
void initMetaStepAttributes(MetaStep *mStep)
Definition engngm.h:677
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
void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
Definition fmelement.C:51
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Definition fmelement.h:63
void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition fmelement.C:44
void registerField(FieldPtr eField, FieldType key)
double computeNorm() const
Definition floatarray.C:861
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 beDifferenceOf(const FloatArray &a, const FloatArray &b)
Definition floatarray.C:403
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void add(const FloatArray &src)
Definition floatarray.C:218
void times(double f)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void zero()
Zeroes all coefficient of receiver.
bool isNotEmpty() const
Tests for empty matrix.
void assemble(const FloatMatrix &src, const IntArray &loc)
FluidModel(int i, EngngModel *master)
Definition fluidmodel.h:49
virtual double evaluateAtTime(double t)
Definition function.C:78
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
virtual void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)=0
virtual void giveLocalPressureDofMap(IntArray &map)
virtual void computeBCLhsPressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
virtual void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep)=0
virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep)=0
virtual void updateElementForNewInterfacePosition(TimeStep *tStep)
Definition supgelement.h:79
virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep)=0
virtual void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeBCLhsTerm_MB(FloatMatrix &answer, TimeStep *tStep)
virtual void computeBCLhsPressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep)=0
virtual void computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void giveLocalVelocityDofMap(IntArray &map)
virtual double computeCriticalTimeStep(TimeStep *tStep)=0
Computes the critical time increment.
virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep)=0
void updateInternalState(TimeStep *tStep) override
MatResponseMode rmode
Definition supg.h:91
int maxiter
Max number of iterations.
Definition supg.h:128
TimeStep * giveSolutionStepWhenIcApply(bool force=false) override
Definition supg.C:412
int deltaTF
Definition supg.h:124
bool equationScalingFlag
Definition supg.h:140
void updateInternalState(TimeStep *tStep)
Definition supg.C:808
std ::unique_ptr< SparseLinearSystemNM > nMethod
Numerical method used to solve the problem.
Definition supg.h:108
int requiresUnknownsDictionaryUpdate() override
Definition supg.h:186
double uscale
Definition supg.h:142
std ::unique_ptr< SparseMtrx > lhs
Definition supg.h:113
void updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm) override
Definition supg.C:378
SparseMtrxType sparseMtrxType
Definition supg.h:111
double Re
Reynolds number.
Definition supg.h:144
int giveUnknownDictHashIndx(ValueModeType mode, TimeStep *tStep) override
Definition supg.C:1127
std ::unique_ptr< MaterialInterface > materialInterface
Definition supg.h:147
int consistentMassFlag
Definition supg.h:138
std ::unique_ptr< PrimaryField > VelocityPressureField
Definition supg.h:114
double lscale
Definition supg.h:142
FloatArray eNorm
Definition supg.h:120
void updateSolutionVectors(FloatArray &solutionVector, FloatArray &accelerationVector, FloatArray &incrementalSolutionVector, TimeStep *tStep)
Definition supg.C:1183
void updateElementsForNewInterfacePosition(TimeStep *tStep)
Definition supg.C:1030
LinSystSolverType solverType
Definition supg.h:110
double alpha
Integration constant.
Definition supg.h:135
bool stopmaxiter
Definition supg.h:133
FloatArray incrementalSolutionVector
Definition supg.h:117
FloatArray accelerationVector
Definition supg.h:116
void updateDofUnknownsDictionary(DofManager *dman, TimeStep *tStep) override
Definition supg.C:1045
double deltaT
Definition supg.h:123
NumericalMethod * giveNumericalMethod(MetaStep *mStep) override
Returns reference to receiver's numerical method.
Definition supg.C:216
FloatArray internalForces
Definition supg.h:119
double atolv
Convergence tolerance.
Definition supg.h:126
void updateDofUnknownsDictionary_predictor(TimeStep *tStep)
Definition supg.C:1059
int initFlag
Definition supg.h:137
void evaluateElementStabilizationCoeffs(TimeStep *tStep)
Definition supg.C:1019
double giveVariableScale(VarScaleType varId) override
Returns the scale factor for given variable type.
Definition supg.C:995
double dscale
Definition supg.h:142
double rtolv
Definition supg.h:126
void updateSolutionVectors_predictor(FloatArray &solutionVector, FloatArray &accelerationVector, TimeStep *tStep)
Definition supg.C:1137
void updateDofUnknownsDictionary_corrector(TimeStep *tStep)
Definition supg.C:1095
void applyIC(TimeStep *tStep)
Definition supg.C:934
virtual void zero()=0
Zeroes the receiver.
void incrementStateCounter()
Updates solution state counter.
Definition timestep.h:213
ConvergedReason convergedReason
Status of solution step (Converged,.
Definition timestep.h:120
double giveTimeIncrement()
Returns solution step associated time increment.
Definition timestep.h:168
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
double getUtime()
Returns total user time elapsed in seconds.
Definition timer.C:105
void startTimer()
Definition timer.C:69
void stopTimer()
Definition timer.C:77
#define THROW_CIOERR(e)
#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_DEBUG(...)
Definition logger.h:144
long ContextMode
Definition contextmode.h:43
FloatArrayF< N > min(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
VarScaleType
Type determining the scale corresponding to particular variable.
@ VST_Length
@ VST_Viscosity
@ VST_Velocity
@ VST_Force
@ VST_Density
@ VST_Pressure
FloatArrayF< N > assemble(const FloatArrayF< M > &x, int const (&c)[M])
Assemble components into zero matrix.
FloatArrayF< N > max(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ CR_DIVERGED_ITS
@ ForceLoadBVT
Definition bcvaltype.h:43
@ VelocityBVT
Definition bcvaltype.h:46
@ PressureBVT
Definition bcvaltype.h:44
ClassFactory & classFactory
@ NonLinearLhs
@ InternalRhs
#define _IFT_SUPG_uscale
Definition supg.h:56
#define _IFT_SUPG_stopmaxiter
Definition supg.h:62
#define _IFT_SUPG_rtolv
Definition supg.h:59
#define _IFT_SUPG_alpha
Definition supg.h:53
#define _IFT_SUPG_cmflag
Definition supg.h:52
#define _IFT_SUPG_atolv
Definition supg.h:60
#define _IFT_SUPG_deltatFunction
Definition supg.h:51
#define _IFT_SUPG_dscale
Definition supg.h:57
#define _IFT_SUPG_miflag
Definition supg.h:58
#define _IFT_SUPG_maxiter
Definition supg.h:61
#define _IFT_SUPG_deltat
Definition supg.h:50
#define _IFT_SUPG_lscale
Definition supg.h:55
#define _IFT_SUPG_scaleflag
Definition supg.h:54

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