OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
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 - 2013 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"
50 
51 
52 namespace oofem {
53 class SparseNonLinearSystemNM;
54 class MeshQualityErrorEstimator;
55 
63 class StokesFlow : public FluidModel
64 {
65 protected:
67  double deltaT;
69  std :: unique_ptr< PrimaryField > velocityPressureField;
73  std :: unique_ptr< SparseNonLinearSystemNM > nMethod;
78 
80  std :: unique_ptr< MeshQualityErrorEstimator > meshqualityee;
82  double maxdef;
83 
84  std :: unique_ptr< SparseMtrx > stiffnessMatrix;
87 
90 
91 public:
92  StokesFlow(int i, EngngModel * _master = NULL);
93  virtual ~StokesFlow();
94 
95  virtual void solveYourselfAt(TimeStep *tStep);
96 
102  virtual void updateYourself(TimeStep *tStep);
103 
104  virtual double giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *domain, Dof *dof);
105 
106  virtual double giveReynoldsNumber();
107 
108  virtual int forceEquationNumbering(int id);
109 
111 
112  virtual int checkConsistency();
113  virtual void doStepOutput(TimeStep *tStep);
114  void updateInternalState(TimeStep *tStep);
115  virtual void updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d);
117  virtual TimeStep *giveNextStep();
118 
119  virtual const char *giveClassName() const { return "StokesFlow"; }
120  virtual const char *giveInputRecordName() const { return _IFT_StokesFlow_Name; }
121 };
122 } // end namespace oofem
123 
124 #endif // stokesflow_h
LinSystSolverType
The values of this type should be related not to specific solvers, but more to specific packages that...
virtual double giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *domain, Dof *dof)
Returns requested unknown.
Definition: stokesflow.C:220
virtual ~StokesFlow()
Definition: stokesflow.C:58
Class and object Domain.
Definition: domain.h:115
Class representing meta step.
Definition: metastep.h:62
virtual double giveReynoldsNumber()
Definition: stokesflow.C:226
virtual const char * giveInputRecordName() const
Definition: stokesflow.h:120
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description in input reader.
Definition: stokesflow.C:62
Base class for fluid problems.
Definition: fluidmodel.h:46
This base class is an abstraction for numerical algorithm.
Definition: nummet.h:80
virtual TimeStep * giveNextStep()
Returns next time step (next to current step) of receiver.
Definition: stokesflow.C:279
std::unique_ptr< SparseMtrx > stiffnessMatrix
Definition: stokesflow.h:84
SparseMtrxType sparseMtrxType
Sparse matrix type.
Definition: stokesflow.h:71
double maxdef
Maximum deformation allowed.
Definition: stokesflow.h:82
virtual void updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d)
Updates components mapped to numerical method if necessary during solution process.
Definition: stokesflow.C:180
StokesFlow(int i, EngngModel *_master=NULL)
Definition: stokesflow.C:53
LinSystSolverType solverType
Linear solver type.
Definition: stokesflow.h:75
NumericalCmpn
Type representing numerical component.
Definition: numericalcmpn.h:46
virtual const char * giveClassName() const
Returns class name of the receiver.
Definition: stokesflow.h:119
virtual int checkConsistency()
Allows programmer to test some receiver&#39;s internal data, before computation begins.
Definition: stokesflow.C:232
FloatArray eNorm
Element norm for nonlinear analysis (squared)
Definition: stokesflow.h:77
SparseMtrxType
Enumerative type used to identify the sparse matrix type.
void updateInternalState(TimeStep *tStep)
Definition: stokesflow.C:247
double deltaT
Time increment read from input record.
Definition: stokesflow.h:67
virtual int forceEquationNumbering()
Forces equation renumbering on all domains associated to engng model.
Definition: fluidmodel.h:52
virtual NumericalMethod * giveNumericalMethod(MetaStep *mStep)
Returns reference to receiver&#39;s numerical method.
Definition: stokesflow.C:271
virtual void doStepOutput(TimeStep *tStep)
Prints the ouput of the solution step (using virtual this->printOutputAtservice) to the stream detemi...
Definition: stokesflow.C:261
Implements the engineering model to solve incompressible Stokes flow.
Definition: stokesflow.h:63
std::unique_ptr< SparseNonLinearSystemNM > nMethod
Numerical method.
Definition: stokesflow.h:73
TopologyState ts
Topology state, most notably used for determining if there is a need to remesh.
Definition: stokesflow.h:89
Class representing vector of real numbers.
Definition: floatarray.h:82
std::unique_ptr< PrimaryField > velocityPressureField
Primary unknowns.
Definition: stokesflow.h:69
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
FloatArray solutionVector
Definition: stokesflow.h:85
Class representing the general Input Record.
Definition: inputrecord.h:101
FloatArray internalForces
Definition: stokesflow.h:86
std::unique_ptr< MeshQualityErrorEstimator > meshqualityee
Used for determining if a new mesh must be created.
Definition: stokesflow.h:80
Abstract base class representing the "problem" under consideration.
Definition: engngm.h:181
#define _IFT_StokesFlow_Name
Definition: stokesflow.h:48
the oofem namespace is to define a context or scope in which all oofem names are defined.
Abstract class Dof represents Degree Of Freedom in finite element mesh.
Definition: dof.h:93
virtual void solveYourselfAt(TimeStep *tStep)
Solves problem for given time step.
Definition: stokesflow.C:90
virtual void updateYourself(TimeStep *tStep)
Updates everything for the problem.
Definition: stokesflow.C:205
TopologyState
Determines the state of the evolving topology.
Class representing solution step.
Definition: timestep.h:80

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:31 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011