OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr1_2d_supg2.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 tr1_2d_supg2_h
36 #define tr1_2d_supg2_h
37 
38 #include "tr1_2d_supg.h"
39 #include "spatiallocalizer.h"
40 #include "zznodalrecoverymodel.h"
42 #include "sprnodalrecoverymodel.h"
43 #include "leplic.h"
44 
45 #define _IFT_TR1_2D_SUPG2_Name "tr1supg2"
46 
47 namespace oofem {
65 class TR1_2D_SUPG2 : public TR1_2D_SUPG
66 {
67 protected:
72  Polygon myPoly [ 2 ];
73  std::vector< FloatArray > vcoords [ 2 ];
74 
80  int mat [ 2 ];
85  std :: unique_ptr< IntegrationRule > defaultIRule;
86 
87 public:
88  TR1_2D_SUPG2(int n, Domain * d);
89  virtual ~TR1_2D_SUPG2();
90 
91  virtual void computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep);
92  virtual void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep);
93  virtual void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep);
94  virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep);
95  virtual void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep);
96  virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep);
97  virtual void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep);
98  virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep);
99  virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep);
100  virtual void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep);
102  answer.resize(3, 6);
103  answer.zero();
104  }
105  virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep) {
106  answer.resize(3);
107  answer.zero();
108  }
109  virtual void computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep);
110  virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep);
111  virtual void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep);
112  virtual void computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep);
113 
114  virtual void updateStabilizationCoeffs(TimeStep *tStep);
116  virtual double computeCriticalTimeStep(TimeStep *tStep);
117 
118  // definition
119  virtual const char *giveClassName() const { return "TR1_2D_SUPG2"; }
120  virtual const char *giveInputRecordName() const { return _IFT_TR1_2D_SUPG2_Name; }
121 
122  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
123  virtual int computeNumberOfDofs();
125  virtual void giveInputRecord(DynamicInputRecord &input);
126  virtual void updateYourself(TimeStep *tStep);
127 
128  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
129  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
130 
132 
134  const FloatArray &coords, IntArray &dofId, ValueModeType mode,
135  TimeStep *tStep);
136 
137  virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag);
138  virtual void formMaterialVolumePoly(Polygon &matvolpoly, LEPlic *matInterface,
139  const FloatArray &normal, const double p, bool updFlag);
140  virtual void formVolumeInterfacePoly(Polygon &matvolpoly, LEPlic *matInterface,
141  const FloatArray &normal, const double p, bool updFlag);
142  virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume);
143  virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool updFlag);
144  virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag);
145  virtual Element *giveElement() { return this; }
146  virtual double computeMyVolume(LEPlic *matInterface, bool updFlag);
147 
149  InternalStateType type, TimeStep *tStep);
150 
152  virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap);
153  virtual int SPRNodalRecoveryMI_giveNumberOfIP();
155 
156  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
157  virtual int giveDefaultIntegrationRule() const { return 0; }
158  virtual IntegrationRule *giveDefaultIntegrationRulePtr() { return defaultIRule.get(); }
159 
160 
161 
162 #ifdef __OOFEG
164  int node, TimeStep *tStep);
165  // Graphics output
166  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
167  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
168  //virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) {}
169 #endif
170 
171  virtual void printOutputAt(FILE *file, TimeStep *tStep);
172 
173 protected:
174  virtual void computeGaussPoints();
175  virtual void postInitialize();
176  virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *);
177  void computeNVector(FloatArray &answer, GaussPoint *gp);
178  virtual void updateVolumePolygons(Polygon &referenceFluidPoly, Polygon &secondFluidPoly, int &rfPoints, int &sfPoints,
179  const FloatArray &normal, const double p, bool updFlag);
180  double computeVolumeAroundID(GaussPoint *gp, integrationDomain id, const std::vector< FloatArray > &idpoly);
181  void updateIntegrationRules();
182  Material *_giveMaterial(int indx) { return domain->giveMaterial(mat [ indx ]); }
183 };
184 } // end namespace oofem
185 #endif // tr1_2d_supg2_h
virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool updFlag)
Computes the receiver center (in updated Lagrangian configuration).
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual Element * giveElement()
Return number of receiver&#39;s element.
Definition: tr1_2d_supg2.h:145
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
virtual const char * giveInputRecordName() const
Definition: tr1_2d_supg2.h:120
Class and object Domain.
Definition: domain.h:115
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: tr1_2d_supg2.C:103
virtual void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
Computes the derivative of diffusion terms for momentum balance equations(s) with respect to nodal ve...
Definition: tr1_2d_supg2.C:511
virtual void computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
Definition: tr1_2d_supg2.C:908
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
Abstract class representing field of primary variables (those, which are unknown and are typically as...
Definition: primaryfield.h:104
Material * _giveMaterial(int indx)
Definition: tr1_2d_supg2.h:182
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes SLIC stabilization term for momentum balance equation(s).
Definition: tr1_2d_supg2.C:646
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
Computes corresponding volume fraction to given interface position.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
Class representing 2d linear triangular element for solving incompressible fluid with SUPG solver...
Definition: tr1_2d_supg2.h:65
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Abstract base class for all finite elements.
Definition: element.h:145
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: tr1_2d_supg2.C:97
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: tr1_2d_supg2.C:134
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Abstract base class representing integration rule.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
int mat[2]
mat[0] reference fluid.
Definition: tr1_2d_supg2.h:80
virtual int giveDefaultIntegrationRule() const
Returns id of default integration rule.
Definition: tr1_2d_supg2.h:157
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
virtual ~TR1_2D_SUPG2()
Definition: tr1_2d_supg2.C:78
Material * giveMaterial(int n)
Service for accessing particular domain material model.
Definition: domain.C:281
virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for momentum balance equations(s).
Definition: tr1_2d_supg2.C:448
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: tr1_2d_supg2.h:158
virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the linear advection term for mass conservation equation.
Definition: tr1_2d_supg2.C:689
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr1_2d_supg2.C:158
TR1_2D_SUPG2(int n, Domain *d)
Definition: tr1_2d_supg2.C:72
virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag)
Assembles receiver volume.
virtual double computeMyVolume(LEPlic *matInterface, bool updFlag)
Computes the volume of receiver.
virtual void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for momentum balance equations(s) with respect to nodal ve...
Definition: tr1_2d_supg2.C:367
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...
#define _IFT_TR1_2D_SUPG2_Name
Definition: tr1_2d_supg2.h:45
Abstract base class for all material models.
Definition: material.h:95
virtual void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for mass conservation equation with respect to nodal veloc...
Definition: tr1_2d_supg2.C:722
virtual void computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms (generalized mass matrix with stabilization terms) for momentum balance e...
Definition: tr1_2d_supg2.C:190
virtual void formVolumeInterfacePoly(Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag)
Assembles receiver material polygon based solely on given interface line.
virtual void computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms for mass conservation equation.
Definition: tr1_2d_supg2.C:763
double computeVolumeAroundID(GaussPoint *gp, integrationDomain id, const std::vector< FloatArray > &idpoly)
Abstract base class representing Lagrangian-Eulerian (moving) material interfaces.
Definition: leplic.h:153
virtual void formMaterialVolumePoly(Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag)
Assembles the true element material polygon (takes receiver vof into accout).
virtual const char * giveClassName() const
Definition: tr1_2d_supg2.h:119
Class representing vector of real numbers.
Definition: floatarray.h:82
Class representing 2D polygon.
Definition: geotoolbox.h:91
virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *)
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 void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for mass conservation equation.
Definition: tr1_2d_supg2.C:769
virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume)
Truncates given material polygon to receiver.
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
virtual void updateElementForNewInterfacePosition(TimeStep *tStep)
Definition: tr1_2d_supg2.h:115
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Class representing the general Input Record.
Definition: inputrecord.h:101
Class Interface.
Definition: interface.h:82
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
Polygon myPoly[2]
myPoly[0] occupied by reference fluid.
Definition: tr1_2d_supg2.h:72
Class representing the a dynamic Input Record.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for momentum balance equations(s).
Definition: tr1_2d_supg2.C:591
std::vector< FloatArray > vcoords[2]
Definition: tr1_2d_supg2.h:73
virtual int EIPrimaryFieldI_evaluateFieldVectorAt(FloatArray &answer, PrimaryField &pf, const FloatArray &coords, IntArray &dofId, ValueModeType mode, TimeStep *tStep)
Evaluates the value of field at given point of interest (should be located inside receiver&#39;s volume) ...
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Class representing 2d linear triangular element for solving incompressible fluid with SUPG solver...
Definition: tr1_2d_supg.h:68
virtual void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes diffusion derivative terms for mass conservation equation.
Definition: tr1_2d_supg2.h:101
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void postInitialize()
Performs post initialization steps.
Definition: tr1_2d_supg2.C:149
double p
Line constant of line segment representing interface.
Definition: leplic.h:65
std::unique_ptr< IntegrationRule > defaultIRule
default integration rule over element volume.
Definition: tr1_2d_supg2.h:85
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
virtual void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
Definition: tr1_2d_supg2.C:809
void computeNVector(FloatArray &answer, GaussPoint *gp)
Definition: tr1_2d_supg2.C:84
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Updates the stabilization coefficients used for CBS and SUPG algorithms.
Definition: tr1_2d_supg2.C:915
Class representing integration point in finite element program.
Definition: gausspoint.h:93
FloatArray normal
Interface segment normal.
Definition: leplic.h:67
Class representing solution step.
Definition: timestep.h:80
InternalStateMode
Determines the mode of internal variable.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr1_2d_supg2.C:91
virtual double computeCriticalTimeStep(TimeStep *tStep)
Computes the critical time increment.
virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes advection terms for mass conservation equation.
Definition: tr1_2d_supg2.C:695
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for mass conservation equation.
Definition: tr1_2d_supg2.h:105
virtual void updateVolumePolygons(Polygon &referenceFluidPoly, Polygon &secondFluidPoly, int &rfPoints, int &sfPoints, const FloatArray &normal, const double p, bool updFlag)
virtual void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes nonlinear advection terms for momentum balance equations(s).
Definition: tr1_2d_supg2.C:286
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

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