OOFEM 3.0
Loading...
Searching...
No Matches
stokesflow.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 "stokesflow.h"
36#include "fmelement.h"
37#include "inputrecord.h"
38#include "timestep.h"
39#include "classfactory.h"
40#include "domain.h"
41#include "nrsolver.h"
44#include "topologydescription.h"
45#include "parallelcontext.h"
46#include "exportmodulemanager.h"
49
50namespace oofem {
52
53StokesFlow :: StokesFlow(int i, EngngModel *_master) : FluidModel(i, _master)
54{
55 this->ndomains = 1;
56}
57
58
59StokesFlow :: ~StokesFlow() { }
60
61
62void StokesFlow :: initializeFrom(InputRecord &ir)
63{
64 FluidModel :: initializeFrom(ir);
65 int val;
66
67 val = ( int ) SMT_PetscMtrx;
69 this->sparseMtrxType = ( SparseMtrxType ) val;
70
71 val = ( int ) ST_Petsc;
73 this->solverType = ( LinSystSolverType ) val;
74
75 this->deltaT = 1.0;
77
78 this->velocityPressureField = std::make_unique<DofDistributedPrimaryField>(this, 1, FT_VelocityPressure, 1);
79 this->stiffnessMatrix = nullptr;
80 this->meshqualityee = nullptr;
81
82 this->ts = TS_OK;
83
84 this->maxdef = 25;
85}
86
87
88void StokesFlow :: solveYourselfAt(TimeStep *tStep)
89{
90 Domain *d = this->giveDomain(1);
91 FloatArray externalForces;
92 FloatArray incrementOfSolution;
93
94 if ( d->giveNumberOfElements() == 0 && d->giveTopology() ) {
96 this->meshqualityee = std::make_unique<MeshQualityErrorEstimator>( 1, d );
97 }
98
99 if ( d->giveTopology() && this->meshqualityee ) {
100 // Check the quality of the deformed mesh.
101 double meshdeformation = this->meshqualityee->giveValue(globalErrorEEV, tStep);
102 if ( this->giveProblemScale() == macroScale ) {
103 OOFEM_LOG_INFO("StokesFlow :: solveYourselfAt - Mesh deformation at %e\n", meshdeformation);
104 }
105
106 if ( this->ts == TS_NeedsRemeshing || meshdeformation > this->maxdef ) {
108 OOFEM_LOG_INFO( "StokesFlow :: solveYourselfAt - New mesh created (%d elements).\n", d->giveNumberOfElements() );
109 /*meshdeformation =*/ this->meshqualityee->giveValue(globalErrorEEV, tStep);
110 this->giveExportModuleManager()->initialize();
111 }
112 }
113
115
116 velocityPressureField->advanceSolution(tStep);
118
119 // Create "stiffness matrix"
120 if ( !this->stiffnessMatrix ) {
121 this->stiffnessMatrix = classFactory.createSparseMtrx(sparseMtrxType);
122 if ( !this->stiffnessMatrix ) {
123 OOFEM_ERROR("Couldn't create requested sparse matrix of type %d", sparseMtrxType);
124 }
125
126 this->stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
127 }
128
129 incrementOfSolution.resize(neq);
130 this->internalForces.resize(neq);
131
132 // Build initial/external load (LoadVector)
133 externalForces.resize(neq);
134 externalForces.zero();
135 this->assembleVector( externalForces, tStep, ExternalForceAssembler(), VM_Total,
138
139 if ( this->giveProblemScale() == macroScale ) {
140 OOFEM_LOG_INFO("StokesFlow :: solveYourselfAt - Solving step %d, metastep %d, (neq = %d)\n", tStep->giveNumber(), tStep->giveMetaStepNumber(), neq);
141 }
142
145#if 1
146 double loadLevel;
147 int currentIterations=0;
148 this->updateInternalRHS( this->internalForces, tStep, d, &this->eNorm );
149 ConvergedReason status = this->nMethod->solve(*this->stiffnessMatrix,
150 externalForces,
151 NULL,
153 incrementOfSolution,
154 this->internalForces,
155 this->eNorm,
156 loadLevel, // Only relevant for incrementalBCLoadVector?
157 SparseNonLinearSystemNM :: rlm_total,
158 currentIterations,
159 tStep);
160#else
161 this->updateInternalRHS( this->internalForces, tStep, d, nullptr );
162 this->updateMatrix( this->stiffnessMatrix, tStep, d );
163 this->internalForces.negated();
164 this->internalForces.add(externalForces);
165 ConvergedReason status = this->nMethod->giveLinearSolver()->solve(this->stiffnessMatrix.get(), & ( this->internalForces ), & solutionVector);
167 this->updateInternalRHS( this->internalForces, tStep, d, nullptr );
168#endif
169
170 if (status != CR_CONVERGED ) {
171 OOFEM_ERROR("No success in solving problem at time step %d", tStep->giveNumber());
172 }
173 tStep->numberOfIterations = currentIterations;
174 tStep->convergedReason = status;
175}
176
177void StokesFlow :: updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d)
178{
180
181 // update element stabilization
182 for ( auto &elem : d->giveElements() ) {
183 static_cast< FMElement * >( elem.get() )->updateStabilizationCoeffs(tStep);
184 }
185
186 if ( cmpn == InternalRhs ) {
187 this->internalForces.zero();
188 this->assembleVector(this->internalForces, tStep, InternalForceAssembler(), VM_Total,
191 return;
192 } else if ( cmpn == NonLinearLhs ) {
193 this->stiffnessMatrix->zero();
194 this->assemble(*stiffnessMatrix, tStep, TangentAssembler(TangentStiffness),
196 return;
197 } else {
198 OOFEM_ERROR("Unknown component");
199 }
200}
201
202void StokesFlow :: updateSolution(FloatArray &solutionVector, TimeStep *tStep, Domain *d)
203{
205 // update element stabilization
206 for ( auto &elem : d->giveElements() ) {
207 static_cast< FMElement * >( elem.get() )->updateStabilizationCoeffs(tStep);
208 }
209}
210
211
212void StokesFlow :: updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *enorm)
213{
214 answer.zero();
215 this->assembleVector(answer, tStep, InternalForceAssembler(), VM_Total,
218}
219
220
221void StokesFlow :: updateMatrix(SparseMtrx &mat, TimeStep *tStep, Domain *d)
222{
223 mat.zero();
224 this->assemble(mat, tStep, TangentAssembler(TangentStiffness),
226}
227
228
229void StokesFlow :: updateYourself(TimeStep *tStep)
230{
231 this->updateInternalState(tStep);
232 FluidModel :: updateYourself(tStep);
233}
234
235
236int StokesFlow :: forceEquationNumbering(int id)
237{
238 int neq = FluidModel :: forceEquationNumbering(id);
239 this->equationNumberingCompleted = false;
240 this->stiffnessMatrix = nullptr;
241 return neq;
242}
243
244
245double StokesFlow :: giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof)
246{
247 return velocityPressureField->giveUnknownValue(dof, mode, tStep);
248}
249
250
251double StokesFlow :: giveReynoldsNumber()
252{
253 return 1.0;
254}
255
256
257int StokesFlow :: checkConsistency()
258{
259 Domain *domain = this->giveDomain(1);
260
261 // check for proper element type
262 for ( auto &elem : domain->giveElements() ) {
263 if ( dynamic_cast< FMElement * >( elem.get() ) == NULL ) {
264 OOFEM_WARNING("Element %d has no FMElement base", elem->giveLabel());
265 return false;
266 }
267 }
268
269 return EngngModel :: checkConsistency();
270}
271
272void StokesFlow :: updateInternalState(TimeStep *tStep)
273{
274 for ( auto &domain: domainList ) {
275 if ( domain->giveTopology() ) {
276 // Must be done before updating nodal positions
277 this->ts = domain->giveTopology()->updateYourself(tStep);
278 }
279
280 for ( auto &elem : domain->giveElements() ) {
281 elem->updateInternalState(tStep);
282 }
283 }
284}
285
286void StokesFlow :: doStepOutput(TimeStep *tStep)
287{
288 TopologyDescription *tp = this->giveDomain(1)->giveTopology();
289 if ( tp ) {
290 tp->doOutput(tStep);
291 }
292
293 FluidModel :: doStepOutput(tStep);
294}
295
296NumericalMethod *StokesFlow :: giveNumericalMethod(MetaStep *mStep)
297{
298 if ( !this->nMethod ) {
299 this->nMethod = std::make_unique<NRSolver>(this->giveDomain(1), this);
300 }
301 return this->nMethod.get();
302}
303
304TimeStep *StokesFlow :: giveNextStep()
305{
306 if ( !currentStep ) {
307 // first step -> generate initial step
308 //currentStep = std::make_unique<TimeStep>(*giveSolutionStepWhenIcApply());
309 currentStep = std::make_unique<TimeStep>(giveNumberOfTimeStepWhenIcApply(), this, 1, 0., this->deltaT, 0);
310 }
311 previousStep = std :: move(currentStep);
312 currentStep = std::make_unique<TimeStep>(*previousStep, this->deltaT);
313
314 return currentStep.get();
315}
316} // end namespace oofem
#define REGISTER_EngngModel(class)
TopologyDescription * giveTopology()
Definition domain.C:1551
int giveNumberOfElements() const
Returns number of elements in domain.
Definition domain.h:463
std ::vector< std ::unique_ptr< Element > > & giveElements()
Definition domain.h:294
int giveNumberOfTimeStepWhenIcApply()
Returns the time step number, when initial conditions should apply.
Definition engngm.h:787
problemScale giveProblemScale() const
Returns scale in multiscale simulation.
Definition engngm.h:447
std ::vector< std ::unique_ptr< Domain > > domainList
List of problem domains.
Definition engngm.h:217
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Definition engngm.C:452
ExportModuleManager * giveExportModuleManager()
Returns receiver's export module manager.
Definition engngm.h:791
std ::unique_ptr< TimeStep > previousStep
Previous time step.
Definition engngm.h:243
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
@ InternalForcesExchangeTag
Definition engngm.h:321
void initMetaStepAttributes(MetaStep *mStep)
Definition engngm.h:677
int updateSharedDofManagers(FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
Definition engngm.C:2162
void assembleVector(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
Definition engngm.C:1106
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Definition fmelement.h:63
void resize(Index s)
Definition floatarray.C:94
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
FluidModel(int i, EngngModel *master)
Definition fluidmodel.h:49
virtual void zero()=0
Zeroes the receiver.
FloatArray eNorm
Element norm for nonlinear analysis (squared).
Definition stokesflow.h:77
void updateSolution(FloatArray &solutionVector, TimeStep *tStep, Domain *d) override
Definition stokesflow.C:202
TopologyState ts
Topology state, most notably used for determining if there is a need to remesh.
Definition stokesflow.h:89
void updateMatrix(SparseMtrx &mat, TimeStep *tStep, Domain *d) override
Definition stokesflow.C:221
std ::unique_ptr< SparseMtrx > stiffnessMatrix
Definition stokesflow.h:84
double deltaT
Time increment read from input record.
Definition stokesflow.h:67
LinSystSolverType solverType
Linear solver type.
Definition stokesflow.h:75
SparseMtrxType sparseMtrxType
Sparse matrix type.
Definition stokesflow.h:71
void updateInternalState(TimeStep *tStep)
Definition stokesflow.C:272
std ::unique_ptr< MeshQualityErrorEstimator > meshqualityee
Used for determining if a new mesh must be created.
Definition stokesflow.h:80
void updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm) override
Definition stokesflow.C:212
NumericalMethod * giveNumericalMethod(MetaStep *mStep) override
Returns reference to receiver's numerical method.
Definition stokesflow.C:296
FloatArray internalForces
Definition stokesflow.h:86
double maxdef
Maximum deformation allowed.
Definition stokesflow.h:82
std ::unique_ptr< SparseNonLinearSystemNM > nMethod
Numerical method.
Definition stokesflow.h:73
std ::unique_ptr< PrimaryField > velocityPressureField
Primary unknowns.
Definition stokesflow.h:69
FloatArray solutionVector
Definition stokesflow.h:85
int giveMetaStepNumber()
Returns receiver's meta step number.
Definition timestep.h:150
ConvergedReason convergedReason
Status of solution step (Converged,.
Definition timestep.h:120
int numberOfIterations
Number of itarations needed to achieve convergence.
Definition timestep.h:116
int giveNumber()
Returns receiver's number.
Definition timestep.h:144
virtual void doOutput(TimeStep *tStep)
#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 OOFEM_LOG_INFO(...)
Definition logger.h:143
FloatArrayF< N > assemble(const FloatArrayF< M > &x, int const (&c)[M])
Assemble components into zero matrix.
@ SMT_PetscMtrx
PETSc library mtrx representation.
ClassFactory & classFactory
@ globalErrorEEV
@ TS_NeedsRemeshing
Indicates that the topology has reached a need for remeshing, as the case with merging surfaces.
@ TS_OK
Indicates that everything is OK with respect to topology.
@ NonLinearLhs
@ InternalRhs
@ macroScale
Definition problemmode.h:46
#define _IFT_StokesFlow_deltat
Definition stokesflow.h:49

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