OOFEM 3.0
Loading...
Searching...
No Matches
stokesflow.h
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#ifndef stokesflow_h
36#define stokesflow_h
37
38#include "fluidmodel.h"
39#include "sparsemtrxtype.h"
40#include "topologydescription.h"
41#include "linsystsolvertype.h"
42#include "floatmatrix.h"
43#include "primaryfield.h"
44#include "floatarray.h"
45
47
48#define _IFT_StokesFlow_Name "stokesflow"
49#define _IFT_StokesFlow_deltat "deltat"
51
52namespace oofem {
53class SparseNonLinearSystemNM;
54class MeshQualityErrorEstimator;
55
63class StokesFlow : public FluidModel
64{
65protected:
67 double deltaT = 1.;
69 std :: unique_ptr< PrimaryField > velocityPressureField;
73 std :: unique_ptr< SparseNonLinearSystemNM > nMethod;
78
80 std :: unique_ptr< MeshQualityErrorEstimator > meshqualityee;
82 double maxdef = 1.;
83
84 std :: unique_ptr< SparseMtrx > stiffnessMatrix;
87
90
91public:
92 StokesFlow(int i, EngngModel * _master = nullptr);
94
95 void solveYourselfAt(TimeStep *tStep) override;
96
102 void updateYourself(TimeStep *tStep) override;
103
104 double giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *domain, Dof *dof) override;
105 bool newDofHandling() override { return true; }
106
107 double giveReynoldsNumber() override;
108
109 int forceEquationNumbering(int id) override;
110
111 void initializeFrom(InputRecord &ir) override;
112
113 int checkConsistency() override;
114 void doStepOutput(TimeStep *tStep) override;
115 void updateInternalState(TimeStep *tStep);
116 void updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d) override;
117 void updateSolution(FloatArray &solutionVector, TimeStep *tStep, Domain *d) override;
118 void updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm) override;
119 void updateMatrix(SparseMtrx &mat, TimeStep *tStep, Domain *d) override;
121 TimeStep *giveNextStep() override;
122
123 const char *giveClassName() const override { return "StokesFlow"; }
124 const char *giveInputRecordName() const { return _IFT_StokesFlow_Name; }
125};
126} // end namespace oofem
127
128#endif // stokesflow_h
int forceEquationNumbering() override
Definition fluidmodel.h:52
FluidModel(int i, EngngModel *master)
Definition fluidmodel.h:49
void initializeFrom(InputRecord &ir) override
Definition stokesflow.C:62
FloatArray eNorm
Element norm for nonlinear analysis (squared).
Definition stokesflow.h:77
void solveYourselfAt(TimeStep *tStep) override
Definition stokesflow.C:88
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
bool newDofHandling() override
Definition stokesflow.h:105
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
StokesFlow(int i, EngngModel *_master=nullptr)
Definition stokesflow.C:53
int checkConsistency() override
Definition stokesflow.C:257
const char * giveInputRecordName() const
Definition stokesflow.h:124
void updateYourself(TimeStep *tStep) override
Definition stokesflow.C:229
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
double giveReynoldsNumber() override
Definition stokesflow.C:251
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
const char * giveClassName() const override
Returns class name of the receiver.
Definition stokesflow.h:123
double giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *domain, Dof *dof) override
Definition stokesflow.C:245
void updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d) override
Definition stokesflow.C:177
void doStepOutput(TimeStep *tStep) override
Definition stokesflow.C:286
TimeStep * giveNextStep() override
Returns next time step (next to current step) of receiver.
Definition stokesflow.C:304
FloatArray solutionVector
Definition stokesflow.h:85
#define _IFT_StokesFlow_Name
Definition stokesflow.h:48

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