OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
quad10_2d_supg.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 quad10_2d_supg_h
36 #define quad10_2d_supg_h
37 
38 #include "supgelement2.h"
39 #include "spatiallocalizer.h"
40 #include "zznodalrecoverymodel.h"
42 #include "leplic.h"
43 #include "levelsetpcs.h"
44 #include "elementinternaldofman.h"
45 
46 #define _IFT_Quad10_2D_SUPG_Name "quad1supg"
47 
48 namespace oofem {
49 class FEI2dQuadLin;
50 class FEI2dQuadConst;
51 
58 {
59 protected:
63 
64 public:
65  Quad10_2D_SUPG(int n, Domain * d);
66  virtual ~Quad10_2D_SUPG();
67 
68  virtual FEInterpolation *giveInterpolation() const;
69  virtual FEInterpolation *giveInterpolation(DofIDItem id) const;
70 
71  // definition
72  virtual const char *giveClassName() const { return "Quad1_2D_SUPG"; }
73  virtual const char *giveInputRecordName() const { return _IFT_Quad10_2D_SUPG_Name; }
74  virtual MaterialMode giveMaterialMode() { return _2dFlow; }
75 
76  virtual void giveInternalDofManDofIDMask(int i, IntArray &answer) const;
77  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
78  virtual int computeNumberOfDofs();
80  virtual void giveInputRecord(DynamicInputRecord &input);
81  virtual void updateYourself(TimeStep *tStep);
82  virtual int checkConsistency();
83 
84  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
85  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
86 
87  virtual double LS_PCS_computeF(LevelSetPCS *ls, TimeStep *tStep);
88  virtual void LS_PCS_computedN(FloatMatrix &answer);
89  virtual double LS_PCS_computeVolume() { return 0.0; }
90  virtual double LS_PCS_computeS(LevelSetPCS *ls, TimeStep *tStep);
91  virtual void LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi);
92 
94 
96  InternalStateType type, TimeStep *tStep);
97 
99 
100  void computeIntersection(int iedge, FloatArray &intcoords, FloatArray &fi);
101 
102  void computeMiddlePointOnParabolicArc(FloatArray &answer, int iedge, FloatArray borderpoints);
103 
104  void computeCenterOf(FloatArray &C, FloatArray c, int dim);
105 
106  void computeQuadraticRoots(FloatArray Coeff, double &r1, double &r2);
107 
108  void computeCoordsOfEdge(FloatArray &answer, int iedge);
109 
110  void computeQuadraticFunct(FloatArray &answer, int iedge);
111 
112  void computeQuadraticFunct(FloatArray &answer, FloatArray line);
114 
115  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
116 
117 #ifdef __OOFEG
119  int node, TimeStep *tStep);
120  // Graphics output
121  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
122  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
123  //virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) {}
124 #endif
125 
126  virtual double computeCriticalTimeStep(TimeStep *tStep);
127 
128  // three terms for computing their norms due to computing t_supg
129  virtual void computeAdvectionTerm(FloatMatrix &answer, TimeStep *tStep);
130  virtual void computeAdvectionDeltaTerm(FloatMatrix &answer, TimeStep *tStep);
131  virtual void computeMassDeltaTerm(FloatMatrix &answer, TimeStep *tStep);
132  virtual void computeLSICTerm(FloatMatrix &answer, TimeStep *tStep);
133  virtual void computeAdvectionEpsilonTerm(FloatMatrix &answer, TimeStep *tStep);
134  virtual void computeMassEpsilonTerm(FloatMatrix &answer, TimeStep *tStep);
135 
136  //virtual int giveNumberOfDofs() { return 1; }
137  virtual int giveNumberOfInternalDofManagers() const { return 1; }
138  virtual DofManager *giveInternalDofManager(int i) const;
139 
140 protected:
141  virtual void giveLocalVelocityDofMap(IntArray &map);
142  virtual void giveLocalPressureDofMap(IntArray &map);
143 
144  virtual void computeGaussPoints();
145  virtual void computeNuMatrix(FloatMatrix &answer, GaussPoint *gp);
146  virtual void computeUDotGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep);
147  virtual void computeBMatrix(FloatMatrix &anwer, GaussPoint *gp);
148  virtual void computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp);
149  virtual void computeNpMatrix(FloatMatrix &answer, GaussPoint *gp);
150  virtual void computeGradPMatrix(FloatMatrix &answer, GaussPoint *gp);
151  virtual void computeDivTauMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep);
152  virtual void computeGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep);
153  virtual int giveNumberOfSpatialDimensions();
154  virtual double computeVolumeAround(GaussPoint *gp);
155 
156  virtual void updateStabilizationCoeffs(TimeStep *tStep);
157 };
158 } // end namespace oofem
159 #endif // quad10_supg_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Class and object Domain.
Definition: domain.h:115
Class representing a 2d quadrilateral with constant interpolation.
virtual void computeMassEpsilonTerm(FloatMatrix &answer, TimeStep *tStep)
virtual double computeCriticalTimeStep(TimeStep *tStep)
Computes the critical time increment.
void computeCoordsOfEdge(FloatArray &answer, int iedge)
virtual void giveInternalDofManDofIDMask(int i, IntArray &answer) const
Returns internal dofmanager dof mask for node.
Class implementing internal element dof manager having some DOFs.
The element interface required by ZZNodalRecoveryModel.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
void computeIntersection(int iedge, FloatArray &intcoords, FloatArray &fi)
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual void computeAdvectionDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
virtual const char * giveInputRecordName() const
virtual void computeNuMatrix(FloatMatrix &answer, GaussPoint *gp)
virtual void computeGradPMatrix(FloatMatrix &answer, GaussPoint *gp)
Element interface for LevelSetPCS class representing level-set like material interface.
Definition: levelsetpcs.h:68
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual int giveNumberOfSpatialDimensions()
void computeQuadraticRoots(FloatArray Coeff, double &r1, double &r2)
Base class for dof managers.
Definition: dofmanager.h:113
virtual void computeNpMatrix(FloatMatrix &answer, GaussPoint *gp)
virtual int checkConsistency()
Performs consistency check.
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
virtual void giveLocalVelocityDofMap(IntArray &map)
Class implementing an array of integers.
Definition: intarray.h:61
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
This abstract class represent a general base element class for fluid dynamic problems.
Definition: supgelement2.h:57
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
virtual void computeAdvectionEpsilonTerm(FloatMatrix &answer, TimeStep *tStep)
virtual double LS_PCS_computeS(LevelSetPCS *ls, TimeStep *tStep)
Evaluates S in level set equation of the form where .
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
virtual double LS_PCS_computeF(LevelSetPCS *ls, TimeStep *tStep)
Evaluates F in level set equation of the form where for interface position driven by flow with speed...
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Quad10_2D_SUPG(int n, Domain *d)
virtual void LS_PCS_computedN(FloatMatrix &answer)
Returns gradient of shape functions.
virtual void computeLSICTerm(FloatMatrix &answer, TimeStep *tStep)
ElementDofManager pressureNode
virtual FEInterpolation * giveInterpolation() const
virtual void computeAdvectionTerm(FloatMatrix &answer, TimeStep *tStep)
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
static FEI2dQuadConst pressureInterpolation
virtual int giveNumberOfInternalDofManagers() const
void computeCenterOf(FloatArray &C, FloatArray c, int dim)
virtual void giveLocalPressureDofMap(IntArray &map)
virtual const char * giveClassName() const
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Updates the stabilization coefficients used for CBS and SUPG algorithms.
Class representing vector of real numbers.
Definition: floatarray.h:82
void computeQuadraticFunct(FloatArray &answer, int iedge)
#define _IFT_Quad10_2D_SUPG_Name
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual Interface * giveInterface(InterfaceType t)
Interface requesting service.
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Abstract base class representing Level Set representation of material interfaces. ...
Definition: levelsetpcs.h:114
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi)
Returns VOF fractions for each material on element according to nodal values of level set function (p...
Class Interface.
Definition: interface.h:82
virtual void computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp)
Class representing the a dynamic Input Record.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual void computeMassDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
Class representing a 2d isoparametric linear interpolation based on natural coordinates for quadrilat...
Definition: fei2dquadlin.h:45
int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep)
Returns internal state variable (like stress,strain) at node of element in Reduced form...
static FEI2dQuadLin velocityInterpolation
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void computeDivTauMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
Class representing 2d quadrilateral element with linear velocity and constant pressure approximations...
virtual void computeBMatrix(FloatMatrix &anwer, GaussPoint *gp)
virtual void computeUDotGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
void computeMiddlePointOnParabolicArc(FloatArray &answer, int iedge, FloatArray borderpoints)
virtual double LS_PCS_computeVolume()
Returns receiver's volume.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
InternalStateMode
Determines the mode of internal variable.
virtual DofManager * giveInternalDofManager(int i) const
Returns i-th internal element dof manager of the receiver.
virtual void computeGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.

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:30 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011