OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr1_2d_cbs.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_cbs_h
36 #define tr1_2d_cbs_h
37 
38 #include "cbselement.h"
39 #include "spatiallocalizer.h"
40 #include "zznodalrecoverymodel.h"
42 #include "sprnodalrecoverymodel.h"
43 #include "primaryfield.h"
44 //<RESTRICTED_SECTION>
45 #include "leplic.h"
46 //</RESTRICTED_SECTION>
47 
49 
50 #define _IFT_TR1_2D_CBS_Name "tr1cbs"
51 #define _IFT_Tr1CBS_vof "vof"
52 #define _IFT_Tr1CBS_pvof "pvof"
53 
54 
55 
56 namespace oofem {
57 class TimeStep;
58 class Node;
59 class Material;
60 class GaussPoint;
61 class FloatMatrix;
62 class FloatArray;
63 class IntArray;
64 class FEI2dTrLin;
65 
72 //<RESTRICTED_SECTION>
74 //</RESTRICTED_SECTION>
75 {
76 protected:
78  //double a[3];
79  double b [ 3 ];
80  double c [ 3 ];
81  double area;
82 
83 public:
84  TR1_2D_CBS(int n, Domain * aDomain);
85  virtual ~TR1_2D_CBS();
86 
87  virtual FEInterpolation *giveInterpolation() const;
88 
89  virtual void computeConsistentMassMtrx(FloatMatrix &answer, TimeStep *);
90  virtual void computeDiagonalMassMtrx(FloatArray &answer, TimeStep *);
91  virtual void computeConvectionTermsI(FloatArray &answer, TimeStep *);
92  virtual void computeDiffusionTermsI(FloatArray &answer, TimeStep *);
93  virtual void computeDensityRhsVelocityTerms(FloatArray &answer, TimeStep *tStep);
94  virtual void computeDensityRhsPressureTerms(FloatArray &answer, TimeStep *tStep);
95  virtual void computePrescribedTractionPressure(FloatArray &answer, TimeStep *tStep);
97  virtual void computePressureLhs(FloatMatrix &answer, TimeStep *tStep);
98  virtual void computeCorrectionRhs(FloatArray &answer, TimeStep *tStep);
99  virtual double computeCriticalTimeStep(TimeStep *tStep);
100 
101  // definition
102  virtual const char *giveClassName() const { return "TR1_2D_CBS"; }
103  virtual const char *giveInputRecordName() const { return _IFT_TR1_2D_CBS_Name; }
104  virtual MaterialMode giveMaterialMode() { return _2dFlow; }
105 
106  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
107  virtual int computeNumberOfDofs();
109  virtual void giveInputRecord(DynamicInputRecord &input);
110  virtual void updateYourself(TimeStep *tStep);
112  virtual int checkConsistency();
113 
114  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
115  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
116 
118 
120  const FloatArray &coords, IntArray &dofId, ValueModeType mode,
121  TimeStep *tStep);
122 
123  //<RESTRICTED_SECTION>
124  virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag);
125  virtual void formMaterialVolumePoly(Polygon &matvolpoly, LEPlic *matInterface,
126  const FloatArray &normal, const double p, bool updFlag);
127  virtual void formVolumeInterfacePoly(Polygon &matvolpoly, LEPlic *matInterface,
128  const FloatArray &normal, const double p, bool updFlag);
129  virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume);
130  virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool upd);
131  virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag);
132  virtual Element *giveElement() { return this; }
133  virtual double computeMyVolume(LEPlic *matInterface, bool updFlag);
134  virtual double computeCriticalLEPlicTimeStep(TimeStep *tStep) { return 1.e6; }
135  //</RESTRICTED_SECTION>
136 
138  InternalStateType type, TimeStep *tStep);
139 
141  virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap);
142  virtual int SPRNodalRecoveryMI_giveNumberOfIP();
144 
145  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
146  virtual void printOutputAt(FILE *file, TimeStep *tStep);
147 
148 #ifdef __OOFEG
150  int node, TimeStep *tStep);
151  // Graphics output
152  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
153  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
154  //virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) {}
155 #endif
156 
157 protected:
158  virtual void computeGaussPoints();
159  virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
160 };
161 } // end namespace oofem
162 #endif // tr1_2d_cbs_h
virtual void computeDiffusionTermsI(FloatArray &answer, TimeStep *)
Calculates contribution from diffusion terms for (*) velocities.
Definition: tr1_2d_cbs.C:274
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: tr1_2d_cbs.C:103
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
virtual int checkConsistency()
Used to check consistency and initialize some element geometry data (area,b,c)
Definition: tr1_2d_cbs.C:692
virtual ~TR1_2D_CBS()
Definition: tr1_2d_cbs.C:83
The element interface required by ZZNodalRecoveryModel.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr1_2d_cbs.C:90
virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag)
Assembles receiver volume.
Definition: tr1_2d_cbs.C:928
Element interface for LEPlic class representing Lagrangian-Eulerian (moving) material interface...
Definition: leplic.h:58
virtual Element * giveElement()
Return number of receiver&#39;s element.
Definition: tr1_2d_cbs.h:132
Class and object Domain.
Definition: domain.h:115
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
Definition: tr1_2d_cbs.C:1100
The element interface required by ZZNodalRecoveryModel.
Abstract class representing field of primary variables (those, which are unknown and are typically as...
Definition: primaryfield.h:104
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: tr1_2d_cbs.h:104
virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume)
Truncates given material polygon to receiver.
Definition: tr1_2d_cbs.C:911
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void computePressureLhs(FloatMatrix &answer, TimeStep *tStep)
Calculates the pressure LHS.
Definition: tr1_2d_cbs.C:589
virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool upd)
Computes the receiver center (in updated Lagrangian configuration).
Definition: tr1_2d_cbs.C:969
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
Definition: tr1_2d_cbs.C:1105
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
Abstract base class for all finite elements.
Definition: element.h:145
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Definition: tr1_2d_cbs.C:1078
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr1_2d_cbs.C:1049
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr1_2d_cbs.C:138
virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Computes deviatoric stress vector in given integration point and solution step from given total strai...
Definition: tr1_2d_cbs.C:678
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: tr1_2d_cbs.C:1125
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual void computeNumberOfNodalPrescribedTractionPressureContributions(FloatArray &answer, TimeStep *tStep)
Computes number of edges/sides with prescribed traction contributing to node with prescribed pressure...
Definition: tr1_2d_cbs.C:538
#define _IFT_TR1_2D_CBS_Name
Definition: tr1_2d_cbs.h:50
Class implementing an array of integers.
Definition: intarray.h:61
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: tr1_2d_cbs.C:126
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) ...
Definition: tr1_2d_cbs.C:994
virtual FEInterpolation * giveInterpolation() const
Definition: tr1_2d_cbs.C:87
virtual void computeDensityRhsVelocityTerms(FloatArray &answer, TimeStep *tStep)
Computes velocity terms on RHS for density equation.
Definition: tr1_2d_cbs.C:380
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Element interface class.
Definition: primaryfield.h:58
virtual const char * giveClassName() const
Definition: tr1_2d_cbs.h:102
Class representing a 2d triangular linear interpolation based on area coordinates.
Definition: fei2dtrlin.h:44
TR1_2D_CBS(int n, Domain *aDomain)
Definition: tr1_2d_cbs.C:71
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
Definition: tr1_2d_cbs.C:1087
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: tr1_2d_cbs.C:1058
This class is the implementation of triangular CFD element with linear (and equal order) interpolatio...
Definition: tr1_2d_cbs.h:70
virtual void computeConsistentMassMtrx(FloatMatrix &answer, TimeStep *)
Calculates consistent mass matrix.
Definition: tr1_2d_cbs.C:149
virtual double computeMyVolume(LEPlic *matInterface, bool updFlag)
Computes the volume of receiver.
Definition: tr1_2d_cbs.C:951
virtual double computeCriticalLEPlicTimeStep(TimeStep *tStep)
Computes critical time step.
Definition: tr1_2d_cbs.h:134
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr1_2d_cbs.C:1225
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
Definition: tr1_2d_cbs.C:654
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).
Definition: tr1_2d_cbs.C:770
virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
Computes corresponding volume fraction to given interface position.
Definition: tr1_2d_cbs.C:756
Abstract base class representing Lagrangian-Eulerian (moving) material interfaces.
Definition: leplic.h:153
virtual void computePrescribedTractionPressure(FloatArray &answer, TimeStep *tStep)
Computes prescribed pressure due to applied tractions.
Definition: tr1_2d_cbs.C:464
Class representing vector of real numbers.
Definition: floatarray.h:82
Class representing 2D polygon.
Definition: geotoolbox.h:91
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: tr1_2d_cbs.C:1113
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: tr1_2d_cbs.C:1070
virtual void computeDensityRhsPressureTerms(FloatArray &answer, TimeStep *tStep)
Computes pressure terms on RHS for density equation.
Definition: tr1_2d_cbs.C:566
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: tr1_2d_cbs.C:1149
virtual void computeDiagonalMassMtrx(FloatArray &answer, TimeStep *)
Calculates diagonal mass matrix.
Definition: tr1_2d_cbs.C:173
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual double computeCriticalTimeStep(TimeStep *tStep)
Calculates critical time step.
Definition: tr1_2d_cbs.C:723
Class Interface.
Definition: interface.h:82
This abstract class represent a general base element class for fluid dynamic problems solved using CB...
Definition: cbselement.h:56
Class representing the a dynamic Input Record.
The spatial localizer element interface associated to spatial localizer.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr1_2d_cbs.C:1195
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: tr1_2d_cbs.C:96
virtual void computeConvectionTermsI(FloatArray &answer, TimeStep *)
Calculates convection component for (*) velocities.
Definition: tr1_2d_cbs.C:187
double p
Line constant of line segment representing interface.
Definition: leplic.h:65
the oofem namespace is to define a context or scope in which all oofem names are defined.
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...
Definition: tr1_2d_cbs.C:1175
Class representing integration point in finite element program.
Definition: gausspoint.h:93
FloatArray normal
Interface segment normal.
Definition: leplic.h:67
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.
Definition: tr1_2d_cbs.C:802
Class representing solution step.
Definition: timestep.h:80
InternalStateMode
Determines the mode of internal variable.
virtual const char * giveInputRecordName() const
Definition: tr1_2d_cbs.h:103
static FEI2dTrLin interp
Definition: tr1_2d_cbs.h:77
virtual void computeCorrectionRhs(FloatArray &answer, TimeStep *tStep)
Calculates the RHS of velocity correction step.
Definition: tr1_2d_cbs.C:606

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