OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
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 - 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 #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"
42 #include "sparsenonlinsystemnm.h"
44 #include "topologydescription.h"
45 #include "parallelcontext.h"
46 #include "exportmodulemanager.h"
48 #include "unknownnumberingscheme.h"
49 
50 namespace oofem {
51 REGISTER_EngngModel(StokesFlow);
52 
53 StokesFlow :: StokesFlow(int i, EngngModel *_master) : FluidModel(i, _master)
54 {
55  this->ndomains = 1;
56 }
57 
59 {
60 }
61 
63 {
64  IRResultType result;
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.reset( new DofDistributedPrimaryField(this, 1, FT_VelocityPressure, 1) );
79  this->stiffnessMatrix.reset( NULL );
80  this->meshqualityee.reset( NULL );
81 
82  this->ts = TS_OK;
83 
84  this->maxdef = 25;
85 
87 }
88 
89 
91 {
92  Domain *d = this->giveDomain(1);
93  FloatArray externalForces;
94  FloatArray incrementOfSolution;
95 
96  if ( d->giveNumberOfElements() == 0 && d->giveTopology() ) {
98  this->meshqualityee.reset( new MeshQualityErrorEstimator( 1, d ) );
99  }
100 
101  if ( d->giveTopology() && this->meshqualityee ) {
102  // Check the quality of the deformed mesh.
103  double meshdeformation = this->meshqualityee->giveValue(globalErrorEEV, tStep);
104  if ( this->giveProblemScale() == macroScale ) {
105  OOFEM_LOG_INFO("StokesFlow :: solveYourselfAt - Mesh deformation at %e\n", meshdeformation);
106  }
107 
108  if ( this->ts == TS_NeedsRemeshing || meshdeformation > this->maxdef ) {
109  d->giveTopology()->replaceFEMesh();
110  OOFEM_LOG_INFO( "StokesFlow :: solveYourselfAt - New mesh created (%d elements).\n", d->giveNumberOfElements() );
111  /*meshdeformation =*/ this->meshqualityee->giveValue(globalErrorEEV, tStep);
113  }
114  }
115 
117 
118  // Move solution space to current time step
119  velocityPressureField->advanceSolution(tStep);
120 
121  // Point pointer SolutionVector to current solution in velocityPressureField
122  velocityPressureField->initialize(VM_Total, tStep, solutionVector, EModelDefaultEquationNumbering() );
123 
124  // Create "stiffness matrix"
125  if ( !this->stiffnessMatrix ) {
127  if ( !this->stiffnessMatrix ) {
128  OOFEM_ERROR("Couldn't create requested sparse matrix of type %d", sparseMtrxType);
129  }
130 
131  this->stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() );
132  }
133 
134  incrementOfSolution.resize(neq);
135  this->internalForces.resize(neq);
136 
137  // Build initial/external load (LoadVector)
138  externalForces.resize(neq);
139  externalForces.zero();
140  this->assembleVector( externalForces, tStep, ExternalForceAssembler(), VM_Total,
143 
144  if ( this->giveProblemScale() == macroScale ) {
145  OOFEM_LOG_INFO("StokesFlow :: solveYourselfAt - Solving step %d, metastep %d, (neq = %d)\n", tStep->giveNumber(), tStep->giveMetaStepNumber(), neq);
146  }
147 
148  this->giveNumericalMethod( this->giveCurrentMetaStep() );
150 #if 1
151  double loadLevel;
152  int currentIterations;
153  this->updateComponent( tStep, InternalRhs, d );
154  NM_Status status = this->nMethod->solve(*this->stiffnessMatrix,
155  externalForces,
156  NULL,
158  incrementOfSolution,
159  this->internalForces,
160 
161  this->eNorm,
162  loadLevel, // Only relevant for incrementalBCLoadVector?
164  currentIterations,
165  tStep);
166 #else
167  this->updateComponent( tStep, InternalRhs, d );
168  this->updateComponent( tStep, NonLinearLhs, d );
169  this->internalForces.negated();
170  this->internalForces.add(externalForces);
171  NM_Status status = this->nMethod->giveLinearSolver()->solve(this->stiffnessMatrix.get(), & ( this->internalForces ), & solutionVector);
172  this->updateComponent( tStep, NonLinearLhs, d );
173 #endif
174 
175  if ( !( status & NM_Success ) ) {
176  OOFEM_ERROR("No success in solving problem at time step", tStep->giveNumber());
177  }
178 }
179 
181 {
183 
184  // update element stabilization
185  for ( auto &elem : d->giveElements() ) {
186  static_cast< FMElement * >( elem.get() )->updateStabilizationCoeffs(tStep);
187  }
188 
189  if ( cmpn == InternalRhs ) {
190  this->internalForces.zero();
191  this->assembleVector(this->internalForces, tStep, InternalForceAssembler(), VM_Total,
192  EModelDefaultEquationNumbering(), d, & this->eNorm);
194  return;
195  } else if ( cmpn == NonLinearLhs ) {
196  this->stiffnessMatrix->zero();
197  this->assemble(*stiffnessMatrix, tStep, TangentAssembler(TangentStiffness),
199  return;
200  } else {
201  OOFEM_ERROR("Unknown component");
202  }
203 }
204 
206 {
207  this->updateInternalState(tStep);
209 }
210 
212 {
214  this->equationNumberingCompleted = false;
215  this->stiffnessMatrix.reset( NULL );
216  return neq;
217 }
218 
219 
221 {
222  return velocityPressureField->giveUnknownValue(dof, mode, tStep);
223 }
224 
225 
227 {
228  return 1.0;
229 }
230 
231 
233 {
234  Domain *domain = this->giveDomain(1);
235 
236  // check for proper element type
237  for ( auto &elem : domain->giveElements() ) {
238  if ( dynamic_cast< FMElement * >( elem.get() ) == NULL ) {
239  OOFEM_WARNING("Element %d has no FMElement base", elem->giveLabel());
240  return false;
241  }
242  }
243 
245 }
246 
248 {
249  for ( auto &domain: domainList ) {
250  if ( domain->giveTopology() ) {
251  // Must be done before updating nodal positions
252  this->ts = domain->giveTopology()->updateYourself(tStep);
253  }
254 
255  for ( auto &elem : domain->giveElements() ) {
256  elem->updateInternalState(tStep);
257  }
258  }
259 }
260 
262 {
263  TopologyDescription *tp = this->giveDomain(1)->giveTopology();
264  if ( tp ) {
265  tp->doOutput(tStep);
266  }
267 
269 }
270 
272 {
273  if ( !this->nMethod ) {
274  this->nMethod.reset( new NRSolver(this->giveDomain(1), this) );
275  }
276  return this->nMethod.get();
277 }
278 
280 {
281  if ( !currentStep ) {
282  // first step -> generate initial step
283  //currentStep.reset( new TimeStep(*giveSolutionStepWhenIcApply()) );
284  currentStep.reset( new TimeStep(giveNumberOfTimeStepWhenIcApply(), this, 1, 0., this->deltaT, 0) );
285  }
286  previousStep = std :: move(currentStep);
287  currentStep.reset( new TimeStep(*previousStep, this->deltaT) );
288 
289  return currentStep.get();
290 }
291 } // end namespace oofem
LinSystSolverType
The values of this type should be related not to specific solvers, but more to specific packages that...
The representation of EngngModel default unknown numbering.
#define NM_Success
Numerical method exited with success.
Definition: nmstatus.h:47
virtual double giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *domain, Dof *dof)
Returns requested unknown.
Definition: stokesflow.C:220
std::unique_ptr< TimeStep > currentStep
Current time step.
Definition: engngm.h:231
virtual ~StokesFlow()
Definition: stokesflow.C:58
void initMetaStepAttributes(MetaStep *mStep)
Update e-model attributes attributes according to step metaStep attributes.
Definition: engngm.C:580
Class and object Domain.
Definition: domain.h:115
Implementation for assembling internal forces vectors in standard monolithic, nonlinear FE-problems...
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Returns number of equations for given domain in active (current time step) time step.
Definition: engngm.C:391
problemScale giveProblemScale()
Returns scale in multiscale simulation.
Definition: engngm.h:418
TopologyDescription * giveTopology()
Returns receiver&#39;s associated topology description.
Definition: domain.C:1443
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description in input reader.
Definition: engngm.C:262
Class representing meta step.
Definition: metastep.h:62
This class implements Newton-Raphson Method, derived from abstract NumericalMethod class for solving ...
Definition: nrsolver.h:95
SparseMtrx * createSparseMtrx(SparseMtrxType type)
Creates new instance of sparse matrix corresponding to given keyword.
Definition: classfactory.C:105
virtual void doOutput(TimeStep *tStep)
File output of the current state of the topology description.
std::vector< std::unique_ptr< Domain > > domainList
List of problem domains.
Definition: engngm.h:207
virtual void doStepOutput(TimeStep *tStep)
Prints the ouput of the solution step (using virtual this->printOutputAtservice) to the stream detemi...
Definition: engngm.C:670
virtual double giveReynoldsNumber()
Definition: stokesflow.C:226
std::unique_ptr< TimeStep > previousStep
Previous time step.
Definition: engngm.h:233
ExportModuleManager * giveExportModuleManager()
Returns receiver&#39;s export module manager.
Definition: engngm.h:758
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
Indicates that everything is OK with respect to topology.
Implementation for assembling tangent matrices in standard monolithic FE-problems.
std::unique_ptr< SparseMtrx > stiffnessMatrix
Definition: stokesflow.h:84
#define _IFT_StokesFlow_deltat
Definition: stokesflow.h:49
int giveNumberOfElements() const
Returns number of elements in domain.
Definition: domain.h:434
Indicates that the topology has reached a need for remeshing, as the case with merging surfaces...
SparseMtrxType sparseMtrxType
Sparse matrix type.
Definition: stokesflow.h:71
REGISTER_EngngModel(ProblemSequence)
unsigned long NM_Status
Mask defining NumMetod Status; which can be asked after finishing computation by Numerical Method...
Definition: nmstatus.h:44
double maxdef
Maximum deformation allowed.
Definition: stokesflow.h:82
Class representing field of primary variables, which are typically allocated on nodes.
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain)
Assembles characteristic matrix of required type into given sparse matrix.
Definition: engngm.C:803
MetaStep * giveCurrentMetaStep()
Returns current meta step.
Definition: engngm.C:1684
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
int giveNumber()
Returns receiver&#39;s number.
Definition: timestep.h:129
#define OOFEM_LOG_INFO(...)
Definition: logger.h:127
LinSystSolverType solverType
Linear solver type.
Definition: stokesflow.h:75
NumericalCmpn
Type representing numerical component.
Definition: numericalcmpn.h:46
#define _IFT_EngngModel_smtype
Definition: engngm.h:82
#define OOFEM_ERROR(...)
Definition: error.h:61
int ndomains
Number of receiver domains.
Definition: engngm.h:205
virtual int checkConsistency()
Allows programmer to test some receiver&#39;s internal data, before computation begins.
Definition: stokesflow.C:232
int updateSharedDofManagers(FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
Exchanges necessary remote DofManagers data.
Definition: engngm.C:1957
FloatArray eNorm
Element norm for nonlinear analysis (squared)
Definition: stokesflow.h:77
SparseMtrxType
Enumerative type used to identify the sparse matrix type.
The reference incremental load vector is defined as totalLoadVector assembled at given time...
Implementation for assembling external forces vectors in standard monolithic FE-problems.
virtual void updateYourself(TimeStep *tStep)
Updates internal state after finishing time step.
Definition: engngm.C:612
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
PETSc library mtrx representation.
This error estimator measures the quality of the elements.
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
virtual int checkConsistency()
Allows programmer to test some receiver&#39;s internal data, before computation begins.
Definition: engngm.h:995
This abstract class represent a general base element class for fluid dynamic problems.
Definition: fmelement.h:54
std::unique_ptr< PrimaryField > velocityPressureField
Primary unknowns.
Definition: stokesflow.h:69
#define _IFT_EngngModel_lstype
Definition: engngm.h:81
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
int giveMetaStepNumber()
Returns receiver&#39;s meta step number.
Definition: timestep.h:135
int giveNumberOfTimeStepWhenIcApply()
Returns the time step number, when initial conditions should apply.
Definition: engngm.h:754
FloatArray solutionVector
Definition: stokesflow.h:85
Class representing the general Input Record.
Definition: inputrecord.h:101
FloatArray internalForces
Definition: stokesflow.h:86
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
Abstract class for topology description.
ClassFactory & classFactory
Definition: classfactory.C:59
std::unique_ptr< MeshQualityErrorEstimator > meshqualityee
Used for determining if a new mesh must be created.
Definition: stokesflow.h:80
void initialize()
Initializes output manager.
std::vector< std::unique_ptr< Element > > & giveElements()
Definition: domain.h:279
Abstract base class representing the "problem" under consideration.
Definition: engngm.h:181
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
the oofem namespace is to define a context or scope in which all oofem names are defined.
Domain * giveDomain(int n)
Service for accessing particular problem domain.
Definition: engngm.C:1720
Abstract class Dof represents Degree Of Freedom in finite element mesh.
Definition: dof.h:93
void negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
virtual void solveYourselfAt(TimeStep *tStep)
Solves problem for given time step.
Definition: stokesflow.C:90
virtual void replaceFEMesh()
Generates the FE components from the bare mesh.
virtual void updateYourself(TimeStep *tStep)
Updates everything for the problem.
Definition: stokesflow.C:205
#define OOFEM_WARNING(...)
Definition: error.h:62
void assembleVector(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
Assembles characteristic vector of required type from dofManagers, element, and active boundary condi...
Definition: engngm.C:986
Class representing solution step.
Definition: timestep.h:80
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
int equationNumberingCompleted
Equation numbering completed flag.
Definition: engngm.h:223

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